| 1 | subroutine improvedCO2clouds(ngrid,nlay,microtimestep, |
|---|
| 2 | & pplay,pplev,pteff,sum_subpdt, |
|---|
| 3 | & pqeff,sum_subpdq,subpdqcloudco2,subpdtcloudco2, |
|---|
| 4 | & nq,tauscaling, |
|---|
| 5 | & memdMMccn,memdMMh2o,memdNNccn) |
|---|
| 6 | USE comcstfi_h, only: pi, g, cpp |
|---|
| 7 | USE updaterad, only: updaterice_micro, updaterice_microco2 |
|---|
| 8 | use tracer_mod, only: igcm_dust_mass, igcm_dust_number, rho_dust, |
|---|
| 9 | & igcm_h2o_ice, igcm_ccn_mass, |
|---|
| 10 | & igcm_ccn_number, nuice_sed, |
|---|
| 11 | & igcm_co2, igcm_co2_ice, igcm_ccnco2_mass, |
|---|
| 12 | & igcm_ccnco2_number, nuiceco2_sed, |
|---|
| 13 | & rho_ice_co2 |
|---|
| 14 | use conc_mod, only: mmean |
|---|
| 15 | use datafile_mod, only: datadir |
|---|
| 16 | |
|---|
| 17 | implicit none |
|---|
| 18 | |
|---|
| 19 | c------------------------------------------------------------------ |
|---|
| 20 | c This routine is used to form CO2 clouds when a parcel of the GCM is |
|---|
| 21 | c saturated. It includes the ability to have supersaturation, a |
|---|
| 22 | c computation of the nucleation rates, growthrates and the |
|---|
| 23 | c scavenging of dust particles by clouds. |
|---|
| 24 | c It is worth noting that the amount of dust is computed using the |
|---|
| 25 | c dust optical depth computed in aeropacity.F. That's why |
|---|
| 26 | c the variable called "tauscaling" is used to convert |
|---|
| 27 | c pq(dust_mass) and pq(dust_number), which are relative |
|---|
| 28 | c quantities, to absolute and realistic quantities stored in zq. |
|---|
| 29 | c This has to be done to convert the inputs into absolute |
|---|
| 30 | c values, but also to convert the outputs back into relative |
|---|
| 31 | c values which are then used by the sedimentation and advection |
|---|
| 32 | c schemes. |
|---|
| 33 | c CO2 ice particles can nucleate on both dust and on water ice particles |
|---|
| 34 | c When CO2 ice is deposited onto a water ice particles, the particle is |
|---|
| 35 | c removed from the water tracers. |
|---|
| 36 | c Memory of the origin of the co2 particles is kept and thus the |
|---|
| 37 | c water cycle shouldn't be modified by this. |
|---|
| 38 | cWARNING: no sedimentation of the water ice origin is performed |
|---|
| 39 | c in the microphysical timestep in co2cloud.F. |
|---|
| 40 | |
|---|
| 41 | c Authors of the water ice clouds microphysics |
|---|
| 42 | c J.-B. Madeleine, based on the work by Franck Montmessin |
|---|
| 43 | c (October 2011) |
|---|
| 44 | c T. Navarro, debug,correction, new scheme (October-April 2011) |
|---|
| 45 | c A. Spiga, optimization (February 2012) |
|---|
| 46 | c Adaptation for CO2 clouds by Joachim Audouard (09/16), based on the work |
|---|
| 47 | c of Constantino Listowski |
|---|
| 48 | c There is an energy limit to how much co2 can sublimate/condensate. It is |
|---|
| 49 | c defined by the difference of the GCM temperature with the co2 condensation |
|---|
| 50 | c temperature. |
|---|
| 51 | c Warning: |
|---|
| 52 | c If meteoritic particles are activated and turn into co2 ice particles, |
|---|
| 53 | c then they will be reversed in the dust tracers if the cloud sublimates |
|---|
| 54 | |
|---|
| 55 | c------------------------------------------------------------------ |
|---|
| 56 | include "callkeys.h" |
|---|
| 57 | include "microphys.h" |
|---|
| 58 | c------------------------------------------------------------------ |
|---|
| 59 | c Arguments: |
|---|
| 60 | |
|---|
| 61 | INTEGER,INTENT(in) :: ngrid,nlay |
|---|
| 62 | integer,intent(in) :: nq ! number of tracers |
|---|
| 63 | REAL,INTENT(in) :: microtimestep ! physics time step (s) |
|---|
| 64 | REAL,INTENT(in) :: pplay(ngrid,nlay) ! mid-layer pressure (Pa) |
|---|
| 65 | REAL,INTENT(in) :: pplev(ngrid,nlay+1) ! inter-layer pressure (Pa) |
|---|
| 66 | REAL,INTENT(in) :: pteff(ngrid,nlay) ! temperature at the middle of the |
|---|
| 67 | ! layers (K) |
|---|
| 68 | REAL,INTENT(in) :: sum_subpdt(ngrid,nlay) ! tendency on temperature from |
|---|
| 69 | ! previous physical parametrizations |
|---|
| 70 | REAL,INTENT(in) :: pqeff(ngrid,nlay,nq) ! tracers (kg/kg) |
|---|
| 71 | REAL,INTENT(in) :: sum_subpdq(ngrid,nlay,nq) ! tendencies on tracers |
|---|
| 72 | ! before condensation (kg/kg.s-1) |
|---|
| 73 | REAL,INTENT(in) :: tauscaling(ngrid) ! Convertion factor for qdust and Ndust |
|---|
| 74 | c Outputs: |
|---|
| 75 | REAL,INTENT(out) :: subpdqcloudco2(ngrid,nlay,nq) ! tendency on tracers |
|---|
| 76 | ! due to CO2 condensation (kg/kg.s-1) |
|---|
| 77 | ! condensation si igcm_co2_ice |
|---|
| 78 | REAL,INTENT(out) :: subpdtcloudco2(ngrid,nlay) ! tendency on temperature due |
|---|
| 79 | ! to latent heat |
|---|
| 80 | |
|---|
| 81 | c------------------------------------------------------------------ |
|---|
| 82 | c Local variables: |
|---|
| 83 | LOGICAL,SAVE :: firstcall=.true. |
|---|
| 84 | REAL*8 derf ! Error function |
|---|
| 85 | INTEGER ig,l,i |
|---|
| 86 | |
|---|
| 87 | real masse (ngrid,nlay) ! Layer mass (kg.m-2) |
|---|
| 88 | REAL rice(ngrid,nlay) ! Water Ice mass mean radius (m) |
|---|
| 89 | ! used for nucleation of CO2 on ice-coated ccns |
|---|
| 90 | REAL rccnh2o(ngrid,nlay) ! Water Ice mass mean radius (m) |
|---|
| 91 | REAL zq(ngrid,nlay,nq) ! local value of tracers |
|---|
| 92 | REAL zq0(ngrid,nlay,nq) ! local initial value of tracers |
|---|
| 93 | REAL zt(ngrid,nlay) ! local value of temperature |
|---|
| 94 | REAL zqsat(ngrid,nlay) ! saturation vapor pressure for CO2 |
|---|
| 95 | real tcond(ngrid,nlay) |
|---|
| 96 | real zqco2(ngrid,nlay) |
|---|
| 97 | REAL lw !Latent heat of sublimation (J.kg-1) |
|---|
| 98 | REAL,save :: l0,l1,l2,l3,l4 |
|---|
| 99 | DOUBLE PRECISION dMice ! mass of condensed ice |
|---|
| 100 | DOUBLE PRECISION sumcheck |
|---|
| 101 | DOUBLE PRECISION facteurmax!for energy limit on mass growth |
|---|
| 102 | DOUBLE PRECISION pco2,psat ! Co2 vapor partial pressure (Pa) |
|---|
| 103 | DOUBLE PRECISION satu ! Co2 vapor saturation ratio over ice |
|---|
| 104 | DOUBLE PRECISION Mo,No,No_dust,Mo_dust |
|---|
| 105 | DOUBLE PRECISION Rn, Rm, dev2,dev3, n_derf, m_derf |
|---|
| 106 | DOUBLE PRECISION memdMMccn(ngrid,nlay) |
|---|
| 107 | DOUBLE PRECISION memdMMh2o(ngrid,nlay) |
|---|
| 108 | DOUBLE PRECISION memdNNccn(ngrid,nlay) |
|---|
| 109 | |
|---|
| 110 | ! Radius used by the microphysical scheme (m) |
|---|
| 111 | DOUBLE PRECISION n_aer(nbinco2_cld) ! number concentration volume-1 of particle/each size bin |
|---|
| 112 | DOUBLE PRECISION m_aer(nbinco2_cld) ! mass mixing ratio of particle/each size bin |
|---|
| 113 | DOUBLE PRECISION m_aer_h2oice2(nbinco2_cld) ! mass mixing ratio of particle/each size bin |
|---|
| 114 | |
|---|
| 115 | DOUBLE PRECISION n_aer_h2oice(nbinco2_cld) ! Same - for CO2 nucleation |
|---|
| 116 | DOUBLE PRECISION m_aer_h2oice(nbinco2_cld) ! Same - for CO2 nucleation |
|---|
| 117 | DOUBLE PRECISION rad_h2oice(nbinco2_cld) |
|---|
| 118 | |
|---|
| 119 | c REAL*8 sigco2 ! Co2-ice/air surface tension (N.m) |
|---|
| 120 | c EXTERNAL sigco2 |
|---|
| 121 | |
|---|
| 122 | DOUBLE PRECISION dN,dM, dNh2o, dMh2o, dNN,dMM,dNNh2o,dMMh2o |
|---|
| 123 | DOUBLE PRECISION dMh2o_ice,dMh2o_ccn |
|---|
| 124 | |
|---|
| 125 | DOUBLE PRECISION rate(nbinco2_cld) ! nucleation rate |
|---|
| 126 | DOUBLE PRECISION rateh2o(nbinco2_cld) ! nucleation rate |
|---|
| 127 | REAL seq |
|---|
| 128 | DOUBLE PRECISION rho_ice_co2T(ngrid,nlay) |
|---|
| 129 | DOUBLE PRECISION riceco2(ngrid,nlay) ! CO2Ice mean radius (m) |
|---|
| 130 | REAL rhocloud(ngrid,nlay) ! Cloud density (kg.m-3) |
|---|
| 131 | |
|---|
| 132 | REAL rhocloudco2(ngrid,nlay) ! Cloud density (kg.m-3) |
|---|
| 133 | REAL rdust(ngrid,nlay) ! Dust geometric mean radius (m) |
|---|
| 134 | |
|---|
| 135 | c REAL res ! Resistance growth |
|---|
| 136 | DOUBLE PRECISION Ic_rice ! Mass transfer rate CO2 ice crystal |
|---|
| 137 | DOUBLE PRECISION ratioh2o_ccn |
|---|
| 138 | DOUBLE PRECISION vo2co2 |
|---|
| 139 | |
|---|
| 140 | c Parameters of the size discretization used by the microphysical scheme |
|---|
| 141 | DOUBLE PRECISION, PARAMETER :: rmin_cld = 1.e-9 ! Minimum radius (m) |
|---|
| 142 | DOUBLE PRECISION, PARAMETER :: rmax_cld = 5.e-6 ! Maximum radius (m) |
|---|
| 143 | DOUBLE PRECISION, PARAMETER :: rbmin_cld =1.e-10 ! Minimum bounary radius (m) |
|---|
| 144 | DOUBLE PRECISION, PARAMETER :: rbmax_cld = 2.e-4 ! Maximum boundary radius (m) |
|---|
| 145 | DOUBLE PRECISION vrat_cld ! Volume ratio |
|---|
| 146 | DOUBLE PRECISION rb_cldco2(nbinco2_cld+1) ! boundary values of each rad_cldco2 bin (m) |
|---|
| 147 | SAVE rb_cldco2 |
|---|
| 148 | DOUBLE PRECISION dr_cld(nbinco2_cld) ! width of each rad_cldco2 bin (m) |
|---|
| 149 | DOUBLE PRECISION vol_cld(nbinco2_cld) ! particle volume for each bin (m3) |
|---|
| 150 | |
|---|
| 151 | DOUBLE PRECISION Proba,Masse_atm,drsurdt,reff,Probah2o |
|---|
| 152 | REAL sigma_iceco2 ! Variance of the co2 ice and CCN distributions |
|---|
| 153 | SAVE sigma_iceco2 |
|---|
| 154 | REAL sigma_ice ! Variance of the h2o ice and CCN distributions |
|---|
| 155 | SAVE sigma_ice |
|---|
| 156 | DOUBLE PRECISION Niceco2,Qccnco2,Nccnco2 |
|---|
| 157 | |
|---|
| 158 | ! Variables for the meteoritic flux: |
|---|
| 159 | integer,parameter :: nbin_meteor=100 |
|---|
| 160 | integer,parameter :: nlev_meteor=130 |
|---|
| 161 | double precision meteor_ccn(ngrid,nlay,100) !100=nbinco2_cld !!! |
|---|
| 162 | double precision,save :: meteor(130,100) |
|---|
| 163 | double precision mtemp(100),pression_meteor(130) |
|---|
| 164 | logical file_ok |
|---|
| 165 | integer read_ok |
|---|
| 166 | integer nelem,lebon1,lebon2 |
|---|
| 167 | double precision :: ltemp1(130),ltemp2(130) |
|---|
| 168 | integer ibin,j |
|---|
| 169 | integer,parameter :: uMeteor=666 |
|---|
| 170 | |
|---|
| 171 | IF (firstcall) THEN |
|---|
| 172 | !============================================================= |
|---|
| 173 | ! 0. Definition of the size grid |
|---|
| 174 | !============================================================= |
|---|
| 175 | c rad_cldco2 is the primary radius grid used for microphysics computation. |
|---|
| 176 | c The grid spacing is computed assuming a constant volume ratio |
|---|
| 177 | c between two consecutive bins; i.e. vrat_cld. |
|---|
| 178 | c vrat_cld is determined from the boundary values of the size grid: |
|---|
| 179 | c rmin_cld and rmax_cld. |
|---|
| 180 | c The rb_cldco2 array contains the boundary values of each rad_cldco2 bin. |
|---|
| 181 | c dr_cld is the width of each rad_cldco2 bin. |
|---|
| 182 | |
|---|
| 183 | vrat_cld = dlog(rmax_cld/rmin_cld) / float(nbinco2_cld-1) *3. |
|---|
| 184 | vrat_cld = dexp(vrat_cld) |
|---|
| 185 | rb_cldco2(1) = rbmin_cld |
|---|
| 186 | rad_cldco2(1) = rmin_cld |
|---|
| 187 | vol_cld(1) = 4./3. * dble(pi) * rmin_cld*rmin_cld*rmin_cld |
|---|
| 188 | do i=1,nbinco2_cld-1 |
|---|
| 189 | rad_cldco2(i+1) = rad_cldco2(i) * vrat_cld**(1./3.) |
|---|
| 190 | vol_cld(i+1) = vol_cld(i) * vrat_cld |
|---|
| 191 | enddo |
|---|
| 192 | do i=1,nbinco2_cld |
|---|
| 193 | rb_cldco2(i+1)= ( (2.*vrat_cld) / (vrat_cld+1.) )**(1./3.) * |
|---|
| 194 | & rad_cldco2(i) |
|---|
| 195 | dr_cld(i) = rb_cldco2(i+1) - rb_cldco2(i) |
|---|
| 196 | enddo |
|---|
| 197 | rb_cldco2(nbinco2_cld+1) = rbmax_cld |
|---|
| 198 | dr_cld(nbinco2_cld) = rb_cldco2(nbinco2_cld+1) - |
|---|
| 199 | & rb_cldco2(nbinco2_cld) |
|---|
| 200 | print*, ' ' |
|---|
| 201 | print*,'Microphysics co2: size bin information:' |
|---|
| 202 | print*,'i,rb_cldco2(i), rad_cldco2(i),dr_cld(i)' |
|---|
| 203 | print*,'-----------------------------------' |
|---|
| 204 | do i=1,nbinco2_cld |
|---|
| 205 | write(*,'(i3,3x,3(e13.6,4x))') i,rb_cldco2(i), rad_cldco2(i), |
|---|
| 206 | & dr_cld(i) |
|---|
| 207 | enddo |
|---|
| 208 | write(*,'(i3,3x,e13.6)') nbinco2_cld+1,rb_cldco2(nbinco2_cld+1) |
|---|
| 209 | print*,'-----------------------------------' |
|---|
| 210 | do i=1,nbinco2_cld+1 |
|---|
| 211 | rb_cldco2(i) = dlog(rb_cldco2(i)) !! we save that so that it is not computed |
|---|
| 212 | !! at each timestep and gridpoint |
|---|
| 213 | enddo |
|---|
| 214 | c Contact parameter of co2 ice on dst ( m=cos(theta) ) |
|---|
| 215 | c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
|---|
| 216 | c mteta = 0.952 |
|---|
| 217 | c mtetaco2 = 0.952 |
|---|
| 218 | c write(*,*) 'co2_param contact parameter:', mtetaco2 |
|---|
| 219 | |
|---|
| 220 | c Volume of a co2 molecule (m3) |
|---|
| 221 | vo1co2 = m0co2 / dble(rho_ice_co2) ! m0co2 et non mco2 |
|---|
| 222 | |
|---|
| 223 | c Variance of the ice and CCN distributions |
|---|
| 224 | sigma_iceco2 = sqrt(log(1.+nuiceco2_sed)) |
|---|
| 225 | sigma_ice = sqrt(log(1.+nuice_sed)) |
|---|
| 226 | |
|---|
| 227 | |
|---|
| 228 | write(*,*) 'Variance of ice & CCN CO2 distribs :', sigma_iceco2 |
|---|
| 229 | write(*,*) 'nuice for co2 ice sedimentation:', nuiceco2_sed |
|---|
| 230 | write(*,*) 'Volume of a co2 molecule:', vo1co2 |
|---|
| 231 | |
|---|
| 232 | if (co2useh2o) then |
|---|
| 233 | write(*,*) |
|---|
| 234 | write(*,*) "co2useh2o=.true. in callphys.def" |
|---|
| 235 | write(*,*) "This means water ice particles can " |
|---|
| 236 | write(*,*) "serve as CCN for CO2 microphysics" |
|---|
| 237 | endif |
|---|
| 238 | |
|---|
| 239 | if (meteo_flux) then |
|---|
| 240 | write(*,*) |
|---|
| 241 | write(*,*) "meteo_flux=.true. in callphys.def" |
|---|
| 242 | write(*,*) "meteoritic dust particles are available" |
|---|
| 243 | write(*,*) "for co2 ice nucleation! " |
|---|
| 244 | write(*,*) "Flux given by J. Plane (pressions,size bins)" |
|---|
| 245 | ! Initialisation of the flux: it is constant and is it saved |
|---|
| 246 | !We must interpolate the table to the GCM pressures |
|---|
| 247 | INQUIRE(FILE=TRIM(datadir)// |
|---|
| 248 | & '/Meteo_flux_Plane.dat', EXIST=file_ok) |
|---|
| 249 | IF (.not. file_ok) THEN |
|---|
| 250 | write(*,*) 'file Meteo_flux_Plane.dat should be in ' |
|---|
| 251 | & ,trim(datadir) |
|---|
| 252 | STOP |
|---|
| 253 | endif |
|---|
| 254 | !used Variables |
|---|
| 255 | open(unit=uMeteor,file=trim(datadir)// |
|---|
| 256 | & '/Meteo_flux_Plane.dat' |
|---|
| 257 | & ,FORM='formatted') |
|---|
| 258 | !13000 records (130 pressions x 100 bin sizes) |
|---|
| 259 | read(uMeteor,*) !skip 1 line |
|---|
| 260 | do i=1,130 |
|---|
| 261 | read(uMeteor,*,iostat=read_ok) pression_meteor(i) |
|---|
| 262 | if (read_ok==0) then |
|---|
| 263 | write(*,*) pression_meteor(i) |
|---|
| 264 | else |
|---|
| 265 | write(*,*) 'Error reading Meteo_flux_Plane.dat' |
|---|
| 266 | call abort_physic("CO2clouds", |
|---|
| 267 | & "Error reading Meteo_flux_Plane.dat",1) |
|---|
| 268 | endif |
|---|
| 269 | enddo |
|---|
| 270 | read(uMeteor,*) !skip 1 line |
|---|
| 271 | do i=1,130 |
|---|
| 272 | do j=1,100 ! les mêmes 100 bins size que la distri nuclea: on touche pas |
|---|
| 273 | read(uMeteor,'(F12.6)',iostat=read_ok) meteor(i,j) |
|---|
| 274 | if (read_ok/=0) then |
|---|
| 275 | write(*,*) 'Error reading Meteo_flux_Plane.dat' |
|---|
| 276 | call abort_physic("CO2clouds", |
|---|
| 277 | & "Error reading Meteo_flux_Plane.dat",1) |
|---|
| 278 | endif |
|---|
| 279 | enddo |
|---|
| 280 | !On doit maintenant réinterpoler le tableau(130,100) sur les pressions du GCM (nlay,100) |
|---|
| 281 | enddo |
|---|
| 282 | close(uMeteor) |
|---|
| 283 | write(*,*) "File meteo_flux read, end of firstcall in co2 micro" |
|---|
| 284 | endif !of if meteo_flux |
|---|
| 285 | |
|---|
| 286 | ! Parameter values for Latent heat computation |
|---|
| 287 | l0=595594d0 |
|---|
| 288 | l1=903.111d0 |
|---|
| 289 | l2=-11.5959d0 |
|---|
| 290 | l3=0.0528288d0 |
|---|
| 291 | l4=-0.000103183d0 |
|---|
| 292 | firstcall=.false. |
|---|
| 293 | END IF |
|---|
| 294 | !============================================================= |
|---|
| 295 | ! 1. Initialisation |
|---|
| 296 | !============================================================= |
|---|
| 297 | |
|---|
| 298 | meteor_ccn(:,:,:)=0. |
|---|
| 299 | rice(:,:) = 1.e-8 |
|---|
| 300 | riceco2(:,:) = 1.e-11 |
|---|
| 301 | |
|---|
| 302 | c Initialize the tendencies |
|---|
| 303 | subpdqcloudco2(1:ngrid,1:nlay,1:nq)=0. |
|---|
| 304 | subpdtcloudco2(1:ngrid,1:nlay)=0. |
|---|
| 305 | |
|---|
| 306 | c pteff temperature layer; sum_subpdt dT.s-1 |
|---|
| 307 | c pqeff traceur kg/kg; sum_subpdq tendance idem .s-1 |
|---|
| 308 | zt(1:ngrid,1:nlay) = |
|---|
| 309 | & pteff(1:ngrid,1:nlay) + |
|---|
| 310 | & sum_subpdt(1:ngrid,1:nlay) * microtimestep |
|---|
| 311 | zq(1:ngrid,1:nlay,1:nq) = |
|---|
| 312 | & pqeff(1:ngrid,1:nlay,1:nq) + |
|---|
| 313 | & sum_subpdq(1:ngrid,1:nlay,1:nq) * microtimestep |
|---|
| 314 | WHERE( zq(1:ngrid,1:nlay,1:nq) < 1.e-30 ) |
|---|
| 315 | & zq(1:ngrid,1:nlay,1:nq) = 1.e-30 |
|---|
| 316 | |
|---|
| 317 | zq0(1:ngrid,1:nlay,1:nq) = zq(1:ngrid,1:nlay,1:nq) |
|---|
| 318 | !============================================================= |
|---|
| 319 | ! 2. Compute saturation |
|---|
| 320 | !============================================================= |
|---|
| 321 | dev2 = 1. / ( sqrt(2.) * sigma_iceco2 ) |
|---|
| 322 | dev3 = 1. / ( sqrt(2.) * sigma_ice ) |
|---|
| 323 | call co2sat(ngrid*nlay,zt,pplay,zqsat) !zqsat is psat(co2) |
|---|
| 324 | zqco2=zq(:,:,igcm_co2)+zq(:,:,igcm_co2_ice) |
|---|
| 325 | CALL tcondco2(ngrid,nlay,pplay,zqco2,tcond) |
|---|
| 326 | !============================================================= |
|---|
| 327 | ! 3. Bonus: additional meteoritic particles for nucleation |
|---|
| 328 | !============================================================= |
|---|
| 329 | if (meteo_flux) then |
|---|
| 330 | !pression_meteo(130) |
|---|
| 331 | !pplev(ngrid,nlay+1) |
|---|
| 332 | !meteo(130,100) |
|---|
| 333 | !resultat: meteo_ccn(ngrid,nlay,100) |
|---|
| 334 | do l=1,nlay |
|---|
| 335 | do ig=1,ngrid |
|---|
| 336 | masse(ig,l)=(pplev(ig,l) - pplev(ig,l+1)) /g |
|---|
| 337 | ltemp1=abs(pression_meteor(:)-pplev(ig,l)) |
|---|
| 338 | ltemp2=abs(pression_meteor(:)-pplev(ig,l+1)) |
|---|
| 339 | lebon1=minloc(ltemp1,DIM=1) |
|---|
| 340 | lebon2=minloc(ltemp2,DIM=1) |
|---|
| 341 | nelem=lebon2-lebon1+1. |
|---|
| 342 | mtemp(:)=0d0 !mtemp(100) : valeurs pour les 100bins |
|---|
| 343 | do ibin=1,100 |
|---|
| 344 | mtemp(ibin)=sum(meteor(lebon1:lebon2,ibin)) |
|---|
| 345 | enddo |
|---|
| 346 | meteor_ccn(ig,l,:)=mtemp(:)/nelem/masse(ig,l) !Par kg air |
|---|
| 347 | csi par m carre, x epaisseur/masse pour par kg/air |
|---|
| 348 | !write(*,*) "masse air ig l=",masse(ig,l) |
|---|
| 349 | !check original unit with J. Plane |
|---|
| 350 | enddo |
|---|
| 351 | enddo |
|---|
| 352 | endif |
|---|
| 353 | c ------------------------------------------------------------------------ |
|---|
| 354 | c --------- Actual microphysics : Main loop over the GCM's grid --------- |
|---|
| 355 | c ------------------------------------------------------------------------ |
|---|
| 356 | DO l=1,nlay |
|---|
| 357 | DO ig=1,ngrid |
|---|
| 358 | c Get the partial pressure of co2 vapor and its saturation ratio |
|---|
| 359 | pco2 = zq(ig,l,igcm_co2) * (mmean(ig,l)/44.01) * pplay(ig,l) |
|---|
| 360 | satu = pco2 / zqsat(ig,l) |
|---|
| 361 | |
|---|
| 362 | rho_ice_co2T(ig,l)=1000.*(1.72391-2.53e-4*zt(ig,l) |
|---|
| 363 | & -2.87e-6*zt(ig,l)*zt(ig,l)) !T-dependant CO2 ice density |
|---|
| 364 | vo2co2 = m0co2 / dble(rho_ice_co2T(ig,l)) |
|---|
| 365 | rho_ice_co2=rho_ice_co2T(ig,l) |
|---|
| 366 | |
|---|
| 367 | !============================================================= |
|---|
| 368 | !4. Nucleation |
|---|
| 369 | !============================================================= |
|---|
| 370 | IF ( satu .ge. 1 ) THEN ! if there is condensation |
|---|
| 371 | |
|---|
| 372 | c We do rdust computation "by hand" because we don't want to |
|---|
| 373 | c change the mininumum rdust value in updaterad... It would |
|---|
| 374 | c mess up various part of the GCM ! |
|---|
| 375 | |
|---|
| 376 | rdust(ig,l)= zq(ig,l,igcm_dust_mass) |
|---|
| 377 | & *0.75/(pi*rho_dust |
|---|
| 378 | & * zq(ig,l,igcm_dust_number)) |
|---|
| 379 | rdust(ig,l)= rdust(ig,l)**(1./3.) |
|---|
| 380 | if (zq(ig,l,igcm_dust_mass)*tauscaling(ig) .le. 1.e-20 |
|---|
| 381 | & .or. zq(ig,l,igcm_dust_number)*tauscaling(ig) .le.1 |
|---|
| 382 | & .or. rdust(ig,l) .le. 1.e-9) then |
|---|
| 383 | rdust(ig,l)=1.e-9 |
|---|
| 384 | endif |
|---|
| 385 | rdust(ig,l)=min(5.e-4,rdust(ig,l)) |
|---|
| 386 | |
|---|
| 387 | c Expand the dust moments into a binned distribution |
|---|
| 388 | |
|---|
| 389 | n_aer(:)=0. |
|---|
| 390 | m_aer(:)=0. |
|---|
| 391 | |
|---|
| 392 | Mo = zq(ig,l,igcm_dust_mass)* tauscaling(ig)+1.e-30 |
|---|
| 393 | No = zq(ig,l,igcm_dust_number)* tauscaling(ig)+1.e-30 |
|---|
| 394 | |
|---|
| 395 | No_dust=No |
|---|
| 396 | Mo_dust=Mo |
|---|
| 397 | |
|---|
| 398 | Rn = rdust(ig,l) |
|---|
| 399 | Rn = -dlog(Rn) |
|---|
| 400 | Rm = Rn - 3. * sigma_iceco2*sigma_iceco2 |
|---|
| 401 | n_derf = derf( (rb_cldco2(1)+Rn) *dev2) |
|---|
| 402 | m_derf = derf( (rb_cldco2(1)+Rm) *dev2) |
|---|
| 403 | |
|---|
| 404 | do i = 1, nbinco2_cld |
|---|
| 405 | n_aer(i) = -0.5 * No * n_derf !! this ith previously computed |
|---|
| 406 | m_aer(i) = -0.5 * Mo * m_derf !! this ith previously computed |
|---|
| 407 | n_derf = derf((rb_cldco2(i+1)+Rn) *dev2) |
|---|
| 408 | m_derf = derf((rb_cldco2(i+1)+Rm) *dev2) |
|---|
| 409 | n_aer(i) = n_aer(i) + 0.5 * No * n_derf |
|---|
| 410 | m_aer(i) = m_aer(i) + 0.5 * Mo * m_derf |
|---|
| 411 | enddo |
|---|
| 412 | |
|---|
| 413 | c Ajout meteor_ccn particles aux particules de poussière background |
|---|
| 414 | if (meteo_flux) then |
|---|
| 415 | do i = 1, nbinco2_cld |
|---|
| 416 | n_aer(i) = n_aer(i) + meteor_ccn(ig,l,i) |
|---|
| 417 | m_aer(i) = m_aer(i) + 4./3.*pi*rho_dust |
|---|
| 418 | & *meteor_ccn(ig,l,i)*rad_cldco2(i)*rad_cldco2(i) |
|---|
| 419 | & *rad_cldco2(i) |
|---|
| 420 | enddo |
|---|
| 421 | endif |
|---|
| 422 | |
|---|
| 423 | c Same but with h2o particles as CCN only if co2useh2o=.true. |
|---|
| 424 | |
|---|
| 425 | n_aer_h2oice(:)=0. |
|---|
| 426 | m_aer_h2oice(:)=0. |
|---|
| 427 | |
|---|
| 428 | if (co2useh2o) then |
|---|
| 429 | call updaterice_micro(zq(ig,l,igcm_h2o_ice), |
|---|
| 430 | & zq(ig,l,igcm_ccn_mass),zq(ig,l,igcm_ccn_number), |
|---|
| 431 | & tauscaling(ig),rice(ig,l),rhocloud(ig,l)) |
|---|
| 432 | Mo = zq(ig,l,igcm_h2o_ice) + |
|---|
| 433 | & zq(ig,l,igcm_ccn_mass)*tauscaling(ig)+1.e-30 |
|---|
| 434 | ! Total mass of H20 crystals,CCN included |
|---|
| 435 | No = zq(ig,l,igcm_ccn_number)* tauscaling(ig) + 1.e-30 |
|---|
| 436 | Rn = rice(ig,l) |
|---|
| 437 | Rn = -dlog(Rn) |
|---|
| 438 | Rm = Rn - 3. * sigma_ice*sigma_ice |
|---|
| 439 | n_derf = derf( (rb_cldco2(1)+Rn) *dev3) |
|---|
| 440 | m_derf = derf( (rb_cldco2(1)+Rm) *dev3) |
|---|
| 441 | do i = 1, nbinco2_cld |
|---|
| 442 | n_aer_h2oice(i) = -0.5 * No * n_derf |
|---|
| 443 | m_aer_h2oice(i) = -0.5 * Mo * m_derf |
|---|
| 444 | n_derf = derf( (rb_cldco2(i+1)+Rn) *dev3) |
|---|
| 445 | m_derf = derf( (rb_cldco2(i+1)+Rm) *dev3) |
|---|
| 446 | n_aer_h2oice(i) = n_aer_h2oice(i) + 0.5 * No * n_derf |
|---|
| 447 | m_aer_h2oice(i) = m_aer_h2oice(i) + 0.5 * Mo * m_derf |
|---|
| 448 | rad_h2oice(i) = rad_cldco2(i) |
|---|
| 449 | enddo |
|---|
| 450 | endif |
|---|
| 451 | |
|---|
| 452 | |
|---|
| 453 | ! Call to nucleation routine |
|---|
| 454 | call nucleaCO2(dble(pco2),zt(ig,l),dble(satu) |
|---|
| 455 | & ,n_aer,rate,n_aer_h2oice |
|---|
| 456 | & ,rad_h2oice,rateh2o,vo2co2) |
|---|
| 457 | dN = 0. |
|---|
| 458 | dM = 0. |
|---|
| 459 | dNh2o = 0. |
|---|
| 460 | dMh2o = 0. |
|---|
| 461 | do i = 1, nbinco2_cld |
|---|
| 462 | Proba =1.0-dexp(-1.*microtimestep*rate(i)) |
|---|
| 463 | dN = dN + n_aer(i) * Proba |
|---|
| 464 | dM = dM + m_aer(i) * Proba |
|---|
| 465 | enddo |
|---|
| 466 | if (co2useh2o) then |
|---|
| 467 | do i = 1, nbinco2_cld |
|---|
| 468 | Probah2o = 1.0-dexp(-1.*microtimestep*rateh2o(i)) |
|---|
| 469 | dNh2o = dNh2o + n_aer_h2oice(i) * Probah2o |
|---|
| 470 | dMh2o = dMh2o + m_aer_h2oice(i) * Probah2o |
|---|
| 471 | enddo |
|---|
| 472 | endif |
|---|
| 473 | |
|---|
| 474 | ! dM masse activée (kg) et dN nb particules par kg d'air |
|---|
| 475 | ! Now increment CCN tracers and update dust tracers |
|---|
| 476 | dNN= dN/tauscaling(ig) |
|---|
| 477 | dMM= dM/tauscaling(ig) |
|---|
| 478 | dNN=min(dNN,zq(ig,l,igcm_dust_number)) |
|---|
| 479 | dMM=min(dMM,zq(ig,l,igcm_dust_mass)) |
|---|
| 480 | zq(ig,l,igcm_ccnco2_mass) = |
|---|
| 481 | & zq(ig,l,igcm_ccnco2_mass) + dMM |
|---|
| 482 | zq(ig,l,igcm_ccnco2_number) = |
|---|
| 483 | & zq(ig,l,igcm_ccnco2_number) + dNN |
|---|
| 484 | zq(ig,l,igcm_dust_mass)= zq(ig,l,igcm_dust_mass)-dMM |
|---|
| 485 | zq(ig,l,igcm_dust_number)=zq(ig,l,igcm_dust_number)-dNN |
|---|
| 486 | |
|---|
| 487 | |
|---|
| 488 | c Update CCN for CO2 nucleating on H2O CCN : |
|---|
| 489 | ! Warning: must keep memory of it |
|---|
| 490 | if (co2useh2o) then |
|---|
| 491 | dNNh2o=dNh2o/tauscaling(ig) |
|---|
| 492 | dNNh2o=min(dNNh2o,zq(ig,l,igcm_ccn_number)) |
|---|
| 493 | ratioh2o_ccn=1./(zq(ig,l,igcm_h2o_ice) |
|---|
| 494 | & +zq(ig,l,igcm_ccn_mass)*tauscaling(ig)) |
|---|
| 495 | dMh2o_ice=dMh2o*zq(ig,l,igcm_h2o_ice)*ratioh2o_ccn |
|---|
| 496 | dMh2o_ccn=dMh2o*zq(ig,l,igcm_ccn_mass)* |
|---|
| 497 | & tauscaling(ig)*ratioh2o_ccn |
|---|
| 498 | dMh2o_ccn=dMh2o_ccn/tauscaling(ig) |
|---|
| 499 | dMh2o_ccn=min(dMh2o_ccn,zq(ig,l,igcm_ccn_mass)) |
|---|
| 500 | dMh2o_ice=min(dMh2o_ice,zq(ig,l,igcm_h2o_ice)) |
|---|
| 501 | zq(ig,l,igcm_ccnco2_mass) = |
|---|
| 502 | & zq(ig,l,igcm_ccnco2_mass) + dMh2o_ice+dMh2o_ccn |
|---|
| 503 | zq(ig,l,igcm_ccnco2_number) = |
|---|
| 504 | & zq(ig,l,igcm_ccnco2_number) + dNNh2o |
|---|
| 505 | zq(ig,l,igcm_ccn_number)=zq(ig,l,igcm_ccn_number)-dNNh2o |
|---|
| 506 | zq(ig,l,igcm_h2o_ice) = zq(ig,l,igcm_h2o_ice)-dMh2o_ice |
|---|
| 507 | zq(ig,l,igcm_ccn_mass)= zq(ig,l,igcm_ccn_mass)-dMh2o_ccn |
|---|
| 508 | memdMMh2o(ig,l)=memdMMh2o(ig,l)+dMh2o_ice |
|---|
| 509 | memdMMccn(ig,l)=memdMMccn(ig,l)+dMh2o_ccn |
|---|
| 510 | memdNNccn(ig,l)=memdNNccn(ig,l)+dNNh2o |
|---|
| 511 | endif ! of if co2useh2o |
|---|
| 512 | ENDIF ! of is satu >1 |
|---|
| 513 | |
|---|
| 514 | !============================================================= |
|---|
| 515 | ! 5. Ice growth: scheme for radius evolution |
|---|
| 516 | !============================================================= |
|---|
| 517 | |
|---|
| 518 | c We trigger crystal growth if and only if there is at least one nuclei (N>1). |
|---|
| 519 | c Indeed, if we are supersaturated and still don't have at least one nuclei, we should better wait |
|---|
| 520 | c to avoid unrealistic value for nuclei radius and so on for cases that remain negligible. |
|---|
| 521 | No = zq(ig,l,igcm_ccnco2_number)* tauscaling(ig)+1.e-30 |
|---|
| 522 | IF (No .ge. 1)THEN ! we trigger crystal growth |
|---|
| 523 | call updaterice_microco2(zq(ig,l,igcm_co2_ice), |
|---|
| 524 | & zq(ig,l,igcm_ccnco2_mass),zq(ig,l,igcm_ccnco2_number), |
|---|
| 525 | & tauscaling(ig),riceco2(ig,l),rhocloudco2(ig,l)) |
|---|
| 526 | |
|---|
| 527 | Ic_rice=0. |
|---|
| 528 | lw = l0 + l1 * zt(ig,l) + l2 * zt(ig,l)**2 + |
|---|
| 529 | & l3 * zt(ig,l)**3 + l4 * zt(ig,l)**4 !J.kg-1 |
|---|
| 530 | facteurmax=abs(Tcond(ig,l)-zt(ig,l))*cpp/lw |
|---|
| 531 | !specific heat of co2 ice = 1000 J.kg-1.K-1 |
|---|
| 532 | !specific heat of atm cpp = 744.5 J.kg-1.K-1 |
|---|
| 533 | |
|---|
| 534 | c call scheme of microphys. mass growth for CO2 |
|---|
| 535 | call massflowrateCO2(pplay(ig,l),zt(ig,l), |
|---|
| 536 | & satu,riceco2(ig,l),mmean(ig,l),Ic_rice) |
|---|
| 537 | c Ic_rice Mass transfer rate (kg/s) for a rice particle >0 si croissance ! |
|---|
| 538 | |
|---|
| 539 | if (isnan(Ic_rice) .or. Ic_rice .eq. 0.) then |
|---|
| 540 | Ic_rice=0. |
|---|
| 541 | subpdtcloudco2(ig,l)=-sum_subpdt(ig,l) |
|---|
| 542 | dMice=0 |
|---|
| 543 | |
|---|
| 544 | else |
|---|
| 545 | dMice=zq(ig,l,igcm_ccnco2_number)*Ic_rice*microtimestep |
|---|
| 546 | & *tauscaling(ig) ! Kg par kg d'air, >0 si croissance ! |
|---|
| 547 | !kg.s-1 par particule * nb particule par kg air*s |
|---|
| 548 | ! = kg par kg air |
|---|
| 549 | |
|---|
| 550 | dMice = max(dMice,max(-facteurmax,-zq(ig,l,igcm_co2_ice))) |
|---|
| 551 | dMice = min(dMice,min(facteurmax,zq(ig,l,igcm_co2))) |
|---|
| 552 | ! facteurmax maximum quantity of CO2 that can sublime/condense according to available thermal energy |
|---|
| 553 | ! latent heat release >0 if growth i.e. if dMice >0 |
|---|
| 554 | subpdtcloudco2(ig,l)=dMice*lw/cpp/microtimestep |
|---|
| 555 | ! kgco2/kgair* J/kgco2 * 1/(J.kgair-1.K-1)/s= K par seconde |
|---|
| 556 | !Now update tracers |
|---|
| 557 | zq(ig,l,igcm_co2_ice) = zq(ig,l,igcm_co2_ice)+dMice |
|---|
| 558 | zq(ig,l,igcm_co2) = zq(ig,l,igcm_co2)-dMice |
|---|
| 559 | endif |
|---|
| 560 | |
|---|
| 561 | !============================================================= |
|---|
| 562 | ! 6. Dust cores releasing if no more co2 ice : |
|---|
| 563 | !============================================================= |
|---|
| 564 | |
|---|
| 565 | if (zq(ig,l,igcm_co2_ice).le. 1.e-25)THEN |
|---|
| 566 | ! On sublime tout |
|---|
| 567 | if (co2useh2o) then |
|---|
| 568 | if (memdMMccn(ig,l) .gt. 0) then |
|---|
| 569 | zq(ig,l,igcm_ccn_mass)=zq(ig,l,igcm_ccn_mass) |
|---|
| 570 | & +memdMMccn(ig,l) |
|---|
| 571 | endif |
|---|
| 572 | if (memdMMh2o(ig,l) .gt. 0) then |
|---|
| 573 | zq(ig,l,igcm_h2o_ice)=zq(ig,l,igcm_h2o_ice) |
|---|
| 574 | & +memdMMh2o(ig,l) |
|---|
| 575 | endif |
|---|
| 576 | |
|---|
| 577 | if (memdNNccn(ig,l) .gt. 0) then |
|---|
| 578 | zq(ig,l,igcm_ccn_number)=zq(ig,l,igcm_ccn_number) |
|---|
| 579 | & +memdNNccn(ig,l) |
|---|
| 580 | endif |
|---|
| 581 | endif |
|---|
| 582 | zq(ig,l,igcm_dust_mass) = |
|---|
| 583 | & zq(ig,l,igcm_dust_mass) |
|---|
| 584 | & + zq(ig,l,igcm_ccnco2_mass)- |
|---|
| 585 | & (memdMMh2o(ig,l)+memdMMccn(ig,l)) |
|---|
| 586 | zq(ig,l,igcm_dust_number) = |
|---|
| 587 | & zq(ig,l,igcm_dust_number) |
|---|
| 588 | & + zq(ig,l,igcm_ccnco2_number)-memdNNccn(ig,l) |
|---|
| 589 | |
|---|
| 590 | zq(ig,l,igcm_co2) = zq(ig,l,igcm_co2) |
|---|
| 591 | & + zq(ig,l,igcm_co2_ice) |
|---|
| 592 | |
|---|
| 593 | zq(ig,l,igcm_ccnco2_mass)=0. |
|---|
| 594 | zq(ig,l,igcm_co2_ice)=0. |
|---|
| 595 | zq(ig,l,igcm_ccnco2_number)=0. |
|---|
| 596 | memdNNccn(ig,l)=0. |
|---|
| 597 | memdMMh2o(ig,l)=0. |
|---|
| 598 | memdMMccn(ig,l)=0. |
|---|
| 599 | riceco2(ig,l)=0. |
|---|
| 600 | |
|---|
| 601 | endif !of if co2_ice <1e-25 |
|---|
| 602 | |
|---|
| 603 | ENDIF ! of if NCCN > 1 |
|---|
| 604 | ENDDO ! of ig loop |
|---|
| 605 | ENDDO ! of nlayer loop |
|---|
| 606 | |
|---|
| 607 | !============================================================= |
|---|
| 608 | ! 7. END: get cloud tendencies |
|---|
| 609 | !============================================================= |
|---|
| 610 | |
|---|
| 611 | ! Get cloud tendencies |
|---|
| 612 | subpdqcloudco2(1:ngrid,1:nlay,igcm_co2) = |
|---|
| 613 | & (zq(1:ngrid,1:nlay,igcm_co2) - |
|---|
| 614 | & zq0(1:ngrid,1:nlay,igcm_co2))/microtimestep |
|---|
| 615 | subpdqcloudco2(1:ngrid,1:nlay,igcm_co2_ice) = |
|---|
| 616 | & (zq(1:ngrid,1:nlay,igcm_co2_ice) - |
|---|
| 617 | & zq0(1:ngrid,1:nlay,igcm_co2_ice))/microtimestep |
|---|
| 618 | |
|---|
| 619 | if (co2useh2o) then |
|---|
| 620 | subpdqcloudco2(1:ngrid,1:nlay,igcm_h2o_ice) = |
|---|
| 621 | & (zq(1:ngrid,1:nlay,igcm_h2o_ice) - |
|---|
| 622 | & zq0(1:ngrid,1:nlay,igcm_h2o_ice))/microtimestep |
|---|
| 623 | subpdqcloudco2(1:ngrid,1:nlay,igcm_ccn_mass) = |
|---|
| 624 | & (zq(1:ngrid,1:nlay,igcm_ccn_mass) - |
|---|
| 625 | & zq0(1:ngrid,1:nlay,igcm_ccn_mass))/microtimestep |
|---|
| 626 | subpdqcloudco2(1:ngrid,1:nlay,igcm_ccn_number) = |
|---|
| 627 | & (zq(1:ngrid,1:nlay,igcm_ccn_number) - |
|---|
| 628 | & zq0(1:ngrid,1:nlay,igcm_ccn_number))/microtimestep |
|---|
| 629 | endif |
|---|
| 630 | |
|---|
| 631 | subpdqcloudco2(1:ngrid,1:nlay,igcm_ccnco2_mass) = |
|---|
| 632 | & (zq(1:ngrid,1:nlay,igcm_ccnco2_mass) - |
|---|
| 633 | & zq0(1:ngrid,1:nlay,igcm_ccnco2_mass))/microtimestep |
|---|
| 634 | subpdqcloudco2(1:ngrid,1:nlay,igcm_ccnco2_number) = |
|---|
| 635 | & (zq(1:ngrid,1:nlay,igcm_ccnco2_number) - |
|---|
| 636 | & zq0(1:ngrid,1:nlay,igcm_ccnco2_number))/microtimestep |
|---|
| 637 | subpdqcloudco2(1:ngrid,1:nlay,igcm_dust_mass) = |
|---|
| 638 | & (zq(1:ngrid,1:nlay,igcm_dust_mass) - |
|---|
| 639 | & zq0(1:ngrid,1:nlay,igcm_dust_mass))/microtimestep |
|---|
| 640 | subpdqcloudco2(1:ngrid,1:nlay,igcm_dust_number) = |
|---|
| 641 | & (zq(1:ngrid,1:nlay,igcm_dust_number) - |
|---|
| 642 | & zq0(1:ngrid,1:nlay,igcm_dust_number))/microtimestep |
|---|
| 643 | |
|---|
| 644 | |
|---|
| 645 | end |
|---|
| 646 | |
|---|
| 647 | |
|---|
| 648 | |
|---|