source: trunk/LMDZ.MARS/libf/phymars/vdifc_mod.F @ 2156

Last change on this file since 2156 was 2080, checked in by mvals, 6 years ago

Mars GCM:

  • typo in vdifc.F (visible in debug mode) in the case of rdstorm=false, condition for defining pdqsdif(stormdust_mass/number)
  • typo in physiq_mod.F (visible in debug mode): wrong indices for tracers in aerosol(:,:,:)

MV

File size: 30.8 KB
RevLine 
[1969]1      MODULE vdifc_mod
2     
3      IMPLICIT NONE
4     
5      CONTAINS
6     
[38]7      SUBROUTINE vdifc(ngrid,nlay,nq,co2ice,ppopsk,
8     $                ptimestep,pcapcal,lecrit,
9     $                pplay,pplev,pzlay,pzlev,pz0,
10     $                pu,pv,ph,pq,ptsrf,pemis,pqsurf,
11     $                pdufi,pdvfi,pdhfi,pdqfi,pfluxsrf,
12     $                pdudif,pdvdif,pdhdif,pdtsrf,pq2,
[660]13     $                pdqdif,pdqsdif,wstar,zcdv_true,zcdh_true,
[1996]14     $                hfmax,pcondicea_co2microp,sensibFlux,
15     $                dustliftday,local_time)
[1236]16
[1036]17      use tracer_mod, only: noms, igcm_dust_mass, igcm_dust_number,
18     &                      igcm_dust_submicron, igcm_h2o_vap,
[1974]19     &                      igcm_h2o_ice, alpha_lift,
20     &                      igcm_stormdust_mass, igcm_stormdust_number
[1224]21      use surfdat_h, only: watercaptag, frost_albedo_threshold, dryness
[1969]22      USE comcstfi_h, ONLY: cpp, r, rcp, g
[1996]23      use watersat_mod, only: watersat
[1242]24      use turb_mod, only: turb_resolved, ustar, tstar
[1974]25      use compute_dtau_mod, only: ti_injection,tf_injection
[38]26      IMPLICIT NONE
27
28c=======================================================================
29c
30c   subject:
31c   --------
32c   Turbulent diffusion (mixing) for potential T, U, V and tracer
33c
34c   Shema implicite
35c   On commence par rajouter au variables x la tendance physique
36c   et on resoult en fait:
37c      x(t+1) =  x(t) + dt * (dx/dt)phys(t)  +  dt * (dx/dt)difv(t+1)
38c
39c   author:
40c   ------
41c      Hourdin/Forget/Fournier
42c=======================================================================
43
44c-----------------------------------------------------------------------
45c   declarations:
46c   -------------
47
[1944]48      include "callkeys.h"
49      include "microphys.h"
[38]50
51c
52c   arguments:
53c   ----------
54
[1036]55      INTEGER,INTENT(IN) :: ngrid,nlay
56      REAL,INTENT(IN) :: ptimestep
57      REAL,INTENT(IN) :: pplay(ngrid,nlay),pplev(ngrid,nlay+1)
58      REAL,INTENT(IN) :: pzlay(ngrid,nlay),pzlev(ngrid,nlay+1)
59      REAL,INTENT(IN) :: pu(ngrid,nlay),pv(ngrid,nlay)
60      REAL,INTENT(IN) :: ph(ngrid,nlay)
61      REAL,INTENT(IN) :: ptsrf(ngrid),pemis(ngrid)
62      REAL,INTENT(IN) :: pdufi(ngrid,nlay),pdvfi(ngrid,nlay)
63      REAL,INTENT(IN) :: pdhfi(ngrid,nlay)
64      REAL,INTENT(IN) :: pfluxsrf(ngrid)
65      REAL,INTENT(OUT) :: pdudif(ngrid,nlay),pdvdif(ngrid,nlay)
66      REAL,INTENT(OUT) :: pdtsrf(ngrid),pdhdif(ngrid,nlay)
[1130]67      REAL,INTENT(IN) :: pcapcal(ngrid)
68      REAL,INTENT(INOUT) :: pq2(ngrid,nlay+1)
[38]69
70c    Argument added for condensation:
[1036]71      REAL,INTENT(IN) :: co2ice (ngrid), ppopsk(ngrid,nlay)
72      logical,INTENT(IN) :: lecrit
[1996]73      REAL,INTENT(IN) :: pcondicea_co2microp(ngrid,nlay)! tendency due to CO2 condensation (kg/kg.s-1)
74     
[1036]75      REAL,INTENT(IN) :: pz0(ngrid) ! surface roughness length (m)
[224]76
[256]77c    Argument added to account for subgrid gustiness :
78
[1944]79      REAL,INTENT(IN) :: wstar(ngrid), hfmax(ngrid)!, zi(ngrid)
[256]80
[38]81c    Traceurs :
[1036]82      integer,intent(in) :: nq
83      REAL,INTENT(IN) :: pqsurf(ngrid,nq)
84      real,intent(in) :: pq(ngrid,nlay,nq), pdqfi(ngrid,nlay,nq)
85      real,intent(out) :: pdqdif(ngrid,nlay,nq)
[1974]86      real,intent(out) :: pdqsdif(ngrid,nq)
87      REAL,INTENT(in) :: dustliftday(ngrid)
88      REAL,INTENT(in) :: local_time(ngrid)
[38]89     
90c   local:
91c   ------
92
[1130]93      REAL :: pt(ngrid,nlay)
[473]94 
[38]95      INTEGER ilev,ig,ilay,nlev
96
[1047]97      REAL z4st,zdplanck(ngrid)
98      REAL zkv(ngrid,nlay+1),zkh(ngrid,nlay+1)
99      REAL zkq(ngrid,nlay+1)
100      REAL zcdv(ngrid),zcdh(ngrid)
101      REAL zcdv_true(ngrid),zcdh_true(ngrid)    ! drag coeff are used by the LES to recompute u* and hfx
102      REAL zu(ngrid,nlay),zv(ngrid,nlay)
103      REAL zh(ngrid,nlay)
104      REAL ztsrf2(ngrid)
105      REAL z1(ngrid),z2(ngrid)
106      REAL za(ngrid,nlay),zb(ngrid,nlay)
107      REAL zb0(ngrid,nlay)
108      REAL zc(ngrid,nlay),zd(ngrid,nlay)
[38]109      REAL zcst1
[1047]110      REAL zu2(ngrid)
[38]111
112      EXTERNAL SSUM,SCOPY
113      REAL SSUM
[1036]114      LOGICAL,SAVE :: firstcall=.true.
[38]115
[626]116
[38]117c     variable added for CO2 condensation:
118c     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[1047]119      REAL hh , zhcond(ngrid,nlay)
120      REAL,PARAMETER :: latcond=5.9e5
121      REAL,PARAMETER :: tcond1mb=136.27
122      REAL,SAVE :: acond,bcond
[669]123     
124c     For latent heat release from ground ice sublimation   
[1047]125!      REAL tsrf_lw(ngrid)
[1036]126!      REAL alpha
[669]127      REAL T1,T2
128      SAVE T1,T2
[734]129      DATA T1,T2/-604.3,1080.7/ ! zeros of latent heat equation for ice
[38]130
131c    Tracers :
132c    ~~~~~~~
133      INTEGER iq
[1047]134      REAL zq(ngrid,nlay,nq)
135      REAL zq1temp(ngrid)
136      REAL rho(ngrid) ! near surface air density
137      REAL qsat(ngrid)
[38]138
139      REAL kmixmin
140
[473]141c    Mass-variation scheme :
142c    ~~~~~~~
143
144      INTEGER j,l
[1047]145      REAL zcondicea(ngrid,nlay)
146      REAL zt(ngrid,nlay),ztcond(ngrid,nlay+1)
147      REAL betam(ngrid,nlay),dmice(ngrid,nlay)
148      REAL pdtc(ngrid,nlay)
149      REAL zhs(ngrid,nlay)
150      REAL,SAVE :: ccond
[473]151
152c     Theta_m formulation for mass-variation scheme :
153c     ~~~~~~~
154
[1047]155      INTEGER,SAVE :: ico2
156      INTEGER llnt(ngrid)
157      REAL,SAVE :: m_co2, m_noco2, A , B
158      REAL vmr_co2(ngrid,nlay)
[473]159      REAL qco2,mmean
160
[1047]161      REAL,INTENT(OUT) :: sensibFlux(ngrid)
[660]162
[38]163c    ** un petit test de coherence
164c       --------------------------
165
[1779]166      ! AS: OK firstcall absolute
[38]167      IF (firstcall) THEN
168c        To compute: Tcond= 1./(bcond-acond*log(.0095*p)) (p in pascal)
169         bcond=1./tcond1mb
170         acond=r/latcond
[473]171         ccond=cpp/(g*latcond)
[38]172         PRINT*,'In vdifc: Tcond(P=1mb)=',tcond1mb,' Lcond=',latcond
[473]173         PRINT*,'          acond,bcond,ccond',acond,bcond,ccond
[38]174
[473]175
176         ico2=0
177
178         if (tracer) then
179c          Prepare Special treatment if one of the tracer is CO2 gas
[1036]180           do iq=1,nq
[473]181             if (noms(iq).eq."co2") then
182                ico2=iq
183                m_co2 = 44.01E-3  ! CO2 molecular mass (kg/mol)
184                m_noco2 = 33.37E-3  ! Non condensible mol mass (kg/mol)
185c               Compute A and B coefficient use to compute
186c               mean molecular mass Mair defined by
187c               1/Mair = q(ico2)/m_co2 + (1-q(ico2))/m_noco2
188c               1/Mair = A*q(ico2) + B
189                A =(1/m_co2 - 1/m_noco2)
190                B=1/m_noco2
191             endif
192           enddo
193         end if
194
[38]195        firstcall=.false.
196      ENDIF
197
198
199c-----------------------------------------------------------------------
200c    1. initialisation
201c    -----------------
202
203      nlev=nlay+1
204
[1035]205      ! initialize output tendencies to zero:
206      pdudif(1:ngrid,1:nlay)=0
207      pdvdif(1:ngrid,1:nlay)=0
208      pdhdif(1:ngrid,1:nlay)=0
209      pdtsrf(1:ngrid)=0
210      pdqdif(1:ngrid,1:nlay,1:nq)=0
211      pdqsdif(1:ngrid,1:nq)=0
212
[38]213c    ** calcul de rho*dz et dt*rho/dz=dt*rho**2 g/dp
214c       avec rho=p/RT=p/ (R Theta) (p/ps)**kappa
215c       ----------------------------------------
216
217      DO ilay=1,nlay
218         DO ig=1,ngrid
219            za(ig,ilay)=(pplev(ig,ilay)-pplev(ig,ilay+1))/g
[473]220! Mass variation scheme:
221            betam(ig,ilay)=-za(ig,ilay)*latcond/(cpp*ppopsk(ig,ilay))
[38]222         ENDDO
223      ENDDO
224
225      zcst1=4.*g*ptimestep/(r*r)
226      DO ilev=2,nlev-1
227         DO ig=1,ngrid
228            zb0(ig,ilev)=pplev(ig,ilev)*
229     s      (pplev(ig,1)/pplev(ig,ilev))**rcp /
230     s      (ph(ig,ilev-1)+ph(ig,ilev))
231            zb0(ig,ilev)=zcst1*zb0(ig,ilev)*zb0(ig,ilev)/
232     s      (pplay(ig,ilev-1)-pplay(ig,ilev))
233         ENDDO
234      ENDDO
235      DO ig=1,ngrid
[473]236         zb0(ig,1)=ptimestep*pplev(ig,1)/(r*ptsrf(ig))
[38]237      ENDDO
238
239c    ** diagnostique pour l'initialisation
240c       ----------------------------------
241
242      IF(lecrit) THEN
243         ig=ngrid/2+1
244         PRINT*,'Pression (mbar) ,altitude (km),u,v,theta, rho dz'
245         DO ilay=1,nlay
246            WRITE(*,'(6f11.5)')
247     s      .01*pplay(ig,ilay),.001*pzlay(ig,ilay),
248     s      pu(ig,ilay),pv(ig,ilay),ph(ig,ilay),za(ig,ilay)
249         ENDDO
250         PRINT*,'Pression (mbar) ,altitude (km),zb'
251         DO ilev=1,nlay
252            WRITE(*,'(3f15.7)')
253     s      .01*pplev(ig,ilev),.001*pzlev(ig,ilev),
254     s      zb0(ig,ilev)
255         ENDDO
256      ENDIF
257
[473]258c     -----------------------------------
[38]259c     Potential Condensation temperature:
260c     -----------------------------------
261
[473]262c     Compute CO2 Volume mixing ratio
263c     -------------------------------
264      if (callcond.and.(ico2.ne.0)) then
265         DO ilev=1,nlay
266            DO ig=1,ngrid
[529]267              qco2=MAX(1.E-30
268     &             ,pq(ig,ilev,ico2)+pdqfi(ig,ilev,ico2)*ptimestep)
[473]269c             Mean air molecular mass = 1/(q(ico2)/m_co2 + (1-q(ico2))/m_noco2)
270              mmean=1/(A*qco2 +B)
271              vmr_co2(ig,ilev) = qco2*mmean/m_co2
272            ENDDO
273         ENDDO
274      else
275         DO ilev=1,nlay
276            DO ig=1,ngrid
277              vmr_co2(ig,ilev)=0.95
278            ENDDO
279         ENDDO
280      end if
[38]281
[473]282c     forecast of atmospheric temperature zt and frost temperature ztcond
283c     --------------------------------------------------------------------
[38]284
[473]285      if (callcond) then
286        DO ilev=1,nlay
287          DO ig=1,ngrid
[884]288              ztcond(ig,ilev)=
289     &      1./(bcond-acond*log(.01*vmr_co2(ig,ilev)*pplay(ig,ilev)))
290            if (pplay(ig,ilev).lt.1e-4) ztcond(ig,ilev)=0.0 !mars Monica
[473]291!            zhcond(ig,ilev) =
292!     &  (1./(bcond-acond*log(.0095*pplay(ig,ilev))))/ppopsk(ig,ilev)
293              zhcond(ig,ilev) = ztcond(ig,ilev)/ppopsk(ig,ilev)
294          END DO
295        END DO
[884]296        ztcond(:,nlay+1)=ztcond(:,nlay)
[473]297      else
[884]298         zhcond(:,:) = 0
299         ztcond(:,:) = 0
[473]300      end if
301
302
[38]303c-----------------------------------------------------------------------
304c   2. ajout des tendances physiques
305c      -----------------------------
306
307      DO ilev=1,nlay
308         DO ig=1,ngrid
309            zu(ig,ilev)=pu(ig,ilev)+pdufi(ig,ilev)*ptimestep
310            zv(ig,ilev)=pv(ig,ilev)+pdvfi(ig,ilev)*ptimestep
311            zh(ig,ilev)=ph(ig,ilev)+pdhfi(ig,ilev)*ptimestep
[473]312!            zh(ig,ilev)=max(zh(ig,ilev),zhcond(ig,ilev))
[38]313         ENDDO
314      ENDDO
315      if(tracer) then
316        DO iq =1, nq
317         DO ilev=1,nlay
318           DO ig=1,ngrid
319              zq(ig,ilev,iq)=pq(ig,ilev,iq)+pdqfi(ig,ilev,iq)*ptimestep
320           ENDDO
321         ENDDO
322        ENDDO
323      end if
324
325c-----------------------------------------------------------------------
326c   3. schema de turbulence
327c      --------------------
328
329c    ** source d'energie cinetique turbulente a la surface
330c       (condition aux limites du schema de diffusion turbulente
331c       dans la couche limite
332c       ---------------------
333
[499]334      CALL vdif_cd(ngrid,nlay,pz0,g,pzlay,pu,pv,wstar,ptsrf,ph
[1236]335     &             ,zcdv_true,zcdh_true)
[256]336
[290]337        zu2(:)=pu(:,1)*pu(:,1)+pv(:,1)*pv(:,1)
[256]338
[291]339        IF (callrichsl) THEN
[545]340          zcdv(:)=zcdv_true(:)*sqrt(zu2(:)+
341     &     (log(1.+0.7*wstar(:) + 2.3*wstar(:)**2))**2)
342          zcdh(:)=zcdh_true(:)*sqrt(zu2(:)+
343     &     (log(1.+0.7*wstar(:) + 2.3*wstar(:)**2))**2)
[496]344
[1242]345           ustar(:)=sqrt(zcdv_true(:))*sqrt(zu2(:)+
[545]346     &     (log(1.+0.7*wstar(:) + 2.3*wstar(:)**2))**2)
347
[1242]348           tstar(:)=0.
[545]349           DO ig=1,ngrid
350              IF (zcdh_true(ig) .ne. 0.) THEN      ! When Cd=Ch=0, u*=t*=0
[1242]351                 tstar(ig)=(ph(ig,1)-ptsrf(ig))*zcdh(ig)/ustar(ig)
[545]352              ENDIF
353           ENDDO
354
[284]355        ELSE
[545]356           zcdv(:)=zcdv_true(:)*sqrt(zu2(:))     ! 1 / bulk aerodynamic momentum conductance
357           zcdh(:)=zcdh_true(:)*sqrt(zu2(:))     ! 1 / bulk aerodynamic heat conductance
[1242]358           ustar(:)=sqrt(zcdv_true(:))*sqrt(zu2(:))
359           tstar(:)=(ph(:,1)-ptsrf(:))*zcdh_true(:)/sqrt(zcdv_true(:))
[290]360        ENDIF
[38]361
[529]362! Some usefull diagnostics for the new surface layer parametrization :
[290]363
[1130]364!        call WRITEDIAGFI(ngrid,'vdifc_zcdv_true',
[256]365!     &              'momentum drag','no units',
366!     &                         2,zcdv_true)
[1130]367!        call WRITEDIAGFI(ngrid,'vdifc_zcdh_true',
[256]368!     &              'heat drag','no units',
369!     &                         2,zcdh_true)
[1130]370!        call WRITEDIAGFI(ngrid,'vdifc_ust',
[473]371!     &              'friction velocity','m/s',2,ust)
[1130]372!       call WRITEDIAGFI(ngrid,'vdifc_tst',
[473]373!     &              'friction temperature','K',2,tst)
[1130]374!        call WRITEDIAGFI(ngrid,'vdifc_zcdv',
[268]375!     &              'aerodyn momentum conductance','m/s',
[256]376!     &                         2,zcdv)
[1130]377!        call WRITEDIAGFI(ngrid,'vdifc_zcdh',
[268]378!     &              'aerodyn heat conductance','m/s',
[256]379!     &                         2,zcdh)
380
[38]381c    ** schema de diffusion turbulente dans la couche limite
382c       ----------------------------------------------------
[1240]383       IF (.not. callyamada4) THEN
[529]384
[1130]385       CALL vdif_kc(ngrid,nlay,nq,ptimestep,g,pzlev,pzlay
[38]386     &              ,pu,pv,ph,zcdv_true
[325]387     &              ,pq2,zkv,zkh,zq)
[529]388
[544]389      ELSE
[38]390
[555]391      pt(:,:)=ph(:,:)*ppopsk(:,:)
[1036]392      CALL yamada4(ngrid,nlay,nq,ptimestep,g,r,pplev,pt
[1242]393     s   ,pzlev,pzlay,pu,pv,ph,pq,zcdv_true,pq2,zkv,zkh,zkq,ustar
[652]394     s   ,9)
[544]395      ENDIF
396
[38]397      if ((doubleq).and.(ngrid.eq.1)) then
398        kmixmin = 80. !80.! minimum eddy mix coeff in 1D
399        do ilev=1,nlay
400          do ig=1,ngrid
401           zkh(ig,ilev) = max(kmixmin,zkh(ig,ilev))
402           zkv(ig,ilev) = max(kmixmin,zkv(ig,ilev))
403          end do
404        end do
405      end if
406
407c    ** diagnostique pour le schema de turbulence
408c       -----------------------------------------
409
410      IF(lecrit) THEN
411         PRINT*
412         PRINT*,'Diagnostic for the vertical turbulent mixing'
413         PRINT*,'Cd for momentum and potential temperature'
414
415         PRINT*,zcdv(ngrid/2+1),zcdh(ngrid/2+1)
416         PRINT*,'Mixing coefficient for momentum and pot.temp.'
417         DO ilev=1,nlay
418            PRINT*,zkv(ngrid/2+1,ilev),zkh(ngrid/2+1,ilev)
419         ENDDO
420      ENDIF
421
422
423
424
425c-----------------------------------------------------------------------
426c   4. inversion pour l'implicite sur u
427c      --------------------------------
428
429c    ** l'equation est
430c       u(t+1) =  u(t) + dt * {(du/dt)phys}(t)  +  dt * {(du/dt)difv}(t+1)
431c       avec
432c       /zu/ = u(t) + dt * {(du/dt)phys}(t)   (voir paragraphe 2.)
433c       et
434c       dt * {(du/dt)difv}(t+1) = dt * {(d/dz)[ Ku (du/dz) ]}(t+1)
435c       donc les entrees sont /zcdv/ pour la condition a la limite sol
436c       et /zkv/ = Ku
437 
438      CALL multipl((nlay-1)*ngrid,zkv(1,2),zb0(1,2),zb(1,2))
439      CALL multipl(ngrid,zcdv,zb0,zb)
440
441      DO ig=1,ngrid
442         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
443         zc(ig,nlay)=za(ig,nlay)*zu(ig,nlay)*z1(ig)
444         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
445      ENDDO
446
447      DO ilay=nlay-1,1,-1
448         DO ig=1,ngrid
449            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
450     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
451            zc(ig,ilay)=(za(ig,ilay)*zu(ig,ilay)+
452     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
453            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
454         ENDDO
455      ENDDO
456
457      DO ig=1,ngrid
458         zu(ig,1)=zc(ig,1)
459      ENDDO
460      DO ilay=2,nlay
461         DO ig=1,ngrid
462            zu(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zu(ig,ilay-1)
463         ENDDO
464      ENDDO
465
466
467
468
469
470c-----------------------------------------------------------------------
471c   5. inversion pour l'implicite sur v
472c      --------------------------------
473
474c    ** l'equation est
475c       v(t+1) =  v(t) + dt * {(dv/dt)phys}(t)  +  dt * {(dv/dt)difv}(t+1)
476c       avec
477c       /zv/ = v(t) + dt * {(dv/dt)phys}(t)   (voir paragraphe 2.)
478c       et
479c       dt * {(dv/dt)difv}(t+1) = dt * {(d/dz)[ Kv (dv/dz) ]}(t+1)
480c       donc les entrees sont /zcdv/ pour la condition a la limite sol
481c       et /zkv/ = Kv
482
483      DO ig=1,ngrid
484         z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
485         zc(ig,nlay)=za(ig,nlay)*zv(ig,nlay)*z1(ig)
486         zd(ig,nlay)=zb(ig,nlay)*z1(ig)
487      ENDDO
488
489      DO ilay=nlay-1,1,-1
490         DO ig=1,ngrid
491            z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
492     $         zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
493            zc(ig,ilay)=(za(ig,ilay)*zv(ig,ilay)+
494     $         zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
495            zd(ig,ilay)=zb(ig,ilay)*z1(ig)
496         ENDDO
497      ENDDO
498
499      DO ig=1,ngrid
500         zv(ig,1)=zc(ig,1)
501      ENDDO
502      DO ilay=2,nlay
503         DO ig=1,ngrid
504            zv(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zv(ig,ilay-1)
505         ENDDO
506      ENDDO
507
508
509
510
511
512c-----------------------------------------------------------------------
513c   6. inversion pour l'implicite sur h sans oublier le couplage
514c      avec le sol (conduction)
515c      ------------------------
516
517c    ** l'equation est
518c       h(t+1) =  h(t) + dt * {(dh/dt)phys}(t)  +  dt * {(dh/dt)difv}(t+1)
519c       avec
520c       /zh/ = h(t) + dt * {(dh/dt)phys}(t)   (voir paragraphe 2.)
521c       et
522c       dt * {(dh/dt)difv}(t+1) = dt * {(d/dz)[ Kh (dh/dz) ]}(t+1)
523c       donc les entrees sont /zcdh/ pour la condition de raccord au sol
524c       et /zkh/ = Kh
525c       -------------
526
[473]527c Mass variation scheme:
[38]528      CALL multipl((nlay-1)*ngrid,zkh(1,2),zb0(1,2),zb(1,2))
529      CALL multipl(ngrid,zcdh,zb0,zb)
530
[473]531c on initialise dm c
532     
533      pdtc(:,:)=0.
534      zt(:,:)=0.
535      dmice(:,:)=0.
[38]536
537c    ** calcul de (d Planck / dT) a la temperature d'interface
538c       ------------------------------------------------------
539
540      z4st=4.*5.67e-8*ptimestep
[544]541      IF (tke_heat_flux .eq. 0.) THEN
[38]542      DO ig=1,ngrid
543         zdplanck(ig)=z4st*pemis(ig)*ptsrf(ig)*ptsrf(ig)*ptsrf(ig)
544      ENDDO
[544]545      ELSE
546         zdplanck(:)=0.
547      ENDIF
[38]548
[473]549! calcul de zc et zd pour la couche top en prenant en compte le terme
[2080]550! de variation de masse (on fait une boucle pour que \E7a converge)
[473]551
552! Identification des points de grilles qui ont besoin de la correction
553
554      llnt(:)=1
[1236]555      IF (.not.turb_resolved) THEN
[884]556      IF (callcond) THEN
557       DO ig=1,ngrid
[473]558         DO l=1,nlay
559            if(zh(ig,l) .lt. zhcond(ig,l)) then
560               llnt(ig)=300 
561! 200 and 100 do not go beyond month 9 with normal dissipation
562               goto 5
563            endif
564         ENDDO
[884]5655      continue
566       ENDDO
567      ENDIF
[473]568
[529]569      ENDIF
570
[473]571      DO ig=1,ngrid
572
573! Initialization of z1 and zd, which do not depend on dmice
574
575      z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
576      zd(ig,nlay)=zb(ig,nlay)*z1(ig)
577
578      DO ilay=nlay-1,1,-1
579          z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
580     $        zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
581          zd(ig,ilay)=zb(ig,ilay)*z1(ig)
582      ENDDO
583
584! Convergence loop
585
586      DO j=1,llnt(ig)
587
588            z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
589            zc(ig,nlay)=za(ig,nlay)*zh(ig,nlay)
590     &      -betam(ig,nlay)*dmice(ig,nlay)
591            zc(ig,nlay)=zc(ig,nlay)*z1(ig)
592!            zd(ig,nlay)=zb(ig,nlay)*z1(ig)
593
594! calcul de zc et zd pour les couches du haut vers le bas
595
596             DO ilay=nlay-1,1,-1
597               z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
598     $            zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
599               zc(ig,ilay)=(za(ig,ilay)*zh(ig,ilay)+
600     $            zb(ig,ilay+1)*zc(ig,ilay+1)-
601     $            betam(ig,ilay)*dmice(ig,ilay))*z1(ig)
602!               zd(ig,ilay)=zb(ig,ilay)*z1(ig)
603            ENDDO
604
[38]605c    ** calcul de la temperature_d'interface et de sa tendance.
606c       on ecrit que la somme des flux est nulle a l'interface
607c       a t + \delta t,
608c       c'est a dire le flux radiatif a {t + \delta t}
609c       + le flux turbulent a {t + \delta t}
610c            qui s'ecrit K (T1-Tsurf) avec T1 = d1 Tsurf + c1
611c            (notation K dt = /cpp*b/)       
612c       + le flux dans le sol a t
613c       + l'evolution du flux dans le sol lorsque la temperature d'interface
614c       passe de sa valeur a t a sa valeur a {t + \delta t}.
615c       ----------------------------------------------------
616
617         z1(ig)=pcapcal(ig)*ptsrf(ig)+cpp*zb(ig,1)*zc(ig,1)
618     s     +zdplanck(ig)*ptsrf(ig)+ pfluxsrf(ig)*ptimestep
619         z2(ig)= pcapcal(ig)+cpp*zb(ig,1)*(1.-zd(ig,1))+zdplanck(ig)
620         ztsrf2(ig)=z1(ig)/z2(ig)
[473]621!         pdtsrf(ig)=(ztsrf2(ig)-ptsrf(ig))/ptimestep  !incremented outside loop
622            zhs(ig,1)=zc(ig,1)+zd(ig,1)*ztsrf2(ig)
[38]623
624c    ** et a partir de la temperature au sol on remonte
625c       -----------------------------------------------
626
[473]627            DO ilay=2,nlay
628               zhs(ig,ilay)=zc(ig,ilay)+zd(ig,ilay)*zhs(ig,ilay-1)
629            ENDDO
630            DO ilay=1,nlay
631               zt(ig,ilay)=zhs(ig,ilay)*ppopsk(ig,ilay)
632            ENDDO
633
634c      Condensation/sublimation in the atmosphere
635c      ------------------------------------------
636c      (computation of zcondicea and dmice)
637
[1996]638      IF (.NOT. co2clouds) then
639       DO l=nlay , 1, -1
[473]640           IF(zt(ig,l).LT.ztcond(ig,l)) THEN
641               pdtc(ig,l)=(ztcond(ig,l) - zt(ig,l))/ptimestep
642               zcondicea(ig,l)=(pplev(ig,l)-pplev(ig,l+1))
643     &                        *ccond*pdtc(ig,l)
644              dmice(ig,l)= dmice(ig,l) + zcondicea(ig,l)*ptimestep
645            END IF
[1996]646       ENDDO
647      ELSE
648       DO l=nlay , 1, -1
649         zcondicea(ig,l)= 0.!pcondicea_co2microp(ig,l)*
650c     &                        (pplev(ig,l) - pplev(ig,l+1))/g
651         dmice(ig,l)= 0.!dmice(ig,l) + zcondicea(ig,l)*ptimestep
652         pdtc(ig,l)=0.
653       ENDDO   
654      ENDIF
655     
656       ENDDO!of Do j=1,XXX
[38]657
[1996]658      ENDDO   !of Do ig=1,ngrid
[38]659
[473]660      pdtsrf(:)=(ztsrf2(:)-ptsrf(:))/ptimestep
[669]661     
[660]662      DO ig=1,ngrid  ! computing sensible heat flux (atm => surface)
663         sensibFlux(ig)=cpp*zb(ig,1)/ptimestep*(zhs(ig,1)-ztsrf2(ig))
664      ENDDO
[473]665
[38]666c-----------------------------------------------------------------------
667c   TRACERS
668c   -------
669
670      if(tracer) then
671
672c     Using the wind modified by friction for lifting and  sublimation
673c     ----------------------------------------------------------------
674
[529]675!     This is computed above and takes into account surface-atmosphere flux
676!     enhancement by subgrid gustiness and atmospheric-stability related
677!     variations of transfer coefficients.
678
679!        DO ig=1,ngrid
680!          zu2(ig)=zu(ig,1)*zu(ig,1)+zv(ig,1)*zv(ig,1)
681!          zcdv(ig)=zcdv_true(ig)*sqrt(zu2(ig))
682!          zcdh(ig)=zcdh_true(ig)*sqrt(zu2(ig))
683!        ENDDO
684
[38]685c       Calcul du flux vertical au bas de la premiere couche (dust) :
686c       -----------------------------------------------------------
[1047]687        do ig=1,ngrid 
[38]688          rho(ig) = zb0(ig,1) /ptimestep
689c          zb(ig,1) = 0.
690        end do
691c       Dust lifting:
692        if (lifting) then
[310]693#ifndef MESOSCALE
[38]694           if (doubleq.AND.submicron) then
695             do ig=1,ngrid
696c              if(co2ice(ig).lt.1) then
697                 pdqsdif(ig,igcm_dust_mass) =
698     &             -alpha_lift(igcm_dust_mass) 
699                 pdqsdif(ig,igcm_dust_number) =
700     &             -alpha_lift(igcm_dust_number) 
701                 pdqsdif(ig,igcm_dust_submicron) =
702     &             -alpha_lift(igcm_dust_submicron)
703c              end if
704             end do
705           else if (doubleq) then
[1974]706             if (dustinjection.eq.0) then !injection scheme 0 (old)
707                                          !or 2 (injection in CL)
708              do ig=1,ngrid
[1455]709               if(co2ice(ig).lt.1) then ! pas de soulevement si glace CO2
[38]710                 pdqsdif(ig,igcm_dust_mass) =
711     &             -alpha_lift(igcm_dust_mass) 
712                 pdqsdif(ig,igcm_dust_number) =
[520]713     &             -alpha_lift(igcm_dust_number)
714               end if
[1974]715              end do
716             elseif(dustinjection.eq.1)then ! dust injection scheme = 1 injection from surface
717              do ig=1,ngrid
718               if(co2ice(ig).lt.1) then ! pas de soulevement si glace CO2
719                IF((ti_injection.LE.local_time(ig)).and.
720     &                  (local_time(ig).LE.tf_injection)) THEN
721                  if (rdstorm) then !Rocket dust storm scheme
722                        pdqsdif(ig,igcm_stormdust_mass) =
723     &                          -alpha_lift(igcm_stormdust_mass)
724     &                          *dustliftday(ig)
725                        pdqsdif(ig,igcm_stormdust_number) =
726     &                          -alpha_lift(igcm_stormdust_number)
727     &                          *dustliftday(ig)
728                        pdqsdif(ig,igcm_dust_mass)= 0.
729                        pdqsdif(ig,igcm_dust_number)= 0.
730                  else
731                        pdqsdif(ig,igcm_dust_mass)=
732     &                          -dustliftday(ig)*
733     &                          alpha_lift(igcm_dust_mass)               
734                        pdqsdif(ig,igcm_dust_number)=
735     &                          -dustliftday(ig)*
736     &                          alpha_lift(igcm_dust_number)
737                  endif
738                  if (submicron) then
739                        pdqsdif(ig,igcm_dust_submicron) = 0.
740                  endif ! if (submicron)
741                ELSE ! outside dust injection time frame
[2080]742                  pdqsdif(ig,igcm_dust_mass)= 0.
743                  pdqsdif(ig,igcm_dust_number)= 0.
744                  if (rdstorm) then     
[1974]745                        pdqsdif(ig,igcm_stormdust_mass)= 0.
746                        pdqsdif(ig,igcm_stormdust_number)= 0.
[2080]747                  end if
[1974]748                ENDIF
749               
750                end if ! of if(co2ice(ig).lt.1)
751              end do
752             endif ! end if dustinjection
[38]753           else if (submicron) then
754             do ig=1,ngrid
755                 pdqsdif(ig,igcm_dust_submicron) =
756     &             -alpha_lift(igcm_dust_submicron)
757             end do
758           else
[1236]759#endif
[38]760            call dustlift(ngrid,nlay,nq,rho,zcdh_true,zcdh,co2ice,
761     &                    pdqsdif)
[1236]762#ifndef MESOSCALE
[38]763           endif !doubleq.AND.submicron
[310]764#endif
[38]765        else
766           pdqsdif(1:ngrid,1:nq) = 0.
767        end if
768
769c       OU calcul de la valeur de q a la surface (water)  :
770c       ----------------------------------------
771        if (water) then
[1047]772            call watersat(ngrid,ptsrf,pplev(1,1),qsat)
[38]773        end if
774
775c      Inversion pour l'implicite sur q
776c       --------------------------------
[1974]777        do iq=1,nq  !for all tracers including stormdust
[38]778          CALL multipl((nlay-1)*ngrid,zkh(1,2),zb0(1,2),zb(1,2))
779
780          if ((water).and.(iq.eq.igcm_h2o_vap)) then
781c            This line is required to account for turbulent transport
782c            from surface (e.g. ice) to mid-layer of atmosphere:
783             CALL multipl(ngrid,zcdv,zb0,zb(1,1))
784             CALL multipl(ngrid,dryness,zb(1,1),zb(1,1))
785          else ! (re)-initialize zb(:,1)
786             zb(1:ngrid,1)=0
787          end if
788
789          DO ig=1,ngrid
790               z1(ig)=1./(za(ig,nlay)+zb(ig,nlay))
791               zc(ig,nlay)=za(ig,nlay)*zq(ig,nlay,iq)*z1(ig)
792               zd(ig,nlay)=zb(ig,nlay)*z1(ig)
793          ENDDO
794 
795          DO ilay=nlay-1,2,-1
796               DO ig=1,ngrid
797                z1(ig)=1./(za(ig,ilay)+zb(ig,ilay)+
798     $           zb(ig,ilay+1)*(1.-zd(ig,ilay+1)))
799                zc(ig,ilay)=(za(ig,ilay)*zq(ig,ilay,iq)+
800     $           zb(ig,ilay+1)*zc(ig,ilay+1))*z1(ig)
801                zd(ig,ilay)=zb(ig,ilay)*z1(ig)
802               ENDDO
803          ENDDO
804
805          if (water.and.(iq.eq.igcm_h2o_ice)) then
806            ! special case for water ice tracer: do not include
807            ! h2o ice tracer from surface (which is set when handling
[473]808           ! h2o vapour case (see further down).
[38]809            DO ig=1,ngrid
810                z1(ig)=1./(za(ig,1)+zb(ig,1)+
811     $           zb(ig,2)*(1.-zd(ig,2)))
812                zc(ig,1)=(za(ig,1)*zq(ig,1,iq)+
813     $         zb(ig,2)*zc(ig,2)) *z1(ig)
814            ENDDO
815          else ! general case
816            DO ig=1,ngrid
817                z1(ig)=1./(za(ig,1)+zb(ig,1)+
818     $           zb(ig,2)*(1.-zd(ig,2)))
819                zc(ig,1)=(za(ig,1)*zq(ig,1,iq)+
820     $         zb(ig,2)*zc(ig,2) +
821     $        (-pdqsdif(ig,iq)) *ptimestep) *z1(ig)  !tracer flux from surface
822            ENDDO
823          endif ! of if (water.and.(iq.eq.igcm_h2o_ice))
824 
825          IF ((water).and.(iq.eq.igcm_h2o_vap)) then
826c           Calculation for turbulent exchange with the surface (for ice)
827            DO ig=1,ngrid
828              zd(ig,1)=zb(ig,1)*z1(ig)
829              zq1temp(ig)=zc(ig,1)+ zd(ig,1)*qsat(ig)
830
831              pdqsdif(ig,igcm_h2o_ice)=rho(ig)*dryness(ig)*zcdv(ig)
832     &                       *(zq1temp(ig)-qsat(ig))
833c             write(*,*)'flux vers le sol=',pdqsdif(ig,nq)
834            END DO
835
836            DO ig=1,ngrid
837              if(.not.watercaptag(ig)) then
838                if ((-pdqsdif(ig,igcm_h2o_ice)*ptimestep)
839     &             .gt.pqsurf(ig,igcm_h2o_ice)) then
840c                 write(*,*)'on sublime plus que qsurf!'
841                  pdqsdif(ig,igcm_h2o_ice)=
842     &                         -pqsurf(ig,igcm_h2o_ice)/ptimestep
843c                 write(*,*)'flux vers le sol=',pdqsdif(ig,nq)
844                  z1(ig)=1./(za(ig,1)+ zb(ig,2)*(1.-zd(ig,2)))
845                  zc(ig,1)=(za(ig,1)*zq(ig,1,igcm_h2o_vap)+
846     $            zb(ig,2)*zc(ig,2) +
847     $            (-pdqsdif(ig,igcm_h2o_ice)) *ptimestep) *z1(ig)
848                  zq1temp(ig)=zc(ig,1)
849                endif   
850              endif ! if (.not.watercaptag(ig))
851c             Starting upward calculations for water :
[669]852              zq(ig,1,igcm_h2o_vap)=zq1temp(ig)
853             
[742]854!c             Take into account H2O latent heat in surface energy budget
855!c             We solve dT/dt = (2834.3-0.28*(T-To)-0.004*(T-To)^2)*1e3*iceflux/cpp
856!              tsrf_lw(ig) = ptsrf(ig) + pdtsrf(ig) *ptimestep
857!             
858!              alpha = exp(-4*abs(T1-T2)*pdqsdif(ig,igcm_h2o_ice)
859!     &            *ptimestep/pcapcal(ig))
860!
861!              tsrf_lw(ig) = (tsrf_lw(ig)*(T2-alpha*T1)+T1*T2*(alpha-1))
862!     &                      /(tsrf_lw(ig)*(1-alpha)+alpha*T2-T1)  ! surface temperature at t+1
863!
864!              pdtsrf(ig) = (tsrf_lw(ig)-ptsrf(ig))/ptimestep
[669]865
866               if(pqsurf(ig,igcm_h2o_ice)
867     &           +pdqsdif(ig,igcm_h2o_ice)*ptimestep
868     &           .gt.frost_albedo_threshold) ! if there is still ice, T cannot exceed To
869     &           pdtsrf(ig) = min(pdtsrf(ig),(To-ptsrf(ig))/ptimestep) ! ice melt case
870     
[38]871            ENDDO ! of DO ig=1,ngrid
872          ELSE
873c           Starting upward calculations for simple mixing of tracer (dust)
874            DO ig=1,ngrid
875               zq(ig,1,iq)=zc(ig,1)
876            ENDDO
877          END IF ! of IF ((water).and.(iq.eq.igcm_h2o_vap))
878
879          DO ilay=2,nlay
880             DO ig=1,ngrid
[626]881               zq(ig,ilay,iq)=zc(ig,ilay)+zd(ig,ilay)*zq(ig,ilay-1,iq)
[38]882             ENDDO
883          ENDDO
884        enddo ! of do iq=1,nq
885      end if ! of if(tracer)
[669]886     
[38]887
888c-----------------------------------------------------------------------
889c   8. calcul final des tendances de la diffusion verticale
890c      ----------------------------------------------------
891
892      DO ilev = 1, nlay
893         DO ig=1,ngrid
894            pdudif(ig,ilev)=(    zu(ig,ilev)-
895     $      (pu(ig,ilev)+pdufi(ig,ilev)*ptimestep)    )/ptimestep
896            pdvdif(ig,ilev)=(    zv(ig,ilev)-
897     $      (pv(ig,ilev)+pdvfi(ig,ilev)*ptimestep)    )/ptimestep
[473]898            hh = ph(ig,ilev)+pdhfi(ig,ilev)*ptimestep
899     $  + (latcond*dmice(ig,ilev)/cpp)/ppopsk(ig,ilev)
900            pdhdif(ig,ilev)=( zhs(ig,ilev)- hh )/ptimestep
[38]901         ENDDO
902      ENDDO
903
904      if (tracer) then
905        DO iq = 1, nq
906          DO ilev = 1, nlay
907            DO ig=1,ngrid
908              pdqdif(ig,ilev,iq)=(zq(ig,ilev,iq)-
909     $      (pq(ig,ilev,iq) + pdqfi(ig,ilev,iq)*ptimestep))/ptimestep
910            ENDDO
911          ENDDO
912        ENDDO
913      end if
914
915c    ** diagnostique final
916c       ------------------
917
918      IF(lecrit) THEN
919         PRINT*,'In vdif'
920         PRINT*,'Ts (t) and Ts (t+st)'
921         WRITE(*,'(a10,3a15)')
922     s   'theta(t)','theta(t+dt)','u(t)','u(t+dt)'
923         PRINT*,ptsrf(ngrid/2+1),ztsrf2(ngrid/2+1)
924         DO ilev=1,nlay
925            WRITE(*,'(4f15.7)')
[473]926     s      ph(ngrid/2+1,ilev),zhs(ngrid/2+1,ilev),
[38]927     s      pu(ngrid/2+1,ilev),zu(ngrid/2+1,ilev)
928
929         ENDDO
930      ENDIF
931
[1036]932      END SUBROUTINE vdifc
[1969]933
934      END MODULE vdifc_mod
Note: See TracBrowser for help on using the repository browser.