source: LMDZ6/trunk/libf/phylmdiso/thermcell_main.F90 @ 3940

Last change on this file since 3940 was 3940, checked in by crisi, 3 years ago

replace files by symbloic liks from phylmdiso towards phylmd.
Many files at once

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