source: trunk/LMDZ.MARS/libf/phymars/vdifc.F @ 272

Last change on this file since 272 was 268, checked in by acolaitis, 14 years ago

--- AC 03/08/2011 ---
M 267 libf/phymars/physiq.F
<> Minor modification to pass Ch from vdifc to meso_inc_les

M 267 libf/phymars/surflayer_interpol.F
<> Major modification to the formulation of integrals

Now stable for most cases. Some cases with highly negative Monin Obukhov length
remain to be explored.

M 267 libf/phymars/vdif_cd.F
<> Added gustiness to the Richardson computation. Gustiness factor is for now of beta=1., after

several comparisons with LES aerodynamic conductances. May be subject to a minor change (+/- 0.1)
in the near future. (almost transparent for the user)

M 267 libf/phymars/vdifc.F
<> Minor modifications relative to variables.

M 267 libf/phymars/calltherm_mars.F90
<> Added a comment on a "sensitive" parameter that should not be changed without knowing the consequence !

M 267 libf/phymars/meso_inc/meso_inc_les.F
<> Changed the definition for HFX computation in the LES (to be discussed with Aymeric). New definition yields

very similar results to old one and follows a strict definition of what HFX should be.


File size: 21.9 KB
RevLine 
[38]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,
[268]7     $                pdqdif,pdqsdif,wmax,zcdv_true,zcdh_true)
[38]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
[224]59      REAL pz0(ngridmx) ! surface roughness length (m)
60
[256]61c    Argument added to account for subgrid gustiness :
62
63      REAL wmax(ngridmx)
64
[38]65c    Traceurs :
66      integer nq
67      REAL pqsurf(ngrid,nq)
68      real pq(ngrid,nlay,nq), pdqfi(ngrid,nlay,nq)
69      real pdqdif(ngrid,nlay,nq)
70      real pdqsdif(ngrid,nq)
71     
72c   local:
73c   ------
74
75      INTEGER ilev,ig,ilay,nlev
76
77      REAL z4st,zdplanck(ngridmx)
78      REAL zkv(ngridmx,nlayermx+1),zkh(ngridmx,nlayermx+1)
79      REAL zcdv(ngridmx),zcdh(ngridmx)
[268]80      REAL zcdv_true(ngridmx),zcdh_true(ngridmx)    ! drag coeff are used by the LES to recompute u* and hfx
[38]81      REAL zu(ngridmx,nlayermx),zv(ngridmx,nlayermx)
82      REAL zh(ngridmx,nlayermx)
83      REAL ztsrf2(ngridmx)
84      REAL z1(ngridmx),z2(ngridmx)
85      REAL za(ngridmx,nlayermx),zb(ngridmx,nlayermx)
86      REAL zb0(ngridmx,nlayermx)
87      REAL zc(ngridmx,nlayermx),zd(ngridmx,nlayermx)
88      REAL zcst1
89      REAL zu2
90
91      EXTERNAL SSUM,SCOPY
92      REAL SSUM
93      LOGICAL firstcall
94      SAVE firstcall
95
96c     variable added for CO2 condensation:
97c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
98      REAL hh , zhcond(ngridmx,nlayermx)
99      REAL latcond,tcond1mb
100      REAL acond,bcond
101      SAVE acond,bcond
102      DATA latcond,tcond1mb/5.9e5,136.27/
103
104c    Tracers :
105c    ~~~~~~~
106      INTEGER iq
107      REAL zq(ngridmx,nlayermx,nqmx)
108      REAL zq1temp(ngridmx)
109      REAL rho(ngridmx) ! near surface air density
110      REAL qsat(ngridmx)
111      DATA firstcall/.true./
112
113      REAL kmixmin
114
115c    ** un petit test de coherence
116c       --------------------------
117
118      IF (firstcall) THEN
119         IF(ngrid.NE.ngridmx) THEN
120            PRINT*,'STOP dans vdifc'
121            PRINT*,'probleme de dimensions :'
122            PRINT*,'ngrid  =',ngrid
123            PRINT*,'ngridmx  =',ngridmx
124            STOP
125         ENDIF
126c        To compute: Tcond= 1./(bcond-acond*log(.0095*p)) (p in pascal)
127         bcond=1./tcond1mb
128         acond=r/latcond
129         PRINT*,'In vdifc: Tcond(P=1mb)=',tcond1mb,' Lcond=',latcond
130         PRINT*,'          acond,bcond',acond,bcond
131
132        firstcall=.false.
133      ENDIF
134
135
136
137
138
139c-----------------------------------------------------------------------
140c    1. initialisation
141c    -----------------
142
143      nlev=nlay+1
144
145c    ** calcul de rho*dz et dt*rho/dz=dt*rho**2 g/dp
146c       avec rho=p/RT=p/ (R Theta) (p/ps)**kappa
147c       ----------------------------------------
148
149      DO ilay=1,nlay
150         DO ig=1,ngrid
151            za(ig,ilay)=(pplev(ig,ilay)-pplev(ig,ilay+1))/g
152         ENDDO
153      ENDDO
154
155      zcst1=4.*g*ptimestep/(r*r)
156      DO ilev=2,nlev-1
157         DO ig=1,ngrid
158            zb0(ig,ilev)=pplev(ig,ilev)*
159     s      (pplev(ig,1)/pplev(ig,ilev))**rcp /
160     s      (ph(ig,ilev-1)+ph(ig,ilev))
161            zb0(ig,ilev)=zcst1*zb0(ig,ilev)*zb0(ig,ilev)/
162     s      (pplay(ig,ilev-1)-pplay(ig,ilev))
163         ENDDO
164      ENDDO
165      DO ig=1,ngrid
166                 zb0(ig,1)=ptimestep*pplev(ig,1)/(r*ptsrf(ig))
167      ENDDO
168
169c    ** diagnostique pour l'initialisation
170c       ----------------------------------
171
172      IF(lecrit) THEN
173         ig=ngrid/2+1
174         PRINT*,'Pression (mbar) ,altitude (km),u,v,theta, rho dz'
175         DO ilay=1,nlay
176            WRITE(*,'(6f11.5)')
177     s      .01*pplay(ig,ilay),.001*pzlay(ig,ilay),
178     s      pu(ig,ilay),pv(ig,ilay),ph(ig,ilay),za(ig,ilay)
179         ENDDO
180         PRINT*,'Pression (mbar) ,altitude (km),zb'
181         DO ilev=1,nlay
182            WRITE(*,'(3f15.7)')
183     s      .01*pplev(ig,ilev),.001*pzlev(ig,ilev),
184     s      zb0(ig,ilev)
185         ENDDO
186      ENDIF
187
188c     Potential Condensation temperature:
189c     -----------------------------------
190
191c     if (callcond) then
192c       DO ilev=1,nlay
193c         DO ig=1,ngrid
194c           zhcond(ig,ilev) =
195c    &  (1./(bcond-acond*log(.0095*pplay(ig,ilev))))/ppopsk(ig,ilev)
196c         END DO
197c       END DO
198c     else
199        call zerophys(ngrid*nlay,zhcond)
200c     end if
201
202
203c-----------------------------------------------------------------------
204c   2. ajout des tendances physiques
205c      -----------------------------
206
207      DO ilev=1,nlay
208         DO ig=1,ngrid
209            zu(ig,ilev)=pu(ig,ilev)+pdufi(ig,ilev)*ptimestep
210            zv(ig,ilev)=pv(ig,ilev)+pdvfi(ig,ilev)*ptimestep
211            zh(ig,ilev)=ph(ig,ilev)+pdhfi(ig,ilev)*ptimestep
212            zh(ig,ilev)=max(zh(ig,ilev),zhcond(ig,ilev))
213         ENDDO
214      ENDDO
215      if(tracer) then
216        DO iq =1, nq
217         DO ilev=1,nlay
218           DO ig=1,ngrid
219              zq(ig,ilev,iq)=pq(ig,ilev,iq)+pdqfi(ig,ilev,iq)*ptimestep
220           ENDDO
221         ENDDO
222        ENDDO
223      end if
224
225c-----------------------------------------------------------------------
226c   3. schema de turbulence
227c      --------------------
228
229c    ** source d'energie cinetique turbulente a la surface
230c       (condition aux limites du schema de diffusion turbulente
231c       dans la couche limite
232c       ---------------------
233
[268]234      CALL vdif_cd(ngrid,nlay,pz0,g,pzlay,pu,pv,wmax,ptsrf,ph
[38]235     &             ,zcdv_true,zcdh_true)
[256]236
[38]237      DO ig=1,ngrid
[256]238
[38]239        zu2=pu(ig,1)*pu(ig,1)+pv(ig,1)*pv(ig,1)
[268]240     &             +(1.*wmax(ig))**2        ! subgrid gustiness added by quadratic interpolation with a coeff beta, here beta=1. (tuned from
241                                            ! LES comparisons. This parameter is linked to the thermals settings)
[267]242       
[268]243        zcdv(ig)=zcdv_true(ig)*sqrt(zu2)     ! 1 / bulk aerodynamic momentum conductance
244        zcdh(ig)=zcdh_true(ig)*sqrt(zu2)     ! 1 / bulk aerodynamic heat conductance
245                                             ! these are the quantities to be looked at when comparing surface layers of different models
[38]246      ENDDO
247
[256]248! Some usefull diagnostics for the new surface layer parametrization :
249
[268]250!         call WRITEDIAGFI(ngridmx,'pcdv',
[256]251!     &              'momentum drag','no units',
252!     &                         2,zcdv_true)
[268]253!        call WRITEDIAGFI(ngridmx,'pcdh',
[256]254!     &              'heat drag','no units',
255!     &                         2,zcdh_true)
[268]256!        call WRITEDIAGFI(ngridmx,'u*',
257!    &              'friction velocity','m/s',
258!     &                         2,sqrt(zcdv_true(:))*sqrt(zu2))
259!        call WRITEDIAGFI(ngridmx,'T*',
260!     &              'friction temperature','K',
261!     &          2,-zcdh_true(:)*(ptsrf(:)-ph(:,1))/sqrt(zcdv_true(:)))
262!        call WRITEDIAGFI(ngridmx,'rm-1',
263!     &              'aerodyn momentum conductance','m/s',
[256]264!     &                         2,zcdv)
[268]265!        call WRITEDIAGFI(ngridmx,'rh-1',
266!     &              'aerodyn heat conductance','m/s',
[256]267!     &                         2,zcdh)
268
[38]269c    ** schema de diffusion turbulente dans la couche limite
270c       ----------------------------------------------------
271
272        CALL vdif_kc(ptimestep,g,pzlev,pzlay
273     &              ,pu,pv,ph,zcdv_true
274     &              ,pq2,zkv,zkh)
275
276      if ((doubleq).and.(ngrid.eq.1)) then
277        kmixmin = 80. !80.! minimum eddy mix coeff in 1D
278        do ilev=1,nlay
279          do ig=1,ngrid
280           zkh(ig,ilev) = max(kmixmin,zkh(ig,ilev))
281           zkv(ig,ilev) = max(kmixmin,zkv(ig,ilev))
282          end do
283        end do
284      end if
285
286c    ** diagnostique pour le schema de turbulence
287c       -----------------------------------------
288
289      IF(lecrit) THEN
290         PRINT*
291         PRINT*,'Diagnostic for the vertical turbulent mixing'
292         PRINT*,'Cd for momentum and potential temperature'
293
294         PRINT*,zcdv(ngrid/2+1),zcdh(ngrid/2+1)
295         PRINT*,'Mixing coefficient for momentum and pot.temp.'
296         DO ilev=1,nlay
297            PRINT*,zkv(ngrid/2+1,ilev),zkh(ngrid/2+1,ilev)
298         ENDDO
299      ENDIF
300
301
302
303
304c-----------------------------------------------------------------------
305c   4. inversion pour l'implicite sur u
306c      --------------------------------
307
308c    ** l'equation est
309c       u(t+1) =  u(t) + dt * {(du/dt)phys}(t)  +  dt * {(du/dt)difv}(t+1)
310c       avec
311c       /zu/ = u(t) + dt * {(du/dt)phys}(t)   (voir paragraphe 2.)
312c       et
313c       dt * {(du/dt)difv}(t+1) = dt * {(d/dz)[ Ku (du/dz) ]}(t+1)
314c       donc les entrees sont /zcdv/ pour la condition a la limite sol
315c       et /zkv/ = Ku
316 
317      CALL multipl((nlay-1)*ngrid,zkv(1,2),zb0(1,2),zb(1,2))
318      CALL multipl(ngrid,zcdv,zb0,zb)
319
320      DO ig=1,ngrid
321         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
322         zc(ig,nlay)=za(ig,nlay)*zu(ig,nlay)*z1(ig)
323         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
324      ENDDO
325
326      DO ilay=nlay-1,1,-1
327         DO ig=1,ngrid
328            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
329     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
330            zc(ig,ilay)=(za(ig,ilay)*zu(ig,ilay)+
331     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
332            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
333         ENDDO
334      ENDDO
335
336      DO ig=1,ngrid
337         zu(ig,1)=zc(ig,1)
338      ENDDO
339      DO ilay=2,nlay
340         DO ig=1,ngrid
341            zu(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zu(ig,ilay-1)
342         ENDDO
343      ENDDO
344
345
346
347
348
349c-----------------------------------------------------------------------
350c   5. inversion pour l'implicite sur v
351c      --------------------------------
352
353c    ** l'equation est
354c       v(t+1) =  v(t) + dt * {(dv/dt)phys}(t)  +  dt * {(dv/dt)difv}(t+1)
355c       avec
356c       /zv/ = v(t) + dt * {(dv/dt)phys}(t)   (voir paragraphe 2.)
357c       et
358c       dt * {(dv/dt)difv}(t+1) = dt * {(d/dz)[ Kv (dv/dz) ]}(t+1)
359c       donc les entrees sont /zcdv/ pour la condition a la limite sol
360c       et /zkv/ = Kv
361
362      DO ig=1,ngrid
363         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
364         zc(ig,nlay)=za(ig,nlay)*zv(ig,nlay)*z1(ig)
365         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
366      ENDDO
367
368      DO ilay=nlay-1,1,-1
369         DO ig=1,ngrid
370            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
371     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
372            zc(ig,ilay)=(za(ig,ilay)*zv(ig,ilay)+
373     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
374            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
375         ENDDO
376      ENDDO
377
378      DO ig=1,ngrid
379         zv(ig,1)=zc(ig,1)
380      ENDDO
381      DO ilay=2,nlay
382         DO ig=1,ngrid
383            zv(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zv(ig,ilay-1)
384         ENDDO
385      ENDDO
386
387
388
389
390
391c-----------------------------------------------------------------------
392c   6. inversion pour l'implicite sur h sans oublier le couplage
393c      avec le sol (conduction)
394c      ------------------------
395
396c    ** l'equation est
397c       h(t+1) =  h(t) + dt * {(dh/dt)phys}(t)  +  dt * {(dh/dt)difv}(t+1)
398c       avec
399c       /zh/ = h(t) + dt * {(dh/dt)phys}(t)   (voir paragraphe 2.)
400c       et
401c       dt * {(dh/dt)difv}(t+1) = dt * {(d/dz)[ Kh (dh/dz) ]}(t+1)
402c       donc les entrees sont /zcdh/ pour la condition de raccord au sol
403c       et /zkh/ = Kh
404c       -------------
405
406      CALL multipl((nlay-1)*ngrid,zkh(1,2),zb0(1,2),zb(1,2))
407      CALL multipl(ngrid,zcdh,zb0,zb)
408
409      DO ig=1,ngrid
410         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
411         zc(ig,nlay)=za(ig,nlay)*zh(ig,nlay)*z1(ig)
412         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
413      ENDDO
414
415      DO ilay=nlay-1,1,-1
416         DO ig=1,ngrid
417            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
418     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
419            zc(ig,ilay)=(za(ig,ilay)*zh(ig,ilay)+
420     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
421            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
422         ENDDO
423      ENDDO
424
425c    ** calcul de (d Planck / dT) a la temperature d'interface
426c       ------------------------------------------------------
427
428      z4st=4.*5.67e-8*ptimestep
429      DO ig=1,ngrid
430         zdplanck(ig)=z4st*pemis(ig)*ptsrf(ig)*ptsrf(ig)*ptsrf(ig)
431      ENDDO
432
433c    ** calcul de la temperature_d'interface et de sa tendance.
434c       on ecrit que la somme des flux est nulle a l'interface
435c       a t + \delta t,
436c       c'est a dire le flux radiatif a {t + \delta t}
437c       + le flux turbulent a {t + \delta t}
438c            qui s'ecrit K (T1-Tsurf) avec T1 = d1 Tsurf + c1
439c            (notation K dt = /cpp*b/)       
440c       + le flux dans le sol a t
441c       + l'evolution du flux dans le sol lorsque la temperature d'interface
442c       passe de sa valeur a t a sa valeur a {t + \delta t}.
443c       ----------------------------------------------------
444
445      DO ig=1,ngrid
446         z1(ig)=pcapcal(ig)*ptsrf(ig)+cpp*zb(ig,1)*zc(ig,1)
447     s     +zdplanck(ig)*ptsrf(ig)+ pfluxsrf(ig)*ptimestep
448         z2(ig)= pcapcal(ig)+cpp*zb(ig,1)*(1.-zd(ig,1))+zdplanck(ig)
449         ztsrf2(ig)=z1(ig)/z2(ig)
450         pdtsrf(ig)=(ztsrf2(ig)-ptsrf(ig))/ptimestep
451
452c        Modif speciale CO2 condensation:
453c        tconds = 1./(bcond-acond*log(.0095*pplev(ig,1)))
454c        if ((callcond).and.
455c    &      ((co2ice(ig).ne.0).or.(ztsrf2(ig).lt.tconds)))then
456c           zh(ig,1)=zc(ig,1)+zd(ig,1)*tconds
457c        else
458            zh(ig,1)=zc(ig,1)+zd(ig,1)*ztsrf2(ig)
459c        end if
460      ENDDO
461
462c    ** et a partir de la temperature au sol on remonte
463c       -----------------------------------------------
464
465      DO ilay=2,nlay
466         DO ig=1,ngrid
467            hh = max( zh(ig,ilay-1) , zhcond(ig,ilay-1) ) ! modif co2cond
468            zh(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*hh
469         ENDDO
470      ENDDO
471
472
473c-----------------------------------------------------------------------
474c   TRACERS
475c   -------
476
477      if(tracer) then
478
479c     Using the wind modified by friction for lifting and  sublimation
480c     ----------------------------------------------------------------
481        DO ig=1,ngrid
482          zu2=zu(ig,1)*zu(ig,1)+zv(ig,1)*zv(ig,1)
483          zcdv(ig)=zcdv_true(ig)*sqrt(zu2)
484          zcdh(ig)=zcdh_true(ig)*sqrt(zu2)
485        ENDDO
486
487c       Calcul du flux vertical au bas de la premiere couche (dust) :
488c       -----------------------------------------------------------
489        do ig=1,ngridmx 
490          rho(ig) = zb0(ig,1) /ptimestep
491c          zb(ig,1) = 0.
492        end do
493c       Dust lifting:
494        if (lifting) then
495           if (doubleq.AND.submicron) then
496             do ig=1,ngrid
497c              if(co2ice(ig).lt.1) then
498                 pdqsdif(ig,igcm_dust_mass) =
499     &             -alpha_lift(igcm_dust_mass) 
500                 pdqsdif(ig,igcm_dust_number) =
501     &             -alpha_lift(igcm_dust_number) 
502                 pdqsdif(ig,igcm_dust_submicron) =
503     &             -alpha_lift(igcm_dust_submicron)
504c              end if
505             end do
506           else if (doubleq) then
507             do ig=1,ngrid
[77]508           !!! soulevement constant
[38]509                 pdqsdif(ig,igcm_dust_mass) =
510     &             -alpha_lift(igcm_dust_mass) 
511                 pdqsdif(ig,igcm_dust_number) =
512     &             -alpha_lift(igcm_dust_number) 
513             end do
514           else if (submicron) then
515             do ig=1,ngrid
516                 pdqsdif(ig,igcm_dust_submicron) =
517     &             -alpha_lift(igcm_dust_submicron)
518             end do
519           else
520            call dustlift(ngrid,nlay,nq,rho,zcdh_true,zcdh,co2ice,
521     &                    pdqsdif)
522           endif !doubleq.AND.submicron
523        else
524           pdqsdif(1:ngrid,1:nq) = 0.
525        end if
526
527c       OU calcul de la valeur de q a la surface (water)  :
528c       ----------------------------------------
529        if (water) then
530            call watersat(ngridmx,ptsrf,pplev(1,1),qsat)
531        end if
532
533c      Inversion pour l'implicite sur q
534c       --------------------------------
535        do iq=1,nq
536          CALL multipl((nlay-1)*ngrid,zkh(1,2),zb0(1,2),zb(1,2))
537
538          if ((water).and.(iq.eq.igcm_h2o_vap)) then
539c            This line is required to account for turbulent transport
540c            from surface (e.g. ice) to mid-layer of atmosphere:
541             CALL multipl(ngrid,zcdv,zb0,zb(1,1))
542             CALL multipl(ngrid,dryness,zb(1,1),zb(1,1))
543          else ! (re)-initialize zb(:,1)
544             zb(1:ngrid,1)=0
545          end if
546
547          DO ig=1,ngrid
548               z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
549               zc(ig,nlay)=za(ig,nlay)*zq(ig,nlay,iq)*z1(ig)
550               zd(ig,nlay)=zb(ig,nlay)*z1(ig)
551          ENDDO
552 
553          DO ilay=nlay-1,2,-1
554               DO ig=1,ngrid
555                z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
556     $           zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
557                zc(ig,ilay)=(za(ig,ilay)*zq(ig,ilay,iq)+
558     $           zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
559                zd(ig,ilay)=zb(ig,ilay)*z1(ig)
560               ENDDO
561          ENDDO
562
563          if (water.and.(iq.eq.igcm_h2o_ice)) then
564            ! special case for water ice tracer: do not include
565            ! h2o ice tracer from surface (which is set when handling
566            ! h2o vapour case (see further down).
567            DO ig=1,ngrid
568                z1(ig)=1./(za(ig,1)+zb(ig,1)+
569     $           zb(ig,2)*(1.-zd(ig,2)))
570                zc(ig,1)=(za(ig,1)*zq(ig,1,iq)+
571     $         zb(ig,2)*zc(ig,2)) *z1(ig)
572            ENDDO
573          else ! general case
574            DO ig=1,ngrid
575                z1(ig)=1./(za(ig,1)+zb(ig,1)+
576     $           zb(ig,2)*(1.-zd(ig,2)))
577                zc(ig,1)=(za(ig,1)*zq(ig,1,iq)+
578     $         zb(ig,2)*zc(ig,2) +
579     $        (-pdqsdif(ig,iq)) *ptimestep) *z1(ig)  !tracer flux from surface
580            ENDDO
581          endif ! of if (water.and.(iq.eq.igcm_h2o_ice))
582 
583          IF ((water).and.(iq.eq.igcm_h2o_vap)) then
584c           Calculation for turbulent exchange with the surface (for ice)
585            DO ig=1,ngrid
586              zd(ig,1)=zb(ig,1)*z1(ig)
587              zq1temp(ig)=zc(ig,1)+ zd(ig,1)*qsat(ig)
588
589              pdqsdif(ig,igcm_h2o_ice)=rho(ig)*dryness(ig)*zcdv(ig)
590     &                       *(zq1temp(ig)-qsat(ig))
591c             write(*,*)'flux vers le sol=',pdqsdif(ig,nq)
592            END DO
593
594            DO ig=1,ngrid
595              if(.not.watercaptag(ig)) then
596                if ((-pdqsdif(ig,igcm_h2o_ice)*ptimestep)
597     &             .gt.pqsurf(ig,igcm_h2o_ice)) then
598c                 write(*,*)'on sublime plus que qsurf!'
599                  pdqsdif(ig,igcm_h2o_ice)=
600     &                         -pqsurf(ig,igcm_h2o_ice)/ptimestep
601c                 write(*,*)'flux vers le sol=',pdqsdif(ig,nq)
602                  z1(ig)=1./(za(ig,1)+ zb(ig,2)*(1.-zd(ig,2)))
603                  zc(ig,1)=(za(ig,1)*zq(ig,1,igcm_h2o_vap)+
604     $            zb(ig,2)*zc(ig,2) +
605     $            (-pdqsdif(ig,igcm_h2o_ice)) *ptimestep) *z1(ig)
606                  zq1temp(ig)=zc(ig,1)
607                endif   
608              endif ! if (.not.watercaptag(ig))
609c             Starting upward calculations for water :
610               zq(ig,1,igcm_h2o_vap)=zq1temp(ig)
611            ENDDO ! of DO ig=1,ngrid
612          ELSE
613c           Starting upward calculations for simple mixing of tracer (dust)
614            DO ig=1,ngrid
615               zq(ig,1,iq)=zc(ig,1)
616            ENDDO
617          END IF ! of IF ((water).and.(iq.eq.igcm_h2o_vap))
618
619          DO ilay=2,nlay
620             DO ig=1,ngrid
621                zq(ig,ilay,iq)=zc(ig,ilay)+zd(ig,ilay)*zq(ig,ilay-1,iq)
622             ENDDO
623          ENDDO
624        enddo ! of do iq=1,nq
625      end if ! of if(tracer)
626
627c-----------------------------------------------------------------------
628c   8. calcul final des tendances de la diffusion verticale
629c      ----------------------------------------------------
630
631      DO ilev = 1, nlay
632         DO ig=1,ngrid
633            pdudif(ig,ilev)=(    zu(ig,ilev)-
634     $      (pu(ig,ilev)+pdufi(ig,ilev)*ptimestep)    )/ptimestep
635            pdvdif(ig,ilev)=(    zv(ig,ilev)-
636     $      (pv(ig,ilev)+pdvfi(ig,ilev)*ptimestep)    )/ptimestep
637            hh = max(ph(ig,ilev)+pdhfi(ig,ilev)*ptimestep ,
638     $           zhcond(ig,ilev))        ! modif co2cond
639            pdhdif(ig,ilev)=( zh(ig,ilev)- hh )/ptimestep
640         ENDDO
641      ENDDO
642
643   
644      if (tracer) then
645        DO iq = 1, nq
646          DO ilev = 1, nlay
647            DO ig=1,ngrid
648              pdqdif(ig,ilev,iq)=(zq(ig,ilev,iq)-
649     $      (pq(ig,ilev,iq) + pdqfi(ig,ilev,iq)*ptimestep))/ptimestep
650            ENDDO
651          ENDDO
652        ENDDO
653      end if
654
655c    ** diagnostique final
656c       ------------------
657
658      IF(lecrit) THEN
659         PRINT*,'In vdif'
660         PRINT*,'Ts (t) and Ts (t+st)'
661         WRITE(*,'(a10,3a15)')
662     s   'theta(t)','theta(t+dt)','u(t)','u(t+dt)'
663         PRINT*,ptsrf(ngrid/2+1),ztsrf2(ngrid/2+1)
664         DO ilev=1,nlay
665            WRITE(*,'(4f15.7)')
666     s      ph(ngrid/2+1,ilev),zh(ngrid/2+1,ilev),
667     s      pu(ngrid/2+1,ilev),zu(ngrid/2+1,ilev)
668
669         ENDDO
670      ENDIF
671
[161]672      if (calltherm .and. outptherm) then
673      if (ngrid .eq. 1) then
674        call WRITEDIAGFI(ngrid,'zh','zh inside vdifc',
675     &                       'SI',1,ph(:,:)+pdhfi(:,:)*ptimestep)
676        call WRITEDIAGFI(ngrid,'zkh','zkh',
677     &                       'SI',1,zkh)
678      endif
679      endif
680
[38]681      RETURN
682      END
Note: See TracBrowser for help on using the repository browser.