Changeset 3006
- Timestamp:
- Jul 21, 2023, 11:47:04 AM (16 months ago)
- Location:
- trunk/LMDZ.MARS
- Files:
-
- 1 deleted
- 3 edited
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/changelog.txt
r3004 r3006 4111 4111 specifically output for Grads in 1D). 4112 4112 Also turned lwi and lwflux into modules while at it. 4113 4114 == 21/07/2023 == EM 4115 More code cleanup. Turn "nirdata.h" common into module "nirdata.F90" and 4116 include "nir_leedat.F" (reading/loading of the data) in the module. 4117 Also turn nirco2abs.F in a module. -
trunk/LMDZ.MARS/libf/phymars/nirco2abs.F
r2616 r3006 1 MODULE nirco2abs_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE nirco2abs(ngrid,nlayer,pplay,dist_sol,nq,pq, 2 8 $ mu0,fract,declin,pdtnirco2) … … 6 12 USE comcstfi_h, ONLY: pi 7 13 USE time_phylmdz_mod, ONLY: daysec 14 use nirdata_mod, only: npres, alfa, corgcm, oco21d, pres1d 8 15 IMPLICIT NONE 9 16 c======================================================================= … … 49 56 c 50 57 include "callkeys.h" 51 include "nirdata.h"52 58 53 59 c----------------------------------------------------------------------- … … 237 243 END IF 238 244 239 end245 END SUBROUTINE nirco2abs 240 246 241 247 … … 247 253 C escout(nlayer) on pressure grid p(nlayer). 248 254 C 249 real escout(nlayer),p(nlayer) 250 real escin(nl),pin(nl),wm,wp 251 integer nl,nlayer,n1,n,nm,np 255 real,intent(out) :: escout(nlayer) 256 real,intent(in) :: p(nlayer) 257 integer,intent(in) :: nlayer 258 real,intent(in) :: escin(nl) 259 real,intent(in) :: pin(nl) 260 integer,intent(in) :: nl 261 262 real :: wm,wp 263 integer :: n1,n,nm,np 264 252 265 do n1=1,nlayer 253 266 if(p(n1) .gt. 1500. .or. p(n1) .lt. 1.0e-13) then … … 265 278 endif 266 279 enddo 267 end 280 281 end subroutine interpnir 282 283 END MODULE nirco2abs_mod -
trunk/LMDZ.MARS/libf/phymars/nirdata.F90
r3004 r3006 1 module nirdata_mod 1 2 2 integer npres ! Number of pressures in NIR correction 3 parameter (npres=42) ! table 3 implicit none 4 4 5 common /NIRdata/ pres1d,corgcm,oco21d,alfa,p1999 6 !$OMP THREADPRIVATE(/NIRdata/) 7 real pres1d(npres) 8 real corgcm(npres) 9 real oco21d(npres),alfa(npres),p1999(npres) 5 ! Number of pressures in NIR correction table 6 integer,parameter :: npres=42 7 8 real,save,protected :: pres1d(npres) 9 real,save,protected :: corgcm(npres) 10 real,save,protected :: oco21d(npres) 11 real,save,protected :: alfa(npres) 12 real,save,protected :: p1999(npres) 13 !$OMP THREADPRIVATE(pres1d,corgcm,oco21d,alfa,p1999) 14 15 contains 16 17 subroutine NIR_leedat 18 19 ! reads parameters for NIR NLTE calculation 20 21 ! nov 2011 fgg+malv first version 22 !*********************************************************************** 23 24 use datafile_mod, only: datadir 25 USE mod_phys_lmdz_para, ONLY: is_master 26 USE mod_phys_lmdz_transfert_para, ONLY: bcast 27 28 implicit none 29 30 ! local variables 31 32 integer :: ind 33 34 if (is_master) then 35 open(43,file=trim(datadir)//'/NIRcorrection_feb2011.dat', & 36 status='old') 37 do ind=1,9 38 read(43,*) 39 enddo 40 41 do ind=1,npres 42 read(43,*)pres1d(ind),corgcm(ind),oco21d(ind),p1999(ind), & 43 alfa(ind) 44 !Tabulated pression to Pa 45 pres1d(ind)=pres1d(ind)*100. 46 enddo 47 close(43) 48 49 endif ! if(is_master) then 50 51 call bcast(pres1d) 52 call bcast(corgcm) 53 call bcast(oco21d) 54 call bcast(p1999) 55 call bcast(alfa) 56 57 end subroutine NIR_leedat 58 59 end module nirdata_mod -
trunk/LMDZ.MARS/libf/phymars/physiq_mod.F
r3004 r3006 56 56 use comsaison_h, only: dist_sol, declin, zls, 57 57 & mu0, fract, local_time 58 use nirdata_mod, only: NIR_leedat 59 use nirco2abs_mod, only: nirco2abs 58 60 use slope_mod, only: theta_sl, psi_sl 59 61 use conc_mod, only: rnew, cpnew, mmean … … 109 111 & major_slope,compute_meshgridavg, 110 112 & ini_comslope_h 111 USE ioipsl_getincom, only: getin112 113 use write_output_mod, only: write_output 113 114 IMPLICIT NONE
Note: See TracChangeset
for help on using the changeset viewer.