Changeset 334 for trunk/LMDZ.MARS/libf/phymars
- Timestamp:
- Oct 28, 2011, 5:00:48 PM (13 years ago)
- Location:
- trunk/LMDZ.MARS/libf/phymars
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/physiq.F
r333 r334 198 198 REAL nuice(ngridmx,nlayermx) ! Estimated effective variance 199 199 ! of the size distribution 200 REAL surfdust(ngridmx,nlayermx) ! dust surface area (micron2/cm3, if photochemistry) 201 REAL surfice(ngridmx,nlayermx) ! ice surface area (micron2/cm3, if photochemistry) 200 202 201 203 c Variables used by the slope model … … 306 308 real rho(ngridmx,nlayermx) ! density 307 309 real vmr(ngridmx,nlayermx) ! volume mixing ratio 310 !real colden(ngridmx,nqmx) ! vertical column !FL 308 311 REAL mtot(ngridmx) ! Total mass of water vapor (kg/m2) 309 312 REAL icetot(ngridmx) ! Total mass of water ice (kg/m2) … … 918 921 & pq,pdq,zdqcloud,zdqscloud,zdtcloud, 919 922 & nq,naerkind,tau, 920 & ccn,rdust,rice,nuice) 923 & ccn,rdust,rice,nuice, 924 & surfdust,surfice) 925 921 926 if (activice) then 922 927 c Temperature variation due to latent heat release … … 955 960 c -------------- 956 961 IF (photochem .or. thermochem) then 957 call calchim(ptimestep,pplay,pplev,pt,pdt,dist_sol,mu0,958 & zzlay,zday,pq,pdq,rice,959 & zdqchim,zdqschim,zdqcloud,zdqscloud)960 !NB: Photochemistry includes condensation of H2O2 962 call calchim(ptimestep,pplay,pplev,pt,pdt,dist_sol,mu0, 963 $ zzlev,zzlay,zday,pq,pdq,zdqchim,zdqschim, 964 $ zdqcloud,zdqscloud,tauref,co2ice, 965 $ pu,pdu,pv,pdv,surfdust,surfice) 961 966 962 967 ! increment values of tracers: … … 1409 1414 . "Volume mixing ratio","mol/mol",3,vmr) 1410 1415 endif 1411 enddo 1412 endif ! of if (thermochem.or.photochem) 1413 1414 endif ! of if (tracer) 1416 ! do ig = 1,ngrid 1417 ! colden(ig,iq) = 0. !FL 1418 ! end do 1419 ! do l=1,nlayer !FL 1420 ! do ig=1,ngrid !FL 1421 ! colden(ig,iq) = colden(ig,iq) + zq(ig,l,iq) !FL 1422 ! $ *(pplev(ig,l)-pplev(ig,l+1)) !FL 1423 ! $ *6.022e22/(mmol(iq)*g) !FL 1424 ! end do !FL 1425 ! end do !FL 1426 ! call wstats(ngrid,"c_"//trim(noms(iq)), !FL 1427 ! $ "column","mol cm-2",2,colden(1,iq)) !FL 1428 end do 1429 end if ! of if (thermochem.or.photochem) 1430 1431 end if ! of if (tracer) 1415 1432 1416 1433 IF(lastcall) THEN … … 1533 1550 call WRITEDIAGFI(ngrid,"co2col","CO2 column","kg.m-2",2, 1534 1551 & co2col) 1552 !!!!! FL 1553 ! do iq = 1,nq 1554 ! if (noms(iq) .ne. "dust_mass" .and. 1555 ! $ noms(iq) .ne. "dust_number") then 1556 ! call writediagfi(ngrid,"c_"//trim(noms(iq)), 1557 ! $ "column","mol cm-2",2,colden(1,iq)) 1558 ! end if 1559 ! end do 1560 !!!!! FL 1535 1561 endif ! of if (tracer.and.(igcm_co2.ne.0)) 1536 1562 -
trunk/LMDZ.MARS/libf/phymars/watercloud.F
r120 r334 3 3 & pq,pdq,pdqcloud,pdqscloud,pdtcloud, 4 4 & nq,naersize,tau, 5 & ccn,rdust,rice,nuice) 5 & ccn,rdust,rice,nuice, 6 & surfdust,surfice) 6 7 IMPLICIT NONE 7 8 … … 70 71 ! of the size distribution 71 72 73 REAL surfdust(ngrid,nlay) ! dust surface area (micron2/cm3, if photochemistry) 74 REAL surfice(ngrid,nlay) ! ice surface area (micron2/cm3, if photochemistry) 75 72 76 c local: 73 77 c ------ … … 85 89 REAL masse (ngridmx,nlayermx) 86 90 REAL epaisseur (ngridmx,nlayermx) 91 REAL npart(ngridmx,nlayermx) ! Cloud condensation nuclei (#.cm-3) 87 92 ! REAL rfinal ! Ice crystal radius after condensation(m) 88 93 REAL*8 seq ! Equilibrium saturation ration (accounting for curvature effect) … … 179 184 enddo ! of dol=1,nlay 180 185 186 if (photochem) then 187 c computation of dust and ice surface area (micron2/cm3) 188 c for heterogeneous chemistry 189 190 do l = 1,nlay 191 do ig = 1,ngrid 192 c 193 c npart: number density of ccn in #/cm3 194 c 195 npart(ig,l) = 1.e-6*ccn(ig,l) 196 $ *masse(ig,l)/epaisseur(ig,l) 197 c 198 c dust and ice surface area 199 c 200 surfdust(ig,l) = npart(ig,l)*4.*pi*1.e12*rdust(ig,l)**2 201 c 202 if (rice(ig,l) .ge. rdust(ig,l)) then 203 surfice(ig,l) = npart(ig,l)*4.*pi*1.e12*rice(ig,l)**2 204 surfdust(ig,l) = 0. 205 else 206 surfice(ig,l) = 0. 207 end if 208 end do ! of do ig=1,ngrid 209 end do ! of do l=1,nlay 210 end if ! of photochem 211 c 181 212 pdqscloud(1:ngrid,1:nq)=0 182 213 pdqcloud(1:ngrid,1:nlay,1:nq)=0
Note: See TracChangeset
for help on using the changeset viewer.