source: trunk/LMDZ.GENERIC/libf/phystd/vdifc.F @ 135

Last change on this file since 135 was 135, checked in by aslmd, 14 years ago

CHANGEMENT ARBORESCENCE ETAPE 2 -- NON COMPLET

File size: 19.5 KB
Line 
1      SUBROUTINE vdifc(ngrid,nlay,nq,ppopsk,
2     $                ptimestep,pcapcal,lecrit,
3     $                pplay,pplev,pzlay,pzlev,pz0,
4     $                pu,pv,ph,pq,ptsrf,pemis,pqsurf,
5     $                pdufi,pdvfi,pdhfi,pdqfi,pfluxsrf,
6     $                pdudif,pdvdif,pdhdif,pdtsrf,pq2,
7     $                pdqdif,pdqsdif)
8
9      use watercommon_h, only : RLVTT
10
11      implicit none
12
13c=======================================================================
14c
15c   subject:
16c   --------
17c   Turbulent diffusion (mixing) for potential T, U, V and tracer
18c
19c   Shema implicite
20c   On commence par rajouter au variables x la tendance physique
21c   et on resoult en fait:
22c      x(t+1) =  x(t) + dt * (dx/dt)phys(t)  +  dt * (dx/dt)difv(t+1)
23c
24c   author:
25c   ------
26c      Hourdin/Forget/Fournier
27c=======================================================================
28
29c-----------------------------------------------------------------------
30c   declarations:
31c   -------------
32
33#include "dimensions.h"
34#include "dimphys.h"
35#include "comcstfi.h"
36#include "callkeys.h"
37#include "surfdat.h"
38#include "comgeomfi.h"
39#include "tracer.h"
40
41#include "watercap.h"
42c
43c   arguments:
44c   ----------
45
46      INTEGER ngrid,nlay
47      REAL ptimestep
48      REAL pplay(ngrid,nlay),pplev(ngrid,nlay+1)
49      REAL pzlay(ngrid,nlay),pzlev(ngrid,nlay+1)
50      REAL pu(ngrid,nlay),pv(ngrid,nlay),ph(ngrid,nlay)
51      REAL ptsrf(ngrid),pemis(ngrid)
52      REAL pdufi(ngrid,nlay),pdvfi(ngrid,nlay),pdhfi(ngrid,nlay)
53      REAL pfluxsrf(ngrid)
54      REAL pdudif(ngrid,nlay),pdvdif(ngrid,nlay),pdhdif(ngrid,nlay)
55      REAL pdtsrf(ngrid),pcapcal(ngrid)
56      REAL pq2(ngrid,nlay+1)
57
58c    Argument added for condensation:
59!      REAL co2ice (ngrid)
60      REAL ppopsk(ngrid,nlay)
61      logical lecrit
62      REAL pz0
63
64c    Traceurs :
65      integer nq
66      REAL pqsurf(ngrid,nq)
67      real pq(ngrid,nlay,nq), pdqfi(ngrid,nlay,nq)
68      real pdqdif(ngrid,nlay,nq)
69      real pdqsdif(ngrid,nq)
70     
71c   local:
72c   ------
73
74      INTEGER ilev,ig,ilay,nlev
75
76      REAL z4st,zdplanck(ngridmx)
77      REAL zkv(ngridmx,nlayermx+1),zkh(ngridmx,nlayermx+1)
78      REAL zcdv(ngridmx),zcdh(ngridmx)
79      REAL zcdv_true(ngridmx),zcdh_true(ngridmx)
80      REAL zu(ngridmx,nlayermx),zv(ngridmx,nlayermx)
81      REAL zh(ngridmx,nlayermx)
82      REAL ztsrf2(ngridmx)
83      REAL z1(ngridmx),z2(ngridmx)
84      REAL za(ngridmx,nlayermx),zb(ngridmx,nlayermx)
85      REAL zb0(ngridmx,nlayermx)
86      REAL zc(ngridmx,nlayermx),zd(ngridmx,nlayermx)
87      REAL zcst1
88      REAL zu2
89
90      EXTERNAL SSUM,SCOPY
91      REAL SSUM
92      LOGICAL firstcall
93      SAVE firstcall
94
95c     variable added for CO2 condensation:
96c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
97      REAL hh , zhcond(ngridmx,nlayermx)
98      REAL latcond,tcond1mb
99      REAL acond,bcond
100      SAVE acond,bcond
101      DATA latcond,tcond1mb/5.9e5,136.27/
102
103c    Tracers :
104c    ~~~~~~~
105      INTEGER iq
106      REAL zq(ngridmx,nlayermx,nqmx)
107      REAL zq1temp(ngridmx)
108      REAL rho(ngridmx) ! near surface air density
109      REAL qsat(ngridmx)
110      DATA firstcall/.true./
111      REAL kmixmin
112
113c    ** un petit test de coherence
114c       --------------------------
115
116      IF (firstcall) THEN
117         IF(ngrid.NE.ngridmx) THEN
118            PRINT*,'STOP dans vdifc'
119            PRINT*,'probleme de dimensions :'
120            PRINT*,'ngrid  =',ngrid
121            PRINT*,'ngridmx  =',ngridmx
122            STOP
123         ENDIF
124c        To compute: Tcond= 1./(bcond-acond*log(.0095*p)) (p in pascal)
125         bcond=1./tcond1mb
126         acond=r/latcond
127         PRINT*,'In vdifc: Tcond(P=1mb)=',tcond1mb,' Lcond=',latcond
128         PRINT*,'          acond,bcond',acond,bcond
129
130        firstcall=.false.
131      ENDIF
132
133
134
135
136
137c-----------------------------------------------------------------------
138c    1. initialisation
139c    -----------------
140
141      nlev=nlay+1
142
143c    ** calcul de rho*dz et dt*rho/dz=dt*rho**2 g/dp
144c       avec rho=p/RT=p/ (R Theta) (p/ps)**kappa
145c       ----------------------------------------
146
147      DO ilay=1,nlay
148         DO ig=1,ngrid
149            za(ig,ilay)=(pplev(ig,ilay)-pplev(ig,ilay+1))/g
150         ENDDO
151      ENDDO
152
153      zcst1=4.*g*ptimestep/(r*r)
154      DO ilev=2,nlev-1
155         DO ig=1,ngrid
156            zb0(ig,ilev)=pplev(ig,ilev)*
157     s      (pplev(ig,1)/pplev(ig,ilev))**rcp /
158     s      (ph(ig,ilev-1)+ph(ig,ilev))
159            zb0(ig,ilev)=zcst1*zb0(ig,ilev)*zb0(ig,ilev)/
160     s      (pplay(ig,ilev-1)-pplay(ig,ilev))
161         ENDDO
162      ENDDO
163      DO ig=1,ngrid
164                 zb0(ig,1)=ptimestep*pplev(ig,1)/(r*ptsrf(ig))
165      ENDDO
166
167c    ** diagnostique pour l'initialisation
168c       ----------------------------------
169
170      IF(lecrit) THEN
171         ig=ngrid/2+1
172         PRINT*,'Pression (mbar) ,altitude (km),u,v,theta, rho dz'
173         DO ilay=1,nlay
174            WRITE(*,'(6f11.5)')
175     s      .01*pplay(ig,ilay),.001*pzlay(ig,ilay),
176     s      pu(ig,ilay),pv(ig,ilay),ph(ig,ilay),za(ig,ilay)
177         ENDDO
178         PRINT*,'Pression (mbar) ,altitude (km),zb'
179         DO ilev=1,nlay
180            WRITE(*,'(3f15.7)')
181     s      .01*pplev(ig,ilev),.001*pzlev(ig,ilev),
182     s      zb0(ig,ilev)
183         ENDDO
184      ENDIF
185
186c     Potential Condensation temperature:
187c     -----------------------------------
188
189c     if (callcond) then
190c       DO ilev=1,nlay
191c         DO ig=1,ngrid
192c           zhcond(ig,ilev) =
193c    &  (1./(bcond-acond*log(.0095*pplay(ig,ilev))))/ppopsk(ig,ilev)
194c         END DO
195c       END DO
196c     else
197        call zerophys(ngrid*nlay,zhcond)
198c     end if
199
200
201c-----------------------------------------------------------------------
202c   2. ajout des tendances physiques
203c      -----------------------------
204
205      DO ilev=1,nlay
206         DO ig=1,ngrid
207            zu(ig,ilev)=pu(ig,ilev)+pdufi(ig,ilev)*ptimestep
208            zv(ig,ilev)=pv(ig,ilev)+pdvfi(ig,ilev)*ptimestep
209            zh(ig,ilev)=ph(ig,ilev)+pdhfi(ig,ilev)*ptimestep
210            zh(ig,ilev)=max(zh(ig,ilev),zhcond(ig,ilev))
211         ENDDO
212      ENDDO
213      if(tracer) then
214        DO iq =1, nq
215         DO ilev=1,nlay
216           DO ig=1,ngrid
217              zq(ig,ilev,iq)=pq(ig,ilev,iq)+pdqfi(ig,ilev,iq)*ptimestep
218           ENDDO
219         ENDDO
220        ENDDO
221      end if
222
223c-----------------------------------------------------------------------
224c   3. schema de turbulence
225c      --------------------
226
227c    ** source d'energie cinetique turbulente a la surface
228c       (condition aux limites du schema de diffusion turbulente
229c       dans la couche limite
230c       ---------------------
231
232      CALL vdif_cd( ngrid,nlay,pz0,g,pzlay,pu,pv,ptsrf,ph
233     &             ,zcdv_true,zcdh_true)
234      DO ig=1,ngrid
235        zu2=pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1)
236        zcdv(ig)=zcdv_true(ig)*sqrt(zu2)
237        zcdh(ig)=zcdh_true(ig)*sqrt(zu2)
238      ENDDO
239
240c    ** schema de diffusion turbulente dans la couche limite
241c       ----------------------------------------------------
242
243        CALL vdif_kc(ptimestep,g,pzlev,pzlay
244     &              ,pu,pv,ph,zcdv_true
245     &              ,pq2,zkv,zkh)
246
247
248c    Adding eddy mixing to mimic 3D general circulation in 1D
249c    RW FF 2010
250      if ((ngrid.eq.1)) then
251        kmixmin = 1.0e-2   ! minimum eddy mix coeff in 1D
252        do ilev=1,nlay
253          do ig=1,ngrid
254           zkh(ig,ilev) = max(kmixmin,zkh(ig,ilev))
255           zkv(ig,ilev) = max(kmixmin,zkv(ig,ilev))
256           !zkh(ig,ilev) = kmixmin
257           !zkv(ig,ilev) = kmixmin
258          end do
259        end do
260      end if
261
262
263c    ** diagnostique pour le schema de turbulence
264c       -----------------------------------------
265
266      IF(lecrit) THEN
267         PRINT*
268         PRINT*,'Diagnostic for the vertical turbulent mixing'
269         PRINT*,'Cd for momentum and potential temperature'
270
271         PRINT*,zcdv(ngrid/2+1),zcdh(ngrid/2+1)
272         PRINT*,'Mixing coefficient for momentum and pot.temp.'
273         DO ilev=1,nlay
274            PRINT*,zkv(ngrid/2+1,ilev),zkh(ngrid/2+1,ilev)
275         ENDDO
276      ENDIF
277
278c-----------------------------------------------------------------------
279c   4. inversion pour l'implicite sur u
280c      --------------------------------
281
282c    ** l'equation est
283c       u(t+1) =  u(t) + dt * {(du/dt)phys}(t)  +  dt * {(du/dt)difv}(t+1)
284c       avec
285c       /zu/ = u(t) + dt * {(du/dt)phys}(t)   (voir paragraphe 2.)
286c       et
287c       dt * {(du/dt)difv}(t+1) = dt * {(d/dz)[ Ku (du/dz) ]}(t+1)
288c       donc les entrees sont /zcdv/ pour la condition a la limite sol
289c       et /zkv/ = Ku
290 
291      CALL multipl((nlay-1)*ngrid,zkv(1,2),zb0(1,2),zb(1,2))
292      CALL multipl(ngrid,zcdv,zb0,zb)
293
294      DO ig=1,ngrid
295         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
296         zc(ig,nlay)=za(ig,nlay)*zu(ig,nlay)*z1(ig)
297         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
298      ENDDO
299
300      DO ilay=nlay-1,1,-1
301         DO ig=1,ngrid
302            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
303     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
304            zc(ig,ilay)=(za(ig,ilay)*zu(ig,ilay)+
305     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
306            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
307         ENDDO
308      ENDDO
309
310      DO ig=1,ngrid
311         zu(ig,1)=zc(ig,1)
312      ENDDO
313      DO ilay=2,nlay
314         DO ig=1,ngrid
315            zu(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zu(ig,ilay-1)
316         ENDDO
317      ENDDO
318
319c-----------------------------------------------------------------------
320c   5. inversion pour l'implicite sur v
321c      --------------------------------
322
323c    ** l'equation est
324c       v(t+1) =  v(t) + dt * {(dv/dt)phys}(t)  +  dt * {(dv/dt)difv}(t+1)
325c       avec
326c       /zv/ = v(t) + dt * {(dv/dt)phys}(t)   (voir paragraphe 2.)
327c       et
328c       dt * {(dv/dt)difv}(t+1) = dt * {(d/dz)[ Kv (dv/dz) ]}(t+1)
329c       donc les entrees sont /zcdv/ pour la condition a la limite sol
330c       et /zkv/ = Kv
331
332      DO ig=1,ngrid
333         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
334         zc(ig,nlay)=za(ig,nlay)*zv(ig,nlay)*z1(ig)
335         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
336      ENDDO
337
338      DO ilay=nlay-1,1,-1
339         DO ig=1,ngrid
340            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
341     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
342            zc(ig,ilay)=(za(ig,ilay)*zv(ig,ilay)+
343     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
344            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
345         ENDDO
346      ENDDO
347
348      DO ig=1,ngrid
349         zv(ig,1)=zc(ig,1)
350      ENDDO
351      DO ilay=2,nlay
352         DO ig=1,ngrid
353            zv(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zv(ig,ilay-1)
354         ENDDO
355      ENDDO
356
357c-----------------------------------------------------------------------
358c   6. inversion pour l'implicite sur h sans oublier le couplage
359c      avec le sol (conduction)
360c      ------------------------
361
362c    ** l'equation est
363c       h(t+1) =  h(t) + dt * {(dh/dt)phys}(t)  +  dt * {(dh/dt)difv}(t+1)
364c       avec
365c       /zh/ = h(t) + dt * {(dh/dt)phys}(t)   (voir paragraphe 2.)
366c       et
367c       dt * {(dh/dt)difv}(t+1) = dt * {(d/dz)[ Kh (dh/dz) ]}(t+1)
368c       donc les entrees sont /zcdh/ pour la condition de raccord au sol
369c       et /zkh/ = Kh
370c       -------------
371
372      CALL multipl((nlay-1)*ngrid,zkh(1,2),zb0(1,2),zb(1,2))
373      CALL multipl(ngrid,zcdh,zb0,zb)
374
375      DO ig=1,ngrid
376         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
377         zc(ig,nlay)=za(ig,nlay)*zh(ig,nlay)*z1(ig)
378         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
379      ENDDO
380
381      DO ilay=nlay-1,1,-1
382         DO ig=1,ngrid
383            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
384     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
385            zc(ig,ilay)=(za(ig,ilay)*zh(ig,ilay)+
386     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
387            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
388         ENDDO
389      ENDDO
390
391c    ** calcul de (d Planck / dT) a la temperature d'interface
392c       ------------------------------------------------------
393
394      z4st=4.*5.67e-8*ptimestep
395      DO ig=1,ngrid
396         zdplanck(ig)=z4st*pemis(ig)*ptsrf(ig)*ptsrf(ig)*ptsrf(ig)
397      ENDDO
398
399c    ** calcul de la temperature_d'interface et de sa tendance.
400c       on ecrit que la somme des flux est nulle a l'interface
401c       a t + \delta t,
402c       c'est a dire le flux radiatif a {t + \delta t}
403c       + le flux turbulent a {t + \delta t}
404c            qui s'ecrit K (T1-Tsurf) avec T1 = d1 Tsurf + c1
405c            (notation K dt = /cpp*b/)       
406c       + le flux dans le sol a t
407c       + l'evolution du flux dans le sol lorsque la temperature d'interface
408c       passe de sa valeur a t a sa valeur a {t + \delta t}.
409c       ----------------------------------------------------
410
411      DO ig=1,ngrid
412         z1(ig)=pcapcal(ig)*ptsrf(ig)+cpp*zb(ig,1)*zc(ig,1)
413     s     +zdplanck(ig)*ptsrf(ig)+ pfluxsrf(ig)*ptimestep
414         z2(ig)= pcapcal(ig)+cpp*zb(ig,1)*(1.-zd(ig,1))+zdplanck(ig)
415         ztsrf2(ig)=z1(ig)/z2(ig)
416         pdtsrf(ig)=(ztsrf2(ig)-ptsrf(ig))/ptimestep
417
418c        Modif speciale CO2 condensation:
419c        tconds = 1./(bcond-acond*log(.0095*pplev(ig,1)))
420c        if ((callcond).and.
421c    &      ((co2ice(ig).ne.0).or.(ztsrf2(ig).lt.tconds)))then
422c           zh(ig,1)=zc(ig,1)+zd(ig,1)*tconds
423c        else
424            zh(ig,1)=zc(ig,1)+zd(ig,1)*ztsrf2(ig)
425c        end if
426      ENDDO
427
428c    ** et a partir de la temperature au sol on remonte
429c       -----------------------------------------------
430
431      DO ilay=2,nlay
432         DO ig=1,ngrid
433            hh = max( zh(ig,ilay-1) , zhcond(ig,ilay-1) ) ! modif co2cond
434            zh(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*hh
435         ENDDO
436      ENDDO
437
438
439c-----------------------------------------------------------------------
440c   TRACERS
441c   -------
442
443      if(tracer) then
444
445c     Using the wind modified by friction for lifting and  sublimation
446c     ----------------------------------------------------------------
447        DO ig=1,ngrid
448          zu2=zu(ig,1)*zu(ig,1)+zv(ig,1)*zv(ig,1)
449          zcdv(ig)=zcdv_true(ig)*sqrt(zu2)
450          zcdh(ig)=zcdh_true(ig)*sqrt(zu2)
451        ENDDO
452
453c       Calcul du flux vertical au bas de la premiere couche (dust) :
454c       -----------------------------------------------------------
455        do ig=1,ngridmx 
456          rho(ig) = zb0(ig,1) /ptimestep
457!          zb(ig,1) = 0.
458        end do
459
460        call zerophys(ngrid*nq,pdqsdif)
461
462c       OU calcul de la valeur de q a la surface (water)  :
463c       ----------------------------------------
464        if (water) then
465           !call watersat(ngridmx,ptsrf,pplev(1,1),qsat)
466           do ig=1,ngridmx 
467              call watersat_2(ptsrf(ig),pplev(ig,1),qsat(ig))
468           end do
469        end if
470
471c      Inversion pour l'implicite sur q
472c       --------------------------------
473        do iq=1,nq
474          CALL multipl((nlay-1)*ngrid,zkh(1,2),zb0(1,2),zb(1,2))
475
476          if ((water).and.(iq.eq.igcm_h2o_vap)) then
477c            This line is required to account for turbulent transport
478c            from surface (e.g. ice) to mid-layer of atmosphere:
479             CALL multipl(ngrid,zcdv,zb0,zb(1,1))
480             CALL multipl(ngrid,dryness,zb(1,1),zb(1,1))
481          else ! (re)-initialize zb(:,1)
482             zb(1:ngrid,1)=0
483          end if
484
485          DO ig=1,ngrid
486               z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
487               zc(ig,nlay)=za(ig,nlay)*zq(ig,nlay,iq)*z1(ig)
488               zd(ig,nlay)=zb(ig,nlay)*z1(ig)
489          ENDDO
490 
491          DO ilay=nlay-1,2,-1
492               DO ig=1,ngrid
493                z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
494     $           zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
495                zc(ig,ilay)=(za(ig,ilay)*zq(ig,ilay,iq)+
496     $           zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
497                zd(ig,ilay)=zb(ig,ilay)*z1(ig)
498               ENDDO
499          ENDDO
500
501          if (water.and.(iq.eq.igcm_h2o_ice)) then
502            ! special case for water ice tracer: do not include
503            ! h2o ice tracer from surface (which is set when handling
504            ! h2o vapour case (see further down).
505            DO ig=1,ngrid
506                z1(ig)=1./(za(ig,1)+zb(ig,1)+
507     $           zb(ig,2)*(1.-zd(ig,2)))
508                zc(ig,1)=(za(ig,1)*zq(ig,1,iq)+
509     $         zb(ig,2)*zc(ig,2)) *z1(ig)
510            ENDDO
511          else ! general case
512            DO ig=1,ngrid
513                z1(ig)=1./(za(ig,1)+zb(ig,1)+
514     $           zb(ig,2)*(1.-zd(ig,2)))
515                zc(ig,1)=(za(ig,1)*zq(ig,1,iq)+
516     $         zb(ig,2)*zc(ig,2) +
517     $        (-pdqsdif(ig,iq)) *ptimestep) *z1(ig)  !tracer flux from surface
518            ENDDO
519          endif ! of if (water.and.(iq.eq.igcm_h2o_ice))
520 
521          IF ((water).and.(iq.eq.igcm_h2o_vap)) then
522c           Calculation for turbulent exchange with the surface (for ice)
523            DO ig=1,ngrid
524              zd(ig,1)=zb(ig,1)*z1(ig)
525              zq1temp(ig)=zc(ig,1)+ zd(ig,1)*qsat(ig)
526
527              pdqsdif(ig,igcm_h2o_ice)=rho(ig)*dryness(ig)*zcdv(ig)
528     &                       *(zq1temp(ig)-qsat(ig))
529c             write(*,*)'flux vers le sol=',pdqsdif(ig,nq)
530            END DO
531
532            DO ig=1,ngrid
533              if(.not.watercaptag(ig)) then
534                if ((-pdqsdif(ig,igcm_h2o_ice)*ptimestep)
535     &             .gt.pqsurf(ig,igcm_h2o_ice)) then
536c                 write(*,*)'on sublime plus que qsurf!'
537                  pdqsdif(ig,igcm_h2o_ice)=
538     &                         -pqsurf(ig,igcm_h2o_ice)/ptimestep
539c                 write(*,*)'flux vers le sol=',pdqsdif(ig,nq)
540                  z1(ig)=1./(za(ig,1)+ zb(ig,2)*(1.-zd(ig,2)))
541                  zc(ig,1)=(za(ig,1)*zq(ig,1,igcm_h2o_vap)+
542     $            zb(ig,2)*zc(ig,2) +
543     $            (-pdqsdif(ig,igcm_h2o_ice)) *ptimestep) *z1(ig)
544                  zq1temp(ig)=zc(ig,1)
545                endif   
546              endif ! if (.not.watercaptag(ig))
547c             Starting upward calculations for water :
548               zq(ig,1,igcm_h2o_vap)=zq1temp(ig)
549            ENDDO ! of DO ig=1,ngrid
550
551
552!            ! ADDED BY RW: Water latent heat effect
553!     this doesnt work; causes instability. What do they do on earth?
554!            DO ig=1,ngrid
555!               pdtsrf(ig)=
556!     &             pdtsrf(ig)+RLVTT*pdqsdif(ig,igcm_h2o_ice)/pcapcal(ig)
557!            ENDDO
558!            print*,'Surface latent heat release in vdifc'
559!            print*,'due to H2O NOT taken into account.'
560
561          ELSE
562c           Starting upward calculations for simple mixing of tracer (dust)
563            DO ig=1,ngrid
564               zq(ig,1,iq)=zc(ig,1)
565            ENDDO
566          END IF ! of IF ((water).and.(iq.eq.igcm_h2o_vap))
567
568          DO ilay=2,nlay
569             DO ig=1,ngrid
570                zq(ig,ilay,iq)=zc(ig,ilay)+zd(ig,ilay)*zq(ig,ilay-1,iq)
571             ENDDO
572          ENDDO
573        enddo ! of do iq=1,nq
574      end if ! of if(tracer)
575
576c-----------------------------------------------------------------------
577c   8. calcul final des tendances de la diffusion verticale
578c      ----------------------------------------------------
579
580      DO ilev = 1, nlay
581         DO ig=1,ngrid
582            pdudif(ig,ilev)=(    zu(ig,ilev)-
583     $      (pu(ig,ilev)+pdufi(ig,ilev)*ptimestep)    )/ptimestep
584            pdvdif(ig,ilev)=(    zv(ig,ilev)-
585     $      (pv(ig,ilev)+pdvfi(ig,ilev)*ptimestep)    )/ptimestep
586            hh = max(ph(ig,ilev)+pdhfi(ig,ilev)*ptimestep ,
587     $           zhcond(ig,ilev))        ! modif co2cond
588            pdhdif(ig,ilev)=( zh(ig,ilev)- hh )/ptimestep
589         ENDDO
590      ENDDO
591
592   
593      if (tracer) then
594        DO iq = 1, nq
595          DO ilev = 1, nlay
596            DO ig=1,ngrid
597              pdqdif(ig,ilev,iq)=(zq(ig,ilev,iq)-
598     $      (pq(ig,ilev,iq) + pdqfi(ig,ilev,iq)*ptimestep))/ptimestep
599            ENDDO
600          ENDDO
601        ENDDO
602      end if
603
604c    ** diagnostique final
605c       ------------------
606
607      IF(lecrit) THEN
608         PRINT*,'In vdif'
609         PRINT*,'Ts (t) and Ts (t+st)'
610         WRITE(*,'(a10,3a15)')
611     s   'theta(t)','theta(t+dt)','u(t)','u(t+dt)'
612         PRINT*,ptsrf(ngrid/2+1),ztsrf2(ngrid/2+1)
613         DO ilev=1,nlay
614            WRITE(*,'(4f15.7)')
615     s      ph(ngrid/2+1,ilev),zh(ngrid/2+1,ilev),
616     s      pu(ngrid/2+1,ilev),zu(ngrid/2+1,ilev)
617
618         ENDDO
619      ENDIF
620
621      RETURN
622      END
Note: See TracBrowser for help on using the repository browser.