Changeset 858 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- Dec 20, 2012, 5:54:26 PM (12 years ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 9 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/aeropacity.F90
r804 r858 5 5 use aerosol_mod 6 6 USE comgeomfi_h 7 USE tracer_h 7 USE tracer_h, only: noms,rho_co2,rho_ice 8 8 9 9 implicit none … … 46 46 #include "comvert.h" 47 47 48 INTEGER ngrid,nlayer,nq 49 50 REAL pplay(ngrid,nlayer) 51 REAL pplev(ngrid,nlayer+1) 52 REAL pq(ngrid,nlayer,nq) 53 REAL aerosol(ngrid,nlayer,naerkind) 54 REAL reffrad(ngrid,nlayer,naerkind) 55 REAL QREFvis3d(ngrid,nlayermx,naerkind) 56 REAL QREFir3d(ngrid,nlayermx,naerkind) 57 58 REAL tau_col(ngrid) 59 ! REAL tauref(ngrid), tau_col(ngrid) 60 61 real cloudfrac(ngrid,nlayermx) 48 INTEGER,INTENT(IN) :: ngrid ! number of atmospheric columns 49 INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers 50 INTEGER,INTENT(IN) :: nq ! number of tracers 51 REAL,INTENT(IN) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa) 52 REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) 53 REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers (.../kg_of_air) 54 REAL,INTENT(OUT) :: aerosol(ngrid,nlayer,naerkind) ! aerosol optical depth 55 REAL,INTENT(IN) :: reffrad(ngrid,nlayer,naerkind) ! aerosol effective radius 56 REAL,INTENT(IN) :: QREFvis3d(ngrid,nlayermx,naerkind) ! extinction coefficient in the visible 57 REAL,INTENT(IN) :: QREFir3d(ngrid,nlayermx,naerkind) 58 REAL,INTENT(OUT):: tau_col(ngrid) !column integrated visible optical depth 59 ! BENJAMIN MODIFS 60 real,intent(in) :: cloudfrac(ngrid,nlayermx) ! cloud fraction 61 real,intent(out) :: totcloudfrac(ngrid) ! total cloud fraction 62 logical,intent(in) :: clearsky 63 62 64 real aerosol0 63 65 64 66 INTEGER l,ig,iq,iaer 65 67 66 LOGICAL firstcall 67 DATA firstcall/.true./ 68 SAVE firstcall 68 LOGICAL,SAVE :: firstcall=.true. 69 69 70 70 REAL CBRT … … 79 79 REAL taudusttmp(ngrid) ! Temporary dust opacity used before scaling 80 80 REAL tauh2so4tmp(ngrid) ! Temporary h2so4 opacity used before scaling 81 ! BENJAMIN MODIFS 81 82 82 real CLFtot 83 real totcloudfrac(ngrid)84 logical clearsky85 83 86 84 ! identify tracers -
trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90
r787 r858 1 1 subroutine callcorrk(ngrid,nlayer,pq,nq,qsurf, & 2 2 albedo,emis,mu0,pplev,pplay,pt, & 3 tsurf,fract,dist_star,aerosol,muvar, &3 tsurf,fract,dist_star,aerosol,muvar, & 4 4 dtlw,dtsw,fluxsurf_lw, & 5 5 fluxsurf_sw,fluxtop_lw,fluxabs_sw,fluxtop_dn, & 6 6 OLR_nu,OSR_nu, & 7 reffrad,nueffrad,tau_col,cloudfrac,totcloudfrac,&7 tau_col,cloudfrac,totcloudfrac, & 8 8 clearsky,firstcall,lastcall) 9 9 … … 15 15 use gases_h 16 16 use radii_mod, only : su_aer_radii,co2_reffrad,h2o_reffrad,dust_reffrad,h2so4_reffrad 17 use aerosol_mod 18 17 use aerosol_mod, only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4 19 18 USE tracer_h 20 19 … … 45 44 ! Layer #nlayermx is the layer at the top. 46 45 47 ! INPUT 48 INTEGER icount 49 INTEGER ngrid,nlayer 50 REAL aerosol(ngrid,nlayermx,naerkind) ! aerosol tau (kg/kg) 51 REAL albedo(ngrid) ! SW albedo 52 REAL emis(ngrid) ! LW emissivity 53 REAL pplay(ngrid,nlayermx) ! pres. level in GCM mid of layer 54 REAL pplev(ngrid,nlayermx+1) ! pres. level at GCM layer boundaries 55 56 REAL pt(ngrid,nlayermx) ! air temperature (K) 57 REAL tsurf(ngrid) ! surface temperature (K) 58 REAL dist_star,mu0(ngrid) ! distance star-planet (AU) 59 REAL fract(ngrid) ! fraction of day 46 INTEGER,INTENT(IN) :: ngrid ! number of atmospheric columns 47 INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers 48 REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers (.../kg_of_air) 49 integer,intent(in) :: nq ! number of tracers 50 REAL,INTENT(IN) :: qsurf(ngrid,nq) ! tracer on surface (kg.m-2) 51 REAL,INTENT(IN) :: albedo(ngrid) ! SW albedo 52 REAL,INTENT(IN) :: emis(ngrid) ! LW emissivity 53 real,intent(in) :: mu0(ngrid) ! cosine of sun incident angle 54 REAL,INTENT(IN) :: pplev(ngrid,nlayermx+1) ! inter-layer pressure (Pa) 55 REAL,INTENT(IN) :: pplay(ngrid,nlayermx) ! mid-layer pressure (Pa) 56 REAL,INTENT(IN) :: pt(ngrid,nlayermx) ! air temperature (K) 57 REAL,INTENT(IN) :: tsurf(ngrid) ! surface temperature (K) 58 REAL,INTENT(IN) :: fract(ngrid) ! fraction of day 59 REAL,INTENT(IN) :: dist_star ! distance star-planet (AU) 60 REAL,INTENT(OUT) :: aerosol(ngrid,nlayermx,naerkind) ! aerosol tau (kg/kg) 61 real,intent(in) :: muvar(ngrid,nlayermx+1) 62 REAL,INTENT(OUT) :: dtlw(ngrid,nlayermx) ! heating rate (K/s) due to LW 63 REAL,INTENT(OUT) :: dtsw(ngrid,nlayermx) ! heating rate (K/s) due to SW 64 REAL,INTENT(OUT) :: fluxsurf_lw(ngrid) ! incident LW flux to surf (W/m2) 65 REAL,INTENT(OUT) :: fluxsurf_sw(ngrid) ! incident SW flux to surf (W/m2) 66 REAL,INTENT(OUT) :: fluxtop_lw(ngrid) ! outgoing LW flux to space (W/m2) 67 REAL,INTENT(OUT) :: fluxabs_sw(ngrid) ! SW flux absorbed by planet (W/m2) 68 REAL,INTENT(OUT) :: fluxtop_dn(ngrid) ! incident top of atmosphere SW flux (W/m2) 69 REAL,INTENT(OUT) :: OLR_nu(ngrid,L_NSPECTI)! Outgoing LW radition in each band (Normalized to the band width (W/m2/cm-1) 70 REAL,INTENT(OUT) :: OSR_nu(ngrid,L_NSPECTV)! Outgoing SW radition in each band (Normalized to the band width (W/m2/cm-1) 71 REAL,INTENT(OUT) :: tau_col(ngrid) ! diagnostic from aeropacity 72 ! for H2O cloud fraction in aeropacity 73 real,intent(in) :: cloudfrac(ngrid,nlayermx) 74 real,intent(out) :: totcloudfrac(ngrid) 75 logical,intent(in) :: clearsky 76 logical,intent(in) :: firstcall ! signals first call to physics 77 logical,intent(in) :: lastcall ! signals last call to physics 60 78 61 79 ! Globally varying aerosol optical properties on GCM grid … … 72 90 ! REAL :: omegaREFir3d(ngrid,nlayermx,naerkind) ! not sure of the point of these... 73 91 74 !!! this is a local instance of a variable saved in physiq.F and being an argument 75 REAL, INTENT(INOUT) :: reffrad(ngrid,nlayer,naerkind) 76 REAL, INTENT(INOUT) :: nueffrad(ngrid,nlayermx,naerkind) 77 78 !! OUTPUT 79 REAL dtsw(ngrid,nlayermx) ! heating rate (K/s) due to SW 80 REAL dtlw(ngrid,nlayermx) ! heating rate (K/s) due to LW 81 REAL fluxsurf_lw(ngrid) ! incident LW flux to surf (W/m2) 82 REAL fluxtop_lw(ngrid) ! outgoing LW flux to space (W/m2) 83 REAL fluxsurf_sw(ngrid) ! incident SW flux to surf (W/m2) 84 REAL fluxabs_sw(ngrid) ! SW flux absorbed by planet (W/m2) 85 REAL fluxtop_dn(ngrid) ! incident top of atmosphere SW flux (W/m2) 86 REAL OLR_nu(ngrid,L_NSPECTI)! Outgoing LW radition in each band (Normalized to the band width (W/m2/cm-1) 87 REAL OSR_nu(ngrid,L_NSPECTV)! Outgoing SW radition in each band (Normalized to the band width (W/m2/cm-1) 92 REAL,ALLOCATABLE,SAVE :: reffrad(:,:,:) ! aerosol effective radius (m) 93 REAL,ALLOCATABLE,SAVE :: nueffrad(:,:,:) ! aerosol effective variance 88 94 89 95 !----------------------------------------------------------------------- … … 117 123 118 124 INTEGER ig,l,k,nw,iaer,irad 125 INTEGER icount 119 126 120 127 real fluxtoplanet … … 126 133 127 134 real*8 qvar(L_LEVELS) ! mixing ratio of variable component (mol/mol) 128 REAL pq(ngrid,nlayer,nq)129 REAL qsurf(ngrid,nq) ! tracer on surface (e.g. kg.m-2)130 integer nq131 135 132 136 ! Local aerosol optical properties for each column on RADIATIVE grid 133 real*8 QXVAER(L_LEVELS+1,L_NSPECTV,naerkind) 134 real*8 QSVAER(L_LEVELS+1,L_NSPECTV,naerkind) 135 real*8 GVAER(L_LEVELS+1,L_NSPECTV,naerkind) 136 real*8 QXIAER(L_LEVELS+1,L_NSPECTI,naerkind) 137 real*8 QSIAER(L_LEVELS+1,L_NSPECTI,naerkind) 138 real*8 GIAER(L_LEVELS+1,L_NSPECTI,naerkind) 139 140 save qxvaer, qsvaer, gvaer 141 save qxiaer, qsiaer, giaer 137 real*8,save :: QXVAER(L_LEVELS+1,L_NSPECTV,naerkind) 138 real*8,save :: QSVAER(L_LEVELS+1,L_NSPECTV,naerkind) 139 real*8,save :: GVAER(L_LEVELS+1,L_NSPECTV,naerkind) 140 real*8,save :: QXIAER(L_LEVELS+1,L_NSPECTI,naerkind) 141 real*8,save :: QSIAER(L_LEVELS+1,L_NSPECTI,naerkind) 142 real*8,save :: GIAER(L_LEVELS+1,L_NSPECTI,naerkind) 142 143 143 144 !REAL :: QREFvis3d(ngrid,nlayermx,naerkind) … … 147 148 real, dimension(:,:,:), save, allocatable :: QREFir3d 148 149 149 REAL tau_col(ngrid) ! diagnostic from aeropacity150 150 151 151 ! Misc. 152 logical firstcall, lastcall,nantest152 logical nantest 153 153 real*8 tempv(L_NSPECTV) 154 154 real*8 tempi(L_NSPECTI) … … 169 169 real*8 NFLUXGNDV_nu(L_NSPECTV) 170 170 171 ! for H2O cloud fraction in aeropacity172 real cloudfrac(ngrid,nlayermx)173 real totcloudfrac(ngrid)174 logical clearsky175 171 176 172 ! for weird cloud test … … 179 175 real maxrad, minrad 180 176 181 real CBRT 182 external CBRT 177 real,external :: CBRT 183 178 184 179 ! included by RW for runaway greenhouse 1D study 185 real muvar(ngrid,nlayermx+1)186 180 real vtmp(nlayermx) 187 181 REAL*8 muvarrad(L_LEVELS) … … 208 202 !-------------------------------------------------- 209 203 ! Effective radius and variance of the aerosols 204 allocate(reffrad(ngrid,nlayer,naerkind)) 205 allocate(nueffrad(ngrid,nlayer,naerkind)) 206 210 207 if(naerkind.gt.4)then 211 208 print*,'Code not general enough to deal with naerkind > 4 yet.' … … 262 259 endif 263 260 264 firstcall=.false. 265 266 end if 261 end if ! of if (firstcall) 267 262 268 263 !======================================================================= … … 274 269 275 270 if ((iaer.eq.iaero_co2).and.tracer.and.(igcm_co2_ice.gt.0)) then ! treat condensed co2 particles. 276 call co2_reffrad(ngrid,nq,pq,reffrad )271 call co2_reffrad(ngrid,nq,pq,reffrad(1,1,iaero_co2)) 277 272 print*,'Max. CO2 ice particle size = ',maxval(reffrad(1:ngrid,1:nlayermx,iaer))/1.e-6,' um' 278 273 print*,'Min. CO2 ice particle size = ',minval(reffrad(1:ngrid,1:nlayermx,iaer))/1.e-6,' um' 279 274 end if 280 275 if ((iaer.eq.iaero_h2o).and.water) then ! treat condensed water particles. to be generalized for other aerosols 281 call h2o_reffrad(ngrid,nq,pq,pt,reffrad,nueffrad) 276 call h2o_reffrad(ngrid,pq(1,1,igcm_h2o_ice),pt, & 277 reffrad(1,1,iaero_h2o),nueffrad(1,1,iaero_h2o)) 282 278 print*,'Max. H2O cloud particle size = ',maxval(reffrad(1:ngrid,1:nlayermx,iaer))/1.e-6,' um' 283 279 print*,'Min. H2O cloud particle size = ',minval(reffrad(1:ngrid,1:nlayermx,iaer))/1.e-6,' um' 284 280 endif 285 281 if(iaer.eq.iaero_dust)then 286 call dust_reffrad(ngrid,reffrad )282 call dust_reffrad(ngrid,reffrad(1,1,iaero_dust)) 287 283 print*,'Dust particle size = ',reffrad(1,1,iaer)/1.e-6,' um' 288 284 endif 289 285 if(iaer.eq.iaero_h2so4)then 290 call h2so4_reffrad(ngrid,reffrad )286 call h2so4_reffrad(ngrid,reffrad(1,1,iaero_h2so4)) 291 287 print*,'H2SO4 particle size =',reffrad(1,1,iaer)/1.e-6,' um' 292 288 endif … … 638 634 print*,'Minimum temperature is outside the radiative' 639 635 print*,'transfer kmatrix bounds, exiting.' 636 print*,"k=",k," tlevrad(k)=",tlevrad(k) 637 print*,"tgasmin=",tgasmin 640 638 !print*,'WARNING, OVERRIDING FOR TEST' 641 639 call abort -
trunk/LMDZ.GENERIC/libf/phystd/callsedim.F
r787 r858 1 1 SUBROUTINE callsedim(ngrid,nlay, ptimestep, 2 $ pplev,zlev, pt,3 & pq, pdqfi, pdqsed,pdqs_sed,nq ,rfall)2 & pplev,zlev, pt, pdt, 3 & pq, pdqfi, pdqsed,pdqs_sed,nq) 4 4 5 USE tracer_h 5 use radinc_h, only : naerkind 6 use radii_mod, only: h2o_reffrad 7 use aerosol_mod, only : iaero_h2o 8 USE tracer_h, only : igcm_co2_ice,igcm_h2o_ice,radius,rho_q 6 9 7 10 IMPLICIT NONE … … 33 36 c ---------- 34 37 35 INTEGER ngrid! number of horizontal grid points36 INTEGER nlay! number of atmospheric layers37 REAL ptimestep! physics time step (s)38 REAL pplev(ngrid,nlay+1)! pressure at inter-layers (Pa)39 REAL pt(ngrid,nlay)! temperature at mid-layer (K)40 REAL zlev(ngrid,nlay+1) ! altitude at layer boundaries41 42 43 c Traceurs : 44 integer nq ! number of tracers45 real pq(ngrid,nlay,nq) ! tracers (kg/kg)46 real pdqfi(ngrid,nlay,nq) ! tendency before sedimentation (kg/kg.s-1)47 real pdqsed(ngrid,nlay,nq) ! tendency due to sedimentation (kg/kg.s-1)48 real pdqs_sed(ngrid,nq) ! flux at surface (kg.m-2.s-1)38 integer,intent(in):: ngrid ! number of horizontal grid points 39 integer,intent(in):: nlay ! number of atmospheric layers 40 real,intent(in):: ptimestep ! physics time step (s) 41 real,intent(in):: pplev(ngrid,nlay+1) ! pressure at inter-layers (Pa) 42 real,intent(in):: pt(ngrid,nlay) ! temperature at mid-layer (K) 43 real,intent(in):: pdt(ngrid,nlay) ! tendency on temperature 44 real,intent(in):: zlev(ngrid,nlay+1) ! altitude at layer boundaries 45 integer,intent(in) :: nq ! number of tracers 46 real,intent(in) :: pq(ngrid,nlay,nq) ! tracers (kg/kg) 47 real,intent(in) :: pdqfi(ngrid,nlay,nq) ! tendency on tracers before 48 ! sedimentation (kg/kg.s-1) 49 50 real,intent(out) :: pdqsed(ngrid,nlay,nq) ! tendency due to sedimentation (kg/kg.s-1) 51 real,intent(out) :: pdqs_sed(ngrid,nq) ! flux at surface (kg.m-2.s-1) 49 52 50 53 c local: … … 53 56 INTEGER l,ig, iq 54 57 55 real zqi(ngrid,nlayermx,nq) ! to locally store tracers 56 real masse (ngrid,nlayermx) ! Layer mass (kg.m-2) 57 real epaisseur (ngrid,nlayermx) ! Layer thickness (m) 58 real wq(ngrid,nlayermx+1) ! displaced tracer mass (kg.m-2) 58 ! for particles with varying radii: 59 real reffrad(ngrid,nlay,naerkind) ! particle radius (m) 60 real nueffrad(ngrid,nlay,naerkind) ! aerosol effective radius variance 61 62 real zqi(ngrid,nlay,nq) ! to locally store tracers 63 real zt(ngrid,nlay) ! to locally store temperature (K) 64 real masse (ngrid,nlay) ! Layer mass (kg.m-2) 65 real epaisseur (ngrid,nlay) ! Layer thickness (m) 66 real wq(ngrid,nlay+1) ! displaced tracer mass (kg.m-2) 59 67 c real dens(ngrid,nlayermx) ! Mean density of the ice part. accounting for dust core 60 real rfall(ngrid,nlayermx)61 68 62 69 63 LOGICAL firstcall 64 SAVE firstcall 65 DATA firstcall/.true./ 70 LOGICAL,SAVE :: firstcall=.true. 66 71 67 72 c ** un petit test de coherence … … 80 85 masse(ig,l)=(pplev(ig,l) - pplev(ig,l+1)) /g 81 86 epaisseur(ig,l)= zlev(ig,l+1) - zlev(ig,l) 87 zt(ig,l)=pt(ig,l)+pdt(ig,l)*ptimestep 82 88 end do 83 89 end do … … 88 94 ! of general tracers.] 89 95 90 zqi(1:ngrid,1:nlay,1:nq) = 0.91 96 do iq=1,nq 92 97 if((radius(iq).gt.1.e-9).and.(iq.ne.igcm_co2_ice)) then 93 ! (no sedim for gases ; co2_ice sedim is done elsewhere)98 ! (no sedim for gases, and co2_ice sedim is done in condense_cloud) 94 99 95 100 ! store locally updated tracers … … 106 111 ! Water 107 112 if (water.and.(iq.eq.igcm_h2o_ice)) then 113 ! compute radii for h2o_ice 114 call h2o_reffrad(ngrid,zqi(1,1,igcm_h2o_ice),zt, 115 & reffrad(1,1,iaero_h2o),nueffrad(1,1,iaero_h2o)) 116 ! call sedimentation for h2o_ice 108 117 call newsedim(ngrid,nlay,ngrid*nlay,ptimestep, 109 & pplev,masse,epaisseur,pt,rfall,rho_q(iq),zqi,wq) 118 & pplev,masse,epaisseur,zt,reffrad(1,1,iaero_h2o), 119 & rho_q(iq),zqi(1,1,igcm_h2o_ice),wq) 110 120 111 121 ! General Case 112 122 else 113 123 call newsedim(ngrid,nlay,1,ptimestep, 114 & pplev,masse,epaisseur, pt,radius(iq),rho_q(iq),115 & zqi ,wq)124 & pplev,masse,epaisseur,zt,radius(iq),rho_q(iq), 125 & zqi(1,1,iq),wq) 116 126 endif 117 127 -
trunk/LMDZ.GENERIC/libf/phystd/condense_cloud.F90
r787 r858 4 4 piceco2,psolaralb,pemisurf, & 5 5 pdtc,pdtsrfc,pdpsrf,pduc,pdvc, & 6 pdqc ,reffrad)6 pdqc) 7 7 8 8 use radinc_h, only : naerkind … … 65 65 ! Arguments 66 66 67 INTEGER ngrid, nlayer, nq 68 69 REAL ptimestep 70 REAL pplay(ngrid,nlayer),pplev(ngrid,nlayer+1) 71 REAL pcapcal(ngrid) 72 REAL pt(ngrid,nlayer) 73 REAL ptsrf(ngrid) 74 REAL pphi(ngrid,nlayer) 75 REAL pdt(ngrid,nlayer),pdtsrf(ngrid),pdtc(ngrid,nlayer) 76 REAL pdtsrfc(ngrid),pdpsrf(ngrid) 77 ! REAL piceco2(ngrid),psolaralb(ngrid,2),pemisurf(ngrid) 78 REAL piceco2(ngrid),psolaralb(ngrid),pemisurf(ngrid) 79 80 REAL pu(ngrid,nlayer) , pv(ngrid,nlayer) 81 REAL pdu(ngrid,nlayer) , pdv(ngrid,nlayer) 82 REAL pduc(ngrid,nlayer) , pdvc(ngrid,nlayer) 83 REAL pq(ngrid,nlayer,nq),pdq(ngrid,nlayer,nq) 84 REAL pdqc(ngrid,nlayer,nq), pdqsc(ngrid,nq) 85 86 REAL reffrad(ngrid,nlayer,naerkind) 87 88 67 INTEGER,INTENT(IN) :: ngrid 68 INTEGER,INTENT(IN) :: nlayer 69 INTEGER,INTENT(IN) :: nq 70 REAL,INTENT(IN) :: ptimestep 71 REAL,INTENT(IN) :: pcapcal(ngrid) 72 REAL,INTENT(IN) :: pplay(ngrid,nlayer) 73 REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) 74 REAL,INTENT(IN) :: ptsrf(ngrid) 75 REAL,INTENT(IN) :: pt(ngrid,nlayer) 76 REAL,INTENT(IN) :: pphi(ngrid,nlayer) 77 REAL,INTENT(IN) :: pdt(ngrid,nlayer) 78 REAL,INTENT(IN) :: pdu(ngrid,nlayer) 79 REAL,INTENT(IN) :: pdv(ngrid,nlayer) 80 REAL,INTENT(IN) :: pdtsrf(ngrid) 81 REAL,INTENT(IN) :: pu(ngrid,nlayer) 82 REAL,INTENT(IN) :: pv(ngrid,nlayer) 83 REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) 84 REAL,INTENT(IN) :: pdq(ngrid,nlayer,nq) 85 REAL,INTENT(INOUT) :: piceco2(ngrid) 86 REAL,INTENT(OUT) :: psolaralb(ngrid) 87 REAL,INTENT(OUT) :: pemisurf(ngrid) 88 REAL,INTENT(OUT) :: pdtc(ngrid,nlayer) 89 REAL,INTENT(OUT) :: pdtsrfc(ngrid) 90 REAL,INTENT(OUT) :: pdpsrf(ngrid) 91 REAL,INTENT(OUT) :: pduc(ngrid,nlayer) 92 REAL,INTENT(OUT) :: pdvc(ngrid,nlayer) 93 REAL,INTENT(OUT) :: pdqc(ngrid,nlayer,nq) 89 94 90 95 !----------------------------------------------------------------------- … … 93 98 INTEGER l,ig,icap,ilay,it,iq 94 99 100 REAL reffrad(ngrid,nlayer) ! radius (m) of the co2 ice particles 95 101 REAL*8 zt(ngrid,nlayermx) 96 102 REAL zq(ngrid,nlayermx,nq) … … 123 129 ! Saved local variables 124 130 125 REAL latcond126 REAL ccond127 REAL cpice131 REAL,SAVE :: latcond=5.9e5 132 REAL,SAVE :: ccond 133 REAL,SAVE :: cpice=1000. 128 134 REAL,SAVE,ALLOCATABLE,DIMENSION(:) :: emisref 129 SAVE cpice 130 SAVE latcond,ccond 131 132 LOGICAL firstcall 133 SAVE firstcall 134 REAL SSUM 135 EXTERNAL SSUM 136 137 DATA latcond /5.9e5/ 138 DATA cpice /1000./ 139 DATA firstcall/.true./ 140 141 REAL CBRT 142 EXTERNAL CBRT 135 136 LOGICAL,SAVE :: firstcall=.true. 137 REAL,EXTERNAL :: SSUM 138 139 REAL,EXTERNAL :: CBRT 143 140 144 141 INTEGER,SAVE :: i_co2ice=0 ! co2 ice … … 154 151 ! Initializations 155 152 156 call zerophys(ngrid*nlayer*nq,pdqc)157 call zerophys(ngrid*nlayer*nq,pdtc)158 call zerophys(ngrid*nlayermx*nq,zq)159 call zerophys(ngrid*nlayermx,zt)153 pdqc(1:ngrid,1:nlayer,1:nq)=0 154 pdtc(1:ngrid,1:nlayer)=0 155 zq(1:ngrid,1:nlayer,1:nq)=0 156 zt(1:ngrid,1:nlayer)=0 160 157 161 158 ! Initialisation … … 320 317 DO ig=1,ngrid 321 318 322 reff = reffrad(ig,ilay ,iaero_co2)319 reff = reffrad(ig,ilay) 323 320 324 321 call stokes & -
trunk/LMDZ.GENERIC/libf/phystd/newsedim.F
r787 r858 23 23 ! --------- 24 24 25 INTEGER ngrid,nlay,naersize26 REAL ptimestep ! pas de temps physique (s)27 REAL pplev(ngrid,nlay+1) ! pression aux inter-couches (Pa)28 REAL pt(ngrid,nlay) ! temperature au centre des couches (K)29 real masse (ngrid,nlay) ! masse d'une couche (kg)30 real epaisseur (ngrid,nlay) ! epaisseur d'une couche (m)31 real rd(naersize) ! particle radius (m)32 real rho ! particle density (kg.m-3)33 34 35 c Traceurs : 36 real pqi(ngrid,nlay) ! traceur (e.g. ?/kg)37 real wq(ngrid,nlay+1) ! flux de traceur duranttimestep (?/m-2)25 integer,intent(in) :: ngrid ! number of atmospheric columns 26 integer,intent(in) :: nlay ! number of atmospheric layers 27 integer,intent(in) :: naersize ! number of particle sizes (1 or number 28 ! of grid boxes) 29 real,intent(in) :: ptimestep ! physics time step (s) 30 real,intent(in) :: pplev(ngrid,nlay+1) ! inter-layer pressures (Pa) 31 real,intent(in) :: pt(ngrid,nlay) ! mid-layer temperatures (K) 32 real,intent(in) :: masse (ngrid,nlay) ! atmospheric mass (kg) 33 real,intent(in) :: epaisseur (ngrid,nlay) ! thickness of atm. layers (m) 34 real,intent(in) :: rd(naersize) ! particle radius (m) 35 real,intent(in) :: rho ! particle density (kg.m-3) 36 real,intent(inout) :: pqi(ngrid,nlay) ! tracer (e.g. ?/kg) 37 real,intent(out) :: wq(ngrid,nlay+1) ! flux of tracer during timestep (?/m-2) 38 38 39 39 c local: … … 43 43 REAL rfall 44 44 45 LOGICAL firstcall 46 SAVE firstcall 45 LOGICAL,SAVE :: firstcall=.true. 47 46 48 47 c Traceurs : … … 53 52 real ptop, dztop, Ep, Stra 54 53 55 DATA firstcall/.true./56 54 57 55 c Physical constant 58 56 c ~~~~~~~~~~~~~~~~~ 59 REAL visc, molrad60 57 c Gas molecular viscosity (N.s.m-2) 61 data visc/1.e-5/! CO258 real,parameter :: visc=1.e-5 ! CO2 62 59 c Effective gas molecular radius (m) 63 data molrad/2.2e-10/! CO260 real,parameter :: molrad=2.2e-10 ! CO2 64 61 65 62 c local and saved variable 66 real a,b 67 save a,b 63 real,save :: a,b 68 64 69 65 c ** un petit test de coherence -
trunk/LMDZ.GENERIC/libf/phystd/physiq.F90
r804 r858 12 12 use gases_h 13 13 use radcommon_h, only: sigma 14 use radii_mod, only: h2o_reffrad, co2_reffrad 14 use radii_mod, only: h2o_reffrad, co2_reffrad, h2o_cloudrad 15 15 use aerosol_mod 16 16 use surfdat_h … … 132 132 ! inputs: 133 133 ! ------- 134 integer ngrid,nlayer,nq 135 real ptimestep 136 real pplev(ngrid,nlayer+1),pplay(ngrid,nlayer) 137 real pphi(ngrid,nlayer) 138 real pu(ngrid,nlayer),pv(ngrid,nlayer) 139 real pt(ngrid,nlayer),pq(ngrid,nlayer,nq) 140 real pw(ngrid,nlayer) ! pvervel transmitted by dyn3d 141 real zh(ngrid,nlayermx) ! potential temperature (K) 142 logical firstcall,lastcall 143 144 real pday 145 real ptime 146 logical tracerdyn 147 148 character*20 nametrac(nq) ! name of the tracer from dynamics 134 integer,intent(in) :: ngrid ! number of atmospheric columns 135 integer,intent(in) :: nlayer ! number of atmospheric layers 136 integer,intent(in) :: nq ! number of tracers 137 character*20,intent(in) :: nametrac(nq) ! name of the tracer from dynamics 138 logical,intent(in) :: firstcall ! signals first call to physics 139 logical,intent(in) :: lastcall ! signals last call to physics 140 real,intent(in) :: pday ! number of elapsed sols since reference Ls=0 141 real,intent(in) :: ptime ! "universal time", given as fraction of sol (e.g.: 0.5 for noon) 142 real,intent(in) :: ptimestep ! physics timestep (s) 143 real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa) 144 real,intent(in) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa) 145 real,intent(in) :: pphi(ngrid,nlayer) ! geopotential at mid-layer (m2s-2) 146 real,intent(in) :: pu(ngrid,nlayer) ! zonal wind component (m/s) 147 real,intent(in) :: pv(ngrid,nlayer) ! meridional wind component (m/s) 148 real,intent(in) :: pt(ngrid,nlayer) ! temperature (K) 149 real,intent(in) :: pq(ngrid,nlayer,nq) ! tracers (.../kg_of_air) 150 real,intent(in) :: pw(ngrid,nlayer) ! vertical velocity (m/s) 151 152 149 153 150 154 ! outputs: 151 155 ! -------- 152 156 ! physical tendencies 153 real pdu(ngrid,nlayer),pdv(ngrid,nlayer) 154 real pdt(ngrid,nlayer),pdq(ngrid,nlayer,nq) 155 real pdpsrf(ngrid) ! surface pressure tendency 156 157 real,intent(out) :: pdu(ngrid,nlayer) ! zonal wind tendency (m/s/s) 158 real,intent(out) :: pdv(ngrid,nlayer) ! meridional wind tendency (m/s/s) 159 real,intent(out) :: pdt(ngrid,nlayer) ! temperature tendency (K/s) 160 real,intent(out) :: pdq(ngrid,nlayer,nq) ! tracer tendencies (../kg/s) 161 real,intent(out) :: pdpsrf(ngrid) ! surface pressure tendency (Pa/s) 162 logical,intent(out) :: tracerdyn ! signal to the dynamics to advect tracers or not 157 163 158 164 ! Local saved variables: … … 190 196 ! for the "naerkind" optically active aerosols: 191 197 real aerosol(ngrid,nlayermx,naerkind) 198 real zh(ngrid,nlayermx) ! potential temperature (K) 192 199 193 200 character*80 fichier … … 390 397 real,allocatable,dimension(:,:,:),save :: reffrad ! aerosol effective radius (m) 391 398 real,allocatable,dimension(:,:,:),save :: nueffrad ! aerosol effective radius variance 392 real :: nueffrad_dummy(ngrid,nlayermx,naerkind) !! AS. This is temporary. Check below why.399 ! real :: nueffrad_dummy(ngrid,nlayermx,naerkind) !! AS. This is temporary. Check below why. 393 400 real :: reffh2oliq(ngrid,nlayermx) ! liquid water particles effective radius (m) 394 401 real :: reffh2oice(ngrid,nlayermx) ! water ice particles effective radius (m) 395 real reffH2O(ngrid,nlayermx)402 ! real reffH2O(ngrid,nlayermx) 396 403 real reffcol(ngrid,naerkind) 397 404 … … 406 413 !======================================================================= 407 414 408 !!! ALLOCATE SAVED ARRAYS409 IF (.not. ALLOCATED(tsurf)) ALLOCATE(tsurf(ngrid))410 IF (.not. ALLOCATED(tsoil)) ALLOCATE(tsoil(ngrid,nsoilmx))411 IF (.not. ALLOCATED(albedo)) ALLOCATE(albedo(ngrid))412 IF (.not. ALLOCATED(albedo0)) ALLOCATE(albedo0(ngrid))413 IF (.not. ALLOCATED(rnat)) ALLOCATE(rnat(ngrid))414 IF (.not. ALLOCATED(emis)) ALLOCATE(emis(ngrid))415 IF (.not. ALLOCATED(dtrad)) ALLOCATE(dtrad(ngrid,nlayermx))416 IF (.not. ALLOCATED(fluxrad_sky)) ALLOCATE(fluxrad_sky(ngrid))417 IF (.not. ALLOCATED(fluxrad)) ALLOCATE(fluxrad(ngrid))418 IF (.not. ALLOCATED(capcal)) ALLOCATE(capcal(ngrid))419 IF (.not. ALLOCATED(fluxgrd)) ALLOCATE(fluxgrd(ngrid))420 IF (.not. ALLOCATED(qsurf)) ALLOCATE(qsurf(ngrid,nq))421 IF (.not. ALLOCATED(q2)) ALLOCATE(q2(ngrid,nlayermx+1))422 IF (.not. ALLOCATED(ztprevious)) ALLOCATE(ztprevious(ngrid,nlayermx))423 IF (.not. ALLOCATED(cloudfrac)) ALLOCATE(cloudfrac(ngrid,nlayermx))424 IF (.not. ALLOCATED(totcloudfrac)) ALLOCATE(totcloudfrac(ngrid))425 IF (.not. ALLOCATED(qsurf_hist)) ALLOCATE(qsurf_hist(ngrid,nq))426 IF (.not. ALLOCATED(reffrad)) ALLOCATE(reffrad(ngrid,nlayermx,naerkind))427 IF (.not. ALLOCATED(nueffrad)) ALLOCATE(nueffrad(ngrid,nlayermx,naerkind))428 IF (.not. ALLOCATED(ice_initial)) ALLOCATE(ice_initial(ngrid))429 IF (.not. ALLOCATED(ice_min)) ALLOCATE(ice_min(ngrid))430 431 IF (.not. ALLOCATED(fluxsurf_lw)) ALLOCATE(fluxsurf_lw(ngrid))432 IF (.not. ALLOCATED(fluxsurf_sw)) ALLOCATE(fluxsurf_sw(ngrid))433 IF (.not. ALLOCATED(fluxtop_lw)) ALLOCATE(fluxtop_lw(ngrid))434 IF (.not. ALLOCATED(fluxabs_sw)) ALLOCATE(fluxabs_sw(ngrid))435 IF (.not. ALLOCATED(fluxtop_dn)) ALLOCATE(fluxtop_dn(ngrid))436 IF (.not. ALLOCATED(fluxdyn)) ALLOCATE(fluxdyn(ngrid))437 IF (.not. ALLOCATED(OLR_nu)) ALLOCATE(OLR_nu(ngrid,L_NSPECTI))438 IF (.not. ALLOCATED(OSR_nu)) ALLOCATE(OSR_nu(ngrid,L_NSPECTV))439 IF (.not. ALLOCATED(sensibFlux)) ALLOCATE(sensibFlux(ngrid))440 IF (.not. ALLOCATED(zdtlw)) ALLOCATE(zdtlw(ngrid,nlayermx))441 IF (.not. ALLOCATED(zdtsw)) ALLOCATE(zdtsw(ngrid,nlayermx))442 IF (.not. ALLOCATED(tau_col)) ALLOCATE(tau_col(ngrid))443 444 !=======================================================================445 446 415 ! 1. Initialisation 447 416 ! ----------------- … … 450 419 ! --------------------------------------- 451 420 if (firstcall) then 421 422 !!! ALLOCATE SAVED ARRAYS 423 ALLOCATE(tsurf(ngrid)) 424 ALLOCATE(tsoil(ngrid,nsoilmx)) 425 ALLOCATE(albedo(ngrid)) 426 ALLOCATE(albedo0(ngrid)) 427 ALLOCATE(rnat(ngrid)) 428 ALLOCATE(emis(ngrid)) 429 ALLOCATE(dtrad(ngrid,nlayermx)) 430 ALLOCATE(fluxrad_sky(ngrid)) 431 ALLOCATE(fluxrad(ngrid)) 432 ALLOCATE(capcal(ngrid)) 433 ALLOCATE(fluxgrd(ngrid)) 434 ALLOCATE(qsurf(ngrid,nq)) 435 ALLOCATE(q2(ngrid,nlayermx+1)) 436 ALLOCATE(ztprevious(ngrid,nlayermx)) 437 ALLOCATE(cloudfrac(ngrid,nlayermx)) 438 ALLOCATE(totcloudfrac(ngrid)) 439 ALLOCATE(qsurf_hist(ngrid,nq)) 440 ALLOCATE(reffrad(ngrid,nlayermx,naerkind)) 441 ALLOCATE(nueffrad(ngrid,nlayermx,naerkind)) 442 ALLOCATE(ice_initial(ngrid)) 443 ALLOCATE(ice_min(ngrid)) 444 ALLOCATE(fluxsurf_lw(ngrid)) 445 ALLOCATE(fluxsurf_sw(ngrid)) 446 ALLOCATE(fluxtop_lw(ngrid)) 447 ALLOCATE(fluxabs_sw(ngrid)) 448 ALLOCATE(fluxtop_dn(ngrid)) 449 ALLOCATE(fluxdyn(ngrid)) 450 ALLOCATE(OLR_nu(ngrid,L_NSPECTI)) 451 ALLOCATE(OSR_nu(ngrid,L_NSPECTV)) 452 ALLOCATE(sensibFlux(ngrid)) 453 ALLOCATE(zdtlw(ngrid,nlayermx)) 454 ALLOCATE(zdtsw(ngrid,nlayermx)) 455 ALLOCATE(tau_col(ngrid)) 452 456 453 457 !! this is defined in comsaison_h … … 710 714 zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw,fluxtop_lw, & 711 715 fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu, & 712 reffrad,nueffrad,tau_col,cloudfrac,totcloudfrac,&716 tau_col,cloudfrac,totcloudfrac, & 713 717 clearsky,firstcall,lastcall) 714 718 … … 724 728 zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_sw1,fluxtop_lw1, & 725 729 fluxabs_sw1,fluxtop_dn,OLR_nu1,OSR_nu1, & 726 reffrad,nueffrad,tau_col1,cloudfrac,totcloudfrac,&730 tau_col1,cloudfrac,totcloudfrac, & 727 731 clearsky,firstcall,lastcall) 728 732 clearsky = .false. ! just in case. … … 1003 1007 qsurf(1,igcm_co2_ice),albedo,emis, & 1004 1008 zdtc,zdtsurfc,pdpsrf,zduc,zdvc, & 1005 zdqc ,reffrad)1009 zdqc) 1006 1010 1007 1011 … … 1189 1193 1190 1194 call callsedim(ngrid,nlayer,ptimestep, & 1191 pplev,zzlev,pt,p q,pdq,zdqsed,zdqssed,nq,reffH2O)1195 pplev,zzlev,pt,pdt,pq,pdq,zdqsed,zdqssed,nq) 1192 1196 1193 1197 !------------------------- … … 1468 1472 reffcol(1:ngrid,1:naerkind)=0.0 1469 1473 if(co2cond.and.(iaero_co2.ne.0))then 1470 call co2_reffrad(ngrid,nq,zq,reffrad )1474 call co2_reffrad(ngrid,nq,zq,reffrad(1,1,iaero_co2)) 1471 1475 do ig=1,ngrid 1472 1476 reffcol(ig,iaero_co2) = SUM(zq(ig,1:nlayermx,igcm_co2_ice)*reffrad(ig,1:nlayermx,iaero_co2)*mass(ig,1:nlayermx)) … … 1474 1478 endif 1475 1479 if(water.and.(iaero_h2o.ne.0))then 1476 1477 !! AS: This modifies reffrad and nueffrad 1478 !! ... and those are saved in physiq 1479 !! To be conservative with previous versions, 1480 !! ... I put nueffrad_dummy so that nueffrad 1481 !! ... is not modified, but shouldn't it be modified actually? 1482 call h2o_reffrad(ngrid,nq,zq,zt,reffrad,nueffrad_dummy) 1480 call h2o_reffrad(ngrid,zq(1,1,igcm_h2o_ice),zt, & 1481 reffrad(1,1,iaero_h2o),nueffrad(1,1,iaero_h2o)) 1483 1482 do ig=1,ngrid 1484 1483 reffcol(ig,iaero_h2o) = SUM(zq(ig,1:nlayermx,igcm_h2o_ice)*reffrad(ig,1:nlayermx,iaero_h2o)*mass(ig,1:nlayermx)) -
trunk/LMDZ.GENERIC/libf/phystd/radii_mod.F90
r787 r858 35 35 use radinc_h, only: naerkind 36 36 use aerosol_mod 37 USE tracer_h38 Implicit none 39 40 #include "callkeys.h" 41 #include "dimensions.h" 42 #include "dimphys.h" 43 44 integer :: ngrid45 46 real, intent( inout) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K)47 real, intent( inout) :: nueffrad(ngrid,nlayermx,naerkind) !variance37 ! USE tracer_h 38 Implicit none 39 40 #include "callkeys.h" 41 #include "dimensions.h" 42 #include "dimphys.h" 43 44 integer,intent(in) :: ngrid 45 46 real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 47 real, intent(out) :: nueffrad(ngrid,nlayermx,naerkind) !variance 48 48 49 49 logical, save :: firstcall=.true. … … 119 119 120 120 !================================================================== 121 subroutine h2o_reffrad(ngrid,nq,pq,pt,reffrad,nueffrad) 121 ! subroutine h2o_reffrad(ngrid,nq,pq,pt,reffrad,nueffrad) 122 subroutine h2o_reffrad(ngrid,pq,pt,reffrad,nueffrad) 122 123 !================================================================== 123 124 ! Purpose … … 130 131 ! 131 132 !================================================================== 132 use radinc_h, only: naerkind133 ! use radinc_h, only: naerkind 133 134 use watercommon_h, Only: T_h2O_ice_liq,T_h2O_ice_clouds,rhowater,rhowaterice 134 use aerosol_mod, only : iaero_h2o135 USE tracer_h 135 ! use aerosol_mod, only : iaero_h2o 136 ! USE tracer_h, only: igcm_h2o_ice 136 137 Implicit none 137 138 … … 141 142 #include "comcstfi.h" 142 143 143 integer :: ngrid,nq 144 145 real, intent(in) :: pq(ngrid,nlayermx,nq) !tracer mixing ratios (kg/kg) 146 real, intent(in) :: pt(ngrid,nlayermx) !temperature (K) 147 real, intent(inout) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 148 real, intent(inout) :: nueffrad(ngrid,nlayermx,naerkind) ! dispersion 144 integer,intent(in) :: ngrid 145 ! intent,integer(in) :: nq 146 147 ! real, intent(in) :: pq(ngrid,nlayermx,nq) !tracer mixing ratios (kg/kg) 148 real, intent(in) :: pq(ngrid,nlayermx) !water ice mixing ratios (kg/kg) 149 real, intent(in) :: pt(ngrid,nlayermx) !temperature (K) 150 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosol radii 151 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii 152 real, intent(out) :: nueffrad(ngrid,nlayermx) ! dispersion 153 ! real, intent(out) :: nueffrad(ngrid,nlayermx,naerkind) ! dispersion 149 154 150 155 integer :: ig,l … … 158 163 zfice = 1.0 - (pt(ig,l)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds) 159 164 zfice = MIN(MAX(zfice,0.0),1.0) 160 reffrad(ig,l,iaero_h2o)= rad_h2o * (1.-zfice) + rad_h2o_ice * zfice 161 nueffrad(ig,l,iaero_h2o) = coef_chaud * (1.-zfice) + coef_froid * zfice 165 ! reffrad(ig,l,iaero_h2o)= rad_h2o * (1.-zfice) + rad_h2o_ice * zfice 166 reffrad(ig,l)= rad_h2o * (1.-zfice) + rad_h2o_ice * zfice 167 ! nueffrad(ig,l,iaero_h2o) = coef_chaud * (1.-zfice) + coef_froid * zfice 168 nueffrad(ig,l) = coef_chaud * (1.-zfice) + coef_froid * zfice 162 169 enddo 163 170 enddo … … 167 174 zfice = 1.0 - (pt(ig,l)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds) 168 175 zfice = MIN(MAX(zfice,0.0),1.0) 169 zrad_liq = CBRT( 3*pq(ig,l,igcm_h2o_ice)/(4*Nmix_h2o*pi*rhowater) ) 170 zrad_ice = CBRT( 3*pq(ig,l,igcm_h2o_ice)/(4*Nmix_h2o_ice*pi*rhowaterice) ) 171 nueffrad(ig,l,iaero_h2o) = coef_chaud * (1.-zfice) + coef_froid * zfice 176 ! zrad_liq = CBRT( 3*pq(ig,l,igcm_h2o_ice)/(4*Nmix_h2o*pi*rhowater) ) 177 zrad_liq = CBRT( 3*pq(ig,l)/(4*Nmix_h2o*pi*rhowater) ) 178 ! zrad_ice = CBRT( 3*pq(ig,l,igcm_h2o_ice)/(4*Nmix_h2o_ice*pi*rhowaterice) ) 179 zrad_ice = CBRT( 3*pq(ig,l)/(4*Nmix_h2o_ice*pi*rhowaterice) ) 180 ! nueffrad(ig,l,iaero_h2o) = coef_chaud * (1.-zfice) + coef_froid * zfice 181 nueffrad(ig,l) = coef_chaud * (1.-zfice) + coef_froid * zfice 172 182 zrad = zrad_liq * (1.-zfice) + zrad_ice * zfice 173 reffrad(ig,l,iaero_h2o) = min(max(zrad,1.e-6),100.e-6) 183 ! reffrad(ig,l,iaero_h2o) = min(max(zrad,1.e-6),100.e-6) 184 reffrad(ig,l) = min(max(zrad,1.e-6),100.e-6) 174 185 enddo 175 186 enddo … … 193 204 !================================================================== 194 205 use watercommon_h, Only: rhowater,rhowaterice 195 USE tracer_h206 ! USE tracer_h 196 207 Implicit none 197 208 … … 201 212 #include "comcstfi.h" 202 213 203 integer :: ngrid214 integer,intent(in) :: ngrid 204 215 205 216 real, intent(in) :: pql(ngrid,nlayermx) !condensed water mixing ratios (kg/kg) … … 236 247 ! 237 248 !================================================================== 238 use radinc_h, only: naerkind239 use aerosol_mod, only : iaero_co2240 USE tracer_h 249 ! use radinc_h, only: naerkind 250 ! use aerosol_mod, only : iaero_co2 251 USE tracer_h, only:igcm_co2_ice,rho_co2 241 252 Implicit none 242 253 … … 246 257 #include "comcstfi.h" 247 258 248 integer :: ngrid,nq259 integer,intent(in) :: ngrid,nq 249 260 250 261 real, intent(in) :: pq(ngrid,nlayermx,nq) !tracer mixing ratios (kg/kg) 251 real, intent(inout) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 262 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 263 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosols radii (K) 252 264 253 265 integer :: ig,l … … 258 270 259 271 if (radfixed) then 260 reffrad(1:ngrid,1:nlayermx,iaero_co2) = 5.e-5 ! CO2 ice 272 ! reffrad(1:ngrid,1:nlayermx,iaero_co2) = 5.e-5 ! CO2 ice 273 reffrad(1:ngrid,1:nlayermx) = 5.e-5 ! CO2 ice 261 274 else 262 275 do l=1,nlayermx 263 276 do ig=1,ngrid 264 277 zrad = CBRT( 3*pq(ig,l,igcm_co2_ice)/(4*Nmix_co2*pi*rho_co2) ) 265 reffrad(ig,l,iaero_co2) = min(max(zrad,1.e-6),100.e-6) 278 ! reffrad(ig,l,iaero_co2) = min(max(zrad,1.e-6),100.e-6) 279 reffrad(ig,l) = min(max(zrad,1.e-6),100.e-6) 266 280 enddo 267 281 enddo … … 285 299 ! 286 300 !================================================================== 287 use radinc_h, only: naerkind 288 use aerosol_mod, only : iaero_dust 289 Implicit none 290 291 #include "callkeys.h" 292 #include "dimensions.h" 293 #include "dimphys.h" 294 295 integer :: ngrid 296 297 real, intent(inout) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 301 ! use radinc_h, only: naerkind 302 ! use aerosol_mod, only : iaero_dust 303 Implicit none 304 305 #include "callkeys.h" 306 #include "dimensions.h" 307 #include "dimphys.h" 308 309 integer,intent(in) :: ngrid 310 311 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 312 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosols radii (K) 298 313 299 reffrad(1:ngrid,1:nlayermx,iaero_dust) = 2.e-6 ! dust 314 ! reffrad(1:ngrid,1:nlayermx,iaero_dust) = 2.e-6 ! dust 315 reffrad(1:ngrid,1:nlayermx) = 2.e-6 ! dust 300 316 301 317 end subroutine dust_reffrad … … 315 331 ! 316 332 !================================================================== 317 use radinc_h, only: naerkind 318 use aerosol_mod, only : iaero_h2so4 319 Implicit none 320 321 #include "callkeys.h" 322 #include "dimensions.h" 323 #include "dimphys.h" 324 325 integer :: ngrid 326 327 real, intent(inout) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 333 ! use radinc_h, only: naerkind 334 ! use aerosol_mod, only : iaero_h2so4 335 Implicit none 336 337 #include "callkeys.h" 338 #include "dimensions.h" 339 #include "dimphys.h" 340 341 integer,intent(in) :: ngrid 342 343 ! real, intent(out) :: reffrad(ngrid,nlayermx,naerkind) !aerosols radii (K) 344 real, intent(out) :: reffrad(ngrid,nlayermx) !aerosols radii (K) 328 345 329 reffrad(1:ngrid,1:nlayermx,iaero_h2so4) = 1.e-6 ! h2so4 346 ! reffrad(1:ngrid,1:nlayermx,iaero_h2so4) = 1.e-6 ! h2so4 347 reffrad(1:ngrid,1:nlayermx) = 1.e-6 ! h2so4 330 348 331 349 end subroutine h2so4_reffrad -
trunk/LMDZ.GENERIC/libf/phystd/rain.F90
r787 r858 28 28 #include "callkeys.h" 29 29 30 integer ngrid,nq31 32 ! Pre-arguments (for universal model)33 real pq(ngrid,nlayermx,nq) ! tracer (kg/kg)34 real qsurf(ngrid,nq) ! tracer at the surface (kg.m-2)35 REAL pdt(ngrid,nlayermx),pdq(ngrid,nlayermx,nq)36 37 real dqrain(ngrid,nlayermx,nq) ! tendency of H2O precipitation (kg/kg.s-1)38 real dqsrain(ngrid) ! rain flux at the surface (kg.m-2.s-1)39 real dqssnow(ngrid) ! snow flux at the surface (kg.m-2.s-1)40 REAL d_t(ngrid,nlayermx) ! temperature increment41 42 30 ! Arguments 43 REAL ptimestep ! time interval 44 REAL pplev(ngrid,nlayermx+1) ! inter-layer pressure 45 REAL pplay(ngrid,nlayermx) ! mid-layer pressure 46 REAL t(ngrid,nlayermx) ! input temperature (K) 31 integer,intent(in) :: ngrid ! number of atmospherci columns 32 integer,intent(in) :: nq ! number of tracers 33 real,intent(in) :: ptimestep ! time interval 34 real,intent(in) :: pplev(ngrid,nlayermx+1) ! inter-layer pressure (Pa) 35 real,intent(in) :: pplay(ngrid,nlayermx) ! mid-layer pressure (Pa) 36 real,intent(in) :: t(ngrid,nlayermx) ! input temperature (K) 37 real,intent(in) :: pdt(ngrid,nlayermx) ! input tendency on temperature (K/s) 38 real,intent(in) :: pq(ngrid,nlayermx,nq) ! tracers (kg/kg) 39 real,intent(in) :: pdq(ngrid,nlayermx,nq) ! input tendency on tracers 40 real,intent(out) :: d_t(ngrid,nlayermx) ! temperature tendency (K/s) 41 real,intent(out) :: dqrain(ngrid,nlayermx,nq) ! tendency of H2O precipitation (kg/kg.s-1) 42 real,intent(out) :: dqsrain(ngrid) ! rain flux at the surface (kg.m-2.s-1) 43 real,intent(out) :: dqssnow(ngrid) ! snow flux at the surface (kg.m-2.s-1) 44 real,intent(in) :: rneb(ngrid,nlayermx) ! cloud fraction 45 47 46 REAL zt(ngrid,nlayermx) ! working temperature (K) 48 47 REAL ql(ngrid,nlayermx) ! liquid water (Kg/Kg) 49 48 REAL q(ngrid,nlayermx) ! specific humidity (Kg/Kg) 50 REAL rneb(ngrid,nlayermx) ! cloud fraction51 49 REAL d_q(ngrid,nlayermx) ! water vapor increment 52 50 REAL d_ql(ngrid,nlayermx) ! liquid water / ice increment 53 51 54 52 ! Subroutine options 55 REAL seuil_neb ! Nebulosity threshold 56 PARAMETER (seuil_neb=0.001) 53 REAL,PARAMETER :: seuil_neb=0.001 ! Nebulosity threshold 57 54 58 55 INTEGER,save :: precip_scheme ! id number for precipitaion scheme … … 67 64 REAL,SAVE :: c1 68 65 69 INTEGER ninter 70 PARAMETER (ninter=5) 71 72 logical evap_prec ! Does the rain evaporate? 73 parameter(evap_prec=.true.) 66 INTEGER,PARAMETER :: ninter=5 67 68 logical,parameter :: evap_prec=.true. ! Does the rain evaporate? 74 69 75 70 ! for simple scheme 76 real t_crit 77 PARAMETER (t_crit=218.0) 71 real,parameter :: t_crit=218.0 78 72 real lconvert 79 73 … … 100 94 INTEGER, SAVE :: i_ice=0 ! water ice 101 95 102 LOGICAL firstcall 103 SAVE firstcall 96 LOGICAL,SAVE :: firstcall=.true. 104 97 105 98 ! Online functions … … 108 101 fall2v (zzz) =10.6 * ((zzz)**0.31) !for use with radii 109 102 110 DATA firstcall /.true./111 103 112 104 IF (firstcall) THEN -
trunk/LMDZ.GENERIC/libf/phystd/totalcloudfrac.F90
r787 r858 4 4 use comdiurn_h 5 5 USE comgeomfi_h 6 USE tracer_h 6 USE tracer_h, only: igcm_h2o_ice 7 7 implicit none 8 8 … … 24 24 #include "callkeys.h" 25 25 26 integer ngrid, nq 26 integer,intent(in) :: ngrid ! number of atmospheric columns 27 integer,intent(in) :: nq ! number of tracers 28 real,intent(in) :: rneb(ngrid,nlayermx) ! cloud fraction 29 real,intent(out) :: totalrneb(ngrid) ! total cloud fraction 30 real,intent(in) :: pplev(ngrid,nlayermx+1) ! inter-layer pressure (Pa) 31 real,intent(in) :: pq(ngrid,nlayermx,nq) ! tracers (.../kg_of_air) 32 real,intent(in) :: tau(ngrid,nlayermx) 27 33 28 real rneb(ngrid,nlayermx) ! cloud fraction29 real pplev(ngrid,nlayermx+1)30 real pq(ngrid,nlayermx,nq)31 real tau(ngrid,nlayermx)32 33 real totalrneb(ngrid) ! total cloud fraction34 34 real, dimension(nlayermx+1) :: masse 35 35 integer, parameter :: recovery=7
Note: See TracChangeset
for help on using the changeset viewer.