source: LMDZ5/branches/IPSLCM5A2.1_ISO/libf/phyiso/thermcell_main.F90 @ 3699

Last change on this file since 3699 was 3331, checked in by acozic, 7 years ago

Add modification for isotopes

  • Property svn:executable set to *
File size: 33.6 KB
Line 
1!
2! $Id: thermcell_main.F90 2387 2015-11-07 09:27:40Z fhourdin $
3!
4      SUBROUTINE thermcell_main(itap,ngrid,nlay,ptimestep  &
5     &                  ,pplay,pplev,pphi,debut  &
6     &                  ,pu,pv,pt,po  &
7     &                  ,pduadj,pdvadj,pdtadj,pdoadj  &
8     &                  ,fm0,entr0,detr0,zqta,zqla,lmax  &
9     &                  ,ratqscth,ratqsdiff,zqsatth  &
10     &                  ,Ale_bl,Alp_bl,lalim_conv,wght_th &
11     &                  ,zmax0, f0,zw2,fraca,ztv &
12     &                  ,zpspsk,ztla,zthl &
13!!! nrlmd le 10/04/2012
14     &                  ,pbl_tke,pctsrf,omega,airephy &
15     &                  ,zlcl,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &
16     &                  ,n2,s2,ale_bl_stat &
17     &                  ,therm_tke_max,env_tke_max &
18     &                  ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &
19     &                  ,alp_bl_conv,alp_bl_stat &
20!!! fin nrlmd le 10/04/2012
21     &                  ,ztva &
22#ifdef ISO         
23     &      ,xtpo,xtpdoadj &
24#ifdef DIAGISO         
25     &      ,zoa,xtzoa &
26#endif
27#endif         
28     &   )
29
30      USE dimphy
31      USE ioipsl
32      USE indice_sol_mod
33      USE print_control_mod, ONLY: lunout,prt_level
34#ifdef ISO
35  USE infotrac_phy, ONLY : ntraciso
36#ifdef ISOVERIF
37  USE isotopes_mod, ONLY : iso_eau,iso_HDO
38  USE isotopes_verif_mod, ONLY: iso_verif_egalite, &
39        iso_verif_aberrant_encadre
40#endif
41#endif
42      IMPLICIT NONE
43
44!=======================================================================
45!   Auteurs: Frederic Hourdin, Catherine Rio, Anne Mathieu
46!   Version du 09.02.07
47!   Calcul du transport vertical dans la couche limite en presence
48!   de "thermiques" explicitement representes avec processus nuageux
49!
50!   Reecriture a partir d'un listing papier a Habas, le 14/02/00
51!
52!   le thermique est suppose homogene et dissipe par melange avec
53!   son environnement. la longueur l_mix controle l'efficacite du
54!   melange
55!
56!   Le calcul du transport des differentes especes se fait en prenant
57!   en compte:
58!     1. un flux de masse montant
59!     2. un flux de masse descendant
60!     3. un entrainement
61!     4. un detrainement
62!
63! Modif 2013/01/04 (FH hourdin@lmd.jussieu.fr)
64!    Introduction of an implicit computation of vertical advection in
65!    the environment of thermal plumes in thermcell_dq
66!    impl =     0 : explicit, 1 : implicit, -1 : old version
67!    controled by iflag_thermals =
68!       15, 16 run with impl=-1 : numerical convergence with NPv3
69!       17, 18 run with impl=1  : more stable
70!    15 and 17 correspond to the activation of the stratocumulus "bidouille"
71!
72!=======================================================================
73
74
75!-----------------------------------------------------------------------
76!   declarations:
77!   -------------
78
79#include "YOMCST.h"
80#include "YOETHF.h"
81#include "FCTTRE.h"
82#include "thermcell.h"
83
84!   arguments:
85!   ----------
86
87!IM 140508
88      INTEGER itap
89
90      INTEGER ngrid,nlay
91      real ptimestep
92      REAL pt(ngrid,nlay),pdtadj(ngrid,nlay)
93      REAL pu(ngrid,nlay),pduadj(ngrid,nlay)
94      REAL pv(ngrid,nlay),pdvadj(ngrid,nlay)
95      REAL po(ngrid,nlay),pdoadj(ngrid,nlay)
96      REAL pplay(ngrid,nlay),pplev(ngrid,nlay+1)
97      real pphi(ngrid,nlay)
98      LOGICAL debut
99
100!   local:
101!   ------
102
103      integer icount
104
105      integer, save :: dvdq=1,dqimpl=-1
106!$OMP THREADPRIVATE(dvdq,dqimpl)
107      data icount/0/
108      save icount
109!$OMP THREADPRIVATE(icount)
110
111      integer,save :: igout=1
112!$OMP THREADPRIVATE(igout)
113      integer,save :: lunout1=6
114!$OMP THREADPRIVATE(lunout1)
115      integer,save :: lev_out=10
116!$OMP THREADPRIVATE(lev_out)
117
118      REAL susqr2pi, Reuler
119
120      INTEGER ig,k,l,ll,ierr
121      real zsortie1d(klon)
122      INTEGER lmax(klon),lmin(klon),lalim(klon)
123      INTEGER lmix(klon)
124      INTEGER lmix_bis(klon)
125      real linter(klon)
126      real zmix(klon)
127      real zmax(klon),zw2(klon,klev+1),ztva(klon,klev),zw_est(klon,klev+1),ztva_est(klon,klev)
128!      real fraca(klon,klev)
129
130      real zmax_sec(klon)
131!on garde le zmax du pas de temps precedent
132      real zmax0(klon)
133!FH/IM     save zmax0
134
135      real lambda
136
137      real zlev(klon,klev+1),zlay(klon,klev)
138      real deltaz(klon,klev)
139      REAL zh(klon,klev)
140      real zthl(klon,klev),zdthladj(klon,klev)
141      REAL ztv(klon,klev)
142      real zu(klon,klev),zv(klon,klev),zo(klon,klev)
143      real zl(klon,klev)
144      real zsortie(klon,klev)
145      real zva(klon,klev)
146      real zua(klon,klev)
147      real zoa(klon,klev)
148
149      real zta(klon,klev)
150      real zha(klon,klev)
151      real fraca(klon,klev+1)
152      real zf,zf2
153      real thetath2(klon,klev),wth2(klon,klev),wth3(klon,klev)
154      real q2(klon,klev)
155! FH probleme de dimensionnement avec l'allocation dynamique
156!     common/comtherm/thetath2,wth2
157      real wq(klon,klev)
158      real wthl(klon,klev)
159      real wthv(klon,klev)
160   
161      real ratqscth(klon,klev)
162      real var
163      real vardiff
164      real ratqsdiff(klon,klev)
165
166      logical sorties
167      real rho(klon,klev),rhobarz(klon,klev),masse(klon,klev)
168      real zpspsk(klon,klev)
169
170      real wmax(klon)
171      real wmax_tmp(klon)
172      real wmax_sec(klon)
173      real fm0(klon,klev+1),entr0(klon,klev),detr0(klon,klev)
174      real fm(klon,klev+1),entr(klon,klev),detr(klon,klev)
175
176      real ztla(klon,klev),zqla(klon,klev),zqta(klon,klev)
177!niveau de condensation
178      integer nivcon(klon)
179      real zcon(klon)
180      REAL CHI
181      real zcon2(klon)
182      real pcon(klon)
183      real zqsat(klon,klev)
184      real zqsatth(klon,klev)
185
186      real f_star(klon,klev+1),entr_star(klon,klev)
187      real detr_star(klon,klev)
188      real alim_star_tot(klon)
189      real alim_star(klon,klev)
190      real alim_star_clos(klon,klev)
191      real f(klon), f0(klon)
192!FH/IM     save f0
193      real zlevinter(klon)
194       real seuil
195      real csc(klon,klev)
196
197!!! nrlmd le 10/04/2012
198
199!------Entrées
200      real pbl_tke(klon,klev+1,nbsrf)
201      real pctsrf(klon,nbsrf)
202      real omega(klon,klev)
203      real airephy(klon)
204!------Sorties
205      real zlcl(klon),fraca0(klon),w0(klon),w_conv(klon)
206      real therm_tke_max0(klon),env_tke_max0(klon)
207      real n2(klon),s2(klon)
208      real ale_bl_stat(klon)
209      real therm_tke_max(klon,klev),env_tke_max(klon,klev)
210      real alp_bl_det(klon),alp_bl_fluct_m(klon),alp_bl_fluct_tke(klon),alp_bl_conv(klon),alp_bl_stat(klon)
211!------Local
212      integer nsrf
213      real rhobarz0(klon)                    ! Densité au LCL
214      logical ok_lcl(klon)                   ! Existence du LCL des thermiques
215      integer klcl(klon)                     ! Niveau du LCL
216      real interp(klon)                      ! Coef d'interpolation pour le LCL
217!--Triggering
218      real Su                                ! Surface unité: celle d'un updraft élémentaire
219      parameter(Su=4e4)
220      real hcoef                             ! Coefficient directeur pour le calcul de s2
221      parameter(hcoef=1)
222      real hmincoef                          ! Coefficient directeur pour l'ordonnée à l'origine pour le calcul de s2
223      parameter(hmincoef=0.3)
224      real eps1                              ! Fraction de surface occupée par la population 1 : eps1=n1*s1/(fraca0*Sd)
225      parameter(eps1=0.3)
226      real hmin(ngrid)                       ! Ordonnée à l'origine pour le calcul de s2
227      real zmax_moy(ngrid)                   ! Hauteur moyenne des thermiques : zmax_moy = zlcl + 0.33 (zmax-zlcl)
228      real zmax_moy_coef
229      parameter(zmax_moy_coef=0.33)
230      real depth(klon)                       ! Epaisseur moyenne du cumulus
231      real w_max(klon)                       ! Vitesse max statistique
232      real s_max(klon)
233!--Closure
234      real pbl_tke_max(klon,klev)            ! Profil de TKE moyenne
235      real pbl_tke_max0(klon)                ! TKE moyenne au LCL
236      real w_ls(klon,klev)                   ! Vitesse verticale grande échelle (m/s)
237      real coef_m                            ! On considère un rendement pour alp_bl_fluct_m
238      parameter(coef_m=1.)
239      real coef_tke                          ! On considère un rendement pour alp_bl_fluct_tke
240      parameter(coef_tke=1.)
241
242!!! fin nrlmd le 10/04/2012
243
244!
245      !nouvelles variables pour la convection
246      real Ale_bl(klon)
247      real Alp_bl(klon)
248      real alp_int(klon),dp_int(klon),zdp
249      real ale_int(klon)
250      integer n_int(klon)
251      real fm_tot(klon)
252      real wght_th(klon,klev)
253      integer lalim_conv(klon)
254!v1d     logical therm
255!v1d     save therm
256
257      character*2 str2
258      character*10 str10
259
260      character (len=20) :: modname='thermcell_main'
261      character (len=80) :: abort_message
262
263      EXTERNAL SCOPY
264
265#ifdef ISO
266      REAL xtpo(ntraciso,ngrid,nlay),xtpdoadj(ntraciso,ngrid,nlay)
267      REAL xtzo(ntraciso,klon,klev)
268      real xtzoa(ntraciso,klon,klev)
269      REAL xtpdoadj_tmp(ngrid,nlay)
270      REAL xtpo_tmp(klon,klev)
271      real xtzoa_tmp(klon,klev)
272      integer ixt
273      real xtzl(ntraciso,klon,klev)
274#endif
275!
276
277!-----------------------------------------------------------------------
278!   initialisation:
279!   ---------------
280!
281
282   seuil=0.25
283
284   if (debut) then
285      if (iflag_thermals==15.or.iflag_thermals==16) then
286         dvdq=0
287         dqimpl=-1
288      else
289         dvdq=1
290         dqimpl=1
291      endif
292
293      fm0=0.
294      entr0=0.
295      detr0=0.
296   endif
297   fm=0. ; entr=0. ; detr=0.
298   icount=icount+1
299
300!IM 090508 beg
301!print*,'====================================================================='
302!print*,'====================================================================='
303!print*,' PAS ',icount,' PAS ',icount,' PAS ',icount,' PAS ',icount
304!print*,'====================================================================='
305!print*,'====================================================================='
306!IM 090508 end
307
308      if (prt_level.ge.1) print*,'thermcell_main V4'
309
310       sorties=.true.
311      IF(ngrid.NE.klon) THEN
312         PRINT*
313         PRINT*,'STOP dans convadj'
314         PRINT*,'ngrid    =',ngrid
315         PRINT*,'klon  =',klon
316      ENDIF
317!
318!     write(lunout,*)'WARNING thermcell_main f0=max(f0,1.e-2)'
319     do ig=1,klon
320         f0(ig)=max(f0(ig),1.e-2)
321         zmax0(ig)=max(zmax0(ig),40.)
322!IMmarche pas ?!       if (f0(ig)<1.e-2) f0(ig)=1.e-2
323     enddo
324
325      if (prt_level.ge.20) then
326       do ig=1,ngrid
327          print*,'th_main ig f0',ig,f0(ig)
328       enddo
329      endif
330!-----------------------------------------------------------------------
331! Calcul de T,q,ql a partir de Tl et qT dans l environnement
332!   --------------------------------------------------------------------
333!
334      CALL thermcell_env(ngrid,nlay,po,pt,pu,pv,pplay,  &
335     &           pplev,zo,zh,zl,ztv,zthl,zu,zv,zpspsk,zqsat,lev_out &
336#ifdef ISO
337     &            ,xtpo,xtzo,xtzl &
338#endif     
339     &   )
340       
341      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_env'
342
343!------------------------------------------------------------------------
344!                       --------------------
345!
346!
347!                       + + + + + + + + + + +
348!
349!
350!  wa, fraca, wd, fracd --------------------   zlev(2), rhobarz
351!  wh,wt,wo ...
352!
353!                       + + + + + + + + + + +  zh,zu,zv,zo,rho
354!
355!
356!                       --------------------   zlev(1)
357!                       \\\\\\\\\\\\\\\\\\\\
358!
359!
360
361!-----------------------------------------------------------------------
362!   Calcul des altitudes des couches
363!-----------------------------------------------------------------------
364
365      do l=2,nlay
366         zlev(:,l)=0.5*(pphi(:,l)+pphi(:,l-1))/RG
367      enddo
368         zlev(:,1)=0.
369         zlev(:,nlay+1)=(2.*pphi(:,klev)-pphi(:,klev-1))/RG
370      do l=1,nlay
371         zlay(:,l)=pphi(:,l)/RG
372      enddo
373!calcul de l epaisseur des couches
374      do l=1,nlay
375         deltaz(:,l)=zlev(:,l+1)-zlev(:,l)
376      enddo
377
378!     print*,'2 OK convect8'
379!-----------------------------------------------------------------------
380!   Calcul des densites
381!-----------------------------------------------------------------------
382
383     rho(:,:)=pplay(:,:)/(zpspsk(:,:)*RD*ztv(:,:))
384
385     if (prt_level.ge.10)write(lunout,*)                                &
386    &    'WARNING thermcell_main rhobarz(:,1)=rho(:,1)'
387      rhobarz(:,1)=rho(:,1)
388
389      do l=2,nlay
390         rhobarz(:,l)=0.5*(rho(:,l)+rho(:,l-1))
391      enddo
392
393!calcul de la masse
394      do l=1,nlay
395         masse(:,l)=(pplev(:,l)-pplev(:,l+1))/RG
396      enddo
397
398      if (prt_level.ge.1) print*,'thermcell_main apres initialisation'
399
400!------------------------------------------------------------------
401!
402!             /|\
403!    --------  |  F_k+1 -------   
404!                              ----> D_k
405!             /|\              <---- E_k , A_k
406!    --------  |  F_k ---------
407!                              ----> D_k-1
408!                              <---- E_k-1 , A_k-1
409!
410!
411!
412!
413!
414!    ---------------------------
415!
416!    ----- F_lmax+1=0 ----------         \
417!            lmax     (zmax)              |
418!    ---------------------------          |
419!                                         |
420!    ---------------------------          |
421!                                         |
422!    ---------------------------          |
423!                                         |
424!    ---------------------------          |
425!                                         |
426!    ---------------------------          |
427!                                         |  E
428!    ---------------------------          |  D
429!                                         |
430!    ---------------------------          |
431!                                         |
432!    ---------------------------  \       |
433!            lalim                 |      |
434!    ---------------------------   |      |
435!                                  |      |
436!    ---------------------------   |      |
437!                                  | A    |
438!    ---------------------------   |      |
439!                                  |      |
440!    ---------------------------   |      |
441!    lmin  (=1 pour le moment)     |      |
442!    ----- F_lmin=0 ------------  /      /
443!
444!    ---------------------------
445!    //////////////////////////
446!
447!
448!=============================================================================
449!  Calculs initiaux ne faisant pas intervenir les changements de phase
450!=============================================================================
451
452!------------------------------------------------------------------
453!  1. alim_star est le profil vertical de l'alimentation a la base du
454!     panache thermique, calcule a partir de la flotabilite de l'air sec
455!  2. lmin et lalim sont les indices inferieurs et superieurs de alim_star
456!------------------------------------------------------------------
457!
458      entr_star=0. ; detr_star=0. ; alim_star=0. ; alim_star_tot=0.
459      lmin=1
460
461!-----------------------------------------------------------------------------
462!  3. wmax_sec et zmax_sec sont les vitesses et altitudes maximum d'un
463!     panache sec conservatif (e=d=0) alimente selon alim_star
464!     Il s'agit d'un calcul de type CAPE
465!     zmax_sec est utilise pour determiner la geometrie du thermique.
466!------------------------------------------------------------------------------
467!---------------------------------------------------------------------------------
468!calcul du melange et des variables dans le thermique
469!--------------------------------------------------------------------------------
470!
471      if (prt_level.ge.1) print*,'avant thermcell_plume ',lev_out
472!IM 140508   CALL thermcell_plume(ngrid,nlay,ptimestep,ztv,zthl,po,zl,rhobarz,  &
473
474! Gestion temporaire de plusieurs appels à thermcell_plume au travers
475! de la variable iflag_thermals
476
477!      print*,'THERM thermcell_main iflag_thermals_ed=',iflag_thermals_ed
478      if (iflag_thermals_ed<=9) then
479!         print*,'THERM NOUVELLE/NOUVELLE Arnaud'
480         CALL thermcell_plume(itap,ngrid,nlay,ptimestep,ztv,zthl,po,zl,rhobarz,&
481     &    zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot,  &
482     &    lalim,f0,detr_star,entr_star,f_star,csc,ztva,  &
483     &    ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter &
484     &    ,lev_out,lunout1,igout)
485
486      elseif (iflag_thermals_ed>9) then
487!        print*,'THERM RIO et al 2010, version d Arnaud'
488         CALL thermcellV1_plume(itap,ngrid,nlay,ptimestep,ztv,zthl,po,zl,rhobarz,&
489     &    zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot,  &
490     &    lalim,f0,detr_star,entr_star,f_star,csc,ztva,  &
491     &    ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter &
492     &    ,lev_out,lunout1,igout)
493
494      endif
495
496      if (prt_level.ge.1) print*,'apres thermcell_plume ',lev_out
497
498      call test_ltherm(ngrid,nlay,pplev,pplay,lalim,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_plum lalim ')
499      call test_ltherm(ngrid,nlay,pplev,pplay,lmix ,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_plum lmix  ')
500
501      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_plume'
502      if (prt_level.ge.10) then
503         write(lunout1,*) 'Dans thermcell_main 2'
504         write(lunout1,*) 'lmin ',lmin(igout)
505         write(lunout1,*) 'lalim ',lalim(igout)
506         write(lunout1,*) ' ig l alim_star entr_star detr_star f_star '
507         write(lunout1,'(i6,i4,4e15.5)') (igout,l,alim_star(igout,l),entr_star(igout,l),detr_star(igout,l) &
508     &    ,f_star(igout,l+1),l=1,nint(linter(igout))+5)
509      endif
510
511!-------------------------------------------------------------------------------
512! Calcul des caracteristiques du thermique:zmax,zmix,wmax
513!-------------------------------------------------------------------------------
514!
515      CALL thermcell_height(ngrid,nlay,lalim,lmin,linter,lmix,zw2,  &
516     &           zlev,lmax,zmax,zmax0,zmix,wmax,lev_out)
517! Attention, w2 est transforme en sa racine carree dans cette routine
518! Le probleme vient du fait que linter et lmix sont souvent égaux à 1.
519      wmax_tmp=0.
520      do  l=1,nlay
521         wmax_tmp(:)=max(wmax_tmp(:),zw2(:,l))
522      enddo
523!     print*,"ZMAX ",lalim,lmin,linter,lmix,lmax,zmax,zmax0,zmix,wmax
524
525
526
527      call test_ltherm(ngrid,nlay,pplev,pplay,lalim,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lalim ')
528      call test_ltherm(ngrid,nlay,pplev,pplay,lmin ,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lmin  ')
529      call test_ltherm(ngrid,nlay,pplev,pplay,lmix ,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lmix  ')
530      call test_ltherm(ngrid,nlay,pplev,pplay,lmax ,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lmax  ')
531
532      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_height'
533
534!-------------------------------------------------------------------------------
535! Fermeture,determination de f
536!-------------------------------------------------------------------------------
537!
538!
539!!      write(lunout,*)'THERM NOUVEAU XXXXX'
540      CALL thermcell_dry(ngrid,nlay,zlev,pphi,ztv,alim_star,  &
541    &                      lalim,lmin,zmax_sec,wmax_sec,lev_out)
542
543 
544call test_ltherm(ngrid,nlay,pplev,pplay,lmin,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_dry  lmin  ')
545call test_ltherm(ngrid,nlay,pplev,pplay,lalim,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_dry  lalim ')
546
547      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_dry'
548      if (prt_level.ge.10) then
549         write(lunout1,*) 'Dans thermcell_main 1b'
550         write(lunout1,*) 'lmin ',lmin(igout)
551         write(lunout1,*) 'lalim ',lalim(igout)
552         write(lunout1,*) ' ig l alim_star entr_star detr_star f_star '
553         write(lunout1,'(i6,i4,e15.5)') (igout,l,alim_star(igout,l) &
554     &    ,l=1,lalim(igout)+4)
555      endif
556
557
558
559
560! Choix de la fonction d'alimentation utilisee pour la fermeture.
561! Apparemment sans importance
562      alim_star_clos(:,:)=alim_star(:,:)
563      alim_star_clos(:,:)=entr_star(:,:)+alim_star(:,:)
564!
565!CR Appel de la fermeture seche
566      if (iflag_thermals_closure.eq.1) then
567
568      CALL thermcell_closure(ngrid,nlay,r_aspect_thermals,ptimestep,rho,  &
569     &   zlev,lalim,alim_star_clos,f_star,zmax_sec,wmax_sec,f,lev_out)
570
571!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
572! Appel avec les zmax et wmax tenant compte de la condensation
573! Semble moins bien marcher
574     else if (iflag_thermals_closure.eq.2) then
575
576     CALL thermcell_closure(ngrid,nlay,r_aspect_thermals,ptimestep,rho,  &
577    &   zlev,lalim,alim_star,f_star,zmax,wmax,f,lev_out)
578
579     endif
580
581!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
582
583      if(prt_level.ge.1)print*,'thermcell_closure apres thermcell_closure'
584
585      if (tau_thermals>1.) then
586         lambda=exp(-ptimestep/tau_thermals)
587         f0=(1.-lambda)*f+lambda*f0
588      else
589         f0=f
590      endif
591
592! Test valable seulement en 1D mais pas genant
593      if (.not. (f0(1).ge.0.) ) then
594              abort_message = '.not. (f0(1).ge.0.)'
595              CALL abort_physic (modname,abort_message,1)
596      endif
597
598!-------------------------------------------------------------------------------
599!deduction des flux
600!-------------------------------------------------------------------------------
601
602      CALL thermcell_flux2(ngrid,nlay,ptimestep,masse, &
603     &       lalim,lmax,alim_star,  &
604     &       entr_star,detr_star,f,rhobarz,zlev,zw2,fm,entr,  &
605     &       detr,zqla,lev_out,lunout1,igout)
606!IM 060508    &       detr,zqla,zmax,lev_out,lunout,igout)
607
608      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_flux'
609      call test_ltherm(ngrid,nlay,pplev,pplay,lalim,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_flux lalim ')
610      call test_ltherm(ngrid,nlay,pplev,pplay,lmax ,seuil,ztv,po,ztva,zqla,f_star,zw2,'thermcell_flux lmax  ')
611
612!------------------------------------------------------------------
613!   On ne prend pas directement les profils issus des calculs precedents
614!   mais on s'autorise genereusement une relaxation vers ceci avec
615!   une constante de temps tau_thermals (typiquement 1800s).
616!------------------------------------------------------------------
617
618      if (tau_thermals>1.) then
619         lambda=exp(-ptimestep/tau_thermals)
620         fm0=(1.-lambda)*fm+lambda*fm0
621         entr0=(1.-lambda)*entr+lambda*entr0
622         detr0=(1.-lambda)*detr+lambda*detr0
623      else
624         fm0=fm
625         entr0=entr
626         detr0=detr
627      endif
628
629!c------------------------------------------------------------------
630!   calcul du transport vertical
631!------------------------------------------------------------------
632
633      call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse,  &
634     &                    zthl,zdthladj,zta,lev_out)
635#ifdef ISO
636      call thermcell_dq_iso(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse,  &
637     &                   po,pdoadj,zoa,lev_out,xtpo,xtpdoadj,xtzoa)
638#else
639      call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse,  &
640     &                   po,pdoadj,zoa,lev_out)
641#endif
642
643#ifdef ISO     
644#ifdef ISOVERIF
645      DO  ll=1,nlay
646        DO ig=1,ngrid
647          if (iso_eau.gt.0) then
648              call iso_verif_egalite(xtpo(iso_eau,ig,ll), &
649     &          po(ig,ll),'thermcell_main 594')
650              call iso_verif_egalite(xtpdoadj(iso_eau,ig,ll), &
651     &          pdoadj(ig,ll),'thermcell_main 596')
652          endif
653          if (iso_HDO.gt.0) then
654              call iso_verif_aberrant_encadre(xtpo(iso_hdo,ig,ll) &
655     &           /po(ig,ll),'thermcell_main 610')
656          endif
657        enddo
658      enddo !DO  ll=1,nlay
659      write(*,*) 'thermcell_main 600 tmp: apres thermcell_dq'
660#endif     
661#endif   
662
663!------------------------------------------------------------------
664! Calcul de la fraction de l'ascendance
665!------------------------------------------------------------------
666      do ig=1,klon
667         fraca(ig,1)=0.
668         fraca(ig,nlay+1)=0.
669      enddo
670      do l=2,nlay
671         do ig=1,klon
672            if (zw2(ig,l).gt.1.e-10) then
673            fraca(ig,l)=fm(ig,l)/(rhobarz(ig,l)*zw2(ig,l))
674            else
675            fraca(ig,l)=0.
676            endif
677         enddo
678      enddo
679     
680!------------------------------------------------------------------
681!  calcul du transport vertical du moment horizontal
682!------------------------------------------------------------------
683
684!IM 090508 
685      if (dvdq == 0 ) then
686
687! Calcul du transport de V tenant compte d'echange par gradient
688! de pression horizontal avec l'environnement
689
690         call thermcell_dv2(ngrid,nlay,ptimestep,fm0,entr0,masse  &
691!    &    ,fraca*dvdq,zmax &
692     &    ,fraca,zmax &
693     &    ,zu,zv,pduadj,pdvadj,zua,zva,lev_out)
694
695      else
696
697! calcul purement conservatif pour le transport de V
698         call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse  &
699     &    ,zu,pduadj,zua,lev_out)
700         call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse  &
701     &    ,zv,pdvadj,zva,lev_out)
702
703      endif
704
705!     print*,'13 OK convect8'
706      do l=1,nlay
707         do ig=1,ngrid
708           pdtadj(ig,l)=zdthladj(ig,l)*zpspsk(ig,l) 
709         enddo
710      enddo
711
712      if (prt_level.ge.1) print*,'14 OK convect8'
713!------------------------------------------------------------------
714!   Calculs de diagnostiques pour les sorties
715!------------------------------------------------------------------
716!calcul de fraca pour les sorties
717     
718      if (sorties) then
719      if (prt_level.ge.1) print*,'14a OK convect8'
720! calcul du niveau de condensation
721! initialisation
722      do ig=1,ngrid
723         nivcon(ig)=0
724         zcon(ig)=0.
725      enddo
726!nouveau calcul
727      do ig=1,ngrid
728      CHI=zh(ig,1)/(1669.0-122.0*zo(ig,1)/zqsat(ig,1)-zh(ig,1))
729      pcon(ig)=pplay(ig,1)*(zo(ig,1)/zqsat(ig,1))**CHI
730      enddo
731!IM   do k=1,nlay
732      do k=1,nlay-1
733         do ig=1,ngrid
734         if ((pcon(ig).le.pplay(ig,k))  &
735     &      .and.(pcon(ig).gt.pplay(ig,k+1))) then
736            zcon2(ig)=zlay(ig,k)-(pcon(ig)-pplay(ig,k))/(RG*rho(ig,k))/100.
737         endif
738         enddo
739      enddo
740!IM
741      ierr=0
742      do ig=1,ngrid
743        if (pcon(ig).le.pplay(ig,nlay)) then
744           zcon2(ig)=zlay(ig,nlay)-(pcon(ig)-pplay(ig,nlay))/(RG*rho(ig,nlay))/100.
745           ierr=1
746        endif
747      enddo
748      if (ierr==1) then
749           abort_message = 'thermcellV0_main: les thermiques vont trop haut '
750           CALL abort_physic (modname,abort_message,1)
751      endif
752
753      if (prt_level.ge.1) print*,'14b OK convect8'
754      do k=nlay,1,-1
755         do ig=1,ngrid
756            if (zqla(ig,k).gt.1e-10) then
757               nivcon(ig)=k
758               zcon(ig)=zlev(ig,k)
759            endif
760         enddo
761      enddo
762      if (prt_level.ge.1) print*,'14c OK convect8'
763!calcul des moments
764!initialisation
765      do l=1,nlay
766         do ig=1,ngrid
767            q2(ig,l)=0.
768            wth2(ig,l)=0.
769            wth3(ig,l)=0.
770            ratqscth(ig,l)=0.
771            ratqsdiff(ig,l)=0.
772         enddo
773      enddo     
774      if (prt_level.ge.1) print*,'14d OK convect8'
775      if (prt_level.ge.10)write(lunout,*)                                &
776    &     'WARNING thermcell_main wth2=0. si zw2 > 1.e-10'
777      do l=1,nlay
778         do ig=1,ngrid
779            zf=fraca(ig,l)
780            zf2=zf/(1.-zf)
781!
782            thetath2(ig,l)=zf2*(ztla(ig,l)-zthl(ig,l))**2
783            if(zw2(ig,l).gt.1.e-10) then
784             wth2(ig,l)=zf2*(zw2(ig,l))**2
785            else
786             wth2(ig,l)=0.
787            endif
788            wth3(ig,l)=zf2*(1-2.*fraca(ig,l))/(1-fraca(ig,l))  &
789     &                *zw2(ig,l)*zw2(ig,l)*zw2(ig,l)
790            q2(ig,l)=zf2*(zqta(ig,l)*1000.-po(ig,l)*1000.)**2
791!test: on calcul q2/po=ratqsc
792            ratqscth(ig,l)=sqrt(max(q2(ig,l),1.e-6)/(po(ig,l)*1000.))
793         enddo
794      enddo
795!calcul des flux: q, thetal et thetav
796      do l=1,nlay
797         do ig=1,ngrid
798      wq(ig,l)=fraca(ig,l)*zw2(ig,l)*(zqta(ig,l)*1000.-po(ig,l)*1000.)
799      wthl(ig,l)=fraca(ig,l)*zw2(ig,l)*(ztla(ig,l)-zthl(ig,l))
800      wthv(ig,l)=fraca(ig,l)*zw2(ig,l)*(ztva(ig,l)-ztv(ig,l))
801         enddo
802      enddo
803!
804! $Id: thermcell_main.F90 2387 2015-11-07 09:27:40Z fhourdin $
805!
806      CALL thermcell_alp(ngrid,nlay,ptimestep  &
807     &                  ,pplay,pplev  &
808     &                  ,fm0,entr0,lmax  &
809     &                  ,Ale_bl,Alp_bl,lalim_conv,wght_th &
810     &                  ,zw2,fraca &
811!!! necessire en plus
812     &                  ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax &
813!!! nrlmd le 10/04/2012
814     &                  ,pbl_tke,pctsrf,omega,airephy &
815     &                  ,zlcl,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &
816     &                  ,n2,s2,ale_bl_stat &
817     &                  ,therm_tke_max,env_tke_max &
818     &                  ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &
819     &                  ,alp_bl_conv,alp_bl_stat &
820!!! fin nrlmd le 10/04/2012
821     &                   )
822
823
824
825!calcul du ratqscdiff
826      if (prt_level.ge.1) print*,'14e OK convect8'
827      var=0.
828      vardiff=0.
829      ratqsdiff(:,:)=0.
830
831      do l=1,klev
832         do ig=1,ngrid
833            if (l<=lalim(ig)) then
834            var=var+alim_star(ig,l)*zqta(ig,l)*1000.
835            endif
836         enddo
837      enddo
838
839      if (prt_level.ge.1) print*,'14f OK convect8'
840
841      do l=1,klev
842         do ig=1,ngrid
843            if (l<=lalim(ig)) then
844               zf=fraca(ig,l)
845               zf2=zf/(1.-zf)
846               vardiff=vardiff+alim_star(ig,l)*(zqta(ig,l)*1000.-var)**2
847            endif
848         enddo
849      enddo
850
851      if (prt_level.ge.1) print*,'14g OK convect8'
852      do l=1,nlay
853         do ig=1,ngrid
854            ratqsdiff(ig,l)=sqrt(vardiff)/(po(ig,l)*1000.)   
855!           write(11,*)'ratqsdiff=',ratqsdiff(ig,l)
856         enddo
857      enddo
858!--------------------------------------------------------------------   
859!
860!ecriture des fichiers sortie
861!     print*,'15 OK convect8 CCCCCCCCCCCCCCCCCCc'
862
863      endif
864
865      if (prt_level.ge.1) print*,'thermcell_main FIN  OK'
866
867      return
868      end
869
870!-----------------------------------------------------------------------------
871
872      subroutine test_ltherm(klon,klev,pplev,pplay,long,seuil,ztv,po,ztva,zqla,f_star,zw2,comment)
873      USE print_control_mod, ONLY: prt_level
874      IMPLICIT NONE
875
876      integer i, k, klon,klev
877      real pplev(klon,klev+1),pplay(klon,klev)
878      real ztv(klon,klev)
879      real po(klon,klev)
880      real ztva(klon,klev)
881      real zqla(klon,klev)
882      real f_star(klon,klev)
883      real zw2(klon,klev)
884      integer long(klon)
885      real seuil
886      character*21 comment
887
888      if (prt_level.ge.1) THEN
889       print*,'WARNING !!! TEST ',comment
890      endif
891      return
892
893!  test sur la hauteur des thermiques ...
894         do i=1,klon
895!IMtemp           if (pplay(i,long(i)).lt.seuil*pplev(i,1)) then
896           if (prt_level.ge.10) then
897               print*,'WARNING ',comment,' au point ',i,' K= ',long(i)
898               print*,'  K  P(MB)  THV(K)     Qenv(g/kg)THVA        QLA(g/kg)   F*        W2'
899               do k=1,klev
900                  write(6,'(i3,7f10.3)') k,pplay(i,k),ztv(i,k),1000*po(i,k),ztva(i,k),1000*zqla(i,k),f_star(i,k),zw2(i,k)
901               enddo
902           endif
903         enddo
904
905
906      return
907      end
908
909!!! nrlmd le 10/04/2012                          Transport de la TKE par le thermique moyen pour la fermeture en ALP
910!                                                         On transporte pbl_tke pour donner therm_tke
911!                                          Copie conforme de la subroutine DTKE dans physiq.F écrite par Frederic Hourdin
912      subroutine thermcell_tke_transport(ngrid,nlay,ptimestep,fm0,entr0,  &
913     &           rg,pplev,therm_tke_max)
914      USE print_control_mod, ONLY: prt_level
915      implicit none
916
917!=======================================================================
918!
919!   Calcul du transport verticale dans la couche limite en presence
920!   de "thermiques" explicitement representes
921!   calcul du dq/dt une fois qu'on connait les ascendances
922!
923!=======================================================================
924
925      integer ngrid,nlay,nsrf
926
927      real ptimestep
928      real masse0(ngrid,nlay),fm0(ngrid,nlay+1),pplev(ngrid,nlay+1)
929      real entr0(ngrid,nlay),rg
930      real therm_tke_max(ngrid,nlay)
931      real detr0(ngrid,nlay)
932
933
934      real masse(ngrid,nlay),fm(ngrid,nlay+1)
935      real entr(ngrid,nlay)
936      real q(ngrid,nlay)
937      integer lev_out                           ! niveau pour les print
938
939      real qa(ngrid,nlay),detr(ngrid,nlay),wqd(ngrid,nlay+1)
940
941      real zzm
942
943      integer ig,k
944      integer isrf
945
946
947      lev_out=0
948
949
950      if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0'
951
952!   calcul du detrainement
953      do k=1,nlay
954         detr0(:,k)=fm0(:,k)-fm0(:,k+1)+entr0(:,k)
955         masse0(:,k)=(pplev(:,k)-pplev(:,k+1))/RG
956      enddo
957
958
959! Decalage vertical des entrainements et detrainements.
960      masse(:,1)=0.5*masse0(:,1)
961      entr(:,1)=0.5*entr0(:,1)
962      detr(:,1)=0.5*detr0(:,1)
963      fm(:,1)=0.
964      do k=1,nlay-1
965         masse(:,k+1)=0.5*(masse0(:,k)+masse0(:,k+1))
966         entr(:,k+1)=0.5*(entr0(:,k)+entr0(:,k+1))
967         detr(:,k+1)=0.5*(detr0(:,k)+detr0(:,k+1))
968         fm(:,k+1)=fm(:,k)+entr(:,k)-detr(:,k)
969      enddo
970      fm(:,nlay+1)=0.
971
972!!! nrlmd le 16/09/2010
973!   calcul de la valeur dans les ascendances
974!       do ig=1,ngrid
975!          qa(ig,1)=q(ig,1)
976!       enddo
977!!!
978
979!do isrf=1,nsrf
980
981!   q(:,:)=therm_tke(:,:,isrf)
982   q(:,:)=therm_tke_max(:,:)
983!!! nrlmd le 16/09/2010
984      do ig=1,ngrid
985         qa(ig,1)=q(ig,1)
986      enddo
987!!!
988
989    if (1==1) then
990      do k=2,nlay
991         do ig=1,ngrid
992            if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt.  &
993     &         1.e-5*masse(ig,k)) then
994         qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k))  &
995     &         /(fm(ig,k+1)+detr(ig,k))
996            else
997               qa(ig,k)=q(ig,k)
998            endif
999            if (qa(ig,k).lt.0.) then
1000!               print*,'qa<0!!!'
1001            endif
1002            if (q(ig,k).lt.0.) then
1003!               print*,'q<0!!!'
1004            endif
1005         enddo
1006      enddo
1007
1008! Calcul du flux subsident
1009
1010      do k=2,nlay
1011         do ig=1,ngrid
1012            wqd(ig,k)=fm(ig,k)*q(ig,k)
1013            if (wqd(ig,k).lt.0.) then
1014!               print*,'wqd<0!!!'
1015            endif
1016         enddo
1017      enddo
1018      do ig=1,ngrid
1019         wqd(ig,1)=0.
1020         wqd(ig,nlay+1)=0.
1021      enddo
1022
1023! Calcul des tendances
1024      do k=1,nlay
1025         do ig=1,ngrid
1026            q(ig,k)=q(ig,k)+(detr(ig,k)*qa(ig,k)-entr(ig,k)*q(ig,k)  &
1027     &               -wqd(ig,k)+wqd(ig,k+1))  &
1028     &               *ptimestep/masse(ig,k)
1029         enddo
1030      enddo
1031
1032 endif
1033
1034   therm_tke_max(:,:)=q(:,:)
1035
1036      return
1037!!! fin nrlmd le 10/04/2012
1038     end
1039
Note: See TracBrowser for help on using the repository browser.