source: trunk/LMDZ.MARS/libf/phymars/thermcell_main_mars.F90 @ 542

Last change on this file since 542 was 532, checked in by acolaitis, 14 years ago

Most up-to-date thermals parameters. Optimized for high resolution runs.

File size: 51.5 KB
RevLine 
[161]1!
2!
[185]3      SUBROUTINE thermcell_main_mars(ptimestep  &
[161]4     &                  ,pplay,pplev,pphi,zlev,zlay  &
5     &                  ,pu,pv,pt,pq,pq2  &
6     &                  ,pduadj,pdvadj,pdtadj,pdqadj,pdq2adj  &
[185]7     &                  ,fm,entr,detr,lmax,zmax  &
[161]8     &                  ,r_aspect &
9     &                  ,zw2,fraca &
10     &                  ,zpopsk,ztla,heatFlux,heatFlux_down &
11     &                  ,buoyancyOut, buoyancyEst)
12
13      IMPLICIT NONE
14
15!=======================================================================
16! Mars version of thermcell_main. Author : A Colaitis
17!=======================================================================
[185]18
19#include "dimensions.h"
20#include "dimphys.h"
21#include "comcstfi.h"
[508]22#include "tracer.h"
23#include "callkeys.h"
[185]24
[161]25! ============== INPUTS ==============
26
27      REAL, INTENT(IN) :: ptimestep,r_aspect
[185]28      REAL, INTENT(IN) :: pt(ngridmx,nlayermx)
29      REAL, INTENT(IN) :: pu(ngridmx,nlayermx)
30      REAL, INTENT(IN) :: pv(ngridmx,nlayermx)
31      REAL, INTENT(IN) :: pq(ngridmx,nlayermx,nqmx)
32      REAL, INTENT(IN) :: pq2(ngridmx,nlayermx)
33      REAL, INTENT(IN) :: pplay(ngridmx,nlayermx)
34      REAL, INTENT(IN) :: pplev(ngridmx,nlayermx+1)
35      REAL, INTENT(IN) :: pphi(ngridmx,nlayermx)
36      REAL, INTENT(IN) :: zlay(ngridmx,nlayermx)
37      REAL, INTENT(IN) :: zlev(ngridmx,nlayermx+1)
[161]38
39! ============== OUTPUTS ==============
40
[185]41      REAL, INTENT(OUT) :: pdtadj(ngridmx,nlayermx)
[337]42      REAL :: pduadj(ngridmx,nlayermx)
43      REAL :: pdvadj(ngridmx,nlayermx)
44      REAL :: pdqadj(ngridmx,nlayermx,nqmx)
[313]45!      REAL, INTENT(OUT) :: pdq2adj(ngridmx,nlayermx)
46      REAL :: pdq2adj(ngridmx,nlayermx)
[185]47      REAL, INTENT(OUT) :: zw2(ngridmx,nlayermx+1)
[161]48
49! Diagnostics
[185]50      REAL, INTENT(OUT) :: heatFlux(ngridmx,nlayermx)   ! interface heatflux
[313]51     REAL, INTENT(OUT) :: heatFlux_down(ngridmx,nlayermx) ! interface heat flux from downdraft
52!      REAL, INTENT(OUT) :: buoyancyOut(ngridmx,nlayermx)  ! interlayer buoyancy term
53!      REAL, INTENT(OUT) :: buoyancyEst(ngridmx,nlayermx)  ! interlayer estimated buoyancy term
[161]54
55! dummy variables when output not needed :
56
[185]57!      REAL :: heatFlux(ngridmx,nlayermx)   ! interface heatflux
58!      REAL :: heatFlux_down(ngridmx,nlayermx) ! interface heat flux from downdraft
[313]59      REAL :: buoyancyOut(ngridmx,nlayermx)  ! interlayer buoyancy term
60      REAL :: buoyancyEst(ngridmx,nlayermx)  ! interlayer estimated buoyancy term
[161]61
62
63! ============== LOCAL ================
64
65      INTEGER ig,k,l,ll,iq
[185]66      INTEGER lmax(ngridmx),lmin(ngridmx),lalim(ngridmx)
67      REAL linter(ngridmx)
68      REAL zmax(ngridmx)
69      REAL ztva(ngridmx,nlayermx),zw_est(ngridmx,nlayermx+1),ztva_est(ngridmx,nlayermx)
70      REAL zmax_sec(ngridmx)
71      REAL zh(ngridmx,nlayermx)
72      REAL zdthladj(ngridmx,nlayermx)
73      REAL zdthladj_down(ngridmx,nlayermx)
74      REAL ztvd(ngridmx,nlayermx)
75      REAL ztv(ngridmx,nlayermx)
76      REAL zu(ngridmx,nlayermx),zv(ngridmx,nlayermx),zo(ngridmx,nlayermx)
77      REAL zva(ngridmx,nlayermx)
78      REAL zua(ngridmx,nlayermx)
[161]79
[185]80      REAL zta(ngridmx,nlayermx)
81      REAL fraca(ngridmx,nlayermx+1)
82      REAL q2(ngridmx,nlayermx)
83      REAL rho(ngridmx,nlayermx),rhobarz(ngridmx,nlayermx),masse(ngridmx,nlayermx)
84      REAL zpopsk(ngridmx,nlayermx)
[161]85
[185]86      REAL wmax(ngridmx)
87      REAL wmax_sec(ngridmx)
88      REAL fm(ngridmx,nlayermx+1),entr(ngridmx,nlayermx),detr(ngridmx,nlayermx)
[161]89
[185]90      REAL fm_down(ngridmx,nlayermx+1)
[161]91
[185]92      REAL ztla(ngridmx,nlayermx)
[161]93
[185]94      REAL f_star(ngridmx,nlayermx+1),entr_star(ngridmx,nlayermx)
95      REAL detr_star(ngridmx,nlayermx)
96      REAL alim_star_tot(ngridmx)
97      REAL alim_star(ngridmx,nlayermx)
98      REAL alim_star_clos(ngridmx,nlayermx)
99      REAL f(ngridmx)
[161]100
[185]101      REAL teta_th_int(ngridmx,nlayermx)
102      REAL teta_env_int(ngridmx,nlayermx)
103      REAL teta_down_int(ngridmx,nlayermx)
[161]104
105      CHARACTER (LEN=20) :: modname
106      CHARACTER (LEN=80) :: abort_message
107
108! ============= PLUME VARIABLES ============
109
[185]110      REAL w_est(ngridmx,nlayermx+1)
111      REAL wa_moy(ngridmx,nlayermx+1)
112      REAL wmaxa(ngridmx)
113      REAL zdz,zbuoy(ngridmx,nlayermx),zw2m
[508]114      LOGICAL activecell(ngridmx),activetmp(ngridmx)
[532]115      REAL a1,b1,ae,be,ad,bd,fdfu,b1inv,a1inv,omega
[290]116      INTEGER tic
[161]117
118! ==========================================
119
120! ============= HEIGHT VARIABLES ===========
121
[185]122      REAL num(ngridmx)
123      REAL denom(ngridmx)
124      REAL zlevinter(ngridmx)
[161]125
126! =========================================
127
128! ============= CLOSURE VARIABLES =========
129
[185]130      REAL zdenom(ngridmx)
131      REAL alim_star2(ngridmx)
132      REAL alim_star_tot_clos(ngridmx)
[161]133      INTEGER llmax
134
135! =========================================
136
137! ============= FLUX2 VARIABLES ===========
138
139      INTEGER ncorecfm1,ncorecfm2,ncorecfm3,ncorecalpha
140      INTEGER ncorecfm4,ncorecfm5,ncorecfm6,ncorecfm7,ncorecfm8
141      REAL zfm
142      REAL f_old,ddd0,eee0,ddd,eee,zzz
143      REAL fomass_max,alphamax
144
145! =========================================
146
147! ============= DTETA VARIABLES ===========
148
149! rien : on prends la divergence du flux turbulent
150
151! =========================================
152
153! ============= DV2 VARIABLES =============
154!               not used for now
155
[185]156      REAL qa(ngridmx,nlayermx),detr_dv2(ngridmx,nlayermx),zf,zf2
157      REAL wvd(ngridmx,nlayermx+1),wud(ngridmx,nlayermx+1)
158      REAL gamma0(ngridmx,nlayermx+1),gamma(ngridmx,nlayermx+1)
159      REAL ue(ngridmx,nlayermx),ve(ngridmx,nlayermx)
160      LOGICAL ltherm(ngridmx,nlayermx)
161      REAL dua(ngridmx,nlayermx),dva(ngridmx,nlayermx)
[161]162      INTEGER iter
163      INTEGER nlarga0
164
165! =========================================
166
[508]167! ============== Theta_M Variables ========
168
169      INTEGER ico2
170      SAVE ico2
171      REAL m_co2, m_noco2, A , B
172      SAVE A, B
173      LOGICAL firstcall
174      save firstcall
175      data firstcall/.true./
176      REAL zhc(ngridmx,nlayermx)
177      REAL zdzfull(ngridmx,nlayermx)
178      REAL ratiom(ngridmx,nlayermx)
179
180! =========================================
181
[161]182!-----------------------------------------------------------------------
183!   initialisation:
184!   ---------------
185
[336]186      entr(:,:)=0.
187      detr(:,:)=0.
188      fm(:,:)=0.
[337]189!      zu(:,:)=pu(:,:)
190!      zv(:,:)=pv(:,:)
[508]191      zhc(:,:)=pt(:,:)/zpopsk(:,:)
[161]192
[508]193! **********************************************************************
194! Taking into account vertical molar mass gradients
195! **********************************************************************
196
197      if(firstcall) then
198        ico2=0
199        if (tracer) then
200!     Prepare Special treatment if one of the tracers is CO2 gas
201           do iq=1,nqmx
202             if (noms(iq).eq."co2") then
203                ico2=iq
204                m_co2 = 44.01E-3  ! CO2 molecular mass (kg/mol)
205                m_noco2 = 33.37E-3  ! Non condensible mol mass (kg/mol)
206!               Compute A and B coefficient use to compute
207!               mean molecular mass Mair defined by
208!               1/Mair = q(ico2)/m_co2 + (1-q(ico2))/m_noco2
209!               1/Mair = A*q(ico2) + B
210                A =(1/m_co2 - 1/m_noco2)
211                B=1/m_noco2
212             end if
213           enddo
214        endif
215
216      firstcall=.false.
217      endif !of if firstcall
218
219!.......................................................................
220!  Special treatment for co2
221!.......................................................................
222
223      if (ico2.ne.0) then
224!     Special case if one of the tracers is CO2 gas
225         DO l=1,nlayermx
226           DO ig=1,ngridmx
227            ztv(ig,l) = zhc(ig,l)*(A*pq(ig,l,ico2)+B)
228           ENDDO
229         ENDDO
230       else
231          ztv(:,:)=zhc(:,:)
232       end if
233
234
[161]235!------------------------------------------------------------------------
236!                       --------------------
237!
238!
239!                       + + + + + + + + + + +
240!
241!
242!  wa, fraca, wd, fracd --------------------   zlev(2), rhobarz
243!  wh,wt,wo ...
244!
245!                       + + + + + + + + + + +  zh,zu,zv,zo,rho
246!
247!
248!                       --------------------   zlev(1)
249!                       \\\\\\\\\\\\\\\\\\\\
250!
251!
252
253!-----------------------------------------------------------------------
254!   Calcul des altitudes des couches
255!-----------------------------------------------------------------------
256
[185]257!      do l=2,nlayermx
258!         zlev(:,l)=0.5*(pphi(:,l)+pphi(:,l-1))/g
[161]259!      enddo
260!         zlev(:,1)=0.
[185]261!         zlev(:,nlayermx+1)=(2.*pphi(:,nlayermx)-pphi(:,nlayermx-1))/g
[161]262
[185]263!         zlay(:,:)=pphi(:,:)/g
[161]264!-----------------------------------------------------------------------
265!   Calcul des densites
266!-----------------------------------------------------------------------
267
[185]268      rho(:,:)=pplay(:,:)/(r*pt(:,:))
[161]269
270      rhobarz(:,1)=rho(:,1)
271
[185]272      do l=2,nlayermx
[161]273         rhobarz(:,l)=0.5*(rho(:,l)+rho(:,l-1))
274      enddo
275
276!calcul de la masse
[185]277      do l=1,nlayermx
278         masse(:,l)=(pplev(:,l)-pplev(:,l+1))/g
[161]279      enddo
280
281
282!------------------------------------------------------------------
283!
284!             /|\
285!    --------  |  F_k+1 -------   
286!                              ----> D_k
287!             /|\              <---- E_k , A_k
288!    --------  |  F_k ---------
289!                              ----> D_k-1
290!                              <---- E_k-1 , A_k-1
291!
292!
293!    ---------------------------
294!
295!    ----- F_lmax+1=0 ----------         \
296!            lmax     (zmax)              |
297!    ---------------------------          |
298!                                         |
299!    ---------------------------          |
300!                                         |
301!    ---------------------------          |
302!                                         |
303!    ---------------------------          |
304!                                         |
305!    ---------------------------          |
306!                                         |  E
307!    ---------------------------          |  D
308!                                         |
309!    ---------------------------          |
310!                                         |
311!    ---------------------------  \       |
312!            lalim                 |      |
313!    ---------------------------   |      |
314!                                  |      |
315!    ---------------------------   |      |
316!                                  | A    |
317!    ---------------------------   |      |
318!                                  |      |
319!    ---------------------------   |      |
320!    lmin  (=1 pour le moment)     |      |
321!    ----- F_lmin=0 ------------  /      /
322!
323!    ---------------------------
324!    //////////////////////////
325!
326
327!=============================================================================
328!  Calculs initiaux ne faisant pas intervenir les changements de phase
329!=============================================================================
330
331!------------------------------------------------------------------
332!  1. alim_star est le profil vertical de l'alimentation a la base du
333!     panache thermique, calcule a partir de la flotabilite de l'air sec
334!  2. lmin et lalim sont les indices inferieurs et superieurs de alim_star
335!------------------------------------------------------------------
336!
337      entr_star=0. ; detr_star=0. ; alim_star=0. ; alim_star_tot=0.
338      lmin=1
339
340!-----------------------------------------------------------------------------
341!  3. wmax_sec et zmax_sec sont les vitesses et altitudes maximum d'un
342!     panache sec conservatif (e=d=0) alimente selon alim_star
343!     Il s'agit d'un calcul de type CAPE
344!     zmax_sec est utilise pour determiner la geometrie du thermique.
345!------------------------------------------------------------------------------
346!---------------------------------------------------------------------------------
347!calcul du melange et des variables dans le thermique
348!--------------------------------------------------------------------------------
349
350! ===========================================================================
351! ===================== PLUME ===============================================
352! ===========================================================================
353
354! Initialisations des variables reeles
355      ztva(:,:)=ztv(:,:)
356      ztva_est(:,:)=ztva(:,:)
357      ztla(:,:)=0.
[173]358      zdz=0.
[161]359      zbuoy(:,:)=0.
360      w_est(:,:)=0.
361      f_star(:,:)=0.
362      wa_moy(:,:)=0.
363      linter(:)=1.
[185]364
[496]365! --------------------------------------------------------------------------
366! -------------- MAIN PARAMETERS FOR THERMALS MODEL ------------------------
367! --------------  see thermiques.pro and getfit.py -------------------------
368
[185]369!      a1=2.5 ; b1=0.0015 ; ae=0.045 ; be = 0.6  ! svn baseline
370
[496]371! Using broad downdraft selection
372!      a1=1.60226 ; b1=0.0006 ; ae=0.0454 ; be = 0.57
373!      ad = 0.0005114  ; bd = -0.662
374!      fdfu = -1.9
[185]375
[496]376! Using conditional sampling downdraft selection
377      a1=1.4716 ; b1=0.0005698 ; ae=0.03683 ; be = 0.57421
378      ad = 0.00048088  ; bd = -0.6697
379      fdfu = -1.3
[532]380      a1inv=a1
381      b1inv=b1
382      omega=0.
[496]383
[532]384! One good config for 34/35 levels
385!      a1inv=a1
386!      b1inv=b1
387!      be=1.1*be
388
[512]389! Best configuration for 222 levels:
[532]390
391      omega=0.06
[512]392      b1=0.
[532]393      a1=1.
[512]394      a1inv=0.25*a1
395      b1inv=0.0002
396
[532]397! Best config for norad 222 levels:
398
399!       omega=0.06
400!       a1=1.
401!       b1=0.
402!       a1inv=a1
403!       be=1.1*be
404!       ad = 0.0004
405!       b1inv=0.00035
406
407! Trying stuff :
408
409!       omega=0.06
410!       a1=1.
411!       b1=0.
412!       be=be
413!       a1inv=0.25*a1
414!       b1inv=0.0002   
415!       ad=1.1*ad
416!       ae=1.*ae
[496]417! --------------------------------------------------------------------------
418! --------------------------------------------------------------------------
419! --------------------------------------------------------------------------
420
[161]421! Initialisation des variables entieres
422      wmaxa(:)=0.
423      lalim(:)=1
424
425!-------------------------------------------------------------------------
426! On ne considere comme actif que les colonnes dont les deux premieres
427! couches sont instables.
428!-------------------------------------------------------------------------
[508]429      activecell(:)=ztv(:,1)>ztv(:,2)
[185]430          do ig=1,ngridmx
[161]431            if (ztv(ig,1)>=(ztv(ig,2))) then
432               alim_star(ig,1)=MAX((ztv(ig,1)-ztv(ig,2)),0.)  &
433     &                       *sqrt(zlev(ig,2))
[290]434!     &                       /sqrt(zlev(ig,2))
[165]435!      &                       *zlev(ig,2)
[161]436               lalim(ig)=2
437               alim_star_tot(ig)=alim_star_tot(ig)+alim_star(ig,1)
438            endif
439         enddo
440
[185]441       do l=2,nlayermx-1
442!        do l=2,4
443         do ig=1,ngridmx
[513]444           if (ztv(ig,l)>(ztv(ig,l+1)+0.) .and. ztv(ig,1)>=ztv(ig,l) .and. (alim_star(ig,l-1) .ne. 0.)) then ! .and. (zlev(ig,l+1) .lt. 1000.)) then
[161]445               alim_star(ig,l)=MAX((ztv(ig,l)-ztv(ig,l+1)),0.)  &
446     &                       *sqrt(zlev(ig,l+1))
[165]447!      &                       *zlev(ig,2)
[161]448                lalim(ig)=l+1
449               alim_star_tot(ig)=alim_star_tot(ig)+alim_star(ig,l)
450           endif
451        enddo
452      enddo
[185]453      do l=1,nlayermx
454         do ig=1,ngridmx
[161]455            if (alim_star_tot(ig) > 1.e-10 ) then
456               alim_star(ig,l)=alim_star(ig,l)/alim_star_tot(ig)
457            endif
458         enddo
459      enddo
460
461      alim_star_tot(:)=1.
[185]462!      if(alim_star(1,1) .ne. 0.) then
463!      print*, 'lalim star' ,lalim(1)
464!      endif
[161]465
466!------------------------------------------------------------------------------
467! Calcul dans la premiere couche
468! On decide dans cette version que le thermique n'est actif que si la premiere
469! couche est instable.
470! Pourrait etre change si on veut que le thermiques puisse se déclencher
471! dans une couche l>1
472!------------------------------------------------------------------------------
473
[185]474      do ig=1,ngridmx
[161]475! Le panache va prendre au debut les caracteristiques de l'air contenu
476! dans cette couche.
[508]477          if (activecell(ig)) then
[161]478          ztla(ig,1)=ztv(ig,1)
479!cr: attention, prise en compte de f*(1)=1 => AC : what ? f*(1) =0. ! (d'ou f*(2)=a*(1)
480! dans un panache conservatif
481          f_star(ig,1)=0.
[290]482         
[161]483          f_star(ig,2)=alim_star(ig,1)
[290]484
[185]485          zw2(ig,2)=2.*g*(ztv(ig,1)-ztv(ig,2))/ztv(ig,2)  &
[161]486      &                     *(zlev(ig,2)-zlev(ig,1))  &
[290]487      &                    *0.4*pphi(ig,1)/(pphi(ig,2)-pphi(ig,1))
[161]488          w_est(ig,2)=zw2(ig,2)
489
490          endif
491      enddo
492
493
494!==============================================================================
495!boucle de calcul de la vitesse verticale dans le thermique
496!==============================================================================
[185]497      do l=2,nlayermx-1
[161]498!==============================================================================
499
500
501! On decide si le thermique est encore actif ou non
502! AFaire : Il faut sans doute ajouter entr_star a alim_star dans ce test
[185]503          do ig=1,ngridmx
[508]504             activecell(ig)=activecell(ig) &
[161]505      &                 .and. zw2(ig,l)>1.e-10 &
506      &                 .and. f_star(ig,l)+alim_star(ig,l)>1.e-10
507          enddo
508
509!---------------------------------------------------------------------------
510! calcul des proprietes thermodynamiques et de la vitesse de la couche l
511! sans tenir compte du detrainement et de l'entrainement dans cette
512! couche
513! C'est a dire qu'on suppose
514! ztla(l)=ztla(l-1)
515! Ici encore, on doit pouvoir ajouter entr_star (qui peut etre calculer
516! avant) a l'alimentation pour avoir un calcul plus propre
517!---------------------------------------------------------------------------
518
[185]519          do ig=1,ngridmx
[508]520             if(activecell(ig)) then
[161]521
[165]522!                if(l .lt. lalim(ig)) then
[290]523!               ztva_est(ig,l)=(f_star(ig,l)*ztla(ig,l-1)+  &
[165]524!     &            alim_star(ig,l)*ztv(ig,l))  &
525!     &            /(f_star(ig,l)+alim_star(ig,l))
526!                else
[161]527                ztva_est(ig,l)=ztla(ig,l-1)
[165]528!                endif
[161]529
530                zdz=zlev(ig,l+1)-zlev(ig,l)
[185]531                zbuoy(ig,l)=g*(ztva_est(ig,l)-ztv(ig,l))/ztv(ig,l)
532
533                if (((a1*zbuoy(ig,l)/w_est(ig,l)-b1) .gt. 0.) .and. (w_est(ig,l) .ne. 0.)) then
534                w_est(ig,l+1)=Max(0.0001,w_est(ig,l)+2.*zdz*a1*zbuoy(ig,l)-2.*zdz*w_est(ig,l)*b1 &
[532]535     & -2.*(1.-omega)*zdz*w_est(ig,l)*ae*(a1*zbuoy(ig,l)/w_est(ig,l)-b1)**be)
[161]536                else
[512]537                w_est(ig,l+1)=Max(0.0001,w_est(ig,l)+2.*zdz*a1inv*zbuoy(ig,l)-2.*zdz*w_est(ig,l)*b1inv)
[161]538                endif
539                if (w_est(ig,l+1).lt.0.) then
540                w_est(ig,l+1)=zw2(ig,l)
541                endif
542             endif
543          enddo
544
545!-------------------------------------------------
546!calcul des taux d'entrainement et de detrainement
547!-------------------------------------------------
548
[185]549      do ig=1,ngridmx
[508]550        if (activecell(ig)) then
[161]551
552          zw2m=w_est(ig,l+1)
[185]553
554          if((a1*(zbuoy(ig,l)/zw2m)-b1) .gt. 0.) then
[161]555          entr_star(ig,l)=f_star(ig,l)*zdz*  &
[512]556        &   MAX(0.,ae*(a1*(zbuoy(ig,l)/zw2m)-b1)**be)
557!         &  MAX(0.,log(1. + 0.03*sqrt(a1*(zbuoy(ig,l)/zw2m)-b1)))
[161]558          else
559          entr_star(ig,l)=0.
560          endif
[185]561
[161]562          if(zbuoy(ig,l) .gt. 0.) then
563             if(l .lt. lalim(ig)) then
[185]564                detr_star(ig,l)=0.
[161]565             else
[185]566
[161]567!                 detr_star(ig,l) = f_star(ig,l)*zdz*              &
568!              &     0.0105*((zbuoy(ig,l)/zw2m)/0.048)**(1./1.7)
569!                 detr_star(ig,l) = f_star(ig,l)*zdz*              &
570!              &     0.0085*((zbuoy(ig,l)/zw2m)/0.05)**(1./1.55)
571
572! last baseline from direct les
573!                 detr_star(ig,l) = f_star(ig,l)*zdz*              &
574!               &     0.065*(2.5*(zbuoy(ig,l)/zw2m))**0.75
575
576! new param from continuity eq with a fit on dfdz
577                 detr_star(ig,l) = f_star(ig,l)*zdz*              &
[512]578            &  ad
579!             &  Max(0., 0.0005 - 0.55*zbuoy(ig,l)/zw2m)
[161]580
[185]581!               &     MAX(0.,-0.38*zbuoy(ig,l)/zw2m+0.0005)   !svn baseline
582!                &     MAX(0.,-0.38*zbuoy(ig,l)/zw2m+0.0008)     
583
[161]584!              &     0.014*((zbuoy(ig,l)/zw2m)/0.05)**(1./1.35)
585!                 detr_star(ig,l) = f_star(ig,l)*zdz*              &
586!              &     ((zbuoy(ig,l)/zw2m)/2.222)! + 0.0002)
587
588             endif
589          else
590          detr_star(ig,l)=f_star(ig,l)*zdz*                        &
[512]591                &     bd*zbuoy(ig,l)/zw2m
592!             &  Max(0., 0.001 - 0.45*zbuoy(ig,l)/zw2m)
593!             &  Max(0., Min(0.001,0.0005 - 0.55*zbuoy(ig,l)/zw2m))
[161]594
[508]595
[185]596!               &     MAX(0.,-0.38*zbuoy(ig,l)/zw2m+0.0005)      !svn baseline
597
[161]598!              &  *5.*(-afact*zbetalpha*zbuoy(ig,l)/zw2m)
599!              &  *5.*(-afact*zbuoy(ig,l)/zw2m)
600
601! last baseline from direct les
602!               &     0.065*(-2.5*(zbuoy(ig,l)/zw2m))**0.75
603
604! new param from continuity eq with a fit on dfdz
605
606
607          endif
608
609! En dessous de lalim, on prend le max de alim_star et entr_star pour
610! alim_star et 0 sinon
611
[512]612       if (l.lt.lalim(ig)) then
[161]613          alim_star(ig,l)=max(alim_star(ig,l),entr_star(ig,l))
[512]614          entr_star(ig,l)=0.
615       endif
[161]616
617! Calcul du flux montant normalise
618
619      f_star(ig,l+1)=f_star(ig,l)+alim_star(ig,l)+entr_star(ig,l)  &
620     &              -detr_star(ig,l)
621
622      endif
623      enddo
624
625
626!----------------------------------------------------------------------------
627!calcul de la vitesse verticale en melangeant Tl et qt du thermique
628!---------------------------------------------------------------------------
629
[532]630      DO tic=0,0  ! internal convergence loop
[508]631      activetmp(:)=activecell(:) .and. f_star(:,l+1)>1.e-10
[185]632      do ig=1,ngridmx
[161]633       if (activetmp(ig)) then
634
635           ztla(ig,l)=(f_star(ig,l)*ztla(ig,l-1)+  &
636     &            (alim_star(ig,l)+entr_star(ig,l))*ztv(ig,l))  &
637     &            /(f_star(ig,l+1)+detr_star(ig,l))
638
639        endif
640      enddo
641
[313]642      activetmp(:)=activetmp(:).and.(abs(ztla(:,l)-ztva(:,l)).gt.0.01)
643
[185]644      do ig=1,ngridmx
[161]645      if (activetmp(ig)) then
646           ztva(ig,l) = ztla(ig,l)
[185]647           zbuoy(ig,l)=g*(ztva(ig,l)-ztv(ig,l))/ztv(ig,l)
[161]648
[185]649           if (((a1*zbuoy(ig,l)/zw2(ig,l)-b1) .gt. 0.) .and. (zw2(ig,l) .ne. 0.) ) then
650           zw2(ig,l+1)=Max(0.,zw2(ig,l)+2.*zdz*a1*zbuoy(ig,l)-         &
[532]651     & 2.*zdz*zw2(ig,l)*b1-2.*(1.-omega)*zdz*zw2(ig,l)*ae*(a1*zbuoy(ig,l)/zw2(ig,l)-b1)**be)
[161]652           else
[512]653           zw2(ig,l+1)=Max(0.,zw2(ig,l)+2.*zdz*a1inv*zbuoy(ig,l)-2.*zdz*zw2(ig,l)*b1inv)
[161]654           endif
655      endif
656      enddo
657
[290]658! ================ RECOMPUTE ENTR, DETR, and F FROM NEW W2 ===================
659
660
661      do ig=1,ngridmx
662        if (activetmp(ig)) then
663
664          zw2m=zw2(ig,l+1)
665          if(zw2m .gt. 0) then
666          if((a1*(zbuoy(ig,l)/zw2m)-b1) .gt. 0.) then
667          entr_star(ig,l)=f_star(ig,l)*zdz*  &
668        &   MAX(0.,ae*(a1*(zbuoy(ig,l)/zw2m)-b1)**be)
[512]669!         &  MAX(0.,log(1. + 0.03*sqrt(a1*(zbuoy(ig,l)/zw2m)-b1)))
[290]670          else
671          entr_star(ig,l)=0.
672          endif
673
674          if(zbuoy(ig,l) .gt. 0.) then
675             if(l .lt. lalim(ig)) then
676                detr_star(ig,l)=0.
677             else
678                 detr_star(ig,l) = f_star(ig,l)*zdz*              &
679            &  ad
[512]680!             &  Max(0., 0.0005 - 0.55*zbuoy(ig,l)/zw2m)
681
[290]682             endif
683          else
684          detr_star(ig,l)=f_star(ig,l)*zdz*                        &
685                &     bd*zbuoy(ig,l)/zw2m
[512]686!             &  Max(0.,Min(0.001,0.0005 - 0.55*zbuoy(ig,l)/zw2m))
687
[290]688          endif
689          else
690          entr_star(ig,l)=0.
691          detr_star(ig,l)=0.
692          endif
693
694! En dessous de lalim, on prend le max de alim_star et entr_star pour
695! alim_star et 0 sinon
696
697        if (l.lt.lalim(ig)) then
698          alim_star(ig,l)=max(alim_star(ig,l),entr_star(ig,l))
699          entr_star(ig,l)=0.
700        endif
701
702! Calcul du flux montant normalise
703
704      f_star(ig,l+1)=f_star(ig,l)+alim_star(ig,l)+entr_star(ig,l)  &
705     &              -detr_star(ig,l)
706
707      endif
708      enddo
709
[313]710      ENDDO   ! of tic
711
[161]712!---------------------------------------------------------------------------
713!initialisations pour le calcul de la hauteur du thermique, de l'inversion et de la vitesse verticale max
714!---------------------------------------------------------------------------
715
[185]716      do ig=1,ngridmx
[161]717            if (zw2(ig,l+1)>0. .and. zw2(ig,l+1).lt.1.e-10) then
718             print*,'On tombe sur le cas particulier de thermcell_plume'
719                zw2(ig,l+1)=0.
720                linter(ig)=l+1
721            endif
722
723        if (zw2(ig,l+1).lt.0.) then
724           linter(ig)=(l*(zw2(ig,l+1)-zw2(ig,l))  &
725     &               -zw2(ig,l))/(zw2(ig,l+1)-zw2(ig,l))
726           zw2(ig,l+1)=0.
727        endif
728           wa_moy(ig,l+1)=sqrt(zw2(ig,l+1))
729
730        if (wa_moy(ig,l+1).gt.wmaxa(ig)) then
731!   lmix est le niveau de la couche ou w (wa_moy) est maximum
732            wmaxa(ig)=wa_moy(ig,l+1)
733        endif
734      enddo
735
736!=========================================================================
737! FIN DE LA BOUCLE VERTICALE
738      enddo
739!=========================================================================
740
741!on recalcule alim_star_tot
[185]742       do ig=1,ngridmx
[161]743          alim_star_tot(ig)=0.
744       enddo
[185]745       do ig=1,ngridmx
[161]746          do l=1,lalim(ig)-1
747          alim_star_tot(ig)=alim_star_tot(ig)+alim_star(ig,l)
748          enddo
749       enddo
750
[185]751      do l=1,nlayermx
752         do ig=1,ngridmx
[165]753            if (alim_star_tot(ig) > 1.e-10 ) then
754               alim_star(ig,l)=alim_star(ig,l)/alim_star_tot(ig)
755            endif
756         enddo
757      enddo
[161]758
759! ===========================================================================
760! ================= FIN PLUME ===============================================
761! ===========================================================================
762
763! ===========================================================================
764! ================= HEIGHT ==================================================
765! ===========================================================================
766
767! Attention, w2 est transforme en sa racine carree dans cette routine
768
769!-------------------------------------------------------------------------------
[356]770! Calcul des caracteristiques du thermique:zmax,wmax
[161]771!-------------------------------------------------------------------------------
772
773!calcul de la hauteur max du thermique
[185]774      do ig=1,ngridmx
[161]775         lmax(ig)=lalim(ig)
776      enddo
[185]777      do ig=1,ngridmx
778         do l=nlayermx,lalim(ig)+1,-1
[161]779            if (zw2(ig,l).le.1.e-10) then
[512]780               lmax(ig)=l-1
[161]781            endif
782         enddo
783      enddo
784
785! On traite le cas particulier qu'il faudrait éviter ou le thermique
786! atteind le haut du modele ...
[185]787      do ig=1,ngridmx
788      if ( zw2(ig,nlayermx) > 1.e-10 ) then
[161]789          print*,'WARNING !!!!! W2 thermiques non nul derniere couche '
[185]790          lmax(ig)=nlayermx
[161]791      endif
792      enddo
793
794! pas de thermique si couche 1 stable
[185]795!      do ig=1,ngridmx
[161]796!         if (lmin(ig).gt.1) then
797!             lmax(ig)=1
798!             lmin(ig)=1
799!             lalim(ig)=1
800!         endif
801!      enddo
802!
803! Determination de zw2 max
[185]804      do ig=1,ngridmx
[161]805         wmax(ig)=0.
806      enddo
807
[185]808      do l=1,nlayermx
809         do ig=1,ngridmx
[161]810            if (l.le.lmax(ig)) then
811                if (zw2(ig,l).lt.0.)then
812                  print*,'pb2 zw2<0'
813                endif
814                zw2(ig,l)=sqrt(zw2(ig,l))
815                wmax(ig)=max(wmax(ig),zw2(ig,l))
816            else
817                 zw2(ig,l)=0.
818            endif
819          enddo
820      enddo
821!   Longueur caracteristique correspondant a la hauteur des thermiques.
[185]822      do  ig=1,ngridmx
[161]823         zmax(ig)=0.
824         zlevinter(ig)=zlev(ig,1)
825      enddo
826
827         num(:)=0.
828         denom(:)=0.
[185]829         do ig=1,ngridmx
830          do l=1,nlayermx
[161]831             num(ig)=num(ig)+zw2(ig,l)*zlev(ig,l)*(zlev(ig,l+1)-zlev(ig,l))
832             denom(ig)=denom(ig)+zw2(ig,l)*(zlev(ig,l+1)-zlev(ig,l))
833          enddo
834       enddo
[185]835       do ig=1,ngridmx
[161]836       if (denom(ig).gt.1.e-10) then
837          zmax(ig)=2.*num(ig)/denom(ig)
838       endif
839       enddo
840
841! Attention, w2 est transforme en sa racine carree dans cette routine
842
843! ===========================================================================
844! ================= FIN HEIGHT ==============================================
845! ===========================================================================
846
847! Choix de la fonction d'alimentation utilisee pour la fermeture.
848
849      alim_star_clos(:,:)=entr_star(:,:)+alim_star(:,:)
850
851! ===========================================================================
852! ============= CLOSURE =====================================================
853! ===========================================================================
854
855!-------------------------------------------------------------------------------
856! Fermeture,determination de f
857!-------------------------------------------------------------------------------
858! Appel avec la version seche
859
860      alim_star2(:)=0.
861      alim_star_tot_clos(:)=0.
862      f(:)=0.
863
864! Indice vertical max (max de lalim) atteint par les thermiques sur le domaine
865      llmax=1
[185]866      do ig=1,ngridmx
[161]867         if (lalim(ig)>llmax) llmax=lalim(ig)
868      enddo
869
870
871! Calcul des integrales sur la verticale de alim_star et de
872!   alim_star^2/(rho dz)
873      do k=1,llmax-1
[185]874         do ig=1,ngridmx
[161]875            if (k<lalim(ig)) then
[185]876         alim_star2(ig)=alim_star2(ig)+alim_star_clos(ig,k)*alim_star_clos(ig,k)  &
[161]877      &                    /(rho(ig,k)*(zlev(ig,k+1)-zlev(ig,k)))
878         alim_star_tot_clos(ig)=alim_star_tot_clos(ig)+alim_star_clos(ig,k)
879      endif
880         enddo
881      enddo
[185]882 
[161]883! WARNING : MARS MODIF : we have added 2. : ratio of wmax/vmoy
884! True ratio is 3.5 but wetake into account the vmoy is the one alimentating
885! the thermal, so there are vs=0 into the vmoy... the true vmoy is lower. (a la louche)
886! And r_aspect has been changed from 2 to 1.5 from observations
[185]887      do ig=1,ngridmx
[161]888         if (alim_star2(ig)>1.e-10) then
[185]889!            f(ig)=wmax_sec(ig)*alim_star_tot_clos(ig)/  &
890!      &     (max(500.,zmax_sec(ig))*r_aspect*alim_star2(ig))
891             f(ig)=wmax(ig)*alim_star_tot_clos(ig)/  &
892      &     (max(500.,zmax(ig))*r_aspect*alim_star2(ig))
893
[161]894         endif
895      enddo
896
897! ===========================================================================
898! ============= FIN CLOSURE =================================================
899! ===========================================================================
900
901
902! ===========================================================================
903! ============= FLUX2 =======================================================
904! ===========================================================================
905
906!-------------------------------------------------------------------------------
907!deduction des flux
908!-------------------------------------------------------------------------------
909
910      fomass_max=0.8
911      alphamax=0.5
912
913      ncorecfm1=0
914      ncorecfm2=0
915      ncorecfm3=0
916      ncorecfm4=0
917      ncorecfm5=0
918      ncorecfm6=0
919      ncorecfm7=0
920      ncorecfm8=0
921      ncorecalpha=0
922
923!-------------------------------------------------------------------------
924! Multiplication par le flux de masse issu de la femreture
925!-------------------------------------------------------------------------
926
[185]927      do l=1,nlayermx
[161]928         entr(:,l)=f(:)*(entr_star(:,l)+alim_star(:,l))
929         detr(:,l)=f(:)*detr_star(:,l)
930      enddo
931
[185]932      do l=1,nlayermx
933         do ig=1,ngridmx
[161]934            if (l.lt.lmax(ig)) then
935               fm(ig,l+1)=fm(ig,l)+entr(ig,l)-detr(ig,l)
936            elseif(l.eq.lmax(ig)) then
937               fm(ig,l+1)=0.
938               detr(ig,l)=fm(ig,l)+entr(ig,l)
939            else
940               fm(ig,l+1)=0.
941            endif
942         enddo
943      enddo
944
945! Test provisoire : pour comprendre pourquoi on corrige plein de fois
946! le cas fm6, on commence par regarder une premiere fois avant les
947! autres corrections.
948
[314]949!      do l=1,nlayermx
950!         do ig=1,ngridmx
951!            if (detr(ig,l).gt.fm(ig,l)) then
952!               ncorecfm8=ncorecfm8+1
953!            endif
954!         enddo
955!      enddo
[161]956
957!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
958! FH Version en cours de test;
959! par rapport a thermcell_flux, on fait une grande boucle sur "l"
960! et on modifie le flux avec tous les contr�les appliques d'affilee
961! pour la meme couche
962! Momentanement, on duplique le calcule du flux pour pouvoir comparer
963! les flux avant et apres modif
964!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
965
[185]966      do l=1,nlayermx
[161]967
[185]968         do ig=1,ngridmx
[161]969            if (l.lt.lmax(ig)) then
970               fm(ig,l+1)=fm(ig,l)+entr(ig,l)-detr(ig,l)
971            elseif(l.eq.lmax(ig)) then
972               fm(ig,l+1)=0.
973               detr(ig,l)=fm(ig,l)+entr(ig,l)
974            else
975               fm(ig,l+1)=0.
976            endif
977         enddo
978
979
980!-------------------------------------------------------------------------
981! Verification de la positivite des flux de masse
982!-------------------------------------------------------------------------
983
[185]984         do ig=1,ngridmx
[190]985
[161]986            if (fm(ig,l+1).lt.0.) then
[190]987               if((l+1) .eq. lmax(ig)) then
988               detr(ig,l)=detr(ig,l)+fm(ig,l+1)
989               fm(ig,l+1)=0.
[336]990               entr(ig,l+1)=0.
[190]991               ncorecfm2=ncorecfm2+1
992               else
993          print*,'fm(l+1)<0 : ig, l+1,lmax :',ig,l+1,lmax(ig),fm(ig,l+1)
994               ncorecfm1=ncorecfm1+1
[161]995               fm(ig,l+1)=fm(ig,l)
996               detr(ig,l)=entr(ig,l)
[190]997               endif
[161]998            endif
[190]999
[161]1000         enddo
1001
[508]1002!  Les "optimisations" du flux sont desactivecelles : moins de bidouilles
[161]1003!  je considere que celles ci ne sont pas justifiees ou trop delicates
1004!  pour MARS, d'apres des observations LES.
1005!-------------------------------------------------------------------------
1006!Test sur fraca croissant
1007!-------------------------------------------------------------------------
1008!      if (iflag_thermals_optflux==0) then
[185]1009!         do ig=1,ngridmx
[161]1010!          if (l.ge.lalim(ig).and.l.le.lmax(ig) &
1011!     &    .and.(zw2(ig,l+1).gt.1.e-10).and.(zw2(ig,l).gt.1.e-10) ) then
1012!!  zzz est le flux en l+1 a frac constant
1013!             zzz=fm(ig,l)*rhobarz(ig,l+1)*zw2(ig,l+1)  &
1014!     &                          /(rhobarz(ig,l)*zw2(ig,l))
1015!             if (fm(ig,l+1).gt.zzz) then
1016!                detr(ig,l)=detr(ig,l)+fm(ig,l+1)-zzz
1017!                fm(ig,l+1)=zzz
1018!                ncorecfm4=ncorecfm4+1
1019!             endif
1020!          endif
1021!        enddo
1022!      endif
1023!
1024!-------------------------------------------------------------------------
1025!test sur flux de masse croissant
1026!-------------------------------------------------------------------------
1027!      if (iflag_thermals_optflux==0) then
[185]1028!         do ig=1,ngridmx
[161]1029!            if ((fm(ig,l+1).gt.fm(ig,l)).and.(l.gt.lalim(ig))) then
1030!              f_old=fm(ig,l+1)
1031!              fm(ig,l+1)=fm(ig,l)
1032!              detr(ig,l)=detr(ig,l)+f_old-fm(ig,l+1)
1033!               ncorecfm5=ncorecfm5+1
1034!            endif
1035!         enddo
1036!      endif
1037!
1038!-------------------------------------------------------------------------
1039!detr ne peut pas etre superieur a fm
1040!-------------------------------------------------------------------------
1041
[185]1042         do ig=1,ngridmx
[161]1043            if (detr(ig,l).gt.fm(ig,l)) then
1044               ncorecfm6=ncorecfm6+1
1045               detr(ig,l)=fm(ig,l)
1046               entr(ig,l)=fm(ig,l+1)
1047
1048! Dans le cas ou on est au dessus de la couche d'alimentation et que le
1049! detrainement est plus fort que le flux de masse, on stope le thermique.
[314]1050!            endif
[161]1051
[336]1052            if(l.gt.lmax(ig)) then
1053!            if(l.gt.lalim(ig)) then
[161]1054               detr(ig,l)=0.
1055               fm(ig,l+1)=0.
1056               entr(ig,l)=0.
1057            endif
[314]1058           
1059            endif
1060
[161]1061         enddo
1062
1063!-------------------------------------------------------------------------
1064!fm ne peut pas etre negatif
1065!-------------------------------------------------------------------------
1066
[185]1067         do ig=1,ngridmx
[161]1068            if (fm(ig,l+1).lt.0.) then
1069               detr(ig,l)=detr(ig,l)+fm(ig,l+1)
1070               fm(ig,l+1)=0.
1071               ncorecfm2=ncorecfm2+1
1072            endif
1073         enddo
1074
1075!-----------------------------------------------------------------------
1076!la fraction couverte ne peut pas etre superieure a 1
1077!-----------------------------------------------------------------------
1078!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1079! FH Partie a revisiter.
1080! Il semble qu'etaient codees ici deux optiques dans le cas
1081! F/ (rho *w) > 1
1082! soit limiter la hauteur du thermique en considerant que c'est
1083! la derniere chouche, soit limiter F a rho w.
1084! Dans le second cas, il faut en fait limiter a un peu moins
1085! que ca parce qu'on a des 1 / ( 1 -alpha) un peu plus loin
1086! dans thermcell_main et qu'il semble de toutes facons deraisonable
1087! d'avoir des fractions de 1..
1088! Ci dessous, et dans l'etat actuel, le premier des  deux if est
1089! sans doute inutile.
1090!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1091
[185]1092        do ig=1,ngridmx
[161]1093           if (zw2(ig,l+1).gt.1.e-10) then
1094           zfm=rhobarz(ig,l+1)*zw2(ig,l+1)*alphamax
1095           if ( fm(ig,l+1) .gt. zfm) then
1096              f_old=fm(ig,l+1)
1097              fm(ig,l+1)=zfm
1098              detr(ig,l)=detr(ig,l)+f_old-fm(ig,l+1)
1099              ncorecalpha=ncorecalpha+1
1100           endif
1101           endif
1102
1103        enddo
1104
1105! Fin de la grande boucle sur les niveaux verticaux
1106      enddo
1107
1108!-----------------------------------------------------------------------
1109! On fait en sorte que la quantite totale d'air entraine dans le
1110! panache ne soit pas trop grande comparee a la masse de la maille
1111!-----------------------------------------------------------------------
1112
[185]1113      do l=1,nlayermx-1
1114         do ig=1,ngridmx
[161]1115            eee0=entr(ig,l)
1116            ddd0=detr(ig,l)
1117            eee=entr(ig,l)-masse(ig,l)*fomass_max/ptimestep
1118            ddd=detr(ig,l)-eee
1119            if (eee.gt.0.) then
1120                ncorecfm3=ncorecfm3+1
1121                entr(ig,l)=entr(ig,l)-eee
1122                if ( ddd.gt.0.) then
1123!   l'entrainement est trop fort mais l'exces peut etre compense par une
1124!   diminution du detrainement)
1125                   detr(ig,l)=ddd
1126                else
1127!   l'entrainement est trop fort mais l'exces doit etre compense en partie
1128!   par un entrainement plus fort dans la couche superieure
1129                   if(l.eq.lmax(ig)) then
1130                      detr(ig,l)=fm(ig,l)+entr(ig,l)
1131                   else
1132                      entr(ig,l+1)=entr(ig,l+1)-ddd
1133                      detr(ig,l)=0.
1134                      fm(ig,l+1)=fm(ig,l)+entr(ig,l)
1135                      detr(ig,l)=0.
1136                   endif
1137                endif
1138            endif
1139         enddo
1140      enddo
1141!
1142!              ddd=detr(ig)-entre
1143!on s assure que tout s annule bien en zmax
[185]1144      do ig=1,ngridmx
[161]1145         fm(ig,lmax(ig)+1)=0.
1146         entr(ig,lmax(ig))=0.
1147         detr(ig,lmax(ig))=fm(ig,lmax(ig))+entr(ig,lmax(ig))
1148      enddo
1149
1150!-----------------------------------------------------------------------
1151! Impression du nombre de bidouilles qui ont ete necessaires
1152!-----------------------------------------------------------------------
1153
1154!IM 090508 beg
[185]1155      if (ncorecfm1+ncorecfm2+ncorecfm3+ncorecfm4+ncorecfm5+ncorecalpha > ngridmx/4. ) then
[161]1156         print*,'thermcell warning : large number of corrections'
1157         print*,'PB thermcell : on a du coriger ',ncorecfm1,'x fm1',&
1158     &     ncorecfm2,'x fm2',ncorecfm3,'x fm3 et', &
1159     &     ncorecfm4,'x fm4',ncorecfm5,'x fm5 et', &
1160     &     ncorecfm6,'x fm6', &
1161     &     ncorecfm7,'x fm7', &
1162     &     ncorecfm8,'x fm8', &
1163     &     ncorecalpha,'x alpha'
1164      endif
1165
1166! ===========================================================================
1167! ============= FIN FLUX2 ===================================================
1168! ===========================================================================
1169
1170
1171! ===========================================================================
1172! ============= TRANSPORT ===================================================
1173! ===========================================================================
1174
1175!------------------------------------------------------------------
1176!   calcul du transport vertical
1177!------------------------------------------------------------------
1178
1179! ------------------------------------------------------------------
1180! Transport de teta dans l'updraft : (remplace thermcell_dq)
1181! ------------------------------------------------------------------
1182
1183      zdthladj(:,:)=0.
1184
[185]1185      do ig=1,ngridmx
[161]1186         if(lmax(ig) .gt. 1) then
1187         do k=1,lmax(ig)
1188            zdthladj(ig,k)=(1./masse(ig,k))*(fm(ig,k+1)*ztv(ig,k+1)-    &
1189     &   fm(ig,k)*ztv(ig,k)+fm(ig,k)*ztva(ig,k)-fm(ig,k+1)*ztva(ig,k+1))
1190            if (ztv(ig,k) + ptimestep*zdthladj(ig,k) .le. 0.) then
[165]1191              print*,'Teta<0 in thermcell_dTeta up: qenv .. dq : ', ztv(ig,k),ptimestep*zdthladj(ig,k)
1192              if(ztv(ig,k) .gt. 0.) then
1193              zdthladj(ig,k)=0.
1194              endif
[161]1195            endif
1196         enddo
1197         endif
1198      enddo
1199
1200! ------------------------------------------------------------------
1201! Prescription des proprietes du downdraft
1202! ------------------------------------------------------------------
1203
1204      ztvd(:,:)=ztv(:,:)
1205      fm_down(:,:)=0.
[185]1206      do ig=1,ngridmx
[161]1207         if (lmax(ig) .gt. 1) then
1208         do l=1,lmax(ig)
[512]1209!              if(zlay(ig,l) .le. 0.8*zmax(ig)) then
1210              if(zlay(ig,l) .le. zmax(ig)) then
[496]1211              fm_down(ig,l) =fm(ig,l)* &
1212     &      max(fdfu,-3*max(0.,(zlay(ig,l)/zmax(ig)))-0.9)
[161]1213              endif
1214
[512]1215!             if(zlay(ig,l) .le. 0.06*zmax(ig)) then
1216!          ztvd(ig,l)=ztv(ig,l)*max(0.,(1.+(sqrt((zlay(ig,l)/zmax(ig))/0.122449) - 1.)*(ztva(ig,l)/ztv(ig,l) - 1.)))
1217!             elseif(zlay(ig,l) .le. 0.4*zmax(ig)) then
1218!          ztvd(ig,l)=ztv(ig,l)*max(0.,1.-0.25*(ztva(ig,l)/ztv(ig,l) - 1.))
1219!             elseif(zlay(ig,l) .le. 0.7*zmax(ig)) then
1220!          ztvd(ig,l)=ztv(ig,l)*max(0.,(1.+(((zlay(ig,l)/zmax(ig))-0.7)/1.)*(ztva(ig,l)/ztv(ig,l) - 1.)))
1221!             else
1222!          ztvd(ig,l)=ztv(ig,l)
1223!            endif
1224
1225!            if(zlay(ig,l) .le. 0.6*zmax(ig)) then
1226!          ztvd(ig,l)=ztv(ig,l)*((zlay(ig,l)/zmax(ig))/179.848 + 0.997832)
1227!            elseif(zlay(ig,l) .le. 0.8*zmax(ig)) then
1228!          ztvd(ig,l)=-ztv(ig,l)*(((zlay(ig,l)/zmax(ig))-171.74)/170.94)
1229!             else
1230!          ztvd(ig,l)=ztv(ig,l)
1231!            endif
1232
1233
1234!            if(zlay(ig,l) .le. 0.8*zmax(ig)) then
1235!          ztvd(ig,l)=ztv(ig,l)*((zlay(ig,l)/zmax(ig))/224.81 + 0.997832)
1236!            elseif(zlay(ig,l) .le. zmax(ig)) then
1237!          ztvd(ig,l)=-ztv(ig,l)*(((zlay(ig,l)/zmax(ig))-144.685)/143.885)
1238!             else
1239!          ztvd(ig,l)=ztv(ig,l)
1240!            endif
1241
1242
1243            if(zlay(ig,l) .le. zmax(ig)) then
1244          ztvd(ig,l)=ztv(ig,l)*((zlay(ig,l)/zmax(ig))/299.7 + 0.997832)
[161]1245             else
1246          ztvd(ig,l)=ztv(ig,l)
1247            endif
1248
1249         enddo
1250         endif
1251      enddo
1252
1253! ------------------------------------------------------------------
1254! Transport de teta dans le downdraft : (remplace thermcell_dq)
1255! ------------------------------------------------------------------
1256
1257       zdthladj_down(:,:)=0.
1258
[185]1259      do ig=1,ngridmx
[161]1260         if(lmax(ig) .gt. 1) then
[290]1261! No downdraft in the very-near surface layer, we begin at k=3
[496]1262 
1263         do k=3,lmax(ig)
[161]1264            zdthladj_down(ig,k)=(1./masse(ig,k))*(fm_down(ig,k+1)*ztv(ig,k+1)- &
1265     & fm_down(ig,k)*ztv(ig,k)+fm_down(ig,k)*ztvd(ig,k)-fm_down(ig,k+1)*ztvd(ig,k+1))
1266            if (ztv(ig,k) + ptimestep*zdthladj_down(ig,k) .le. 0.) then
1267              print*,'q<0 in thermcell_dTeta down: qenv .. dq : ', ztv(ig,k),ptimestep*zdthladj_down(ig,k)
[165]1268              if(ztv(ig,k) .gt. 0.) then
1269              zdthladj(ig,k)=0.
1270              endif
[161]1271            endif
1272         enddo
1273         endif
1274      enddo
1275
1276!------------------------------------------------------------------
1277! Calcul de la fraction de l'ascendance
1278!------------------------------------------------------------------
[185]1279      do ig=1,ngridmx
[161]1280         fraca(ig,1)=0.
[185]1281         fraca(ig,nlayermx+1)=0.
[161]1282      enddo
[185]1283      do l=2,nlayermx
1284         do ig=1,ngridmx
[161]1285            if (zw2(ig,l).gt.1.e-10) then
1286            fraca(ig,l)=fm(ig,l)/(rhobarz(ig,l)*zw2(ig,l))
1287            else
1288            fraca(ig,l)=0.
1289            endif
1290         enddo
1291      enddo
1292
1293
1294
1295! ===========================================================================
1296! ============= DV2 =========================================================
1297! ===========================================================================
1298! ROUTINE OVERIDE : ne prends pas en compte le downdraft
1299! de plus, le gradient de pression horizontal semble tout deregler... A VOIR
1300
1301      if (0 .eq. 1) then
1302
1303!------------------------------------------------------------------
1304!  calcul du transport vertical du moment horizontal
1305!------------------------------------------------------------------
1306
1307! Calcul du transport de V tenant compte d'echange par gradient
1308! de pression horizontal avec l'environnement
1309
1310!   calcul du detrainement
1311!---------------------------
1312
1313      nlarga0=0.
1314
[185]1315      do k=1,nlayermx
1316         do ig=1,ngridmx
[161]1317            detr_dv2(ig,k)=fm(ig,k)-fm(ig,k+1)+entr(ig,k)
1318         enddo
1319      enddo
1320
1321!   calcul de la valeur dans les ascendances
[185]1322      do ig=1,ngridmx
[161]1323         zua(ig,1)=zu(ig,1)
1324         zva(ig,1)=zv(ig,1)
1325         ue(ig,1)=zu(ig,1)
1326         ve(ig,1)=zv(ig,1)
1327      enddo
1328
[185]1329      gamma(1:ngridmx,1)=0.
1330      do k=2,nlayermx
1331         do ig=1,ngridmx
[161]1332            ltherm(ig,k)=(fm(ig,k+1)+detr_dv2(ig,k))*ptimestep > 1.e-5*masse(ig,k)
1333            if(ltherm(ig,k).and.zmax(ig)>0.) then
1334               gamma0(ig,k)=masse(ig,k)  &
1335     &         *sqrt( 0.5*(fraca(ig,k+1)+fraca(ig,k)) )  &
1336     &         *0.5/zmax(ig)  &
1337     &         *1.
1338            else
1339               gamma0(ig,k)=0.
1340            endif
1341            if (ltherm(ig,k).and.zmax(ig)<=0.) nlarga0=nlarga0+1
1342          enddo
1343      enddo
1344
1345      gamma(:,:)=0.
1346
[185]1347      do k=2,nlayermx
[161]1348
[185]1349         do ig=1,ngridmx
[161]1350
1351            if (ltherm(ig,k)) then
1352               dua(ig,k)=zua(ig,k-1)-zu(ig,k-1)
1353               dva(ig,k)=zva(ig,k-1)-zv(ig,k-1)
1354            else
1355               zua(ig,k)=zu(ig,k)
1356               zva(ig,k)=zv(ig,k)
1357               ue(ig,k)=zu(ig,k)
1358               ve(ig,k)=zv(ig,k)
1359            endif
1360         enddo
1361
1362
1363! Debut des iterations
1364!----------------------
1365
1366! AC WARNING : SALE !
1367
1368      do iter=1,5
[185]1369         do ig=1,ngridmx
[161]1370! Pour memoire : calcul prenant en compte la fraction reelle
1371!              zf=0.5*(fraca(ig,k)+fraca(ig,k+1))
1372!              zf2=1./(1.-zf)
1373! Calcul avec fraction infiniement petite
1374               zf=0.
1375               zf2=1.
1376
1377!  la première fois on multiplie le coefficient de freinage
1378!  par le module du vent dans la couche en dessous.
1379!  Mais pourquoi donc ???
1380               if (ltherm(ig,k)) then
1381!   On choisit une relaxation lineaire.
1382!                 gamma(ig,k)=gamma0(ig,k)
1383!   On choisit une relaxation quadratique.
[337]1384                gamma(ig,k)=gamma0(ig,k)*sqrt(dua(ig,k)**2+dva(ig,k)**2)
[161]1385                  zua(ig,k)=(fm(ig,k)*zua(ig,k-1)  &
1386     &               +(zf2*entr(ig,k)+gamma(ig,k))*zu(ig,k))  &
1387     &               /(fm(ig,k+1)+detr_dv2(ig,k)+entr(ig,k)*zf*zf2  &
1388     &                 +gamma(ig,k))
1389                  zva(ig,k)=(fm(ig,k)*zva(ig,k-1)  &
1390     &               +(zf2*entr(ig,k)+gamma(ig,k))*zv(ig,k))  &
1391     &               /(fm(ig,k+1)+detr_dv2(ig,k)+entr(ig,k)*zf*zf2  &
1392     &                 +gamma(ig,k))
1393
1394!                  print*,' OUTPUT DV2 !!!!!!!!!!!!!!!!!!!!',k,zua(ig,k),zva(ig,k),zu(ig,k),zv(ig,k),dua(ig,k),dva(ig,k)
1395                  dua(ig,k)=zua(ig,k)-zu(ig,k)
1396                  dva(ig,k)=zva(ig,k)-zv(ig,k)
1397                  ue(ig,k)=(zu(ig,k)-zf*zua(ig,k))*zf2
1398                  ve(ig,k)=(zv(ig,k)-zf*zva(ig,k))*zf2
1399               endif
1400         enddo
1401! Fin des iterations
1402!--------------------
1403      enddo
1404
[185]1405      enddo ! k=2,nlayermx
[161]1406
1407! Calcul du flux vertical de moment dans l'environnement.
1408!---------------------------------------------------------
[185]1409      do k=2,nlayermx
1410         do ig=1,ngridmx
[161]1411            wud(ig,k)=fm(ig,k)*ue(ig,k)
1412            wvd(ig,k)=fm(ig,k)*ve(ig,k)
1413         enddo
1414      enddo
[185]1415      do ig=1,ngridmx
[161]1416         wud(ig,1)=0.
[185]1417         wud(ig,nlayermx+1)=0.
[161]1418         wvd(ig,1)=0.
[185]1419         wvd(ig,nlayermx+1)=0.
[161]1420      enddo
1421
1422! calcul des tendances.
1423!-----------------------
[185]1424      do k=1,nlayermx
1425         do ig=1,ngridmx
[161]1426            pduadj(ig,k)=((detr_dv2(ig,k)+gamma(ig,k))*zua(ig,k)  &
1427     &               -(entr(ig,k)+gamma(ig,k))*ue(ig,k)  &
1428     &               -wud(ig,k)+wud(ig,k+1))  &
1429     &               /masse(ig,k)
1430            pdvadj(ig,k)=((detr_dv2(ig,k)+gamma(ig,k))*zva(ig,k)  &
1431     &               -(entr(ig,k)+gamma(ig,k))*ve(ig,k)  &
1432     &               -wvd(ig,k)+wvd(ig,k+1))  &
1433     &               /masse(ig,k)
1434         enddo
1435      enddo
1436
1437
1438! Sorties eventuelles.
1439!----------------------
1440
1441!      if (nlarga0>0) then
1442!          print*,'WARNING !!!!!! DANS THERMCELL_DV2 '
1443!      print*,nlarga0,' points pour lesquels laraga=0. dans un thermique'
1444!          print*,'Il faudrait decortiquer ces points'
1445!      endif
1446
1447! ===========================================================================
1448! ============= FIN DV2 =====================================================
1449! ===========================================================================
1450
1451      else
1452
[337]1453!      modname='momentum'
1454!      call thermcell_dqupdown(ngridmx,nlayermx,ptimestep,fm,entr,detr,  &
1455!     &     masse,zu,pduadj,ztvd,fm_down,ztv,modname,lmax)
1456!
1457!      call thermcell_dqupdown(ngridmx,nlayermx,ptimestep,fm,entr,detr,  &
1458!     &     masse,zv,pdvadj,ztvd,fm_down,ztv,modname,lmax)
[161]1459
1460      endif
1461
1462!------------------------------------------------------------------
[508]1463!  calcul du transport vertical de traceurs
[161]1464!------------------------------------------------------------------
1465
[508]1466! We only transport co2 tracer because it is coupled to the scheme through theta_m
1467! The rest is transported outside the sub-timestep loop
1468
1469      if (ico2.ne.0) then
1470!      if (nqmx .ne. 0.) then
[185]1471      do l=1,nlayermx
[508]1472         zdzfull(:,l)=zlev(:,l+1)-zlev(:,l)
[161]1473      enddo
1474
[508]1475      modname='tracer'
1476      call thermcell_dqup(ngridmx,nlayermx,ptimestep     &
1477     &     ,fm,entr,detr,  &
1478     &    masse,pq(:,:,ico2),pdqadj(:,:,ico2),modname,zdzfull)
1479!      endif
1480
1481! Compute the ratio between theta and theta_m
1482     
1483       do l=1,nlayermx
1484          do ig=1,ngridmx
1485             ratiom(ig,l)=1./(A*(pq(ig,l,ico2)+pdqadj(ig,l,ico2)*ptimestep)+B)
1486          enddo
1487       enddo
1488       else
1489          ratiom(:,:)=1.
1490       endif
1491
1492
[161]1493!------------------------------------------------------------------
[508]1494!  incrementation dt
[161]1495!------------------------------------------------------------------
1496
[508]1497      do l=1,nlayermx
1498         do ig=1,ngridmx
1499         pdtadj(ig,l)=(zdthladj(ig,l)+zdthladj_down(ig,l))*zpopsk(ig,l)*ratiom(ig,l)
[512]1500!         pdtadj(ig,l)=zdthladj(ig,l)*zpopsk(ig,l)*ratiom(ig,l)
[508]1501         enddo
1502      enddo
[161]1503
1504!------------------------------------------------------------------
1505!  calcul du transport vertical de la tke
1506!------------------------------------------------------------------
1507
[313]1508!      modname='tke'
1509!      call thermcell_dqupdown(ngridmx,nlayermx,ptimestep,fm,entr,detr,  &
1510!      &      masse,pq2,pdq2adj,ztvd,fm_down,ztv,modname,lmax)
[161]1511
1512! ===========================================================================
1513! ============= FIN TRANSPORT ===============================================
1514! ===========================================================================
1515
1516
1517!------------------------------------------------------------------
1518!   Calculs de diagnostiques pour les sorties
1519!------------------------------------------------------------------
1520! DIAGNOSTIQUE
1521! We compute interface values for teta env and th. The last interface
1522! value does not matter, as the mass flux is 0 there.
1523
1524   
[185]1525      do l=1,nlayermx-1
1526       do ig=1,ngridmx
[508]1527        teta_th_int(ig,l)=0.5*(ztva(ig,l+1)+ztva(ig,l))*ratiom(ig,l)
1528        teta_down_int(ig,l) = 0.5*(ztvd(ig,l+1)+ztvd(ig,l))*ratiom(ig,l)
1529        teta_env_int(ig,l)=0.5*(ztv(ig,l+1)+ztv(ig,l))*ratiom(ig,l)
[161]1530       enddo
1531      enddo
[185]1532      do ig=1,ngridmx
1533       teta_th_int(ig,nlayermx)=teta_th_int(ig,nlayermx-1)
1534       teta_env_int(ig,nlayermx)=teta_env_int(ig,nlayermx-1)
1535       teta_down_int(ig,nlayermx)=teta_down_int(ig,nlayermx-1)
[161]1536      enddo
[185]1537      do l=1,nlayermx
1538       do ig=1,ngridmx
[161]1539        heatFlux(ig,l)=fm(ig,l)*(teta_th_int(ig,l)-teta_env_int(ig,l))/(rhobarz(ig,l))
[508]1540        buoyancyOut(ig,l)=g*(ztva(ig,l)-ztv(ig,l))/ztv(ig,l)
1541        buoyancyEst(ig,l)=g*(ztva_est(ig,l)-ztv(ig,l))/ztv(ig,l)
[161]1542        heatFlux_down(ig,l)=fm_down(ig,l)*(teta_down_int(ig,l)-teta_env_int(ig,l))/rhobarz(ig,l)
1543       enddo
1544      enddo
1545
1546      return
1547      end
Note: See TracBrowser for help on using the repository browser.