Changeset 2032
- Timestamp:
- Oct 30, 2018, 12:22:31 PM (6 years ago)
- Location:
- trunk/LMDZ.GENERIC
- Files:
-
- 6 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/README
r2004 r2032 1397 1397 == 02/10/2018 == JL 1398 1398 - correct a bug introduced in commit 1987 in optcv. 1399 1400 == 30/10/2018 == EM 1401 - correct a bug introduced in revision 2026; now that L_NGAUSS is a parameter 1402 read in via sugas_corrk (called at first call by callcorrk), automatic arrays 1403 of size L_NGAUSS cannot be declared in callcorrk and must be allocated 1404 once the value of L_NGAUSS has been set. 1405 - turned optci, optcv and callcorrk into modules in the process. -
trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90
r2026 r2032 1 MODULE callcorrk_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 subroutine callcorrk(ngrid,nlayer,pq,nq,qsurf, & 2 8 albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt, & … … 10 16 11 17 use mod_phys_lmdz_para, only : is_master 12 use radinc_h 13 use radcommon_h 14 use watercommon_h 18 use radinc_h, only: L_NSPECTV, L_NSPECTI, naerkind, banddir, corrkdir,& 19 L_LEVELS, L_NGAUSS, L_NLEVRAD, L_NLAYRAD, L_REFVAR 20 use radcommon_h, only: wrefvar, Cmk, fzeroi, fzerov, gasi, gasv, & 21 glat_ig, gweight, pfgasref, pgasmax, pgasmin, & 22 pgasref, tgasmax, tgasmin, tgasref, scalep, & 23 ubari, wnoi, stellarf, glat, dwnv, dwni, tauray 24 use watercommon_h, only: psat_water, epsi 15 25 use datafile_mod, only: datadir 16 26 use ioipsl_getin_p_mod, only: getin_p 17 use gases_h 27 use gases_h, only: ngasmx 18 28 use radii_mod, only : su_aer_radii,co2_reffrad,h2o_reffrad,dust_reffrad,h2so4_reffrad,back2lay_reffrad 19 29 use aerosol_mod, only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4, iaero_back2lay, iaero_nh3, iaero_aurora 20 USE tracer_h30 use tracer_h, only: igcm_h2o_vap, igcm_h2o_ice, igcm_co2_ice 21 31 use comcstfi_mod, only: pi, mugaz, cpp 22 32 use callkeys_mod, only: varactive,diurnal,tracer,water,varfixed,satval, & 23 33 kastprof,strictboundcorrk,specOLR,CLFvarying 24 34 use optcv_mod, only: optcv 35 use optci_mod, only: optci 25 36 implicit none 26 37 … … 113 124 ! Optical values for the optci/cv subroutines 114 125 REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV) 115 REAL*8 dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS) 116 REAL*8 dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS) 117 REAL*8 cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS) 118 REAL*8 cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS) 119 REAL*8 wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS) 120 REAL*8 wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS) 121 REAL*8 tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS) 122 REAL*8 taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS) 123 REAL*8 taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS) 126 ! NB: Arrays below are "save" to avoid reallocating them at every call 127 ! not because their content needs be reused from call to the next 128 REAL*8,allocatable,save :: dtaui(:,:,:) 129 REAL*8,allocatable,save :: dtauv(:,:,:) 130 REAL*8,allocatable,save :: cosbv(:,:,:) 131 REAL*8,allocatable,save :: cosbi(:,:,:) 132 REAL*8,allocatable,save :: wbari(:,:,:) 133 REAL*8,allocatable,save :: wbarv(:,:,:) 134 REAL*8,allocatable,save :: tauv(:,:,:) 135 REAL*8,allocatable,save :: taucumv(:,:,:) 136 REAL*8,allocatable,save :: taucumi(:,:,:) 124 137 125 138 REAL*8 tauaero(L_LEVELS,naerkind) … … 136 149 INTEGER ig,l,k,nw,iaer 137 150 138 real szangle 139 logical global1d 140 save szangle,global1d 151 real,save :: szangle 152 logical,save :: global1d 141 153 !$OMP THREADPRIVATE(szangle,global1d) 142 real*8 taugsurf(L_NSPECTV,L_NGAUSS-1) 143 real*8 taugsurfi(L_NSPECTI,L_NGAUSS-1) 144 real*8 qvar(L_LEVELS) ! Mixing ratio of variable component (mol/mol). 154 155 real*8,allocatable,save :: taugsurf(:,:) 156 real*8,allocatable,save :: taugsurfi(:,:) 157 real*8 qvar(L_LEVELS) ! Mixing ratio of variable component (mol/mol). 145 158 146 159 ! Local aerosol optical properties for each column on RADIATIVE grid. … … 232 245 call suaer_corrk ! Set up aerosol optical properties. 233 246 247 248 ! now that L_NGAUSS has been initialized (by sugas_corrk) 249 ! allocate related arrays 250 ALLOCATE(dtaui(L_NLAYRAD,L_NSPECTI,L_NGAUSS)) 251 ALLOCATE(dtauv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)) 252 ALLOCATE(cosbv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)) 253 ALLOCATE(cosbi(L_NLAYRAD,L_NSPECTI,L_NGAUSS)) 254 ALLOCATE(wbari(L_NLAYRAD,L_NSPECTI,L_NGAUSS)) 255 ALLOCATE(wbarv(L_NLAYRAD,L_NSPECTV,L_NGAUSS)) 256 ALLOCATE(tauv(L_NLEVRAD,L_NSPECTV,L_NGAUSS)) 257 ALLOCATE(taucumv(L_LEVELS,L_NSPECTV,L_NGAUSS)) 258 ALLOCATE(taucumi(L_LEVELS,L_NSPECTI,L_NGAUSS)) 259 ALLOCATE(taugsurf(L_NSPECTV,L_NGAUSS-1)) 260 ALLOCATE(taugsurfi(L_NSPECTI,L_NGAUSS-1)) 234 261 235 262 if((igcm_h2o_vap.eq.0) .and. varactive)then … … 795 822 end do 796 823 endif 797 824 798 825 call optcv(dtauv,tauv,taucumv,plevrad, & 799 826 qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero, & … … 973 1000 974 1001 end subroutine callcorrk 1002 1003 END MODULE callcorrk_mod -
trunk/LMDZ.GENERIC/libf/phystd/optci.F90
r1725 r2032 1 MODULE optci_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 subroutine optci(PLEV,TLEV,DTAUI,TAUCUMI, & 2 8 QXIAER,QSIAER,GIAER,COSBI,WBARI,TAUAERO, & 3 9 TMID,PMID,TAUGSURF,QVAR,MUVAR) 4 10 5 use radinc_h 11 use radinc_h, only: L_LEVELS, L_NLAYRAD, L_NSPECTI, L_NGAUSS, & 12 L_NLEVRAD, L_REFVAR, naerkind 6 13 use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi,glat_ig 7 use gases_h 14 use gases_h, only: gfrac, ngasmx, igas_N2, igas_He, igas_H2O, igas_H2 8 15 use comcstfi_mod, only: g, r, mugaz 9 16 use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple … … 444 451 445 452 446 return447 448 449 453 end subroutine optci 450 454 451 452 455 END MODULE optci_mod 456 -
trunk/LMDZ.GENERIC/libf/phystd/optcv.F90
r2004 r2032 1 MODULE optcv_mod 2 3 IMPLICIT NONE 4 5 CONTAINS 6 1 7 SUBROUTINE OPTCV(DTAUV,TAUV,TAUCUMV,PLEV, & 2 8 QXVAER,QSVAER,GVAER,WBARV,COSBV, & 3 9 TAURAY,TAUAERO,TMID,PMID,TAUGSURF,QVAR,MUVAR) 4 10 5 use radinc_h 11 use radinc_h, only: L_NLAYRAD, L_NLEVRAD, L_LEVELS, L_NSPECTV, L_NGAUSS, L_REFVAR, NAERKIND 6 12 use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig 7 use gases_h 13 use gases_h, only: gfrac, ngasmx, igas_H2, igas_H2O, igas_He, igas_N2 8 14 use comcstfi_mod, only: g, r, mugaz 9 15 use callkeys_mod, only: kastprof,continuum,graybody,H2Ocont_simple,callgasvis … … 369 375 370 376 371 return372 373 374 377 end subroutine optcv 378 379 END MODULE optcv_mod -
trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90
r2019 r2032 49 49 use conc_mod 50 50 use phys_state_var_mod 51 use callcorrk_mod, only: callcorrk 51 52 use turb_mod, only : q2,sensibFlux,turb_resolved 52 53 #ifndef MESOSCALE -
trunk/LMDZ.GENERIC/libf/phystd/sugas_corrk.F90
r2026 r2032 20 20 !================================================================== 21 21 22 use radinc_h 22 use radinc_h, only: corrkdir, banddir, L_NSPECTI, L_NSPECTV, & 23 L_NGAUSS, L_NPREF, L_NTREF, L_REFVAR, L_PINT 23 24 use radcommon_h, only : pgasref,pfgasref,pgasmin,pgasmax 24 25 use radcommon_h, only : tgasref,tgasmin,tgasmax … … 27 28 use datafile_mod, only: datadir 28 29 use comcstfi_mod, only: mugaz 29 use gases_h 30 use gases_h, only: gnom, ngasmx, & 31 igas_H2, igas_H2O, igas_He, igas_N2 30 32 use ioipsl_getin_p_mod, only: getin_p 31 33 use callkeys_mod, only: varactive,varfixed,graybody,callgasvis,&
Note: See TracChangeset
for help on using the changeset viewer.