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

Last change on this file since 544 was 544, checked in by acolaitis, 13 years ago

27/02/12 == AC

Continuation of thermals setting, comparisons with mesoscale results on Case C
Added possibility to call gcm (or 1d) with constant prescribed sensible heat flux, in the spirit of direct comparisons with LES

... This is directly comparable to the variable tke_heat_flux in namelist.input
... Requires the use of yamada4.F from terrestrial GCM (mixes more, seem more numerically stable)
... Usually requires high timesteps (>1000) to avoids crashes. Best approach is to compare
height of first model level z1 and teta_1 between LES and 1D, and increase the timesteps until results
between the two models are comparable (might require a slitghly different tke_heat_flux between the two models
due to difference in vertical diffusion schemes and subgrid effects)
... Syntax for use is to add "tke_heat_flux = ..." in callphys.def

Corrected some stuff with tke transport in thermals (which is anyway deactivated but one day maybe...)

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