source: trunk/LMDZ.MARS/libf/phymars/co2condens_mod.F @ 2551

Last change on this file since 2551 was 2524, checked in by cmathe, 4 years ago

minor fix in co2condens4micro; fix co2 sedimentation rate outputs

File size: 36.3 KB
Line 
1      MODULE co2condens_mod
2
3      IMPLICIT NONE
4
5      logical, save :: scavco2cond = .false. ! flag for using scavenging_by_co2
6     
7      CONTAINS
8
9      SUBROUTINE co2condens(ngrid,nlayer,nq,ptimestep,
10     $                  pcapcal,pplay,pplev,ptsrf,pt,
11     $                  pphi,pdt,pdu,pdv,pdtsrf,pu,pv,pq,pdq,
12     $                  piceco2,psolaralb,pemisurf,rdust,
13     $                  pdtc,pdtsrfc,pdpsrf,pduc,pdvc,pdqc,
14     $                  fluxsurf_sw,zls,
15     $                  zdqssed_co2,pcondicea_co2microp,
16     $                  pdqsc)
17                                                   
18       use tracer_mod, only: noms, igcm_h2o_ice, igcm_h2o_vap,
19     &                      igcm_dust_mass, igcm_dust_number,
20     &                      igcm_ccn_mass, igcm_ccn_number,
21     &                      igcm_hdo_ice, igcm_hdo_vap,
22     &                      nqperes,nqfils, ! MVals: variables isotopes
23     &                      qperemin,masseqmin
24       use surfdat_h, only: emissiv, phisfi
25       use geometry_mod, only: latitude, ! grid point latitudes (rad)
26     &                         longitude_deg, latitude_deg
27       use planete_h, only: obliquit
28       use comcstfi_h, only: cpp, g, r, pi
29       use dust_param_mod, only: freedust
30
31#ifndef MESOSCALE
32       USE vertical_layers_mod, ONLY: ap, bp
33#endif
34       IMPLICIT NONE
35c=======================================================================
36c   subject:
37c   --------
38c   Condensation/sublimation of CO2 ice on the ground and in the
39c   atmosphere
40c   (Scheme described in Forget et al., Icarus, 1998)
41c
42c   author:   Francois Forget     1994-1996 ; updated 1996 -- 2018
43c   ------
44c             adapted to external CO2 ice clouds scheme by Deborah Bardet (2018) '
45c
46c=======================================================================
47c
48c    0.  Declarations :
49c    ------------------
50c
51      include "callkeys.h"
52
53c-----------------------------------------------------------------------
54c    Arguments :
55c    ---------
56      INTEGER,INTENT(IN) :: ngrid  ! number of atmospheric columns
57      INTEGER,INTENT(IN) :: nlayer ! number of vertical layers
58      INTEGER,INTENT(IN) :: nq  ! number of tracers
59
60      REAL,INTENT(IN) :: ptimestep ! physics timestep (s)
61      REAL,INTENT(IN) :: pcapcal(ngrid)
62      REAL,INTENT(IN) :: pplay(ngrid,nlayer) !mid-layer pressure (Pa)
63      REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa)
64      REAL,INTENT(IN) :: ptsrf(ngrid) ! surface temperature (K)
65      REAL,INTENT(IN) :: pt(ngrid,nlayer) ! atmospheric temperature (K)
66      REAL,INTENT(IN) :: pphi(ngrid,nlayer) ! geopotential (m2.s-2)
67      REAL,INTENT(IN) :: pdt(ngrid,nlayer) ! tendency on temperature from
68                                           ! previous physical processes (K/s)
69      REAL,INTENT(IN) :: pdu(ngrid,nlayer) ! tendency on zonal wind (m/s2)
70                                           ! from previous physical processes
71      REAL,INTENT(IN) :: pdv(ngrid,nlayer) ! tendency on meridional wind (m/s2)
72                                           ! from previous physical processes
73      REAL,INTENT(IN) :: pdtsrf(ngrid) ! tendency on surface temperature from
74                                       ! previous physical processes (K/s)
75      REAL,INTENT(IN) :: pu(ngrid,nlayer) ! zonal wind (m/s)
76      REAL,INTENT(IN) :: pv(ngrid,nlayer) ! meridional wind (m/s)
77      REAL,INTENT(IN) :: pq(ngrid,nlayer,nq) ! tracers (../kg_air)
78      REAL,INTENT(IN) :: pdq(ngrid,nlayer,nq) ! tendency on tracers from
79                                              ! previous physical processes
80
81      REAL,INTENT(IN) :: zdqssed_co2(ngrid)  ! CO2 flux at the surface  (kg.m-2.s-1)
82      REAL,INTENT(IN) :: pcondicea_co2microp(ngrid,nlayer)! tendency due to CO2 condensation (kg/kg.s-1)
83
84      REAL,INTENT(INOUT) :: piceco2(ngrid) ! CO2 ice on the surface (kg.m-2)
85      REAL,INTENT(INOUT) :: psolaralb(ngrid,2) ! albedo of the surface
86      REAL,INTENT(INOUT) :: pemisurf(ngrid) ! emissivity of the surface
87      REAL,INTENT(IN) :: rdust(ngrid,nlayer) ! dust effective radius
88     
89      ! tendencies due to CO2 condensation/sublimation:
90      REAL,INTENT(OUT) :: pdtc(ngrid,nlayer) ! tendency on temperature (K/s)
91      REAL,INTENT(OUT) :: pdtsrfc(ngrid) ! tendency on surface temperature (K/s)
92      REAL,INTENT(OUT) :: pdpsrf(ngrid) ! tendency on surface pressure (Pa/s)
93      REAL,INTENT(OUT) :: pduc(ngrid,nlayer) ! tendency on zonal wind (m.s-2)
94      REAL,INTENT(OUT) :: pdvc(ngrid,nlayer) ! tendency on meridional wind (m.s-2)
95      REAL,INTENT(OUT) :: pdqc(ngrid,nlayer,nq) ! tendency on tracers
96      REAL,INTENT(OUT) :: pdqsc(ngrid,nq) ! tendency on surface tracers
97     
98      ! added to calculate flux dependent albedo:
99      REAL,intent(in) :: fluxsurf_sw(ngrid,2)
100      real,intent(in) :: zls ! solar longitude (rad)
101
102c
103c    Local variables :
104c    -----------------
105
106      INTEGER i,j
107      INTEGER l,ig,iq,icap
108      REAL zt(ngrid,nlayer)
109      REAL zcpi
110      REAL ztcond (ngrid,nlayer+1) ! CO2 condensation temperature (atm)
111      REAL ztcondsol(ngrid) ! CO2 condensation temperature (surface)
112      REAL zdiceco2(ngrid)
113      REAL zcondicea(ngrid,nlayer) ! condensation rate in layer  l  (kg/m2/s)
114      REAL zcondices(ngrid) ! condensation rate on the ground (kg/m2/s)
115      REAL zfallice(ngrid,nlayer+1) ! amount of ice falling from layer l (kg/m2/s)
116      REAL zfallheat
117      REAL zmflux(nlayer+1)
118      REAL zu(nlayer),zv(nlayer)
119      REAL zqc(nlayer,nq),zq1(nlayer)
120      REAL ztsrf(ngrid)
121      REAL ztc(nlayer), ztm(nlayer+1)
122      REAL zum(nlayer+1) , zvm(nlayer+1)
123      REAL zqm(nlayer+1,nq),zqm1(nlayer+1)
124      REAL masse(nlayer),w(nlayer+1)
125      REAL Sm(nlayer),Smq(nlayer,nq),mixmas,qmix
126      REAL availco2
127      LOGICAL condsub(ngrid)
128
129      real :: emisref(ngrid)
130     
131      REAL zdq_scav(ngrid,nlayer,nq) ! tendancy due to scavenging by co2
132      REAL zq(ngrid,nlayer,nq)
133
134c variable speciale diagnostique
135      real tconda1(ngrid,nlayer)
136      real tconda2(ngrid,nlayer)
137c     REAL zdiceco2a(ngrid) ! for diagnostic only
138      real zdtsig (ngrid,nlayer)
139      real zdt (ngrid,nlayer)
140      real vmr_co2(ngrid,nlayer) ! co2 volume mixing ratio
141! improved_ztcond flag: If set to .true. (AND running with a 'co2' tracer)
142! then condensation temperature is computed using partial pressure of CO2
143      logical,parameter :: improved_ztcond=.true.
144
145c   local saved variables
146      integer,save :: ico2 ! index of CO2 tracer
147      real,save :: qco2,mmean
148      real,parameter :: latcond=5.9e5 ! (J/kg) Latent heat of solid CO2 ice
149      real,parameter :: tcond1mb=136.27 ! condensation temperature (K) at 1 mbar
150      real,parameter :: cpice=1000. ! (J.kg-1.K-1) specific heat of CO2 ice
151      REAL,SAVE :: acond,bcond,ccond
152      real,save :: m_co2, m_noco2, A , B
153
154      LOGICAL,SAVE :: firstcall = .true. !,firstcall2=.true.
155
156c D.BARDET: for debug
157      real ztc3D(ngrid,nlayer)
158      REAL ztm3D(ngrid,nlayer)
159      REAL zmflux3D(ngrid,nlayer)
160
161c MVals: variables isotopes
162      REAL Ratio1(nlayer),Ratiom1(nlayer+1)
163      REAL masseq(nlayer),wq(nlayer+1)
164      INTEGER ifils,iq2
165c----------------------------------------------------------------------
166
167c   Initialisation
168c   --------------
169c
170      ! AS: firstcall OK absolute
171      IF (firstcall) THEN
172         
173         bcond=1./tcond1mb
174         ccond=cpp/(g*latcond)
175         acond=r/latcond
176
177         firstcall=.false.
178         write(*,*) 'CO2condens: improved_ztcond=',improved_ztcond
179         PRINT*,'In co2condens:Tcond(P=1mb)=',tcond1mb,' Lcond=',latcond
180         PRINT*,'acond,bcond,ccond',acond,bcond,ccond
181
182         ico2=0
183
184         if (tracer) then
185c          Prepare Special treatment if one of the tracer is CO2 gas
186           do iq=1,nq
187             if (noms(iq).eq."co2") then
188                ico2=iq
189                m_co2 = 44.01E-3  ! CO2 molecular mass (kg/mol)   
190                m_noco2 = 33.37E-3  ! Non condensible mol mass (kg/mol)   
191c               Compute A and B coefficient use to compute
192c               mean molecular mass Mair defined by
193c               1/Mair = q(ico2)/m_co2 + (1-q(ico2))/m_noco2
194c               1/Mair = A*q(ico2) + B
195                A =(1/m_co2 - 1/m_noco2)
196                B=1/m_noco2
197             endif
198           enddo
199         end if
200      ENDIF ! of IF (firstcall)
201      zcpi=1./cpp
202
203c
204c======================================================================
205c    Calcul of CO2 condensation sublimation
206c    ============================================================
207
208c    Used variable :
209c       piceco2(ngrid)   :  amount of co2 ice on the ground (kg/m2)
210c       zcondicea(ngrid,l):  condensation rate in layer  l  (kg/m2/s)
211c       zcondices(ngrid):  condensation rate on the ground (kg/m2/s)
212c       zfallice(ngrid,l):amount of ice falling from layer l (kg/m2/s)
213c                           
214c       pdtc(ngrid,nlayer)  : dT/dt due to cond/sub
215c
216c
217c     Tendencies set to 0
218c     -------------------------------------
219      zcondicea(1:ngrid,1:nlayer) = 0.
220      zfallice(1:ngrid,1:nlayer+1) = 0.
221      pduc(1:ngrid,1:nlayer)  = 0
222      pdvc(1:ngrid,1:nlayer)  = 0
223      pdtc(1:ngrid,1:nlayer) = 0.
224      pdqsc(1:ngrid,1:nq) = 0
225         
226      pdqc(1:ngrid,1:nlayer,1:nq)  = 0
227
228      zcondices(1:ngrid) = 0.
229      pdtsrfc(1:ngrid) = 0.
230      pdpsrf(1:ngrid) = 0.
231      condsub(1:ngrid) = .false.
232      zdiceco2(1:ngrid) = 0.
233
234      zfallheat=0
235     
236      zdq_scav(:,:,:)=0.
237
238c     Update tendencies from previous processes
239c     -------------------------------------
240      DO l=1,nlayer
241         DO ig=1,ngrid
242            zt(ig,l)=pt(ig,l)+ pdt(ig,l)*ptimestep
243            do iq=1,nq
244             zq(ig,l,iq)=pq(ig,l,iq)+pdq(ig,l,iq)*ptimestep
245            enddo
246         ENDDO
247      ENDDO
248     
249c     *************************
250c     ATMOSPHERIC CONDENSATION
251c     *************************
252
253c     Compute CO2 Volume mixing ratio
254c     -------------------------------
255      if (improved_ztcond.and.(ico2.ne.0)) then
256         DO l=1,nlayer
257            DO ig=1,ngrid
258              qco2=pq(ig,l,ico2)+pdq(ig,l,ico2)*ptimestep
259c             Mean air molecular mass = 1/(q(ico2)/m_co2 + (1-q(ico2))/m_noco2)
260              mmean=1/(A*qco2 +B)
261              vmr_co2(ig,l) = qco2*mmean/m_co2
262            ENDDO
263         ENDDO
264      else
265         DO l=1,nlayer
266            DO ig=1,ngrid
267              vmr_co2(ig,l)=0.95
268            ENDDO
269         ENDDO
270      endif
271
272      IF (.NOT. co2clouds) then
273c     forecast of atmospheric temperature zt and frost temperature ztcond
274c     --------------------------------------------------------------------
275
276      DO l=1,nlayer
277         DO ig=1,ngrid
278!            ztcond(ig,l)=1./(bcond-acond*log(.0095*pplay(ig,l)))
279            if (pplay(ig,l).ge.1e-4) then
280              ztcond(ig,l)=
281     &         1./(bcond-acond*log(.01*vmr_co2(ig,l)*pplay(ig,l)))
282            else
283              ztcond(ig,l)=0.0 !mars Monica
284            endif
285         ENDDO
286      ENDDO
287
288      ztcond(:,nlayer+1)=ztcond(:,nlayer)
289 
290c      Condensation/sublimation in the atmosphere
291c      ------------------------------------------
292c      (calcul of zcondicea , zfallice and pdtc)
293c
294      DO l=nlayer , 1, -1
295         DO ig=1,ngrid
296           pdtc(ig,l)=0.
297           IF((zt(ig,l).LT.ztcond(ig,l)).or.(zfallice(ig,l+1).gt.0))THEN
298               condsub(ig)=.true.
299               IF (zfallice(ig,l+1).gt.0) then 
300                 zfallheat=zfallice(ig,l+1)*
301     &           (pphi(ig,l+1)-pphi(ig,l) +
302     &          cpice*(ztcond(ig,l+1)-ztcond(ig,l)))/latcond
303               ELSE
304                   zfallheat=0.
305               ENDIF
306               pdtc(ig,l)=(ztcond(ig,l) - zt(ig,l))/ptimestep
307               zcondicea(ig,l)=(pplev(ig,l)-pplev(ig,l+1))
308     &                        *ccond*pdtc(ig,l)- zfallheat
309c              Case when the ice from above sublimes entirely
310c              """""""""""""""""""""""""""""""""""""""""""""""
311               IF (zfallice(ig,l+1).lt.- zcondicea(ig,l)) then
312                  pdtc(ig,l)=(-zfallice(ig,l+1)+zfallheat)/
313     &              (ccond*(pplev(ig,l)-pplev(ig,l+1)))
314                  zcondicea(ig,l)= -zfallice(ig,l+1)
315               END IF
316
317               zfallice(ig,l) = zcondicea(ig,l)+zfallice(ig,l+1)
318            END IF
319         ENDDO
320      ENDDO
321     
322      if (scavco2cond) then
323        call scavenging_by_co2(ngrid,nlayer,nq,ptimestep,pplev,zq,
324     &                       rdust,zcondicea,zfallice,zdq_scav,pdqsc)
325      endif
326     
327      ELSE ! if co2 clouds
328        DO ig=1,ngrid
329         zfallice(ig,1) = zdqssed_co2(ig)
330        ENDDO
331        DO l=nlayer , 1, -1
332            DO ig=1,ngrid
333         zcondicea(ig,l) = pcondicea_co2microp(ig,l)*
334     &                        (pplev(ig,l) - pplev(ig,l+1))/g
335            ENDDO
336        ENDDO
337 
338      ENDIF ! end of if co2clouds
339
340      call WRITEdiagfi(ngrid,"pdtc_atm",
341     &         "temperature tendency due to CO2 condensation",
342     &         " ",3,pdtc)
343     
344       call WRITEdiagfi(ngrid,"zcondicea",
345     &         "",
346     &         " ",3,zcondicea)         
347     
348       call WRITEdiagfi(ngrid,"zfallice",
349     &         "",
350     &         " ",2,zfallice(ngrid,1))
351     
352
353c     *************************
354c       SURFACE  CONDENSATION
355c     *************************
356
357c     forecast of ground temperature ztsrf and frost temperature ztcondsol
358c     --------------------------------------------------------------------
359      DO ig=1,ngrid
360         ztcondsol(ig)=
361     &          1./(bcond-acond*log(.01*vmr_co2(ig,1)*pplev(ig,1)))
362         ztsrf(ig) = ptsrf(ig) + pdtsrf(ig)*ptimestep
363      ENDDO
364
365c
366c      Condensation/sublimation on the ground
367c      --------------------------------------
368c      (compute zcondices and pdtsrfc)
369c
370      DO ig=1,ngrid
371         IF(latitude(ig).lt.0) THEN
372           ! Southern hemisphere
373            icap=2
374         ELSE
375           ! Northern hemisphere
376            icap=1
377         ENDIF
378       
379c
380c        Loop on where we have  condensation/ sublimation
381         IF ((ztsrf(ig) .LT. ztcondsol(ig)) .OR.   ! ground cond
382     $       (zfallice(ig,1).NE.0.) .OR.           ! falling snow
383     $      ((ztsrf(ig) .GT. ztcondsol(ig)) .AND.  ! ground sublim.
384     $      ((piceco2(ig)+zfallice(ig,1)*ptimestep) .NE. 0.))) THEN
385            condsub(ig) = .true.
386
387            IF (zfallice(ig,1).gt.0) then 
388                 zfallheat=zfallice(ig,1)*
389     &           (pphi(ig,1)- phisfi(ig) +
390     &           cpice*(ztcond(ig,1)-ztcondsol(ig)))/latcond
391            ELSE
392                 zfallheat=0.
393            ENDIF
394
395c           condensation or partial sublimation of CO2 ice
396c           """""""""""""""""""""""""""""""""""""""""""""""
397            zcondices(ig)=pcapcal(ig)*(ztcondsol(ig)-ztsrf(ig))
398     &      /(latcond*ptimestep)         - zfallheat
399            pdtsrfc(ig) = (ztcondsol(ig) - ztsrf(ig))/ptimestep
400            zdiceco2(ig) = zcondices(ig) + zfallice(ig,1)
401
402#ifdef MESOSCALE
403      print*, "not enough CO2 tracer in 1st layer to condense"
404      print*, ">>> to be implemented in the mesoscale case"
405      print*, "because this uses ap levels..."
406#else
407c           If there is not enough CO2 tracer in 1st layer to condense
408c           """"""""""""""""""""""""""""""""""""""""""""""""""""""
409            IF(ico2.ne.0) then
410c              Available CO2 tracer in layer 1 at end of timestep (kg/m2)
411                availco2= pq(ig,1,ico2)*((ap(1)-ap(2))+
412     &          (bp(1)-bp(2))*(pplev(ig,1)/g-zdiceco2(ig)*ptimestep))
413
414               IF ((zcondices(ig) + zcondicea(ig,1))*ptimestep
415     &           .gt.availco2) then
416                   zcondices(ig) = availco2/ptimestep -zcondicea(ig,1)
417                   zdiceco2(ig) = zcondices(ig) + zfallice(ig,1)
418                   pdtsrfc(ig)=(latcond/pcapcal(ig))*
419     &                          (zcondices(ig)+zfallheat)
420               ENDIF
421            ENDIF
422#endif
423
424c           If the entire CO2 ice layer sublimes
425c           """"""""""""""""""""""""""""""""""""""""""""""""""""
426c           (including what has just condensed in the atmosphere)
427
428            IF((piceco2(ig)/ptimestep+zfallice(ig,1)).LE.
429     &          -zcondices(ig))THEN
430               zcondices(ig) = -piceco2(ig)/ptimestep - zfallice(ig,1)
431               pdtsrfc(ig)=(latcond/pcapcal(ig))*
432     &         (zcondices(ig)+zfallheat)
433               zdiceco2(ig) = zcondices(ig) + zfallice(ig,1)
434            END IF
435
436c           Changing CO2 ice amount and pressure :
437c           """"""""""""""""""""""""""""""""""""
438
439            piceco2(ig) = piceco2(ig) + zdiceco2(ig)*ptimestep
440            pdpsrf(ig) = -zdiceco2(ig)*g
441
442            IF(ABS(pdpsrf(ig)*ptimestep).GT.pplev(ig,1)) THEN
443               PRINT*,'STOP in condens'
444               PRINT*,'condensing more than total mass'
445               PRINT*,'Grid point ',ig
446               PRINT*,'Longitude(degrees): ',longitude_deg(ig)
447               PRINT*,'Latitude(degrees): ',latitude_deg(ig)
448               PRINT*,'Ps = ',pplev(ig,1)
449               PRINT*,'d Ps = ',pdpsrf(ig)
450               CALL abort_physic('co2condens',
451     &              'condensing more than total mass', 1)
452            ENDIF
453         END IF ! if there is condensation/sublimation
454      ENDDO ! of DO ig=1,ngrid
455
456c  ********************************************************************
457c  Surface albedo and emissivity of the surface below the snow (emisref)
458c  ********************************************************************
459
460!     Check that amont of CO2 ice is not problematic
461      DO ig=1,ngrid
462           if(.not.piceco2(ig).ge.0.) THEN
463              if(piceco2(ig).le.-5.e-8) print*,
464     $         'WARNING co2condens piceco2(',ig,')=', piceco2(ig)
465               piceco2(ig)=0.
466           endif
467      ENDDO
468     
469!     Set albedo and emissivity of the surface
470!     ----------------------------------------
471      CALL albedocaps(zls,ngrid,piceco2,psolaralb,emisref)
472
473! set pemisurf() to emissiv when there is bare surface (needed for co2snow)
474      DO ig=1,ngrid
475        if (piceco2(ig).eq.0) then
476          pemisurf(ig)=emissiv
477        endif
478      ENDDO
479
480!         firstcall2=.false.
481c ***************************************************************
482c  Correction to account for redistribution between sigma or hybrid
483c  layers when changing surface pressure (and warming/cooling
484c  of the CO2 currently changing phase).
485c *************************************************************
486
487      DO ig=1,ngrid
488        if (condsub(ig)) then
489           do l=1,nlayer
490             ztc(l)  =zt(ig,l)   +pdtc(ig,l)  *ptimestep
491             zu(l)   =pu(ig,l)   +pdu( ig,l)  *ptimestep
492             zv(l)   =pv(ig,l)   +pdv( ig,l)  *ptimestep
493            do iq=1,nq
494             zqc(l,iq)=zq(ig,l,iq)+zdq_scav(ig,l,iq)*ptimestep ! zdq_scav=0 if watercloud=false
495            enddo
496           enddo
497
498c  Mass fluxes through the sigma levels (kg.m-2.s-1)  (>0 when up)
499c  """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
500
501            zmflux(1) = -zcondices(ig)
502            DO l=1,nlayer
503             zmflux(l+1) = zmflux(l) -zcondicea(ig,l)
504#ifndef MESOSCALE
505     &        + (bp(l)-bp(l+1))*(zfallice(ig,1)-zmflux(1))
506c zmflux set to 0 if very low to avoid: top layer is disappearing in v1ld 
507          if (abs(zmflux(l+1)).lt.1E-13.OR.bp(l+1).eq.0.) zmflux(l+1)=0.
508#else
509          if (abs(zmflux(l+1)).lt.1E-13) zmflux(l+1)=0.
510#endif
511            END DO
512
513#ifdef MESOSCALE
514         print*,"absurd mass set because bp not available"
515         print*,"TO BE FIXED"
516#else
517c Mass of each layer at the end of timestep
518c -----------------------------------------
519            DO l=1,nlayer
520              masse(l)=( pplev(ig,l) - pplev(ig,l+1) + 
521     &                 (bp(l)-bp(l+1))*pdpsrf(ig)*ptimestep)/g
522            END DO
523#endif
524
525c  Corresponding fluxes for T,U,V,Q
526c  """"""""""""""""""""""""""""""""
527
528c           averaging operator for TRANSPORT 
529c           """"""""""""""""""""""""""""""""
530c           Value transfert at the surface interface when condensation
531c           sublimation:
532            ztm(1) = ztsrf(ig) + pdtsrfc(ig)*ptimestep
533            zum(1) = 0 
534            zvm(1) = 0 
535            do iq=1,nq
536              zqm(1,iq)=0. ! most tracer do not condense !
537            enddo
538c           Special case if one of the tracer is CO2 gas
539            if (ico2.ne.0) zqm(1,ico2)=1. ! flux is 100% CO2
540
541c           Van Leer scheme:
542            DO l=1,nlayer+1
543                w(l)=-zmflux(l)*ptimestep
544            END DO
545            call vl1d(nlayer,ztc,2.,masse,w,ztm)
546            call vl1d(nlayer,zu ,2.,masse,w,zum)
547            call vl1d(nlayer,zv ,2.,masse,w,zvm)
548            ! MVals: loop over the fathers ("peres")
549            do iq=1,nqperes
550             do l=1,nlayer
551              zq1(l)=zqc(l,iq)
552             enddo
553             zqm1(1)=zqm(1,iq)
554             call vl1d(nlayer,zq1,2.,masse,w,zqm1)
555             do l=2,nlayer
556              zqc(l,iq)=zq1(l)
557              zqm(l,iq)=zqm1(l)
558             enddo
559             ! MVals: loop over the sons ("fils")
560             if (nqfils(iq).gt.0) then
561              if (iq.eq.igcm_h2o_ice) then
562                 iq2=igcm_hdo_ice
563              else if (iq.eq.igcm_h2o_vap) then
564                 iq2=igcm_hdo_vap
565              else
566                 call abort_physic("co2condens_mod","invalid isotope",1)
567              endif
568              do l=1,nlayer
569               if (zqc(l,iq).gt.qperemin) then
570                 Ratio1(l)=zqc(l,iq2)/zqc(l,iq)
571               else
572                 Ratio1(l)=0.
573               endif
574               masseq(l)=max(masse(l)*zqc(l,iq),masseqmin)
575               wq(l)=w(l)*zqm(l,iq)
576              enddo
577              Ratiom1(1)=zqm(1,iq2)
578              call vl1d(nlayer,Ratio1,2.,masseq,wq,Ratiom1)
579              zqm(1,iq2)=Ratiom1(1)*zqc(1,iq)
580              do l=2,nlayer
581               zqm(l,iq2)=Ratiom1(l)*zqm(l,iq)
582              enddo
583             endif !if (nqfils(iq).gt.0) then
584            enddo !iq=1,nqperes
585
586c           Surface condensation affects low winds
587            if (zmflux(1).lt.0) then
588                zum(1)= zu(1) *  (w(1)/masse(1))
589                zvm(1)= zv(1) *  (w(1)/masse(1))
590                if (w(1).gt.masse(1)) then ! ensure numerical stability
591                  zum(1)= (zu(1)-zum(2))*masse(1)/w(1) + zum(2)
592                  zvm(1)= (zv(1)-zvm(2))*masse(1)/w(1) + zvm(2)
593                end if
594            end if
595                   
596            ztm(nlayer+1)= ztc(nlayer) ! should not be used, but...
597            zum(nlayer+1)= zu(nlayer)  ! should not be used, but...
598            zvm(nlayer+1)= zv(nlayer)  ! should not be used, but...
599            do iq=1,nq
600             zqm(nlayer+1,iq)= zqc(nlayer,iq)
601            enddo
602
603#ifdef MESOSCALE
604!!!! AS: This part must be commented in the mesoscale model
605!!!! AS: ... to avoid instabilities.
606!!!! AS: you have to compile with -DMESOSCALE to do so
607#else 
608c           Tendencies on T, U, V, Q
609c           """"""""""""""""""""""""
610            DO l=1,nlayer
611               IF(.not. co2clouds) THEN
612c             Tendencies on T
613                zdtsig(ig,l) = (1/masse(l)) *
614     &        ( zmflux(l)*(ztm(l) - ztc(l))
615     &        - zmflux(l+1)*(ztm(l+1) - ztc(l))
616     &        + zcondicea(ig,l)*(ztcond(ig,l)-ztc(l))  )
617               ELSE
618                zdtsig(ig,l) = (1/masse(l)) *
619     &        ( zmflux(l)*(ztm(l) - ztc(l))
620     &        - zmflux(l+1)*(ztm(l+1) - ztc(l)))               
621               ENDIF
622c D.BARDET: for diagnotics
623                zmflux3D(ig,l)=zmflux(l)
624                ztm3D(ig,l)=ztm(l)
625                ztc3D(ig,l)=ztc(l)
626               
627                pdtc(ig,l) =  pdtc(ig,l) + zdtsig(ig,l)
628
629c             Tendencies on U
630                pduc(ig,l)   = (1/masse(l)) *
631     &        ( zmflux(l)*(zum(l) - zu(l))
632     &        - zmflux(l+1)*(zum(l+1) - zu(l)) )
633
634
635c             Tendencies on V
636                pdvc(ig,l)   = (1/masse(l)) *
637     &        ( zmflux(l)*(zvm(l) - zv(l))
638     &        - zmflux(l+1)*(zvm(l+1) - zv(l)) )
639
640            END DO
641
642#endif
643
644c           Tendencies on Q
645            do iq=1,nq
646!              if (noms(iq).eq.'co2') then
647              if (iq.eq.ico2) then
648c               SPECIAL Case when the tracer IS CO2 :
649                DO l=1,nlayer
650                  pdqc(ig,l,iq)= (1/masse(l)) *
651     &           ( zmflux(l)*(zqm(l,iq) - zqc(l,iq))
652     &           - zmflux(l+1)*(zqm(l+1,iq) - zqc(l,iq))
653     &           + zcondicea(ig,l)*(zqc(l,iq)-1.) )
654                END DO
655              else
656                DO l=1,nlayer
657                  pdqc(ig,l,iq)= (1/masse(l)) *
658     &           ( zmflux(l)*(zqm(l,iq) - zqc(l,iq))
659     &           - zmflux(l+1)*(zqm(l+1,iq) - zqc(l,iq))
660     &           + zcondicea(ig,l)*zqc(l,iq) )
661     
662                  pdqc(ig,l,iq)=pdqc(ig,l,iq)+zdq_scav(ig,l,iq) ! zdq_scav=0 if watercloud=false
663                END DO
664              end if
665            enddo
666
667       end if ! if (condsub)
668      END DO  ! loop on ig
669
670c ***************************************************************
671c   CO2 snow / clouds scheme
672c ***************************************************************
673
674      call co2snow(ngrid,nlayer,ptimestep,emisref,condsub,pplev,
675     &        zcondicea,zcondices,zfallice,pemisurf)
676
677c ***************************************************************
678c Ecriture des diagnostiques
679c ***************************************************************
680
681c     DO l=1,nlayer
682c        DO ig=1,ngrid
683c         Taux de cond en kg.m-2.pa-1.s-1
684c          tconda1(ig,l)=zcondicea(ig,l)/(pplev(ig,l)-pplev(ig,l+1))
685c          Taux de cond en kg.m-3.s-1
686c          tconda2(ig,l)=tconda1(ig,l)*pplay(ig,l)*g/(r*pt(ig,l))
687c        END DO
688c     END DO
689c     call WRITEDIAGFI(ngrid,'tconda1',
690c    &'Taux de condensation CO2 atmospherique /Pa',
691c    & 'kg.m-2.Pa-1.s-1',3,tconda1)
692c     call WRITEDIAGFI(ngrid,'tconda2',
693c    &'Taux de condensation CO2 atmospherique /m',
694c    & 'kg.m-3.s-1',3,tconda2)
695
696! output falling co2 ice in 1st layer:
697!      call WRITEDIAGFI(ngrid,'fallice',
698!     &'Precipitation of co2 ice',
699!     & 'kg.m-2.s-1',2,zfallice(1,1))
700
701#ifndef MESOSCALE
702! Extra special case for surface temperature tendency pdtsrfc:
703! we want to fix the south pole temperature to CO2 condensation temperature
704         if (caps.and.(obliquit.lt.27.)) then
705           ! check if last grid point is the south pole
706           if (abs(latitude(ngrid)-(-pi/2.)).lt.1.e-5) then
707           ! NB: Updated surface pressure, at grid point 'ngrid', is
708           !     ps(ngrid)=pplev(ngrid,1)+pdpsrf(ngrid)*ptimestep
709!             write(*,*) "co2condens: South pole: latitude(ngrid)=",
710!     &                                           latitude(ngrid)
711             ztcondsol(ngrid)=
712     &          1./(bcond-acond*log(.01*vmr_co2(ngrid,1)*
713     &                    (pplev(ngrid,1)+pdpsrf(ngrid)*ptimestep)))
714             pdtsrfc(ngrid)=(ztcondsol(ngrid)-ztsrf(ngrid))/ptimestep
715           endif
716         endif
717#endif
718
719      END SUBROUTINE co2condens
720
721
722
723c *****************************************************************
724      SUBROUTINE vl1d(nlayer,q,pente_max,masse,w,qm)
725c
726c   
727c     Operateur de moyenne inter-couche pour calcul de transport type
728c     Van-Leer " pseudo amont " dans la verticale
729c    q,w sont des arguments d'entree  pour le s-pg ....
730c    masse : masse de la couche Dp/g
731c    w : masse d'atm ``transferee'' a chaque pas de temps (kg.m-2)
732c    pente_max = 2 conseillee
733c
734c
735c   --------------------------------------------------------------------
736      IMPLICIT NONE
737
738c
739c
740c
741c   Arguments:
742c   ----------
743      integer nlayer
744      real masse(nlayer),pente_max
745      REAL q(nlayer),qm(nlayer+1)
746      REAL w(nlayer+1)
747c
748c      Local
749c   ---------
750c
751      INTEGER l
752c
753      real dzq(nlayer),dzqw(nlayer),adzqw(nlayer),dzqmax
754      real sigw, Mtot, MQtot
755      integer m
756c     integer ismax,ismin
757
758
759c    On oriente tout dans le sens de la pression
760c     W > 0 WHEN DOWN !!!!!!!!!!!!!
761
762      do l=2,nlayer
763            dzqw(l)=q(l-1)-q(l)
764            adzqw(l)=abs(dzqw(l))
765      enddo
766
767      do l=2,nlayer-1
768            if(dzqw(l)*dzqw(l+1).gt.0.) then
769                dzq(l)=0.5*(dzqw(l)+dzqw(l+1))
770            else
771                dzq(l)=0.
772            endif
773            dzqmax=pente_max*min(adzqw(l),adzqw(l+1))
774            dzq(l)=sign(min(abs(dzq(l)),dzqmax),dzq(l))
775      enddo
776
777         dzq(1)=0.
778         dzq(nlayer)=0.
779
780       do l = 1,nlayer-1
781
782c         Regular scheme (transfered mass < layer mass)
783c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
784          if(w(l+1).gt.0. .and. w(l+1).le.masse(l+1)) then
785             sigw=w(l+1)/masse(l+1)
786             qm(l+1)=(q(l+1)+0.5*(1.-sigw)*dzq(l+1))
787          else if(w(l+1).le.0. .and. -w(l+1).le.masse(l)) then
788             sigw=w(l+1)/masse(l)
789             qm(l+1)=(q(l)-0.5*(1.+sigw)*dzq(l))
790
791c         Extended scheme (transfered mass > layer mass)
792c         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
793          else if(w(l+1).gt.0.) then
794             m=l+1
795             Mtot = masse(m)
796             MQtot = masse(m)*q(m)
797             do while ((m.lt.nlayer).and.(w(l+1).gt.(Mtot+masse(m+1))))
798                m=m+1
799                Mtot = Mtot + masse(m)
800                MQtot = MQtot + masse(m)*q(m)
801             end do
802             if (m.lt.nlayer) then
803                sigw=(w(l+1)-Mtot)/masse(m+1)
804                qm(l+1)= (1/w(l+1))*(MQtot + (w(l+1)-Mtot)*
805     &          (q(m+1)+0.5*(1.-sigw)*dzq(m+1)) )
806             else
807                w(l+1) = Mtot
808                qm(l+1) = Mqtot / Mtot
809                CALL abort_physic('co2condens',
810     &               'top layer is disapearing !', 1)
811             end if
812          else      ! if(w(l+1).lt.0)
813             m = l-1
814             Mtot = masse(m+1)
815             MQtot = masse(m+1)*q(m+1)
816             if (m.gt.0) then ! because some compilers will have problems
817                              ! evaluating masse(0)
818              do while ((m.gt.0).and.(-w(l+1).gt.(Mtot+masse(m))))
819                m=m-1
820                Mtot = Mtot + masse(m+1)
821                MQtot = MQtot + masse(m+1)*q(m+1)
822                if (m.eq.0) exit
823              end do
824             endif
825             if (m.gt.0) then
826                sigw=(w(l+1)+Mtot)/masse(m)
827                qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*
828     &          (q(m)-0.5*(1.+sigw)*dzq(m))  )
829             else
830                qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*qm(1))
831             end if   
832          end if
833       enddo
834
835c     boundary conditions (not used in co2condens !!)
836c         qm(nlayer+1)=0.
837c         if(w(1).gt.0.) then
838c            qm(1)=q(1)
839c         else
840c           qm(1)=0.
841c         end if
842
843      END SUBROUTINE vl1d
844         
845c *****************************************************************
846      SUBROUTINE scavenging_by_co2(ngrid,nlayer,nq,ptimestep,pplev,pq,
847     &                          rdust,pcondicea,pfallice,pdq_scav,pdqsc)
848                     
849c
850c   
851c     Calcul de la quantite de poussiere lessivee par les nuages de CO2
852c     
853c   --------------------------------------------------------------------
854      use tracer_mod, only: nqmx, igcm_h2o_vap, igcm_h2o_ice,
855     &                      igcm_dust_mass, igcm_dust_number,
856     &                      igcm_ccn_mass, igcm_ccn_number,
857     &                      rho_dust, nuice_sed, nuice_ref,r3n_q
858      use comcstfi_h, only: g
859      use dust_param_mod, only: freedust
860      IMPLICIT NONE
861      include "callkeys.h" ! for the flags water and microphys
862c
863c
864c   Arguments:
865      INTEGER,INTENT(IN) :: ngrid  ! number of atmospheric columns
866      INTEGER,INTENT(IN) :: nlayer ! number of vertical layers
867      INTEGER,INTENT(IN) :: nq  ! number of tracers
868      REAL,INTENT(IN) :: ptimestep ! physics timestep (s)
869      REAL,INTENT(IN) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa)
870      REAL,INTENT(IN) :: pq(ngrid,nlayer,nq)
871      REAL,INTENT(IN) :: rdust(ngrid,nlayer) ! dust effective radius
872      REAL,INTENT(IN) :: pcondicea(ngrid,nlayer) ! condensation rate in layer  l  (kg/m2/s)
873      REAL,INTENT(IN) :: pfallice(ngrid,nlayer+1) ! amount of ice falling from layer l (kg/m2/s)
874     
875      REAL,INTENT(OUT) :: pdq_scav(ngrid,nlayer,nq) ! tendancy due to scavenging by co2
876      REAL,INTENT(OUT) :: pdqsc(ngrid,nq) ! tendency on surface tracers
877     
878c   Locals:
879      INTEGER l,ig
880      REAL scav_ratio_dust, scav_ratio_wice ! ratio of the dust/water ice mass mixing ratios in condensing CO2 ice and in air
881      REAL scav_dust_mass(nlayer+1) ! dust flux (mass) scavenged towards the lower layer (kg/m2/s) (POSITIVE WHEN DOWNWARD)
882      REAL scav_dust_number(nlayer+1) ! dust flux (number) scavenged towards the lower layer (kg/m2/s) (POSITIVE WHEN DOWNWARD)
883      REAL scav_ccn_mass(nlayer+1) ! ccn flux (mass) scavenged towards the lower layer
884      REAL scav_ccn_number(nlayer+1) ! ccn flux (number) scavenged towards the lower layer
885      REAL scav_h2o_ice(nlayer+1) ! water ice flux (mass) scavenged towards the lower layer
886      REAL massl ! mass of the layer l at point ig (kg/m2)
887     
888c   Initialization:
889      scav_ratio_dust = 100 !1 !10 !100 !1000
890      scav_ratio_wice = scav_ratio_dust
891      pdq_scav(:,:,:)=0.
892     
893      DO ig=1,ngrid
894        scav_dust_mass(nlayer+1)=0.
895        scav_dust_number(nlayer+1)=0.
896        scav_ccn_mass(nlayer+1)=0.
897        scav_ccn_number(nlayer+1)=0.
898        scav_h2o_ice(nlayer+1)=0.
899       
900        DO l=nlayer , 1, -1
901         massl=(pplev(ig,l)-pplev(ig,l+1))/g
902         IF(pcondicea(ig,l).GT.0.)THEN ! if CO2 condenses and traps dust/water ice
903           ! Calculation of the tendencies
904           if (freedust) then
905             pdq_scav(ig,l,igcm_dust_mass)=-pq(ig,l,igcm_dust_mass)
906     &                                     /ptimestep*(1-exp(
907     &              -scav_ratio_dust*pcondicea(ig,l)*ptimestep/massl))
908             
909             pdq_scav(ig,l,igcm_dust_number)=pdq_scav(ig,l,igcm_dust_mass)
910     &                                    *r3n_q/rdust(ig,l)
911           endif
912           if (freedust.AND.microphys) then
913             pdq_scav(ig,l,igcm_ccn_mass)=-pq(ig,l,igcm_ccn_mass)
914     &                                     /ptimestep*(1-exp(
915     &              -scav_ratio_wice*pcondicea(ig,l)*ptimestep/massl))
916             pdq_scav(ig,l,igcm_ccn_number)=pdq_scav(ig,l,igcm_ccn_mass)
917     &                                    *r3n_q/rdust(ig,l)
918           endif
919           if (water) then
920             pdq_scav(ig,l,igcm_h2o_ice)=-pq(ig,l,igcm_h2o_ice)
921     &                                     /ptimestep*(1-exp(
922     &              -scav_ratio_wice*pcondicea(ig,l)*ptimestep/massl))
923           endif
924     
925         ELSE IF(pcondicea(ig,l).LT.0.)THEN ! if CO2 sublimates and releases dust/water ice
926           ! Calculation of the tendencies
927           if (freedust) then
928             pdq_scav(ig,l,igcm_dust_mass)=-pcondicea(ig,l)/massl*
929     &                              scav_dust_mass(l+1)/pfallice(ig,l+1)
930           
931             pdq_scav(ig,l,igcm_dust_number)=-pcondicea(ig,l)/massl*
932     &                            scav_dust_number(l+1)/pfallice(ig,l+1)
933           endif
934           if (freedust.AND.microphys) then
935             pdq_scav(ig,l,igcm_ccn_mass)=-pcondicea(ig,l)/massl*
936     &                              scav_ccn_mass(l+1)/pfallice(ig,l+1)
937           
938             pdq_scav(ig,l,igcm_ccn_number)=-pcondicea(ig,l)/massl*
939     &                            scav_ccn_number(l+1)/pfallice(ig,l+1)
940           endif
941           if (water) then
942             pdq_scav(ig,l,igcm_h2o_ice)=-pcondicea(ig,l)/massl*
943     &                              scav_h2o_ice(l+1)/pfallice(ig,l+1)
944           endif
945 
946         END IF
947         ! Calculation of the scavenged dust/wice flux towards the lower layers
948         if (freedust) then
949           scav_dust_mass(l)=-pdq_scav(ig,l,igcm_dust_mass)*massl
950     &                     +scav_dust_mass(l+1)
951         
952           scav_dust_number(l)=-pdq_scav(ig,l,igcm_dust_number)*massl
953     &                     +scav_dust_number(l+1)
954         endif
955         if (freedust.AND.microphys) then
956           scav_ccn_mass(l)=-pdq_scav(ig,l,igcm_ccn_mass)*massl
957     &                     +scav_ccn_mass(l+1)
958         
959           scav_ccn_number(l)=-pdq_scav(ig,l,igcm_ccn_number)*massl
960     &                     +scav_dust_number(l+1)
961         endif
962         if (water) then
963           scav_h2o_ice(l)=-pdq_scav(ig,l,igcm_h2o_ice)*massl
964     &                     +scav_h2o_ice(l+1)
965         endif
966         
967       ENDDO
968       ! Calculation of the surface tendencies
969       pdqsc(ig,igcm_dust_mass)=0.
970       pdqsc(ig,igcm_dust_number)=0.
971       
972       if (freedust) then
973         pdqsc(ig,igcm_dust_mass)=pdqsc(ig,igcm_dust_mass)
974     &                           +scav_dust_mass(1)
975         pdqsc(ig,igcm_dust_number)=pdqsc(ig,igcm_dust_number)
976     &                             +scav_dust_number(1)
977       endif
978       if (freedust.AND.microphys) then
979         pdqsc(ig,igcm_dust_mass)=pdqsc(ig,igcm_dust_mass)
980     &                           +scav_ccn_mass(1)
981         pdqsc(ig,igcm_dust_number)=pdqsc(ig,igcm_dust_number)
982     &                             +scav_ccn_number(1)
983       endif
984       if (water) then
985         pdqsc(ig,igcm_h2o_ice)=scav_h2o_ice(1)
986       endif
987      ENDDO
988     
989      END SUBROUTINE scavenging_by_co2
990     
991      END MODULE co2condens_mod
Note: See TracBrowser for help on using the repository browser.