source: trunk/LMDZ.MARS/libf/phymars/co2cloud.F90 @ 2616

Last change on this file since 2616 was 2616, checked in by romain.vande, 3 years ago

LMDZ_MARS RV : Open_MP;
Put all the "save" variables as "!$OMP THREADPRIVATE" in phymars.
The code can now be tested, see README for more info

  • Property svn:executable set to *
File size: 60.7 KB
Line 
1!======================================================================================================================!
2! Module: CO2 clouds formation ========================================================================================!
3!----------------------------------------------------------------------------------------------------------------------!
4! Authors: Joachim Audouard, Constantino Listowski, Anni Määttänen
5! Date: 09/2016
6!----------------------------------------------------------------------------------------------------------------------!
7! Contains subroutines:
8!      - co2cloud: of co2 cloud microphysics
9!======================================================================================================================!
10module co2cloud_mod
11
12implicit none
13
14contains
15!======================================================================================================================!
16! SUBROUTINE: co2cloud ================================================================================================!
17!======================================================================================================================!
18! Subject:
19!---------
20!   Main subroutine of co2 cloud microphysics
21!----------------------------------------------------------------------------------------------------------------------!
22! Comments:
23!----------
24!     - Adaptation of the water ice clouds scheme (with specific microphysics) of Montmessin, Navarro et al.
25!
26!     - Microphysics subroutine is improvedCO2clouds.F
27!
28!     - There is a time loop specific to cloud formation due to timescales smaller than the GCM integration timestep
29!
30!     - The microphysics time step is a fraction of the physical one
31!
32!     - The CO2 clouds tracers (co2_ice, ccn mass and concentration) are sedimented at each microtimestep. pdqs_sedco2
33!       keeps track of the CO2 flux at the surface
34!
35!     - The subgrid Temperature distribution is modulated (0 or 1) by Spiga et al. (GRL 2012)
36!
37!     - Saturation Index to account for GW propagation or dissipation upwards
38!
39!     - 4D and column opacities are computed using Qext values at 1µm
40!----------------------------------------------------------------------------------------------------------------------!
41! Papers:
42!--------
43!   "Near-pure vapor condensation in the Martian atmosphere: CO2 ice crystal growth", Listowski et al. (2013), JGRE
44!   "Modeling the microphysics of CO2 ice clouds within wave-induced cold pockets in the martian mesosphere", Listowski
45!     et al. (2014), Icarus
46!   "Global climate modeling of the Martian water cycle with improved microphysics and radiatively active water ice
47!     clouds", Navarro et al. (2014), JGRE
48!   "Martian GCM with complete CO2 clouds microphysics", Audouard et al. (2017), EPSC abstract
49!----------------------------------------------------------------------------------------------------------------------!
50! Algorithm:
51!-----------
52!   0. Firstcall
53!     0.1. Initialization of microtimestep from imicroco2
54!     0.2. Compute the radius grid of CO2 ice particles (rb_cldco2)
55!     0.3. Read file 'optprop_co2ice_1mic.dat' to extract optical properties of CO2 ice at 1 micron (Qext)
56!     0.4. Interpole the radius grid (rb_cldco2) to get the corresponding exctinction coefficients (Qext1bins)
57!     0.5. Save the radius grid of CO2 particles (rb_cldco2)
58!   1. Initialization
59!   2. Compute mass and thickness layers
60!   3. Define the sub-grid cloud (CLFvaryingCO2)
61!     3.1. Representation of sub-grid CO2 ice clouds (CLFvaryingCO2 = True)
62!       3.1.a. Saturation index CO2
63!       3.1.b. Compute tcond
64!       3.1.c. Compute cloud fraction in cells
65!     3.2. No sub-grid cloud representation (CLFvaryingCO2 = False)
66!   4. Microphysics of CO2 cloud formation
67!     4.1. Stepped entry for tendancies: At each micro timestep we add pdt in order to have a stepped entry
68!     4.2. Effective tracers quantities in the cloud
69!     4.3. Gravitational sedimentation
70!       4.3.a. Compute cloud density
71!       4.3.b. Save actualized tracer values to compute sedimentation tendancies
72!       4.3.c. Sedimentation of co2 ice
73!       4.3.d. Sedimentation for other tracers
74!       4.3.e. Compute tendencies due to the sedimation process
75!     4.4. Main call to the cloud scheme
76!     4.5. Updating tendencies after cloud scheme
77!   5. Compute final tendencies after time loop
78!   6. Update clouds physical values in the cloud (for output)
79!     6.1. Update density of co2 ice, riceco2 and opacity
80!     6.2. Update rice and rdust
81!   7. Correction if a lot of subliming CO2 fills the 1st layer
82!   8. Compute water cloud sedimentation radius
83!   9. CO2 saturated quantities
84!     9.1 Compute CO2 saturation in layers
85!     9.2 Compute CO2 saturated quantities in layers
86!   10. Everything modified by CO2 microphysics must be wrt co2cloudfrac
87!   11. Compute opacity at 1 micron
88!   12. Write outputs in diagfi.nc
89!======================================================================================================================!
90  subroutine co2cloud(ngrid, nlay, ptimestep, pplev, pplay, pdpsrf, pzlay, pt, pdt, pq, pdq, pdqcloudco2, pdtcloudco2, &
91                      nq, tau, tauscaling, rdust, rice, riceco2, nuice, rhocloud, rsedcloudco2, rhocloudco2, pzlev,&
92                      pdqs_sedco2, pdqs_sedccn, pdu, pu, pcondicea, co2ice)
93
94  use ioipsl_getincom, only: getin
95  use dimradmars_mod,  only: naerkind
96  use comcstfi_h,      only: pi, g, cpp
97  use updaterad,       only: updaterice_microco2, updaterice_micro, updaterdust
98  use conc_mod,        only: mmean, rnew
99  use tracer_mod,      only: igcm_co2, igcm_co2_ice, igcm_dust_mass, igcm_dust_number, igcm_h2o_ice, &
100                             igcm_ccn_mass, igcm_ccn_number, igcm_ccnco2_mass, igcm_ccnco2_number, &
101                             igcm_ccnco2_h2o_number, igcm_ccnco2_h2o_mass_ice, igcm_ccnco2_h2o_mass_ccn, rho_dust, &
102                             nuiceco2_sed, nuiceco2_ref, r3n_q, rho_ice, nuice_sed, igcm_ccnco2_meteor_mass, &
103                             igcm_ccnco2_meteor_number
104
105  use newsedim_mod,    only: newsedim
106
107  use datafile_mod,    only: datadir
108
109  use improvedCO2clouds_mod, only: improvedCO2clouds
110
111#ifndef MESOSCALE
112  use vertical_layers_mod, only: ap, bp
113#endif
114
115  implicit none
116
117  include "callkeys.h"
118  include "microphys.h"
119!----------------------------------------------------------------------------------------------------------------------!
120! VARIABLES DECLARATION
121!----------------------------------------------------------------------------------------------------------------------!
122! Input arguments:
123!-----------------
124  integer, intent(in) ::&
125     ngrid, &! Number of grid points
126     nlay,  &! Number of layers
127     nq      ! Number of tracers
128
129  real, intent(in) :: &
130     ptimestep,           &! Physical time step (s)
131     pplev(ngrid,nlay+1), &! Inter-layer pressures (Pa)
132     pplay(ngrid,nlay),   &! Mid-layer pressures (Pa)
133     pdpsrf(ngrid),       &! Tendency on surface pressure
134     pzlay(ngrid,nlay),   &! Altitude at the middle of the layers
135     pt(ngrid,nlay),      &! Temperature at the middle of the layers (K)
136     pdt(ngrid,nlay),     &! Tendency on temperature from other parametrizations
137     pq(ngrid,nlay,nq),   &! Tracers (kg/kg)
138     pdq(ngrid,nlay,nq),  &! Tendencies before condensation (kg/kg.s-1)
139     tau(ngrid,naerkind), &! Column dust optical depth at each point
140     tauscaling(ngrid),   &! Convertion factor for dust amount
141     pu(ngrid,nlay),      &! Zonal Wind: zu = pu + (pdu * ptimestep)
142     pdu(ngrid,nlay),     &! Tendency of zonal wind before condensation
143     pzlev(ngrid,nlay+1), &! Altitude at the boundaries of the layers
144     nuice(ngrid,nlay),   &! Estimated effective variance of the size distribution
145     co2ice(ngrid)         ! Amount of co2 ice at the surface
146!----------------------------------------------------------------------------------------------------------------------!
147! Output arguments:
148!------------------
149  real, intent(out) :: &
150     rice(ngrid,nlay),          & ! Water Ice mass mean radius (m)
151!     rsedcloud(ngrid,nlay),     & ! Water Cloud sedimentation radius
152     rhocloud(ngrid,nlay),      & ! Water Cloud density (kg.m-3)
153     pdqs_sedco2(ngrid),        & ! CO2 flux at the surface
154     pdqs_sedccn(ngrid,nq),      & ! CCN flux at the surface
155     pdqcloudco2(ngrid,nlay,nq),& ! Tendency due to CO2 condensation (kg/kg.s-1)
156     pcondicea(ngrid,nlay),     & ! Rate of condensation/sublimation of co2 ice in layers
157     pdtcloudco2(ngrid,nlay),   & ! Tendency on temperature due to latent heat
158     rsedcloudco2(ngrid,nlay)     ! Cloud sedimentation radius
159
160  real, intent(inout) :: &
161     rdust(ngrid,nlay) ! Dust geometric mean radius (m)
162
163  double precision, intent(out) :: &
164     riceco2(ngrid,nlay)          ! Ice mass mean radius (m) r_c in Montmessin et al. (2004)
165!----------------------------------------------------------------------------------------------------------------------!
166! Local:
167!-------
168!-----1) Parameters:
169!-------------------
170  integer, parameter :: &
171     uQext = 555,         &! file_qext unit ID
172     var_dim_qext = 10000  ! Exact dimension of radv and qextv1mic from file_qext
173
174  real, parameter :: &
175     mincloud = 0.1,  &! Minimum cloud fraction
176     beta = 0.85,     &! correction for the shape of the particles (see Murphy et al. JGR 1990 vol.95):
177                       !   beta = 1    for spheres
178                       !   beta = 0.85 for irregular particles
179                       !   beta = 0.5  for disk shaped particles
180     threshold = 1e-30 ! limit value
181
182  double precision, parameter :: &
183     rmin_cld = 1.e-9,  &! Minimum radius (m)
184     rmax_cld = 5.e-6,  &! Maximum radius (m)
185     rbmin_cld = 1.e-10,&! Minimum boundary radius (m)
186     rbmax_cld = 2.e-4, &! Maximum boundary radius (m)
187     Fo = 7.5e-7,       &! for sat index  (J.m-3)
188     lambdaH = 150.e3    ! for sat index (km)
189
190  character(len=23), parameter :: &
191     file_qext = 'optprop_co2ice_1mic.dat' ! File extinction coefficients of CO2 particles
192!----------------------------------------------------------------------------------------------------------------------!
193!-----2) Saved:
194!--------------
195  integer, save :: &
196     imicroco2 ! Time subsampling for coupled water microphysics sedimentation microtimestep timeloop for microphysics:
197               !   if imicroco2 = 1, subpdt is the same as pdt
198  real, save :: &
199     sigma_iceco2, &! Variance of the ice and CCN distributions
200     microtimestep  ! Integration timestep for coupled water microphysics & sedimentation
201
202  double precision, save :: &
203     dev2,                   &! 1. / ( sqrt(2.) * sigma_iceco2 )
204     Qext1bins(nbinco2_cld), &! Extinction coefficients for rb_cldco2 radius of CO2 ice particles
205     Qextv1mic(var_dim_qext), &
206     radv(var_dim_qext),      &    ! radius of CO2 ice at 1 µm (read from file_qext)
207     rb_cldco2(nbinco2_cld+1) ! boundary values of each rad_cldco2 bin (m)
208  logical, save :: &
209     firstcall = .true. ! Used to compute saved variables
210     
211!$OMP THREADPRIVATE(imicroco2,sigma_iceco2,microtimestep)
212!$OMP THREADPRIVATE(dev2,Qext1bins,Qextv1mic,radv,rb_cldco2)
213!$OMP THREADPRIVATE(firstcall)     
214
215!----------------------------------------------------------------------------------------------------------------------!
216!-----3) Variables:
217!------------------
218  integer :: &
219     iq,       &! loop on tracers
220     ig,       &! loop on grid points
221     l,        &! loop on layers
222     i,        &! loop on nbinco2_cld
223     nelem,    &! number of point between lebon1 and lebon2 => interpolation
224     lebon1,   &! bound limit for the interpolation
225     lebon2,   &! bound limit for the interpolation
226     microstep  ! Time subsampling step variable
227
228  real :: &
229! ---Tendency given by clouds inside the micro loop
230     subpdqcloudco2(ngrid,nlay,nq), &! On tracers, cf. pdqcloud
231     subpdtcloudco2(ngrid,nlay),    &! On temperature, cf. pdtcloud
232!  ---Global tendency (clouds+physics)
233     sum_subpdq(ngrid,nlay,nq),     &! On tracers, cf. pdqcloud
234     sum_subpdt(ngrid,nlay),        &! On temperature, cf. pdtcloud
235!  ---Sedimentation
236     ztsed(ngrid,nlay),             &! Temperature with real-time value in microtimeloop
237     zqsed(ngrid,nlay,nq),          &! Tracers with real-time value in µloop
238     zqsed0(ngrid,nlay,nq),         &! For sedimentation tendancy
239     subpdqsed(ngrid,nlay,nq),      &! Tendancy due to sedimentation
240     sum_subpdqs_sedco2(ngrid),     &! CO2 flux at the surface
241     sum_subpdqs_sedccn(ngrid,nq),  &!
242      !  ---For sub grid T distribution
243     zt(ngrid,nlay),                &! Local value of temperature
244     zq_co2vap(ngrid, nlay),        &! Local value of CO2 vap
245     rhocloudco2t(ngrid, nlay),     &! Cloud density (kg.m-3)
246!  ---For Saturation Index computation
247     zdelt,                         &! Delta T for the temperature distribution
248     co2cloudfrac(ngrid,nlay),      &! Cloud fraction used only with CLFvarying is true
249!  ---Misc
250     rhocloudco2(ngrid, nlay),      &! Cloud density (kg.m-3)
251     Nccnco2,                       &! buffer: number of ccn used for co2 condensation
252     Qccnco2,                       &! buffer: mass of ccn used for co2 condensation
253     Nccnco2_h2o, &
254     Qccnco2_h2o, &
255     Niceco2,                       &! buffer: mmr co2 ice
256     epaisseur(ngrid,nlay),         &! Layer thickness (m)
257     masse(ngrid,nlay),             &! Layer  mass (kg.m-2)
258     pteff(ngrid, nlay),            &! Effective temperature in the cloud
259     pqeff(ngrid, nlay, nq),        &! Effective tracers quantities in the cloud
260     wq(ngrid,nlay+1),              &! Displaced tracer mass (kg.m-2) during microtimestep
261     satuco2(ngrid,nlay),           &! CO2 satu ratio for output diagfi
262     zqsatco2(ngrid,nlay),          &! Saturation co2
263     availco2,&
264     masslayer, &
265     tmp, a,b, &
266     new_pdq(ngrid,nlay)
267     
268  double precision :: &
269! ---Extinction coefficients at 1 micron of CO2 particles
270     vrat_cld,                      &! Volume ratio
271     n_derf,                        &! derf( (rb_cldco2(1)-log(riceco2(ig,l))) *dev2)
272     Qtemp,                         &! mean value in the interval during the interpolation
273     ltemp1(var_dim_qext),          &! abs(radv(:)-rb_cldco2(i))
274     ltemp2(var_dim_qext),          &! abs(radv(:)-rb_cldco2(i+1))
275     n_aer(nbinco2_cld),            &! -0.5 * Nccnco2*tauscaling(ig) * n_derf
276     tau1mic(ngrid),                &! CO2 ice column opacity at 1µm
277     Qext1bins2(ngrid,nlay),        &! CO2 ice opacities
278! ---For Saturation Index computation
279     rho,                           &! background density
280     zu,                            &! absolute value of zonal wind field
281     NN,                            &! N^2 static stability
282     gradT,                         &! thermal gradient
283     SatIndex(ngrid,nlay),          &! Saturation index S in Spiga 2012 paper, assuming like in the paper GW phase speed
284                                     !  (stationary waves): c = 0 m.s-1, lambdaH = 150 km, Fo = 7.5e-7 J.m-3
285     SatIndexmap(ngrid),            &! maxval(SatIndex(ig,12:26))
286! ---Misc
287     myT,                           &! Temperature scalar for co2 density computation
288     tcond(ngrid,nlay)               ! CO2 condensation temperature
289
290  logical :: &
291     file_qext_ok ! Check if file_qext exists
292!======================================================================================================================!
293! BEGIN ===============================================================================================================!
294!======================================================================================================================!
295! 0. Firstcall
296!----------------------------------------------------------------------------------------------------------------------!
297  if (firstcall) then
298    firstcall=.false.
299!----------------------------------------------------------------------------------------------------------------------!
300! 0.1. Initialization of microtimestep from imicroco2
301!----------------------------------------------------------------------------------------------------------------------!
302#ifdef MESOSCALE
303    imicroco2 = 2
304#else
305    imicroco2 = 30
306#endif
307    call getin("imicroco2", imicroco2)
308    microtimestep = ptimestep/real(imicroco2)
309    sigma_iceco2 = sqrt(log(1.+nuiceco2_sed))
310    dev2 = 1. / ( sqrt(2.) * sigma_iceco2 )
311!----------------------------------------------------------------------------------------------------------------------!
312! 0.2. Compute the radius grid of CO2 ice particles (rb_cldco2)
313!        > the grid spacing is computed assuming a constant volume ratio between two consecutive bins; i.e. vrat_cld.
314!   - rad_cldco2 is the primary radius grid used for microphysics computation.
315!   - The grid spacing is computed assuming a constant volume ratio between two consecutive bins; i.e. vrat_cld.
316!   - vrat_cld is determined from the boundary values of the size grid: rmin_cld and rmax_cld.
317!   - The rb_cldco2 array contains the boundary values of each rad_cldco2 bin.
318!----------------------------------------------------------------------------------------------------------------------!
319    ! vrat_cld is determined from the boundary values of the size grid: rmin_cld and rmax_cld.
320    vrat_cld = exp(log(rmax_cld/rmin_cld) / float(nbinco2_cld-1) * 3.)
321
322    ! rad_cldco2 is the primary radius grid used for microphysics computation.
323    rad_cldco2(1) = rmin_cld
324    do i = 1, nbinco2_cld-1
325      rad_cldco2(i+1) = rad_cldco2(i) * vrat_cld**(1./3.)
326    end do
327
328    ! rb_cldco2 array contains the boundary values of each rad_cldco2 bin.
329    rb_cldco2(1) = rbmin_cld
330    do i = 1, nbinco2_cld
331      rb_cldco2(i+1) = ( (2.*vrat_cld) / (vrat_cld+1.) )**(1./3.) * rad_cldco2(i)
332    end do
333    rb_cldco2(nbinco2_cld+1) = rbmax_cld
334!----------------------------------------------------------------------------------------------------------------------!
335! 0.3. Read file 'optprop_co2ice_1mic.dat' to extract optical properties of CO2 ice at 1 micron (Qext)
336!----------------------------------------------------------------------------------------------------------------------!
337    ! get information about file_qext
338    inquire(file=trim(datadir)//'/'//file_qext, exist=file_qext_ok)
339
340    ! if file_qext is missing then stop
341    if (.not. file_qext_ok) then
342      write(*,*)'file'//file_qext//'should be in ', trim(datadir)
343      call abort_physic('co2cloud', 'file missing', 1)
344    end if
345
346    ! read file_qext
347    open(unit=uQext,file=trim(datadir)//'/'//file_qext, form='formatted')
348
349    ! skip 1 line
350    read(uQext,*)
351
352    ! extract radv
353    do i = 1, var_dim_qext
354      read(uQext,'(E12.5)')radv(i)
355    end do
356
357    ! skip 1 line
358    read(uQext,*)
359
360    ! Qextv1mic
361    do i = 1 , var_dim_qext
362      read(uQext,'(E12.5)')Qextv1mic(i)
363    end do
364
365    ! close file_qext
366    close(uQext)
367!----------------------------------------------------------------------------------------------------------------------!
368! 0.4. Interpole the radius grid (rb_cldco2) to get the corresponding exctinction coefficients (Qext1bins), using
369!       file_qext values (radv, Qextv1mic)
370!----------------------------------------------------------------------------------------------------------------------!
371    do i = 1, nbinco2_cld
372      ltemp1 = abs(radv(:)-rb_cldco2(i))
373      ltemp2 = abs(radv(:)-rb_cldco2(i+1))
374      lebon1 = minloc(ltemp1,DIM=1)
375      lebon2 = min(minloc(ltemp2,DIM=1), var_dim_qext)
376      nelem = lebon2 - lebon1 + 1.
377
378      ! mean value in the interval
379      Qtemp = 0d0
380      do l = 0, nelem
381        Qtemp = Qtemp + Qextv1mic(min(lebon1+l, var_dim_qext))
382      end do
383
384      Qext1bins(i) = Qtemp / nelem
385    end do
386
387    Qext1bins(:) = Qext1bins(:) * pi * (rad_cldco2(:)**2)
388
389    ! print result of the interpolation
390    write(*,*)'--------------------------------------------'
391    write(*,*)'Microphysics co2: size bin-Qext information:'
392    write(*,*)'   i, rad_cldco2(i), Qext1bins(i)'
393    do i = 1, nbinco2_cld
394      write(*,'(i3,3x,3(e13.6,4x))')i, rad_cldco2(i), Qext1bins(i)
395    end do
396    write(*,*)'--------------------------------------------'
397!----------------------------------------------------------------------------------------------------------------------!
398! 0.5. Save the radius grid of CO2 particles (rb_cldco2)
399!----------------------------------------------------------------------------------------------------------------------!
400    do i = 1, nbinco2_cld+1
401      rb_cldco2(i) = log(rb_cldco2(i))
402    end do
403  end if ! of IF (firstcall)
404!----------------------------------------------------------------------------------------------------------------------!
405! 1. Initialization
406!----------------------------------------------------------------------------------------------------------------------!
407  sum_subpdq(1:ngrid,1:nlay,1:nq) = 0.
408  sum_subpdt(1:ngrid,1:nlay) = 0.
409
410  subpdqcloudco2(1:ngrid,1:nlay,1:nq) = 0.
411  subpdtcloudco2(1:ngrid,1:nlay) = 0.
412
413  pdqcloudco2(1:ngrid,1:nlay,1:nq) = 0.
414  pdtcloudco2(1:ngrid,1:nlay) = 0.
415
416  ! default value if no ice
417  rhocloudco2(1:ngrid,1:nlay) = rho_dust
418  rhocloudco2t(1:ngrid,1:nlay) = rho_dust
419  epaisseur(1:ngrid,1:nlay) = 0.
420  masse(1:ngrid,1:nlay) = 0.
421  riceco2(1:ngrid, 1:nlay) = 0.
422  zqsed0(1:ngrid,1:nlay,1:nq) = 0.
423  sum_subpdqs_sedco2(1:ngrid) = 0.
424  sum_subpdqs_sedccn(1:ngrid,1:nq) = 0.
425  subpdqsed(1:ngrid,1:nlay,1:nq) = 0.
426!----------------------------------------------------------------------------------------------------------------------!
427! 2. Compute mass and thickness layers
428!----------------------------------------------------------------------------------------------------------------------!
429  do l = 1, nlay
430    do ig = 1, ngrid
431#ifndef MESOSCALE
432      masse(ig,l) = (pplev(ig,l) - pplev(ig,l+1) + (bp(l)-bp(l+1)) ) / g
433#else
434      masse(ig,l) = (pplev(ig,l) - pplev(ig,l+1)) / g
435#endif
436      epaisseur(ig,l) = pzlev(ig,l+1) - pzlev(ig,l)
437    end do
438  end do
439!----------------------------------------------------------------------------------------------------------------------!
440! 3. Define the sub-grid cloud (CLFvaryingCO2)
441!----------------------------------------------------------------------------------------------------------------------!
442! 3.1. Representation of sub-grid CO2 ice clouds (CLFvaryingCO2 = True)
443!----------------------------------------------------------------------------------------------------------------------!
444  if (CLFvaryingCO2) then
445    ! effective temperature
446    pteff(:,:) = pt(:,:)
447
448    ! min co2cloudfrac when there is ice
449    co2cloudfrac(:,:) = mincloud
450
451    ! temperature
452    do l=1,nlay
453      do ig=1,ngrid
454        zt(ig,l) = pt(ig,l) + pdt(ig,l)*ptimestep
455      end do
456    end do
457
458    ! Quantities of traceurs
459    if (igcm_co2 /= 0) then
460      do l = 1, nlay
461        do ig = 1, ngrid
462          zq_co2vap(ig,l) = pq(ig,l,igcm_co2) + pdq(ig,l,igcm_co2)*ptimestep
463        end do
464      end do
465    end if
466!----------------------------------------------------------------------------------------------------------------------!
467! 3.1.a. Saturation index CO2
468!----------------------------------------------------------------------------------------------------------------------!
469    ! if saturation index co2 is true
470    if (satindexco2) then
471      ! layers 12 --> 26 ~ 12->85 km
472      do l = 12, 26
473        do ig = 1, ngrid
474          ! compute N^2 static stability
475          gradT = (zt(ig,l+1)-zt(ig,l))/(pzlev(ig,l+1)-pzlev(ig,l))
476          NN = sqrt(g/zt(iq,l) * (max(gradT,-g/cpp) + g/cpp))
477
478          ! compute absolute value of zonal wind field
479          zu = abs(pu(ig,l) + pdu(ig,l)*ptimestep)
480
481          ! compute background density
482          rho = pplay(ig,l) / (rnew(ig,l)*zt(ig,l))
483
484          ! saturation index: Modulate the DeltaT by GW propagation index:
485          ! --------------------------------------------------------------
486          SatIndex(ig,l) = sqrt(Fo*lambdaH/(2.*pi)*NN / (rho*zu**3) )
487        end do
488      end do
489
490      ! Then compute Satindex map in layers 12 --> 26 ~ 12->85 km
491      do ig = 1, ngrid
492        SatIndexmap(ig) = maxval(SatIndex(ig,12:26))
493      end do
494
495      ! Write outputs in diagfi.nc
496      call WRITEDIAGFI(ngrid, "SatIndex", "SatIndex", " ", 3, SatIndex)
497
498      call WRITEDIAGFI(ngrid, "SatIndexmap", "SatIndexmap", "km", 2, SatIndexmap)
499    !------------------------------------------------------------------------------------------------------------------!
500    ! if saturation index co2 is false, set saturation index to 0.05
501    !------------------------------------------------------------------------------------------------------------------!
502    else
503      do ig = 1, ngrid
504        SatIndexmap(ig)=0.05
505      end do
506    end if ! of if (satindexco2)
507!----------------------------------------------------------------------------------------------------------------------!
508! 3.1.b. Compute tcond
509!----------------------------------------------------------------------------------------------------------------------!
510    call tcondco2(ngrid,nlay,pplay,zq_co2vap,tcond)
511!----------------------------------------------------------------------------------------------------------------------!
512! 3.1.c. Compute cloud fraction in cells
513!----------------------------------------------------------------------------------------------------------------------!
514    do ig = 1, ngrid
515      if (SatIndexmap(ig) <= 0.1) then
516        do l = 1, nlay-1
517
518          ! The entire fraction is saturated
519          if (tcond(ig,l) >= (zt(ig,l)+zdelt) .or. tcond(ig,l) <= 0.) then
520            pteff(ig,l) = zt(ig,l)
521            co2cloudfrac(ig,l) = 1.
522
523          ! No saturation at all
524          else if (tcond(ig,l) <= (zt(ig,l)-zdelt)) then
525            pteff(ig,l) = zt(ig,l) - zdelt
526            co2cloudfrac(ig,l) = mincloud
527
528          ! Mean temperature of the cloud fraction
529          else
530            pteff(ig,l) = (tcond(ig,l)+zt(ig,l)-zdelt) / 2.
531            co2cloudfrac(ig,l) = (tcond(ig,l)-zt(ig,l)+zdelt) / (2.0*zdelt)
532          end if
533
534          pteff(ig,l) = pteff(ig,l) - pdt(ig,l)*ptimestep
535
536          ! check boundary values of co2cloudfrac
537          if (co2cloudfrac(ig,l) <= mincloud) then
538            co2cloudfrac(ig,l) = mincloud
539          else if (co2cloudfrac(ig,l)> 1) then
540            co2cloudfrac(ig,l) = 1.
541          end if
542        end do
543
544      ! SatIndex not favorable for GW: leave pt untouched
545      else
546        pteff(ig,l) = pt(ig,l)
547        co2cloudfrac(ig,l) = mincloud
548      end if ! of if (SatIndexmap <= 0.1)
549    end do ! of ngrid
550    ! TODO: Totalcloud frac of the column missing here
551!----------------------------------------------------------------------------------------------------------------------!
552! 3.2. No sub-grid cloud representation (CLFvarying = False)
553!----------------------------------------------------------------------------------------------------------------------!
554  else
555    do l = 1, nlay
556      do ig = 1, ngrid
557        pteff(ig,l) = pt(ig,l)
558      end do
559    end do
560  end if ! end if (CLFvaryingco2)
561!----------------------------------------------------------------------------------------------------------------------!
562! 4. Microphysics of CO2 cloud formation
563!----------------------------------------------------------------------------------------------------------------------!
564  pqeff(:,:,:) = pq(:,:,:)
565  pteff(:,:) = pt(:,:)
566!----------------------------------------------------------------------------------------------------------------------!
567! 4.2. Effective tracers quantities in the cloud
568!----------------------------------------------------------------------------------------------------------------------!
569  if (CLFvaryingCO2) then
570    pqeff(:,:,igcm_ccnco2_mass) = pq(:,:,igcm_ccnco2_mass) / co2cloudfrac(:,:)
571
572    pqeff(:,:,igcm_ccnco2_number) = pq(:,:,igcm_ccnco2_number) / co2cloudfrac(:,:)
573
574    pqeff(:,:,igcm_co2_ice) = pq(:,:,igcm_co2_ice) / co2cloudfrac(:,:)
575  end if
576!----------------------------------------------------------------------------------------------------------------------!
577  do microstep = 1, imicroco2
578!----------------------------------------------------------------------------------------------------------------------!
579! 4.1. Stepped entry for tendancies: At each micro timestep we add pdt in order to have a stepped entry
580!----------------------------------------------------------------------------------------------------------------------!
581   do l = 1, nlay
582     do ig = 1, ngrid
583       ! on temperature
584       sum_subpdt(ig,l) = sum_subpdt(ig,l) + pdt(ig,l)
585       
586       ! on tracers
587       sum_subpdq(ig,l,igcm_dust_mass) = sum_subpdq(ig,l,igcm_dust_mass) + pdq(ig,l,igcm_dust_mass)
588
589       sum_subpdq(ig,l,igcm_dust_number) = sum_subpdq(ig,l,igcm_dust_number) + pdq(ig,l,igcm_dust_number)
590
591       sum_subpdq(ig,l,igcm_ccnco2_mass) = sum_subpdq(ig,l,igcm_ccnco2_mass) + pdq(ig,l,igcm_ccnco2_mass)
592
593       sum_subpdq(ig,l,igcm_ccnco2_number) = sum_subpdq(ig,l,igcm_ccnco2_number) + pdq(ig,l,igcm_ccnco2_number)
594
595       sum_subpdq(ig,l,igcm_co2_ice) = sum_subpdq(ig,l,igcm_co2_ice) + pdq(ig,l,igcm_co2_ice)
596
597       sum_subpdq(ig,l,igcm_co2) = sum_subpdq(ig,l,igcm_co2) + pdq(ig,l,igcm_co2)
598
599       if (meteo_flux) then
600        sum_subpdq(ig,l,igcm_ccnco2_meteor_number) = sum_subpdq(ig,l,igcm_ccnco2_meteor_number) + &
601                                                  pdq(ig,l,igcm_ccnco2_meteor_number)
602
603        sum_subpdq(ig,l,igcm_ccnco2_meteor_mass) = sum_subpdq(ig,l,igcm_ccnco2_meteor_mass) + &
604                                                    pdq(ig,l,igcm_ccnco2_meteor_mass)
605       end if
606       if (co2useh2o) then
607         sum_subpdq(ig,l,igcm_h2o_ice) = sum_subpdq(ig,l,igcm_h2o_ice) + pdq(ig,l,igcm_h2o_ice)
608
609         sum_subpdq(ig,l,igcm_ccn_mass) = sum_subpdq(ig,l,igcm_ccn_mass) + pdq(ig,l,igcm_ccn_mass)
610
611         sum_subpdq(ig,l,igcm_ccn_number) = sum_subpdq(ig,l,igcm_ccn_number) + pdq(ig,l,igcm_ccn_number)
612
613         sum_subpdq(ig,l,igcm_ccnco2_h2o_number) = sum_subpdq(ig,l,igcm_ccnco2_h2o_number) + &
614                                                   pdq(ig,l,igcm_ccnco2_h2o_number)
615
616         sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) = sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) + &
617                                                     pdq(ig,l,igcm_ccnco2_h2o_mass_ice)
618
619         sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) = sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) + &
620                                                     pdq(ig,l,igcm_ccnco2_h2o_mass_ccn)
621       end if
622     end do ! ngrid
623   end do ! nlay
624!----------------------------------------------------------------------------------------------------------------------!
625! 4.4. Main call to the cloud scheme
626!----------------------------------------------------------------------------------------------------------------------!
627   call improvedco2clouds(ngrid, nlay, microtimestep, pplay, pplev, pteff, sum_subpdt, pqeff, sum_subpdq, &
628                          subpdqcloudco2, subpdtcloudco2, nq, tauscaling, rb_cldco2, sigma_iceco2, dev2)
629
630   do l = 1, nlay
631     do ig = 1, ngrid
632       if(pq(ig,l,igcm_co2_ice) + microtimestep*(sum_subpdq(ig,l,igcm_co2_ice)+subpdqcloudco2(ig,l,igcm_co2_ice)) &
633         <= 1.e-12) then
634         subpdqcloudco2(ig,l,igcm_co2_ice) = -pq(ig,l,igcm_co2_ice)/microtimestep - sum_subpdq(ig,l,igcm_co2_ice)
635         subpdqcloudco2(ig,l,igcm_co2) = -subpdqcloudco2(ig,l,igcm_co2_ice)
636       end if
637
638       if(pq(ig,l,igcm_co2) + microtimestep*(sum_subpdq(ig,l,igcm_co2)+subpdqcloudco2(ig,l,igcm_co2)) <= 1.e-12) then
639         subpdqcloudco2(ig,l,igcm_co2) =  - pq(ig,l,igcm_co2)/microtimestep - sum_subpdq(ig,l,igcm_co2)
640         subpdqcloudco2(ig,l,igcm_co2_ice) = -subpdqcloudco2(ig,l,igcm_co2)
641       end if
642
643  ! ccnco2_number and ccnco2_mass
644       if (((pq(ig,l,igcm_ccnco2_number)+(sum_subpdq(ig,l,igcm_ccnco2_number)+subpdqcloudco2(ig,l,igcm_ccnco2_number)) &
645          *microtimestep)<=1.) .or. &
646         (pq(ig,l,igcm_ccnco2_mass)+(sum_subpdq(ig,l,igcm_ccnco2_mass)+subpdqcloudco2(ig,l,igcm_ccnco2_mass)) &
647          *microtimestep<=1e-20)) then
648         subpdqcloudco2(ig,l,igcm_ccnco2_number) = - pq(ig,l,igcm_ccnco2_number)/microtimestep + 1. &
649                                                  - sum_subpdq(ig,l,igcm_ccnco2_number)
650         subpdqcloudco2(ig,l,igcm_dust_number) = - subpdqcloudco2(ig,l,igcm_ccnco2_number)
651
652         subpdqcloudco2(ig,l,igcm_ccnco2_mass) = - pq(ig,l,igcm_ccnco2_mass)/microtimestep + 1e-20 &
653                                                - sum_subpdq(ig,l,igcm_ccnco2_mass)
654         subpdqcloudco2(ig,l,igcm_dust_mass) = - subpdqcloudco2(ig,l,igcm_ccnco2_mass)
655       end if
656
657       ! ccnco2_meteor_number and ccnco2_meteor_mass
658       if (meteo_flux) then
659         if (((pq(ig,l,igcm_ccnco2_meteor_number)+(sum_subpdq(ig,l,igcm_ccnco2_meteor_number)+ &
660               subpdqcloudco2(ig,l,igcm_ccnco2_meteor_number))*microtimestep)<=1.) .or. &
661             (pq(ig,l,igcm_ccnco2_meteor_mass)+(sum_subpdq(ig,l,igcm_ccnco2_meteor_mass)+ &
662               subpdqcloudco2(ig,l,igcm_ccnco2_meteor_mass))*microtimestep<=1e-20)) then
663           subpdqcloudco2(ig,l,igcm_ccnco2_meteor_number) = - pq(ig,l,igcm_ccnco2_meteor_number)/microtimestep + 1. &
664                                                    - sum_subpdq(ig,l,igcm_ccnco2_meteor_number)
665             subpdqcloudco2(ig,l,igcm_ccnco2_meteor_mass) = - pq(ig,l,igcm_ccnco2_meteor_mass)/microtimestep + 1e-20 &
666                                                  - sum_subpdq(ig,l,igcm_ccnco2_meteor_mass)
667         end if
668       end if
669
670       ! ccnco2_h2o_number and masses
671       if (co2useh2o) then
672         if (((pq(ig,l,igcm_ccnco2_h2o_number) + (sum_subpdq(ig,l,igcm_ccnco2_h2o_number) + &
673           subpdqcloudco2(ig,l,igcm_ccnco2_h2o_number)) &
674            *microtimestep)<=1.) .or. &
675           (pq(ig,l,igcm_ccnco2_h2o_mass_ice)+pq(ig,l,igcm_ccnco2_h2o_mass_ccn) +&
676           (sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice)+subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) +&
677           sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn)+subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn)) &
678            *microtimestep<=1e-20)) then
679
680           subpdqcloudco2(ig,l,igcm_ccnco2_h2o_number) = - pq(ig,l,igcm_ccnco2_h2o_number)/microtimestep + 1. &
681                                                         - sum_subpdq(ig,l,igcm_ccnco2_h2o_number)
682           subpdqcloudco2(ig,l,igcm_ccn_number) = - subpdqcloudco2(ig,l,igcm_ccnco2_h2o_number)
683
684           subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) = - pq(ig,l,igcm_ccnco2_h2o_mass_ice)/microtimestep + 1e-20 &
685                                                           - sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice)
686
687           subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn) = - pq(ig,l,igcm_ccnco2_h2o_mass_ccn)/microtimestep + 1e-20 &
688                                                           - sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn)
689
690           subpdqcloudco2(ig,l,igcm_ccn_mass) = - subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) &
691                                                - subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn)
692         end if
693       end if
694     end do
695   end do
696!----------------------------------------------------------------------------------------------------------------------!
697! 4.5. Updating tendencies after cloud scheme
698!----------------------------------------------------------------------------------------------------------------------!
699    do l = 1, nlay
700      do ig = 1, ngrid
701        sum_subpdt(ig,l) = sum_subpdt(ig,l) + subpdtcloudco2(ig,l)
702
703        sum_subpdq(ig,l,igcm_dust_mass) = sum_subpdq(ig,l,igcm_dust_mass) + subpdqcloudco2(ig,l,igcm_dust_mass)
704
705        sum_subpdq(ig,l,igcm_dust_number) = sum_subpdq(ig,l,igcm_dust_number) + subpdqcloudco2(ig,l,igcm_dust_number)
706
707        sum_subpdq(ig,l,igcm_ccnco2_mass) = sum_subpdq(ig,l,igcm_ccnco2_mass) + subpdqcloudco2(ig,l,igcm_ccnco2_mass)
708
709        sum_subpdq(ig,l,igcm_ccnco2_number) = sum_subpdq(ig,l,igcm_ccnco2_number) + &
710                                              subpdqcloudco2(ig,l,igcm_ccnco2_number)
711
712        sum_subpdq(ig,l,igcm_co2_ice) = sum_subpdq(ig,l,igcm_co2_ice) + subpdqcloudco2(ig,l,igcm_co2_ice)
713
714        sum_subpdq(ig,l,igcm_co2) = sum_subpdq(ig,l,igcm_co2) + subpdqcloudco2(ig,l,igcm_co2)
715
716        if (meteo_flux) then
717          sum_subpdq(ig,l,igcm_ccnco2_meteor_mass) = sum_subpdq(ig,l,igcm_ccnco2_meteor_mass) + &
718                                                     subpdqcloudco2(ig,l,igcm_ccnco2_meteor_mass)
719
720          sum_subpdq(ig,l,igcm_ccnco2_meteor_number) = sum_subpdq(ig,l,igcm_ccnco2_meteor_number) + &
721                                                       subpdqcloudco2(ig,l,igcm_ccnco2_meteor_number)
722        end if
723        if (co2useh2o) then
724          sum_subpdq(ig,l,igcm_h2o_ice) = sum_subpdq(ig,l,igcm_h2o_ice) + subpdqcloudco2(ig,l,igcm_h2o_ice)
725
726          sum_subpdq(ig,l,igcm_ccn_mass) = sum_subpdq(ig,l,igcm_ccn_mass) + subpdqcloudco2(ig,l,igcm_ccn_mass)
727
728          sum_subpdq(ig,l,igcm_ccn_number) = sum_subpdq(ig,l,igcm_ccn_number) + subpdqcloudco2(ig,l,igcm_ccn_number)
729
730          sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) = sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) + &
731                                                      subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice)
732
733          sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) = sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) + &
734                                                      subpdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn)
735
736          sum_subpdq(ig,l,igcm_ccnco2_h2o_number) = sum_subpdq(ig,l,igcm_ccnco2_h2o_number) + &
737                                                    subpdqcloudco2(ig,l,igcm_ccnco2_number)
738
739        end if
740      end do ! ngrid
741    end do ! nlay
742!----------------------------------------------------------------------------------------------------------------------!
743! 4.3. Gravitational sedimentation
744!----------------------------------------------------------------------------------------------------------------------!
745    if (sedimentation) then
746!----------------------------------------------------------------------------------------------------------------------!
747! 4.3.a. Compute cloud density
748!----------------------------------------------------------------------------------------------------------------------!
749      do l = 1, nlay
750        do ig = 1, ngrid
751          ! temperature during the sedimentation process
752          ztsed(ig,l) = pteff(ig,l) + sum_subpdt(ig,l) * microtimestep
753
754          ! quantities tracers during the sedimentation process
755          zqsed(ig,l,:) = pqeff(ig,l,:) + sum_subpdq(ig,l,:) * microtimestep
756
757          ! assure positive value of co2_ice mmr, ccnco2 number, ccnco2 mass
758          ! meteoritic particle are considered like dust, rho_dust
759          Niceco2 = max(zqsed(ig,l,igcm_co2_ice), threshold)
760          Nccnco2 = max(zqsed(ig,l,igcm_ccnco2_number), threshold)
761          Qccnco2 = max(zqsed(ig,l,igcm_ccnco2_mass), threshold)
762
763          Nccnco2_h2o = 0.
764          Qccnco2_h2o = 0.
765          if (co2useh2o) then
766            Nccnco2_h2o = max(zqsed(ig,l,igcm_ccnco2_h2o_number), threshold)
767            Qccnco2_h2o = max(zqsed(ig,l,igcm_ccnco2_h2o_mass_ice)+zqsed(ig,l,igcm_ccnco2_h2o_mass_ccn), threshold)
768            Nccnco2 = Nccnco2 - Nccnco2_h2o
769            Qccnco2 = Qccnco2 - Qccnco2_h2o
770          end if
771
772          ! Get density cloud and co2 ice particle radius
773          if (Niceco2/=0d0) then
774            call updaterice_microco2(dble(Niceco2), dble(Qccnco2), dble(Nccnco2), dble(Qccnco2_h2o), dble(Nccnco2_h2o),&
775                                     ztsed(ig,l), tauscaling(ig), riceco2(ig,l), rhocloudco2t(ig,l))
776          else
777            riceco2(ig,l) = 0.
778            rhocloudco2t(ig,l) = 0.
779          end if
780        end do ! ngrid
781      end do ! nlay
782!----------------------------------------------------------------------------------------------------------------------!
783! 4.3.b. Save actualized tracer values to compute sedimentation tendancies
784!----------------------------------------------------------------------------------------------------------------------!
785      zqsed0(:,:,igcm_co2_ice) = zqsed(:,:,igcm_co2_ice)
786      zqsed0(:,:,igcm_ccnco2_mass) = zqsed(:,:,igcm_ccnco2_mass)
787      zqsed0(:,:,igcm_ccnco2_number) = zqsed(:,:,igcm_ccnco2_number)
788
789      if (meteo_flux) then
790        zqsed0(:,:,igcm_ccnco2_meteor_mass) = zqsed(:,:,igcm_ccnco2_meteor_mass)
791        zqsed0(:,:,igcm_ccnco2_meteor_number) = zqsed(:,:,igcm_ccnco2_meteor_number)
792      end if
793
794      if (co2useh2o) then
795        zqsed0(:,:,igcm_ccnco2_h2o_number) = zqsed(:,:,igcm_ccnco2_h2o_number)
796        zqsed0(:,:,igcm_ccnco2_h2o_mass_ice) = zqsed(:,:,igcm_ccnco2_h2o_mass_ice)
797        zqsed0(:,:,igcm_ccnco2_h2o_mass_ccn) = zqsed(:,:,igcm_ccnco2_h2o_mass_ccn)
798      end if
799!----------------------------------------------------------------------------------------------------------------------!
800! 4.3.c. Sedimentation of co2 ice
801!----------------------------------------------------------------------------------------------------------------------!
802      do ig = 1, ngrid
803        do l = 1, nlay
804          rsedcloudco2(ig,l) = max( riceco2(ig,l)*(1.+sigma_iceco2)*(1.+sigma_iceco2)*(1.+sigma_iceco2), &
805                                    rdust(ig,l) )
806        end do
807      end do
808
809      wq(:,:) = 0.
810      call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay, microtimestep, pplev, masse, epaisseur, ztsed, &
811                   rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_co2_ice), wq, beta)
812
813      do ig = 1, ngrid
814        sum_subpdqs_sedco2(ig) = sum_subpdqs_sedco2(ig) + wq(ig,1) / microtimestep !wq est en kg.m-2
815      end do
816!----------------------------------------------------------------------------------------------------------------------!
817! 4.3.d. Sedimentation for other tracers
818!----------------------------------------------------------------------------------------------------------------------!
819      wq(:,:) = 0.
820      ! for ccnco2_mass
821      call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay, microtimestep, pplev, masse, epaisseur, ztsed, &
822                   rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_mass), wq, beta)
823      do ig = 1, ngrid
824        sum_subpdqs_sedccn(ig,igcm_ccnco2_mass) = sum_subpdqs_sedccn(ig,igcm_ccnco2_mass) + wq(ig,1) / microtimestep
825      end do
826
827      wq(:,:) = 0.
828      ! for ccnco2_number
829      call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay,microtimestep, pplev, masse, epaisseur, ztsed, &
830                   rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_number), wq, beta)
831      do ig = 1, ngrid
832        sum_subpdqs_sedccn(ig,igcm_ccnco2_number) = sum_subpdqs_sedccn(ig,igcm_ccnco2_number) + wq(ig,1) / microtimestep
833      end do
834
835      if (meteo_flux) then
836        wq(:,:) = 0.
837        ! for ccnco2_meteor_mass
838        call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay, microtimestep, pplev, masse, epaisseur, ztsed, &
839                     rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_meteor_mass), wq, beta)
840        do ig = 1, ngrid
841          sum_subpdqs_sedccn(ig,igcm_ccnco2_meteor_mass) = sum_subpdqs_sedccn(ig,igcm_ccnco2_meteor_mass) + &
842                                                           wq(ig,1)/ microtimestep
843        end do
844
845        wq(:,:) = 0.
846        ! for ccnco2_meteor_number
847        call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay,microtimestep, pplev, masse, epaisseur, ztsed, &
848                     rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_meteor_number), wq, beta)
849        do ig = 1, ngrid
850          sum_subpdqs_sedccn(ig,igcm_ccnco2_meteor_number) = sum_subpdqs_sedccn(ig,igcm_ccnco2_meteor_number) + &
851                                                             wq(ig,1) / microtimestep
852        end do
853      end if
854      ! for ccnco2_h2o_mass_ice
855      if (co2useh2o) then
856        wq(:,:) = 0.
857        call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay,microtimestep, pplev, masse, epaisseur, ztsed, &
858                     rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_h2o_mass_ice), wq, beta)
859        do ig = 1, ngrid
860          sum_subpdqs_sedccn(ig,igcm_ccnco2_h2o_mass_ice) = sum_subpdqs_sedccn(ig,igcm_ccnco2_h2o_mass_ice) + &
861                                                            wq(ig,1) / microtimestep
862        end do
863
864        wq(:,:) = 0.
865        ! for ccnco2_h2o_mass_ccn
866        call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay,microtimestep, pplev, masse, epaisseur, ztsed, &
867                     rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_h2o_mass_ccn), wq, beta)
868        do ig = 1, ngrid
869          sum_subpdqs_sedccn(ig,igcm_ccnco2_h2o_mass_ccn) = sum_subpdqs_sedccn(ig,igcm_ccnco2_h2o_mass_ccn) + &
870                                                            wq(ig,1) / microtimestep
871        end do
872
873        wq(:,:) = 0.
874        ! for ccnco2_h2o_number
875        call newsedim(ngrid, nlay, ngrid*nlay, ngrid*nlay,microtimestep, pplev, masse, epaisseur, ztsed, &
876                     rsedcloudco2, rhocloudco2t, zqsed(:,:,igcm_ccnco2_h2o_number), wq, beta)
877        do ig = 1, ngrid
878          sum_subpdqs_sedccn(ig,igcm_ccnco2_h2o_number) = sum_subpdqs_sedccn(ig,igcm_ccnco2_h2o_number) + &
879                                                            wq(ig,1) / microtimestep
880        end do
881      end if
882!----------------------------------------------------------------------------------------------------------------------!
883! 4.3.e. Compute tendencies due to the sedimation process
884!----------------------------------------------------------------------------------------------------------------------!
885      do l = 1, nlay
886        do ig = 1, ngrid
887          subpdqsed(ig,l,igcm_ccnco2_mass) = ( zqsed(ig,l,igcm_ccnco2_mass) - zqsed0(ig,l,igcm_ccnco2_mass)  ) &
888                                             / microtimestep
889
890          subpdqsed(ig,l,igcm_ccnco2_number) = ( zqsed(ig,l,igcm_ccnco2_number) - zqsed0(ig,l,igcm_ccnco2_number) )&
891                                               / microtimestep
892
893          subpdqsed(ig,l,igcm_co2_ice) = ( zqsed(ig,l,igcm_co2_ice) - zqsed0(ig,l,igcm_co2_ice) ) / microtimestep
894
895          if (meteo_flux) then
896            subpdqsed(ig,l,igcm_ccnco2_meteor_mass) = ( zqsed(ig,l,igcm_ccnco2_meteor_mass) - &
897                                                        zqsed0(ig,l,igcm_ccnco2_meteor_mass) ) / microtimestep
898
899            subpdqsed(ig,l,igcm_ccnco2_meteor_number) = ( zqsed(ig,l,igcm_ccnco2_meteor_number) - &
900                                                          zqsed0(ig,l,igcm_ccnco2_meteor_number) ) / microtimestep
901          end if
902          if (co2useh2o) then
903          subpdqsed(ig,l,igcm_ccnco2_h2o_number) = ( zqsed(ig,l,igcm_ccnco2_h2o_number) - &
904                                                     zqsed0(ig,l,igcm_ccnco2_h2o_number) ) / microtimestep
905
906          subpdqsed(ig,l,igcm_ccnco2_h2o_mass_ice) = ( zqsed(ig,l,igcm_ccnco2_h2o_mass_ice) - &
907                                                     zqsed0(ig,l,igcm_ccnco2_h2o_mass_ice) ) / microtimestep
908
909          subpdqsed(ig,l,igcm_ccnco2_h2o_mass_ccn) = ( zqsed(ig,l,igcm_ccnco2_h2o_mass_ccn) - &
910                                                     zqsed0(ig,l,igcm_ccnco2_h2o_mass_ccn) ) / microtimestep
911          end if
912        end do
913      end do
914      ! update subtimestep tendencies with sedimentation input
915      do l = 1, nlay
916        do ig = 1, ngrid
917          sum_subpdq(ig,l,igcm_ccnco2_mass) = sum_subpdq(ig,l,igcm_ccnco2_mass) + subpdqsed(ig,l,igcm_ccnco2_mass)
918
919          sum_subpdq(ig,l,igcm_ccnco2_number) = sum_subpdq(ig,l,igcm_ccnco2_number) + subpdqsed(ig,l,igcm_ccnco2_number)
920
921          sum_subpdq(ig,l,igcm_co2_ice) = sum_subpdq(ig,l,igcm_co2_ice) + subpdqsed(ig,l,igcm_co2_ice)
922          if (meteo_flux) then
923            sum_subpdq(ig,l,igcm_ccnco2_meteor_mass) = sum_subpdq(ig,l,igcm_ccnco2_meteor_mass) + &
924                                                       subpdqsed(ig,l,igcm_ccnco2_meteor_mass)
925
926            sum_subpdq(ig,l,igcm_ccnco2_meteor_number) = sum_subpdq(ig,l,igcm_ccnco2_meteor_number) + &
927                                                         subpdqsed(ig,l,igcm_ccnco2_meteor_number)
928          end if
929          if (co2useh2o) then
930            sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) = sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) + &
931                                                        subpdqsed(ig,l,igcm_ccnco2_h2o_mass_ice)
932
933            sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) = sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) + &
934                                                        subpdqsed(ig,l,igcm_ccnco2_h2o_mass_ccn)
935
936            sum_subpdq(ig,l,igcm_ccnco2_h2o_number) = sum_subpdq(ig,l,igcm_ccnco2_h2o_number) + &
937                                                        subpdqsed(ig,l,igcm_ccnco2_h2o_number)
938          end if
939        end do
940      end do
941    end if !(end if sedimentation)
942  end do ! of do microstep = 1, imicroco2
943!----------------------------------------------------------------------------------------------------------------------!
944! 5. Compute final tendencies after time loop
945!----------------------------------------------------------------------------------------------------------------------!
946! condensation/sublimation rate in the atmosphere (kg.kg-1.s-1)
947  do l = nlay, 1, -1
948    do ig = 1, ngrid
949      pcondicea(ig,l) = sum_subpdq(ig,l,igcm_co2_ice) / real(imicroco2)
950    end do
951  end do
952
953  ! CO2 flux at surface (kg.m-2.s-1)
954  do ig = 1, ngrid
955    pdqs_sedco2(ig) = sum_subpdqs_sedco2(ig) / real(imicroco2)
956    pdqs_sedccn(ig,:) = sum_subpdqs_sedccn(ig,:) / real(imicroco2)
957  end do
958
959  ! temperature tendency (T.s-1)
960  do l = 1, nlay
961    do ig = 1, ngrid
962      pdtcloudco2(ig,l) = ( sum_subpdt(ig,l)/real(imicroco2) ) - pdt(ig,l)
963    end do
964  end do
965
966  ! tracers tendencies
967  do l = 1, nlay
968    do ig = 1, ngrid
969      pdqcloudco2(ig,l,igcm_co2) = 0. ! here is the trick, this tendency is computed in co2condens_mod4micro
970
971      pdqcloudco2(ig,l,igcm_co2_ice) = ( sum_subpdq(ig,l,igcm_co2_ice) / real(imicroco2) ) - pdq(ig,l,igcm_co2_ice)
972
973      pdqcloudco2(ig,l,igcm_ccnco2_mass) = ( sum_subpdq(ig,l,igcm_ccnco2_mass)/real(imicroco2) ) - &
974                                           pdq(ig,l,igcm_ccnco2_mass)
975
976      pdqcloudco2(ig,l,igcm_ccnco2_number) = ( sum_subpdq(ig,l,igcm_ccnco2_number) / real(imicroco2) ) - &
977                                             pdq(ig,l,igcm_ccnco2_number)
978
979      pdqcloudco2(ig,l,igcm_dust_mass) = ( sum_subpdq(ig,l,igcm_dust_mass) / real(imicroco2) ) - &
980                                          pdq(ig,l,igcm_dust_mass)
981
982      pdqcloudco2(ig,l,igcm_dust_number) = ( sum_subpdq(ig,l,igcm_dust_number)/real(imicroco2) ) - &
983                                           pdq(ig,l,igcm_dust_number)
984
985      if (meteo_flux) then
986        pdqcloudco2(ig,l,igcm_ccnco2_meteor_mass) = ( sum_subpdq(ig,l,igcm_ccnco2_meteor_mass)/real(imicroco2) ) &
987                                                    - pdq(ig,l,igcm_ccnco2_meteor_mass)
988
989        pdqcloudco2(ig,l,igcm_ccnco2_meteor_number) = ( sum_subpdq(ig,l,igcm_ccnco2_meteor_number) / real(imicroco2) ) &
990                                                       - pdq(ig,l,igcm_ccnco2_meteor_number)
991      end if
992      if (co2useh2o) then
993        pdqcloudco2(ig,l,igcm_h2o_ice) = ( sum_subpdq(ig,l,igcm_h2o_ice) / real(imicroco2) ) - &
994                                          pdq(ig,l,igcm_h2o_ice)
995
996        pdqcloudco2(ig,l,igcm_ccn_mass) = ( sum_subpdq(ig,l,igcm_ccn_mass) / real(imicroco2) ) - &
997                                           pdq(ig,l,igcm_ccn_mass)
998
999        pdqcloudco2(ig,l,igcm_ccn_number) = ( sum_subpdq(ig,l,igcm_ccn_number) / real(imicroco2) ) - &
1000                                            pdq(ig,l,igcm_ccn_number)
1001
1002        pdqcloudco2(ig,l,igcm_ccnco2_h2o_number) = ( sum_subpdq(ig,l,igcm_ccnco2_h2o_number) / real(imicroco2) ) - &
1003                                                   pdq(ig,l,igcm_ccnco2_h2o_number)
1004
1005        pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) = ( sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ice) / real(imicroco2) ) - &
1006                                                   pdq(ig,l,igcm_ccnco2_h2o_mass_ice)
1007
1008        pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn) = ( sum_subpdq(ig,l,igcm_ccnco2_h2o_mass_ccn) / real(imicroco2) )- &
1009                                                        pdq(ig,l,igcm_ccnco2_h2o_mass_ccn)
1010
1011      end if
1012    end do ! ngrid
1013  end do ! nlay
1014!----------------------------------------------------------------------------------------------------------------------!
1015! 6. Update clouds physical values in the cloud (for output)
1016!----------------------------------------------------------------------------------------------------------------------!
1017! 6.1. Update density of co2 ice, riceco2 and opacity
1018!----------------------------------------------------------------------------------------------------------------------!
1019  do l = 1, nlay
1020    do ig = 1, ngrid
1021      Niceco2 = pqeff(ig,l,igcm_co2_ice) + ( pdq(ig,l,igcm_co2_ice) + pdqcloudco2(ig,l,igcm_co2_ice) ) * ptimestep
1022      Niceco2 = max (Niceco2, threshold)
1023
1024      ! meteoritic particles are considered like dust, => rho_dust
1025      Nccnco2 = max( (pqeff(ig,l,igcm_ccnco2_number) + (pdq(ig,l,igcm_ccnco2_number) + &
1026                       pdqcloudco2(ig,l, igcm_ccnco2_number)) * ptimestep) &
1027                    , threshold)
1028
1029      Qccnco2 = max( (pqeff(ig,l,igcm_ccnco2_mass) + (pdq(ig,l,igcm_ccnco2_mass) + &
1030                       pdqcloudco2(ig,l, igcm_ccnco2_mass)) * ptimestep)&
1031                    , threshold)
1032
1033      myT = pteff(ig,l) + (pdt(ig,l)+pdtcloudco2(ig,l))*ptimestep
1034
1035      Nccnco2_h2o = 0.
1036      Qccnco2_h2o = 0.
1037      if (co2useh2o) then
1038      Nccnco2_h2o = max( (pqeff(ig,l,igcm_ccnco2_h2o_number) + (pdq(ig,l,igcm_ccnco2_h2o_number) + &
1039                       pdqcloudco2(ig,l, igcm_ccnco2_h2o_number)) * ptimestep) &
1040                    , threshold)
1041
1042      Qccnco2_h2o = max( (pqeff(ig,l,igcm_ccnco2_h2o_mass_ice) + pqeff(ig,l,igcm_ccnco2_h2o_mass_ccn) + &
1043                     (pdq(ig,l,igcm_ccnco2_h2o_mass_ice) + pdq(ig,l,igcm_ccnco2_h2o_mass_ccn) + &
1044                      pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) + pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn)) * &
1045                      ptimestep)&
1046                    , threshold)
1047      Nccnco2 = Nccnco2 - Nccnco2_h2o
1048      Qccnco2 = Qccnco2 - Qccnco2_h2o
1049      end if
1050      ! Compute particle size
1051      call updaterice_microco2(dble(Niceco2), dble(Qccnco2), dble(Nccnco2), dble(Qccnco2_h2o), dble(Nccnco2_h2o), myT, &
1052                               tauscaling(ig), riceco2(ig,l), rhocloudco2(ig,l))
1053
1054     ! Compute opacities
1055      if ( (Niceco2 <= threshold .or. (Nccnco2)*tauscaling(ig) <= 1.) ) then ! Nccnco2 inclut Nccnco2_h2o
1056        riceco2(ig,l) = 0.
1057        Qext1bins2(ig,l) = 0.
1058      else
1059        n_derf = derf( (rb_cldco2(1)-log(riceco2(ig,l))) *dev2)
1060        Qext1bins2(ig,l) = 0.
1061
1062        do i = 1, nbinco2_cld
1063          n_aer(i) = -0.5 * (Nccnco2)*tauscaling(ig) * n_derf
1064
1065          n_derf = derf((rb_cldco2(i+1)-log(riceco2(ig,l))) *dev2)
1066          n_aer(i) = n_aer(i) + (0.5 * (Nccnco2)*tauscaling(ig) * n_derf)
1067
1068          Qext1bins2(ig,l) = Qext1bins2(ig,l) + Qext1bins(i) * n_aer(i)
1069        end do
1070      end if
1071!----------------------------------------------------------------------------------------------------------------------!
1072! 6.2. Update rice and rdust
1073!----------------------------------------------------------------------------------------------------------------------!
1074      ! update rice water only if co2 use h2o ice as CCN
1075      if (co2useh2o) then
1076        call updaterice_micro( &
1077              pqeff(ig,l,igcm_h2o_ice) + (pdq(ig,l,igcm_h2o_ice) + pdqcloudco2(ig,l,igcm_h2o_ice))*ptimestep, &
1078              pqeff(ig,l,igcm_ccn_mass) + (pdq(ig,l,igcm_ccn_mass) + pdqcloudco2(ig,l,igcm_ccn_mass))*ptimestep, &
1079              pqeff(ig,l,igcm_ccn_number) + (pdq(ig,l,igcm_ccn_number) + pdqcloudco2(ig,l,igcm_ccn_number))*ptimestep, &
1080              tauscaling(ig),rice(ig,l),rhocloud(ig,l))
1081      end if
1082
1083      ! update rdust
1084      call updaterdust( &
1085           pqeff(ig,l,igcm_dust_mass) + (pdq(ig,l,igcm_dust_mass) + pdqcloudco2(ig,l,igcm_dust_mass))*ptimestep, &
1086           pqeff(ig,l,igcm_dust_number) + (pdq(ig,l,igcm_dust_number) + pdqcloudco2(ig,l,igcm_dust_number))*ptimestep, &
1087           rdust(ig,l))
1088    end do ! ngrid
1089  end do ! nlay
1090!----------------------------------------------------------------------------------------------------------------------!
1091! 7. A correction if a lot of subliming CO2 fills the 1st layer FF (04/2005). Then that should not affect the ice
1092!      particle radius
1093!----------------------------------------------------------------------------------------------------------------------!
1094  do ig = 1, ngrid
1095    if ( pdpsrf(ig)*ptimestep > 0.9*(pplev(ig,1)-pplev(ig,2))) then
1096
1097      if ( pdpsrf(ig)*ptimestep > 0.9*(pplev(ig,1)-pplev(ig,3)) ) then
1098        riceco2(ig,2) = riceco2(ig,3)
1099      end if
1100
1101      riceco2(ig,1) = riceco2(ig,2)
1102    end if
1103  end do
1104!----------------------------------------------------------------------------------------------------------------------!
1105! 9. CO2 saturated quantities
1106!----------------------------------------------------------------------------------------------------------------------!
1107! 9.1 Compute CO2 saturation in layers
1108!----------------------------------------------------------------------------------------------------------------------!
1109  call co2sat(ngrid*nlay, pteff+(pdt+pdtcloudco2)*ptimestep, zqsatco2)
1110!----------------------------------------------------------------------------------------------------------------------!
1111! 9.2 Compute CO2 saturated quantities in layers
1112!----------------------------------------------------------------------------------------------------------------------!
1113  do l = 1, nlay
1114    do ig = 1, ngrid
1115      satuco2(ig,l) = ( pqeff(ig,l,igcm_co2)  + (pdq(ig,l,igcm_co2) + pdqcloudco2(ig,l,igcm_co2))*ptimestep ) *  &
1116                      (mmean(ig,l)/(mco2*1e3)) * pplay(ig,l) / zqsatco2(ig,l)
1117    end do
1118  end do
1119!----------------------------------------------------------------------------------------------------------------------!
1120! 10. Everything modified by CO2 microphysics must be wrt co2cloudfrac
1121!----------------------------------------------------------------------------------------------------------------------!
1122  if (CLFvaryingCO2) then
1123    do l = 1, nlay
1124      do ig = 1, ngrid
1125        pdqcloudco2(ig,l,igcm_ccnco2_mass) = pdqcloudco2(ig,l,igcm_ccnco2_mass) * co2cloudfrac(ig,l)
1126
1127        pdqcloudco2(ig,l,igcm_ccnco2_number) = pdqcloudco2(ig,l,igcm_ccnco2_number) * co2cloudfrac(ig,l)
1128
1129        pdqcloudco2(ig,l,igcm_dust_mass) = pdqcloudco2(ig,l,igcm_dust_mass) * co2cloudfrac(ig,l)
1130
1131        pdqcloudco2(ig,l,igcm_dust_number) = pdqcloudco2(ig,l,igcm_dust_number) * co2cloudfrac(ig,l)
1132
1133        pdqcloudco2(ig,l,igcm_co2_ice) = pdqcloudco2(ig,l,igcm_co2_ice) * co2cloudfrac(ig,l)
1134
1135        pdqcloudco2(ig,l,igcm_co2) = pdqcloudco2(ig,l,igcm_co2) * co2cloudfrac(ig,l)
1136
1137        pdtcloudco2(ig,l) = pdtcloudco2(ig,l) * co2cloudfrac(ig,l)
1138
1139        Qext1bins2(ig,l) = Qext1bins2(ig,l) * co2cloudfrac(ig,l)
1140
1141        if (meteo_flux) then
1142          pdqcloudco2(ig,l,igcm_ccnco2_meteor_mass) = pdqcloudco2(ig,l,igcm_ccnco2_meteor_mass) * co2cloudfrac(ig,l)
1143
1144          pdqcloudco2(ig,l,igcm_ccnco2_meteor_number) = pdqcloudco2(ig,l,igcm_ccnco2_meteor_number) * co2cloudfrac(ig,l)
1145        end if
1146        if (co2useh2o) then
1147          pdqcloudco2(ig,l,igcm_h2o_ice) = pdqcloudco2(ig,l,igcm_h2o_ice) * co2cloudfrac(ig,l)
1148
1149          pdqcloudco2(ig,l,igcm_ccn_mass) = pdqcloudco2(ig,l,igcm_ccn_mass) * co2cloudfrac(ig,l)
1150
1151          pdqcloudco2(ig,l,igcm_ccn_number) = pdqcloudco2(ig,l,igcm_ccn_number) * co2cloudfrac(ig,l)
1152
1153          pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) = pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ice) * co2cloudfrac(ig,l)
1154
1155          pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn) = pdqcloudco2(ig,l,igcm_ccnco2_h2o_mass_ccn) * co2cloudfrac(ig,l)
1156
1157          pdqcloudco2(ig,l,igcm_ccnco2_h2o_number) = pdqcloudco2(ig,l,igcm_ccnco2_h2o_number) * co2cloudfrac(ig,l)
1158        end if
1159      end do ! ngrid
1160    end do ! nlay
1161  end if ! if CLFvaryingCO2 is true
1162!----------------------------------------------------------------------------------------------------------------------!
1163! 11. Compute opacity at 1 micron: Opacity in mesh ig is the sum over l of Qext1bins2. Is this true ?
1164!----------------------------------------------------------------------------------------------------------------------!
1165  tau1mic(:)=0.
1166  do l = 1, nlay
1167    do ig = 1, ngrid
1168      tau1mic(ig) = tau1mic(ig) + Qext1bins2(ig,l)
1169    end do
1170  end do
1171!----------------------------------------------------------------------------------------------------------------------!
1172! 12. Write outputs in diagfi.nc
1173!----------------------------------------------------------------------------------------------------------------------!
1174  call WRITEDIAGFI(ngrid, "satuco2", "vap in satu", " ", 3, satuco2)
1175
1176  call WRITEDIAGFI(ngrid, "precip_co2_ice_rate", "surface deposition rate of co2 ice", "kg.m-2.s-1", 2, pdqs_sedco2(:))
1177
1178  call WRITEDIAGFI(ngrid, "co2ice_cond_rate", "CO2 condensation rate in atm layers", "kg.kg-1.s-1", 3, pcondicea)
1179
1180  call WRITEDIAGFI(ngrid, "pdtcloudco2", "temperature variation of CO2 latent heat", "K.s-1", 3, pdtcloudco2)
1181
1182  call writediagfi(ngrid, "riceco2", "ice radius", "m", 3, riceco2)
1183
1184  call WRITEDIAGFI(ngrid, "Tau3D1mic", " co2 ice opacities", " ", 3, Qext1bins2)
1185
1186  call WRITEDIAGFI(ngrid, "tau1mic", "co2 ice opacity 1 micron", " ", 2, tau1mic)
1187
1188  if (CLFvaryingCO2) then
1189    call WRITEDIAGFI(ngrid, "co2cloudfrac", "co2 cloud fraction", " ", 3, co2cloudfrac)
1190  end if
1191!======================================================================================================================!
1192! END =================================================================================================================!
1193!======================================================================================================================!
1194  end subroutine co2cloud
1195
1196end module co2cloud_mod
Note: See TracBrowser for help on using the repository browser.