source: trunk/LMDZ.PLUTO/libf/phypluto/condense_n2.F90 @ 3438

Last change on this file since 3438 was 3438, checked in by afalco, 2 months ago

Pluto PCM: Minor fixes.
AF

  • Property svn:executable set to *
File size: 29.0 KB
Line 
1subroutine condense_n2(klon,klev,nq,ptimestep, &
2      pcapcal,pplay,pplev,ptsrf,pt, &
3      pphi,pdt,pdu,pdv,pdtsrf,pu,pv,pq,pdq, &
4      picen2,psolaralb,pemisurf, &
5      pdtc,pdtsrfc,pdpsrf,pduc,pdvc, &
6      pdqc,pdicen2)
7
8  use radinc_h, only : naerkind
9  use comgeomfi_h
10  use comcstfi_mod, only: g, r, cpp, pi
11  USE surfdat_h, only: phisfi,kp,p00
12  USE tracer_h, only: noms, igcm_n2, lw_n2
13  USE callkeys_mod, only: fast,ch4lag,latlag,nbsub,no_n2frost,tsurfmax,kmixmin,source_haze,vmrlag
14  USE comvert_mod, ONLY: ap,bp
15  use geometry_mod, only: latitude
16
17
18  implicit none
19
20!==================================================================
21!     Purpose
22!     -------
23!     Condense and/or sublime N2 ice on the ground and in the
24!     atmosphere, and sediment the ice.
25!
26!     Inputs
27!     ------
28!     klon                 Number of vertical columns
29!     klev                Number of layers
30!     pplay(klon,klev)   Pressure layers
31!     pplev(klon,klev+1) Pressure levels
32!     pt(klon,klev)      Temperature (in K)
33!     ptsrf(klon)          Surface temperature
34!
35!     pdt(klon,klev)   Time derivative before condensation/sublimation of pt
36!     pdtsrf(klon)         Time derivative before condensation/sublimation of ptsrf
37!     picen2(klon)         n2 ice at the surface (kg/m2)
38!
39!     Outputs
40!     -------
41!     pdpsrf(klon)         \  Contribution of condensation/sublimation
42!     pdtc(klon,klev)  /  to the time derivatives of Ps, pt, and ptsrf
43!     pdtsrfc(klon)       /
44!     pdicen2(klon)         Tendancy n2 ice at the surface (kg/m2)
45!
46!     Both
47!     ----
48!     psolaralb(klon)      Albedo at the surface
49!     pemisurf(klon)       Emissivity of the surface
50!
51!     Authors
52!     -------
53!     Francois Forget (1996,2013)
54!     Converted to Fortran 90 and slightly modified by R. Wordsworth (2009)
55!
56!
57!==================================================================
58
59!-----------------------------------------------------------------------
60!     Arguments
61
62  INTEGER klon, klev, nq
63
64  REAL ptimestep
65  REAL pplay(klon,klev),pplev(klon,klev+1)
66  REAL pcapcal(klon)
67  REAL pt(klon,klev)
68  REAL ptsrf(klon),flu1(klon),flu2(klon),flu3(klon)
69  REAL pphi(klon,klev)
70  REAL pdt(klon,klev),pdtsrf(klon),pdtc(klon,klev)
71  REAL pdtsrfc(klon),pdpsrf(klon)
72  REAL picen2(klon),psolaralb(klon),pemisurf(klon)
73
74
75  REAL pu(klon,klev) ,  pv(klon,klev)
76  REAL pdu(klon,klev) , pdv(klon,klev)
77  REAL pduc(klon,klev) , pdvc(klon,klev)
78  REAL pq(klon,klev,nq),pdq(klon,klev,nq)
79  REAL pdqc(klon,klev,nq)
80
81!-----------------------------------------------------------------------
82!     Local variables
83
84  INTEGER l,ig,ilay,it,iq,ich4_gas
85
86  REAL*8 zt(klon,klev)
87  REAL tcond_n2
88  REAL pcond_n2
89  REAL glob_average2d         ! function
90  REAL zqn2(klon,klev) ! N2 MMR used to compute Tcond/zqn2
91  REAL ztcond (klon,klev)
92  REAL ztcondsol(klon),zfallheat
93  REAL pdicen2(klon)
94  REAL zcondicea(klon,klev), zcondices(klon)
95  REAL zfallice(klon,klev+1)
96  REAL zmflux(klev+1)
97  REAL zu(klev),zv(klev)
98  REAL zq(klev,nq),zq1(klev)
99  REAL ztsrf(klon)
100  REAL ztc(klev), ztm(klev+1)
101  REAL zum(klev+1) , zvm(klev+1)
102  REAL zqm(klev+1,nq),zqm1(klev+1)
103  LOGICAL condsub(klon)
104  REAL subptimestep
105  Integer Ntime
106  real masse (klev),w(klev+1)
107  real wq(klon,klev+1)
108  real vstokes,reff
109  real dWtotsn2
110  real condnconsn2(klon)
111  real nconsMAXn2
112!     Special diagnostic variables
113  real tconda1(klon,klev)
114  real tconda2(klon,klev)
115  real zdtsig (klon,klev)
116  real zdtlatent (klon,klev)
117  real zdt (klon,klev)
118!  REAL albediceF(klon)
119!   SAVE albediceF
120  INTEGER nsubtimestep,itsub    !number of subtimestep when calling vl1d
121  REAL subtimestep        !ptimestep/nsubtimestep
122  REAL dtmax
123
124  REAL zplevhist(klon)
125  REAL zplevnew(klon)
126  REAL zplev(klon)
127  REAL zpicen2(klon)
128  REAL ztsrfhist(klon)
129  REAL zdtsrf(klon)
130  REAL globzplevnew
131
132  real,dimension(:),save,allocatable :: vmrn2
133!$OMP THREADPRIVATE(vmrn2)
134  REAL stephan
135  DATA stephan/5.67e-08/  ! Stephan Boltzman constant
136  SAVE stephan
137!-----------------------------------------------------------------------
138!     Saved local variables
139
140!      REAL latcond
141  REAL ccond
142  REAL cpice  ! for atm condensation
143  SAVE cpice
144!      SAVE latcond,ccond
145  SAVE ccond
146
147  LOGICAL firstcall
148  SAVE firstcall
149  REAL SSUM
150  EXTERNAL SSUM
151
152!      DATA latcond /2.5e5/
153!      DATA latcond /1.98e5/
154  DATA cpice /1300./
155  DATA firstcall/.true./
156
157  INTEGER,SAVE :: i_n2ice=0       ! n2 ice
158  CHARACTER(LEN=20) :: tracername ! to temporarily store text
159  logical olkin  ! option to prevent N2 ice effect in the south
160  DATA olkin/.false./
161  save olkin
162
163  CHARACTER(len=10) :: tname
164
165!-----------------------------------------------------------------------
166
167!     Initialisation
168  IF (firstcall) THEN
169     ccond=cpp/(g*lw_n2)
170     print*,'In condense_n2cloud: ccond=',ccond,' latcond=',lw_n2
171
172     ! calculate global mean surface pressure for the fast mode
173     IF (.not. ALLOCATED(kp)) ALLOCATE(kp(klon))
174     DO ig=1,klon
175        kp(ig) = exp(-phisfi(ig)/(r*38.))
176     ENDDO
177     IF (fast) THEN
178        p00=glob_average2d(kp) ! mean pres at ref level
179     ENDIF
180
181     ALLOCATE(vmrn2(klon))
182     vmrn2(:) = 1.
183     !IF (ch4lag) then
184     !   DO ig=1,klon
185     !      if (latitude(ig)*180./pi.ge.latlag) then
186     !         vmrn2(ig) = vmrlag
187     !      endif
188     !   ENDDO
189     !ENDIF
190     IF (no_n2frost) then
191        DO ig=1,klon
192           if (picen2(ig).eq.0.) then
193              vmrn2(ig) = 1.e-15
194           endif
195        ENDDO
196     ENDIF
197     firstcall=.false.
198  ENDIF
199
200!-----------------------------------------------------------------------
201!     Calculation of n2 condensation / sublimation
202
203!     Variables used:
204!     picen2(klon)            : amount of n2 ice on the ground     (kg/m2)
205!     zcondicea(klon,klev): condensation rate in layer l       (kg/m2/s)
206!     zcondices(klon)         : condensation rate on the ground    (kg/m2/s)
207!     zfallice(klon,klev) : amount of ice falling from layer l (kg/m2/s)
208!     zdtlatent(klon,klev): dT/dt due to phase changes         (K/s)
209
210!     Tendencies initially set to 0
211  zcondices(1:klon) = 0.
212  pdtsrfc(1:klon) = 0.
213  pdpsrf(1:klon) = 0.
214  ztsrfhist(1:klon) = 0.
215  condsub(1:klon) = .false.
216  pdicen2(1:klon) = 0.
217  zfallheat=0
218  pdqc(1:klon,1:klev,1:nq)=0
219  pdtc(1:klon,1:klev)=0
220  pduc(1:klon,1:klev)=0
221  pdvc(1:klon,1:klev)=0
222  zfallice(1:klon,1:klev+1)=0
223  zcondicea(1:klon,1:klev)=0
224  zdtlatent(1:klon,1:klev)=0
225  zt(1:klon,1:klev)=0.
226
227!-----------------------------------------------------------------------
228!     Atmospheric condensation
229
230!     Condensation / sublimation in the atmosphere
231!     --------------------------------------------
232!      (calcul of zcondicea , zfallice and pdtc)
233
234  zt(1:klon,1:klev)=pt(1:klon,1:klev)+ pdt(1:klon,1:klev)*ptimestep
235  if (igcm_n2.ne.0) then
236     zqn2(1:klon,1:klev) = 1. ! & temporaire
237!        zqn2(1:klon,1:klev)= pq(1:klon,1:klev,igcm_n2) + pdq(1:klon,1:klev,igcm_n2)*ptimestep
238  else
239     zqn2(1:klon,1:klev) = 1.
240  end if
241
242  if (.not.fast) then
243!     Forecast the atmospheric frost temperature 'ztcond' with function tcond_n2
244    DO l=1,klev
245      DO ig=1,klon
246          ztcond (ig,l) = tcond_n2(pplay(ig,l),zqn2(ig,l))
247      ENDDO
248    ENDDO
249
250    DO l=klev,1,-1
251      DO ig=1,klon
252       pdtc(ig,l)=0.  ! final tendancy temperature set to 0
253
254       IF((zt(ig,l).LT.ztcond(ig,l)).or.(zfallice(ig,l+1).gt.0))THEN
255           condsub(ig)=.true.  !condensation at level l
256           IF (zfallice(ig,l+1).gt.0) then
257             zfallheat=zfallice(ig,l+1)*&
258            (pphi(ig,l+1)-pphi(ig,l) +&
259           cpice*(ztcond(ig,l+1)-ztcond(ig,l)))/lw_n2
260           ELSE
261               zfallheat=0.
262           ENDIF
263           zdtlatent(ig,l)=(ztcond(ig,l) - zt(ig,l))/ptimestep
264           zcondicea(ig,l)=(pplev(ig,l)-pplev(ig,l+1))&
265                         *ccond*zdtlatent(ig,l)- zfallheat
266!              Case when the ice from above sublimes entirely
267!              """""""""""""""""""""""""""""""""""""""""""""""
268           IF ((zfallice(ig,l+1).lt.-zcondicea(ig,l)) &
269                .AND. (zfallice(ig,l+1).gt.0)) THEN
270
271              zdtlatent(ig,l)=(-zfallice(ig,l+1)+zfallheat)/&
272               (ccond*(pplev(ig,l)-pplev(ig,l+1)))
273              zcondicea(ig,l)= -zfallice(ig,l+1)
274           END IF
275
276           zfallice(ig,l) = zcondicea(ig,l)+zfallice(ig,l+1)
277
278        END IF
279
280      ENDDO
281    ENDDO
282  endif  ! not fast
283
284!-----------------------------------------------------------------------
285!     Condensation/sublimation on the ground
286!     (calculation of zcondices and pdtsrfc)
287
288!     Adding subtimesteps : in fast version, pressures too low lead to
289!     instabilities.
290  IF (fast) THEN
291     IF (pplev(1,1).gt.0.3) THEN
292         nsubtimestep= 1
293     ELSE
294         nsubtimestep= nbsub !max(nint(ptimestep/dtmax),1)
295     ENDIF
296  ELSE
297     nsubtimestep= 1 ! if more then kp and p00 have to be calculated
298                     ! for nofast mode
299  ENDIF
300  subtimestep=ptimestep/float(nsubtimestep)
301
302  do itsub=1,nsubtimestep
303     ! first loop : getting zplev, ztsurf
304     IF (itsub.eq.1) then
305       DO ig=1,klon
306          zplev(ig)=pplev(ig,1)
307          ztsrfhist(ig)=ptsrf(ig) + pdtsrf(ig)*ptimestep
308          ztsrf(ig)=ptsrf(ig) + pdtsrf(ig)*subtimestep    !!
309          zpicen2(ig)=picen2(ig)
310       ENDDO
311     ELSE
312     ! additional loop :
313     ! 1)  pressure updated
314     ! 2)  direct redistribution of pressure over the globe
315     ! 3)  modification pressure for unstable cases
316     ! 4)  pressure update to conserve mass
317     ! 5)  temperature updated with radiative tendancies
318       DO ig=1,klon
319          zplevhist(ig)=zplev(ig)
320          zplevnew(ig)=zplev(ig)+pdpsrf(ig)*subtimestep   ! 1)
321          !IF (zplevnew(ig).le.0.0001) then
322          !   zplevnew(ig)=0.0001*kp(ig)/p00
323          !ENDIF
324       ENDDO
325       ! intermediaire de calcul: valeur moyenne de zplevnew (called twice in the code)
326       globzplevnew=glob_average2d(zplevnew)
327       DO ig=1,klon
328         zplev(ig)=kp(ig)*globzplevnew/p00  ! 2)
329       ENDDO
330       DO ig=1,klon     ! 3) unstable case condensation and sublimation: zplev=zplevhist
331          IF (((pdpsrf(ig).lt.0.).and.(tcond_n2(zplev(ig),zqn2(ig,1)).le.ztsrf(ig))).or.  &
332          ((pdpsrf(ig).gt.0.).and.(tcond_n2(zplev(ig),zqn2(ig,1)).ge.ztsrf(ig)))) then
333             zplev(ig)=zplevhist(ig)
334          ENDIF
335          zplevhist(ig)=zplev(ig)
336       ENDDO
337       zplev=zplev*globzplevnew/glob_average2d(zplevhist)   ! 4)
338       DO ig=1,klon    ! 5)
339         zdtsrf(ig)=pdtsrf(ig) + (stephan/pcapcal(ig))*(ptsrf(ig)**4-ztsrf(ig)**4)
340         ztsrf(ig)=ztsrf(ig)+pdtsrfc(ig)*subtimestep+zdtsrf(ig)*subtimestep
341         zpicen2(ig)=zpicen2(ig)+pdicen2(ig)*subtimestep
342       ENDDO
343     ENDIF  ! (itsub=1)
344
345   DO ig=1,klon
346     ! forecast of frost temperature ztcondsol
347     !ztcondsol(ig) = tcond_n2(zplev(ig),zqn2(ig,1))
348     ztcondsol(ig) = tcond_n2(zplev(ig),vmrn2(ig))
349
350!     Loop over where we have condensation / sublimation
351     IF ((ztsrf(ig) .LT. ztcondsol(ig)) .OR.           &    ! ground cond
352                ((ztsrf(ig) .GT. ztcondsol(ig)) .AND.  &    ! ground sublim
353               (zpicen2(ig) .GT. 0.))) THEN
354        condsub(ig) = .true.    ! condensation or sublimation
355
356!     Condensation or partial sublimation of N2 ice
357        if (ztsrf(ig) .LT. ztcondsol(ig)) then  ! condensation
358!             Include a correction to account for the cooling of air near
359!             the surface before condensing:
360         zcondices(ig)=pcapcal(ig)*(ztcondsol(ig)-ztsrf(ig)) &
361         /((lw_n2+cpp*(zt(ig,1)-ztcondsol(ig)))*subtimestep)
362        else                                    ! sublimation
363          zcondices(ig)=pcapcal(ig)*(ztcondsol(ig)-ztsrf(ig)) &
364          /(lw_n2*subtimestep)
365        end if
366
367        pdtsrfc(ig) = (ztcondsol(ig) - ztsrf(ig))/subtimestep
368
369!     partial sublimation of N2 ice
370!     If the entire N_2 ice layer sublimes
371!     (including what has just condensed in the atmosphere)
372        IF((zpicen2(ig)/subtimestep).LE. &
373            -zcondices(ig))THEN
374           zcondices(ig) = -zpicen2(ig)/subtimestep
375           pdtsrfc(ig)=(lw_n2/pcapcal(ig))*       &
376               (zcondices(ig))
377        END IF
378
379!     Changing N2 ice amount and pressure
380
381        pdicen2(ig) = zcondices(ig)
382        pdpsrf(ig)   = -pdicen2(ig)*g
383    !    pdpsrf(ig)   = 0. ! OPTION to check impact N2 sub/cond
384        IF (fast.and.(zplev(ig)+pdpsrf(ig)*subtimestep.le.0.0000001)) then
385            pdpsrf(ig)=(0.0000001*kp(ig)/p00-zplev(ig))/subtimestep
386            pdicen2(ig)=-pdpsrf(ig)/g
387        ENDIF
388
389     ELSE  ! no condsub
390        pdpsrf(ig)=0.
391        pdicen2(ig)=0.
392        pdtsrfc(ig)=0.
393     ENDIF
394   ENDDO ! ig
395  enddo ! subtimestep
396
397!     Updating pressure, temperature and ice reservoir
398  DO ig=1,klon
399     pdpsrf(ig)=(zplev(ig)+pdpsrf(ig)*subtimestep-pplev(ig,1))/ptimestep
400     ! Two options here : 1 ok, 2 is wrong
401     pdicen2(ig)=(zpicen2(ig)+pdicen2(ig)*subtimestep-picen2(ig))/ptimestep
402     !pdicen2(ig)=-pdpsrf(ig)/g
403
404     pdtsrfc(ig)=((ztsrf(ig)+pdtsrfc(ig)*subtimestep)-(ztsrfhist(ig)))/ptimestep
405
406! security
407     if (picen2(ig) + pdicen2(ig)*ptimestep.lt.0.) then
408       write(*,*) 'WARNING in condense_n2:'
409       write(*,*) picen2(ig),pdicen2(ig)*ptimestep
410       pdicen2(ig)= -picen2(ig)/ptimestep
411       pdpsrf(ig)=-pdicen2(ig)*g
412     endif
413
414     if(.not.picen2(ig).ge.0.) THEN
415!           if(picen2(ig) + pdicen2(ig)*ptimestep.le.-1.e-8) then
416          print*, 'WARNING NEG RESERVOIR in condense_n2: picen2(',ig,')=', picen2(ig) + pdicen2(ig)*ptimestep
417!              pdicen2(ig)= -picen2(ig)/ptimestep
418!           else
419          picen2(ig)=0.0
420!           endif
421     endif
422  ENDDO
423
424! ***************************************************************
425!  Correction to account for redistribution between sigma or hybrid
426!  layers when changing surface pressure (and warming/cooling
427!  of the n2 currently changing phase).
428! *************************************************************
429  if (.not.fast) then
430   DO ig=1,klon
431    if (condsub(ig)) then
432
433!  Mass fluxes through the sigma levels (kg.m-2.s-1)  (>0 when up)
434!  """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
435        zmflux(1) = -zcondices(ig)
436        DO l=1,klev
437         zmflux(l+1) = zmflux(l) -zcondicea(ig,l)  &
438         + (bp(l)-bp(l+1))*(zfallice(ig,1)-zmflux(1))
439! zmflux set to 0 if very low to avoid: top layer is disappearing in v1ld
440      if (abs(zmflux(l+1)).lt.1E-13.OR.bp(l+1).eq.0.) zmflux(l+1)=0.
441        END DO
442
443! Mass of each layer
444! ------------------
445       DO l=1,klev
446         masse(l)=(pplev(ig,l) - pplev(ig,l+1))/g
447       END DO
448
449
450!  Corresponding fluxes for T,U,V,Q
451!  """"""""""""""""""""""""""""""""
452!           averaging operator for TRANSPORT
453!           """"""""""""""""""""""""""""""""
454
455!           Subtimestep loop to perform the redistribution gently and simultaneously with
456!           the other tendencies
457!           Estimation of subtimestep (using only  the first layer, the most critical)
458        dtmax=ptimestep
459        if (zmflux(1).gt.1.e-20) then
460            dtmax=min(dtmax,masse(1)*zqn2(ig,1)/abs(zmflux(1)))
461        endif
462        nsubtimestep= max(nint(ptimestep/dtmax),nint(2.))
463        subtimestep=ptimestep/float(nsubtimestep)
464
465!          New flux for each subtimestep
466       do l=1,klev+1
467            w(l)=-zmflux(l)*subtimestep
468       enddo
469!          initializing variables that will vary during subtimestep:
470       do l=1,klev
471           ztc(l)  =pt(ig,l)
472           zu(l)   =pu(ig,l)
473           zv(l)   =pv(ig,l)
474           do iq=1,nq
475              zq(l,iq) = pq(ig,l,iq)
476           enddo
477       end do
478
479!          loop over nsubtimestep
480!          """"""""""""""""""""""
481       do itsub=1,nsubtimestep
482!             Progressively adding tendancies from other processes.
483          do l=1,klev
484             ztc(l)  =ztc(l)   +(pdt(ig,l) + zdtlatent(ig,l))*subtimestep
485             zu(l)   =zu(l)   +pdu( ig,l) * subtimestep
486             zv(l)   =zv(l)   +pdv( ig,l) * subtimestep
487             do iq=1,nq
488                zq(l,iq) = zq(l,iq) + pdq(ig,l,iq)* subtimestep
489             enddo
490           end do
491
492!             Change of mass in each layer
493          do l=1,klev
494             masse(l)=masse(l)+pdpsrf(ig)*subtimestep*(pplev(ig,l) - pplev(ig,l+1))&
495                                             /(g*pplev(ig,1))
496          end do
497
498!             Value transfert at the surface interface when condensation sublimation:
499
500          if (zmflux(1).lt.0) then
501!               Surface condensation
502            zum(1)= zu(1)
503            zvm(1)= zv(1)
504            ztm(1) = ztc(1)
505          else
506!               Surface sublimation:
507            ztm(1) = ztsrf(ig) + pdtsrfc(ig)*ptimestep
508            zum(1) = 0
509            zvm(1) = 0
510          end if
511          do iq=1,nq
512              zqm(1,iq)=0. ! most tracer do not condense !
513          enddo
514!             Special case if the tracer is n2 gas
515          if (igcm_n2.ne.0) zqm(1,igcm_n2)=1.
516
517           ztm(2:klev+1)=0.
518           zum(2:klev+1)=0.
519           zvm(2:klev+1)=0.
520           zqm1(1:klev+1)=0.
521
522!             Van Leer scheme:
523          call vl1d(klev,ztc,2.,masse,w,ztm)
524          call vl1d(klev,zu ,2.,masse,w,zum)
525          call vl1d(klev,zv ,2.,masse,w,zvm)
526         do iq=1,nq
527           do l=1,klev
528            zq1(l)=zq(l,iq)
529           enddo
530           zqm1(1)=zqm(1,iq)
531           call vl1d(klev,zq1,2.,masse,w,zqm1)
532           do l=2,klev
533            zqm(l,iq)=zqm1(l)
534           enddo
535          enddo
536
537!             Correction to prevent negative value for qn2
538          if (igcm_n2.ne.0) then
539            zqm(1,igcm_n2)=1.
540            do l=1,klev-1
541              if (w(l)*zqm(l,igcm_n2).gt.zq(l,igcm_n2)*masse(l)) then
542                 zqm(l+1,igcm_n2)=max(zqm(l+1,igcm_n2),    &
543                 (zqm(l,igcm_n2)*w(l) -zq(l,igcm_n2)*masse(l))/w(l+1) )
544              else
545                exit
546              endif
547            end do
548          end if
549
550!         Value transfert at the surface interface when condensation sublimation:
551          if (zmflux(1).lt.0) then
552!               Surface condensation
553                zum(1)= zu(1)
554                zvm(1)= zv(1)
555                ztm(1) = ztc(1)
556          else
557!               Surface sublimation:
558                ztm(1) = ztsrf(ig) + pdtsrfc(ig)*ptimestep
559                zum(1) = 0 
560                zvm(1) = 0 
561          end if
562          do iq=1,nq
563                 zqm(1,iq)=0. ! most tracer do not condense !
564          enddo
565!         Special case if the tracer is n2 gas
566          if (igcm_n2.ne.0) zqm(1,igcm_n2)=1.
567
568           !!! Source haze: 0.02 pourcent when n2 sublimes
569          IF (source_haze) THEN
570               IF (pdicen2(ig).lt.0) THEN
571                DO iq=1,nq
572                 tname=noms(iq)
573                 if (tname(1:4).eq."haze") then
574                       !zqm(1,iq)=0.02
575                       !zqm(1,iq)=-pdicen2(ig)*0.02
576                       zqm(1,iq)=-pdicen2(ig)*ptimestep*0.02
577                       !zqm(10,iq)=-pdicen2(ig)*ptimestep*100.
578                       !zqm(1,iq)=-pdicen2(ig)*1000000.
579
580                 endif
581                ENDDO
582               ENDIF     
583          ENDIF     
584          ztm(klev+1)= ztc(klev) ! should not be used, but...
585          zum(klev+1)= zu(klev)  ! should not be used, but...
586          zvm(klev+1)= zv(klev)  ! should not be used, but...
587          do iq=1,nq
588             zqm(klev+1,iq)= zq(klev,iq)
589          enddo
590
591!         Tendencies on T, U, V, Q
592!         """""""""""""""""""""""
593          DO l=1,klev
594
595!               Tendencies on T
596            zdtsig(ig,l) = (1/masse(l)) *   &
597           ( zmflux(l)*(ztm(l) - ztc(l))    &
598           - zmflux(l+1)*(ztm(l+1) - ztc(l)) &
599           + zcondicea(ig,l)*(ztcond(ig,l)-ztc(l))  )
600
601!               Tendencies on U
602              pduc(ig,l)   = (1/masse(l)) *   &
603           ( zmflux(l)*(zum(l) - zu(l))&
604           - zmflux(l+1)*(zum(l+1) - zu(l)) )
605
606!               Tendencies on V
607              pdvc(ig,l)   = (1/masse(l)) *   &
608           ( zmflux(l)*(zvm(l) - zv(l))   &
609           - zmflux(l+1)*(zvm(l+1) - zv(l)) )
610
611          END DO
612
613!             Tendencies on Q
614          do iq=1,nq
615            if (iq.eq.igcm_n2) then
616!                 SPECIAL Case when the tracer IS N2 :
617              DO l=1,klev
618              pdqc(ig,l,iq)= (1/masse(l)) *        &
619               ( zmflux(l)*(zqm(l,iq) - zq(l,iq))    &
620               - zmflux(l+1)*(zqm(l+1,iq) - zq(l,iq))&
621               + zcondicea(ig,l)*(zq(l,iq)-1.) )
622              END DO
623            else
624              DO l=1,klev
625                pdqc(ig,l,iq)= (1/masse(l)) *         &
626               ( zmflux(l)*(zqm(l,iq) - zq(l,iq))     &
627               - zmflux(l+1)*(zqm(l+1,iq) - zq(l,iq)) &
628               + zcondicea(ig,l)*zq(l,iq) )
629              END DO
630            end if
631          enddo
632!             Update variables at the end of each subtimestep.
633          do l=1,klev
634            ztc(l)  =ztc(l)  + zdtsig(ig,l)  *subtimestep
635            zu(l)   =zu(l)   + pduc(ig,l)  *subtimestep
636            zv(l)   =zv(l)   + pdvc(ig,l)  *subtimestep
637            do iq=1,nq
638              zq(l,iq) = zq(l,iq) + pdqc(ig,l,iq)*subtimestep
639            enddo
640          end do
641       enddo   ! loop on nsubtimestep
642!          Recomputing Total tendencies
643       do l=1,klev
644         pdtc(ig,l)   = (ztc(l) - zt(ig,l) )/ptimestep
645         pduc(ig,l)   = (zu(l) - (pu(ig,l) + pdu(ig,l)*ptimestep))/ptimestep
646         pdvc(ig,l)   = (zv(l) - (pv(ig,l) + pdv(ig,l)*ptimestep))/ptimestep
647         do iq=1,nq
648           pdqc(ig,l,iq) = (zq(l,iq) - (pq(ig,l,iq) + pdq(ig,l,iq)*ptimestep))/ptimestep
649
650
651!           Correction temporaire
652            if (iq.eq.igcm_n2) then
653             if((pq(ig,l,iq) +(pdqc(ig,l,iq)+ pdq(ig,l,iq))*ptimestep) &
654                     .lt.0.01) then ! if n2 < 1 %  !
655                     write(*,*) 'Warning: n2 < 1%'
656                     pdqc(ig,l,iq)=(0.01-pq(ig,l,iq))/ptimestep-pdq(ig,l,iq)
657             end if
658            end if
659
660         enddo
661       end do
662! *******************************TEMPORAIRE ******************
663       if (klon.eq.1) then
664              write(*,*) 'nsubtimestep=' ,nsubtimestep
665           write(*,*) 'masse avant' , (pplev(ig,1) - pplev(ig,2))/g
666              write(*,*) 'masse apres' , masse(1)
667              write(*,*) 'zmflux*DT, l=1' ,  zmflux(1)*ptimestep
668              write(*,*) 'zmflux*DT, l=2' ,  zmflux(2)*ptimestep
669          write(*,*) 'pq, l=1,2,3' ,  pq(1,1,1), pq(1,2,1),pq(1,3,1)
670          write(*,*) 'zq, l=1,2,3' ,  zq(1,1), zq(2,1),zq(3,1)
671              write(*,*) 'dq*Dt l=1' ,  pdq(1,1,1)*ptimestep
672              write(*,*) 'dqcond*Dt l=1' ,  pdqc(1,1,1)*ptimestep
673       end if
674
675! ***********************************************************
676    end if ! if (condsub)
677   END DO  ! loop on ig
678  endif ! not fast
679
680! ************ Option Olkin to prevent N2 effect in the south********
681112   continue
682  if (olkin) then
683  DO ig=1,klon
684     if (latitude(ig).lt.0) then
685       pdtsrfc(ig) = max(0.,pdtsrfc(ig))
686       pdpsrf(ig) = 0.
687       pdicen2(ig)  = 0.
688       do l=1,klev
689         pdtc(ig,l)  =  max(0.,zdtlatent(ig,l))
690         pduc(ig,l)  = 0.
691         pdvc(ig,l)  = 0.
692         do iq=1,nq
693           pdqc(ig,l,iq) = 0
694         enddo
695       end do
696     end if
697  END DO
698  end if
699! *******************************************************************
700
701! ***************************************************************
702! Ecriture des diagnostiques
703! ***************************************************************
704
705!     DO l=1,klev
706!        DO ig=1,klon
707!         Taux de cond en kg.m-2.pa-1.s-1
708!          tconda1(ig,l)=zcondicea(ig,l)/(pplev(ig,l)-pplev(ig,l+1))
709!          Taux de cond en kg.m-3.s-1
710!          tconda2(ig,l)=tconda1(ig,l)*pplay(ig,l)*g/(r*pt(ig,l))
711!        END DO
712!     END DO
713!     call WRITEDIAGFI(klon,'tconda1',              &
714!     'Taux de condensation N2 atmospherique /Pa',    &
715!      'kg.m-2.Pa-1.s-1',3,tconda1)
716!     call WRITEDIAGFI(klon,'tconda2',              &
717!     'Taux de condensation N2 atmospherique /m',     &
718!      'kg.m-3.s-1',3,tconda2)
719
720
721  return
722  end subroutine condense_n2
723
724!-------------------------------------------------------------------------
725
726  real function tcond_n2(p,vmr)
727!   Calculates the condensation temperature for N2 at pressure P and vmr
728  implicit none
729  real, intent(in):: p,vmr
730
731!       tcond_n2 = (1.)/(0.026315-0.0011877*log(.7143*p*vmr))
732  IF (p.ge.0.529995) then
733!     tcond Fray and Schmitt for N2 phase beta (T>35.6 K) FIT TB
734  ! tcond_n2 = (1.)/(1./63.1470-296.925/(2.5e5*0.98)*log(1./(0.125570*1.e5)*p*vmr))
735    tcond_n2 = (1.)/(0.01583606505-1.211938776e-3*log(7.963685594e-5*p*vmr))
736  ELSE
737!     tcond Fray and Schmitt for N2 phase alpha(T<35.6 K) FIT BT
738  ! tcond_n2 = (1.)/(1./35.6-296.925/(2.5e5*1.09)*log(1./(0.508059)*p*vmr))
739    tcond_n2 = (1.)/(1./35.6-1.089633028e-3*log(1.968275338*p*vmr))
740  ENDIF
741  return
742  end function tcond_n2
743
744!-------------------------------------------------------------------------
745
746  real function pcond_n2(t,vmr)
747!   Calculates the condensation pressure for N2 at temperature T and vmr
748  implicit none
749  real, intent(in):: t,vmr
750
751!       tcond_n2 = (1.)/(0.026315-0.0011877*log(.7143*p*vmr))
752  IF (t.ge.35.6) then
753!     tcond Fray and Schmitt for N2 phase beta (T>35.6 K) FIT TB
754  ! pcond_n2 = 0.125570*1.e5/vmr*exp((2.5e5*0.98)/296.925*(1./63.1470-1./t))
755    pcond_n2 = 0.125570e5/vmr*exp(825.1241896*(1./63.147-1./t))
756  ELSE
757!     tcond Fray and Schmitt for N2 phase alpha(T<35.6 K) FIT TB
758  ! pcond_n2 = 0.508059/vmr*exp((2.5e5*1.09)/296.925*(1./35.6-1./t))
759    pcond_n2 = 0.508059/vmr*exp(917.7401701*(1./35.6-1./t))
760  ENDIF
761  return
762  end function pcond_n2
763
764!-------------------------------------------------------------------------
765
766  real function glob_average2d(var)
767!   Calculates the global average of variable var
768  use comgeomfi_h
769  use dimphy, only: klon
770  USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat
771  use geometry_mod, only: cell_area, latitude
772
773  implicit none
774
775!      INTEGER klon
776  REAL var(klon)
777  INTEGER ig
778
779  glob_average2d = 0.
780  DO ig=2,klon-1
781       glob_average2d = glob_average2d + var(ig)*cell_area(ig)
782  END DO
783  glob_average2d = glob_average2d + var(1)*cell_area(1)*nbp_lon
784  glob_average2d = glob_average2d + var(klon)*cell_area(klon)*nbp_lon
785  glob_average2d = glob_average2d/(totarea+(cell_area(1)+cell_area(klon))*(nbp_lon-1))
786
787  end function glob_average2d
788
789! *****************************************************************
790
791  subroutine vl1d(klev,q,pente_max,masse,w,qm)
792!
793!     Operateur de moyenne inter-couche pour calcul de transport type
794!     Van-Leer " pseudo amont " dans la verticale
795!    q,w sont des arguments d'entree  pour le s-pg ....
796!    masse : masse de la couche Dp/g
797!    w : masse d'atm ``transferee'' a chaque pas de temps (kg.m-2)
798!    pente_max = 2 conseillee
799!   --------------------------------------------------------------------
800  IMPLICIT NONE
801
802!   Arguments:
803!   ----------
804  integer klev
805  real masse(klev),pente_max
806  REAL q(klev),qm(klev+1)
807  REAL w(klev+1)
808!
809!      Local
810!   ---------
811!
812  INTEGER l
813!
814  real dzq(klev),dzqw(klev),adzqw(klev),dzqmax
815  real sigw, Mtot, MQtot
816  integer m
817
818
819!    On oriente tout dans le sens de la pression
820!     W > 0 WHEN DOWN !!!!!!!!!!!!!
821
822  do l=2,klev
823        dzqw(l)=q(l-1)-q(l)
824        adzqw(l)=abs(dzqw(l))
825  enddo
826
827  do l=2,klev-1
828        if(dzqw(l)*dzqw(l+1).gt.0.) then
829            dzq(l)=0.5*(dzqw(l)+dzqw(l+1))
830        else
831            dzq(l)=0.
832        endif
833        dzqmax=pente_max*min(adzqw(l),adzqw(l+1))
834        dzq(l)=sign(min(abs(dzq(l)),dzqmax),dzq(l))
835  enddo
836
837     dzq(1)=0.
838     dzq(klev)=0.
839
840   do l = 1,klev-1
841
842!         Regular scheme (transfered mass < layer mass)
843!         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
844      if(w(l+1).gt.0. .and. w(l+1).le.masse(l+1)) then
845         sigw=w(l+1)/masse(l+1)
846         qm(l+1)=(q(l+1)+0.5*(1.-sigw)*dzq(l+1))
847      else if(w(l+1).le.0. .and. -w(l+1).le.masse(l)) then
848         sigw=w(l+1)/masse(l)
849         qm(l+1)=(q(l)-0.5*(1.+sigw)*dzq(l))
850
851!         Extended scheme (transfered mass > layer mass)
852!         ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
853      else if(w(l+1).gt.0.) then
854         m=l+1
855         Mtot = masse(m)
856         MQtot = masse(m)*q(m)
857         if (m.lt.klev) then ! because some compilers will have problems
858                              ! evaluating masse(klev+1)
859         do while ((m.lt.klev).and.(w(l+1).gt.(Mtot+masse(m+1))))
860            m=m+1
861            Mtot = Mtot + masse(m)
862            MQtot = MQtot + masse(m)*q(m)
863           if (m.eq.klev) exit
864         end do
865         endif
866         if (m.lt.klev) then
867            sigw=(w(l+1)-Mtot)/masse(m+1)
868            qm(l+1)= (1/w(l+1))*(MQtot + (w(l+1)-Mtot)* &
869           (q(m+1)+0.5*(1.-sigw)*dzq(m+1)) )
870         else
871            w(l+1) = Mtot
872            qm(l+1) = Mqtot / Mtot
873            write(*,*) 'top layer is disapearing !'
874            stop
875         end if
876      else      ! if(w(l+1).lt.0)
877         m = l-1
878         Mtot = masse(m+1)
879         MQtot = masse(m+1)*q(m+1)
880         if (m.gt.0) then ! because some compilers will have problems
881                          ! evaluating masse(0)
882          do while ((m.gt.0).and.(-w(l+1).gt.(Mtot+masse(m))))
883            m=m-1
884            Mtot = Mtot + masse(m+1)
885            MQtot = MQtot + masse(m+1)*q(m+1)
886            if (m.eq.0) exit
887          end do
888         endif
889         if (m.gt.0) then
890            sigw=(w(l+1)+Mtot)/masse(m)
891            qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)* &
892           (q(m)-0.5*(1.+sigw)*dzq(m))  )
893         else
894            qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*qm(1))
895         end if
896      end if
897   enddo
898
899   return
900   end  subroutine vl1d
Note: See TracBrowser for help on using the repository browser.