source: trunk/WRF.COMMON/WRFV2/mars_lmd/libf/phymars/vdifc.F @ 3026

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

spiga@svn-planeto:ajoute le modele meso-echelle martien

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