source: LMDZ6/trunk/libf/phylmd/lmdz_thermcell_main.F90 @ 4660

Last change on this file since 4660 was 4590, checked in by fhourdin, 18 months ago

Passage des thermiques a la nouvelle norme.

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 31.0 KB
Line 
1MODULE lmdz_thermcell_main
2! $Id: lmdz_thermcell_main.F90 4590 2023-06-29 01:03:15Z fairhead $
3!
4CONTAINS
5
6      subroutine thermcell_main(itap,ngrid,nlay,ptimestep  &
7     &                  ,pplay,pplev,pphi,debut  &
8     &                  ,pu,pv,pt,po  &
9     &                  ,pduadj,pdvadj,pdtadj,pdoadj  &
10     &                  ,fm0,entr0,detr0,zqta,zqla,lmax  &
11     &                  ,ratqscth,ratqsdiff,zqsatth  &
12     &                  ,zmax0, f0,zw2,fraca,ztv &
13     &                  ,zpspsk,ztla,zthl,ztva &
14     &                  ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax &
15#ifdef ISO         
16     &      ,xtpo,xtpdoadj &
17#endif         
18     &   )
19
20
21      USE lmdz_thermcell_ini, ONLY: thermcell_ini,dqimpl,dvdq,prt_level,lunout,prt_level
22      USE lmdz_thermcell_ini, ONLY: iflag_thermals_closure,iflag_thermals_ed,tau_thermals,r_aspect_thermals
23      USE lmdz_thermcell_ini, ONLY: iflag_thermals_down,fact_thermals_down
24      USE lmdz_thermcell_ini, ONLY: RD,RG
25
26      USE lmdz_thermcell_down, ONLY: thermcell_updown_dq
27      USE lmdz_thermcell_closure, ONLY: thermcell_closure
28      USE lmdz_thermcell_dq, ONLY: thermcell_dq
29      USE lmdz_thermcell_dry, ONLY: thermcell_dry
30      USE lmdz_thermcell_dv2, ONLY: thermcell_dv2
31      USE lmdz_thermcell_env, ONLY: thermcell_env
32      USE lmdz_thermcell_flux2, ONLY: thermcell_flux2
33      USE lmdz_thermcell_height, ONLY: thermcell_height
34      USE lmdz_thermcell_plume, ONLY: thermcell_plume
35      USE lmdz_thermcell_plume_6A, ONLY: thermcell_plume_6A,thermcell_plume_5B
36
37#ifdef ISO
38  USE infotrac_phy, ONLY : ntiso
39#ifdef ISOVERIF
40  USE isotopes_mod, ONLY : iso_eau,iso_HDO
41  USE isotopes_verif_mod, ONLY: iso_verif_egalite, &
42        iso_verif_aberrant_encadre
43#endif
44#endif
45
46
47      IMPLICIT NONE
48
49!=======================================================================
50!   Auteurs: Frederic Hourdin, Catherine Rio, Anne Mathieu
51!   Version du 09.02.07
52!   Calcul du transport vertical dans la couche limite en presence
53!   de "thermiques" explicitement representes avec processus nuageux
54!
55!   Reecriture a partir d'un listing papier a Habas, le 14/02/00
56!
57!   le thermique est suppose homogene et dissipe par melange avec
58!   son environnement. la longueur l_mix controle l'efficacite du
59!   melange
60!
61!   Le calcul du transport des differentes especes se fait en prenant
62!   en compte:
63!     1. un flux de masse montant
64!     2. un flux de masse descendant
65!     3. un entrainement
66!     4. un detrainement
67!
68! Modif 2013/01/04 (FH hourdin@lmd.jussieu.fr)
69!    Introduction of an implicit computation of vertical advection in
70!    the environment of thermal plumes in thermcell_dq
71!    impl =     0 : explicit, 1 : implicit, -1 : old version
72!    controled by iflag_thermals =
73!       15, 16 run with impl=-1 : numerical convergence with NPv3
74!       17, 18 run with impl=1  : more stable
75!    15 and 17 correspond to the activation of the stratocumulus "bidouille"
76!
77! Using
78!    abort_physic
79!    iso_verif_aberrant_encadre
80!    iso_verif_egalite
81!    test_ltherm
82!    thermcell_closure
83!    thermcell_dq
84!    thermcell_dry
85!    thermcell_dv2
86!    thermcell_env
87!    thermcell_flux2
88!    thermcell_height
89!    thermcell_plume
90!    thermcell_plume_5B
91!    thermcell_plume_6A
92!
93!=======================================================================
94
95
96!-----------------------------------------------------------------------
97!   declarations:
98!   -------------
99
100
101!   arguments:
102!   ----------
103      integer, intent(in) :: itap,ngrid,nlay
104      real, intent(in) ::  ptimestep
105      real, intent(in), dimension(ngrid,nlay)    :: pt,pu,pv,pplay,pphi
106! ATTENTION : po et zpspsk sont inout et out mais c'est pas forcement pour de bonnes raisons (FH, 2023)
107      real, intent(inout), dimension(ngrid,nlay)    :: po
108      real, intent(out), dimension(ngrid,nlay)    :: zpspsk
109      real, intent(in), dimension(ngrid,nlay+1)  :: pplev
110      integer, intent(out), dimension(ngrid) :: lmax
111      real, intent(out), dimension(ngrid,nlay)   :: pdtadj,pduadj,pdvadj,pdoadj,entr0,detr0
112      real, intent(out), dimension(ngrid,nlay)   :: ztla,zqla,zqta,zqsatth,zthl
113      real, intent(out), dimension(ngrid,nlay+1) :: fm0,zw2,fraca
114      real, intent(inout), dimension(ngrid) :: zmax0,f0
115      real, intent(out), dimension(ngrid,nlay) :: ztva,ztv
116      logical, intent(in) :: debut
117      real,intent(out), dimension(ngrid,nlay) :: ratqscth,ratqsdiff
118
119      real, intent(out), dimension(ngrid) :: pcon
120      real, intent(out), dimension(ngrid,nlay) :: rhobarz,wth3
121      real, intent(out), dimension(ngrid) :: wmax_sec
122      integer,intent(out), dimension(ngrid) :: lalim
123      real, intent(out), dimension(ngrid,nlay+1) :: fm
124      real, intent(out), dimension(ngrid,nlay) :: alim_star
125      real, intent(out), dimension(ngrid) :: zmax
126
127!   local:
128!   ------
129
130
131      integer,save :: igout=1
132!$OMP THREADPRIVATE(igout)
133      integer,save :: lunout1=6
134!$OMP THREADPRIVATE(lunout1)
135      integer,save :: lev_out=10
136!$OMP THREADPRIVATE(lev_out)
137
138      real lambda, zf,zf2,var,vardiff,CHI
139      integer ig,k,l,ierr,ll
140      logical sorties
141      real, dimension(ngrid) :: linter,zmix, zmax_sec
142      integer,dimension(ngrid) :: lmin,lmix,lmix_bis,nivcon
143      real, dimension(ngrid,nlay) :: ztva_est
144      real, dimension(ngrid,nlay) :: deltaz,zlay,zh,zdthladj,zu,zv,zo,zl,zva,zua,zoa
145      real, dimension(ngrid,nlay) :: zta,zha,q2,wq,wthl,wthv,thetath2,wth2
146      real, dimension(ngrid,nlay) :: rho,masse
147      real, dimension(ngrid,nlay+1) :: zw_est,zlev
148      real, dimension(ngrid) :: wmax,wmax_tmp
149      real, dimension(ngrid,nlay+1) :: f_star
150      real, dimension(ngrid,nlay) :: entr,detr,entr_star,detr_star,alim_star_clos
151      real, dimension(ngrid,nlay) :: zqsat,csc
152      real, dimension(ngrid) :: zcon,zcon2,alim_star_tot,f
153      real, dimension(ngrid,nlay) :: entrdn,detrdn
154
155      character (len=20) :: modname='thermcell_main'
156      character (len=80) :: abort_message
157
158
159#ifdef ISO
160      REAL xtpo(ntiso,ngrid,nlay),xtpdoadj(ntiso,ngrid,nlay)
161      REAL xtzo(ntiso,ngrid,nlay)
162      REAL xtpdoadj_tmp(ngrid,nlay)
163      REAL xtpo_tmp(ngrid,nlay)
164      REAL xtzo_tmp(ngrid,nlay)
165      integer ixt
166#endif
167
168!
169
170!-----------------------------------------------------------------------
171!   initialisation:
172!   ---------------
173!
174   fm=0. ; entr=0. ; detr=0.
175
176      if (prt_level.ge.1) print*,'thermcell_main V4'
177
178       sorties=.true.
179      IF(ngrid.NE.ngrid) THEN
180         PRINT*
181         PRINT*,'STOP dans convadj'
182         PRINT*,'ngrid    =',ngrid
183         PRINT*,'ngrid  =',ngrid
184      ENDIF
185!
186!     write(lunout,*)'WARNING thermcell_main f0=max(f0,1.e-2)'
187     do ig=1,ngrid
188         f0(ig)=max(f0(ig),1.e-2)
189         zmax0(ig)=max(zmax0(ig),40.)
190!IMmarche pas ?!       if (f0(ig)<1.e-2) f0(ig)=1.e-2
191     enddo
192
193      if (prt_level.ge.20) then
194       do ig=1,ngrid
195          print*,'th_main ig f0',ig,f0(ig)
196       enddo
197      endif
198!-----------------------------------------------------------------------
199! Calcul de T,q,ql a partir de Tl et qT dans l environnement
200!   --------------------------------------------------------------------
201!
202      CALL thermcell_env(ngrid,nlay,po,pt,pu,pv,pplay,  &
203     &           pplev,zo,zh,zl,ztv,zthl,zu,zv,zpspsk,zqsat,lev_out)
204       
205      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_env'
206
207!------------------------------------------------------------------------
208!                       --------------------
209!
210!
211!                       + + + + + + + + + + +
212!
213!
214!  wa, fraca, wd, fracd --------------------   zlev(2), rhobarz
215!  wh,wt,wo ...
216!
217!                       + + + + + + + + + + +  zh,zu,zv,zo,rho
218!
219!
220!                       --------------------   zlev(1)
221!                       \\\\\\\\\\\\\\\\\\\\
222!
223!
224
225!-----------------------------------------------------------------------
226!   Calcul des altitudes des couches
227!-----------------------------------------------------------------------
228
229      do l=2,nlay
230         zlev(:,l)=0.5*(pphi(:,l)+pphi(:,l-1))/RG
231      enddo
232      zlev(:,1)=0.
233      zlev(:,nlay+1)=(2.*pphi(:,nlay)-pphi(:,nlay-1))/RG
234      do l=1,nlay
235         zlay(:,l)=pphi(:,l)/RG
236      enddo
237      do l=1,nlay
238         deltaz(:,l)=zlev(:,l+1)-zlev(:,l)
239      enddo
240
241!-----------------------------------------------------------------------
242!   Calcul des densites et masses
243!-----------------------------------------------------------------------
244
245      rho(:,:)=pplay(:,:)/(zpspsk(:,:)*RD*ztv(:,:))
246      if (prt_level.ge.10) write(lunout,*) 'WARNING thermcell_main rhobarz(:,1)=rho(:,1)'
247      rhobarz(:,1)=rho(:,1)
248      do l=2,nlay
249         rhobarz(:,l)=0.5*(rho(:,l)+rho(:,l-1))
250      enddo
251      do l=1,nlay
252         masse(:,l)=(pplev(:,l)-pplev(:,l+1))/RG
253      enddo
254      if (prt_level.ge.1) print*,'thermcell_main apres initialisation'
255
256!------------------------------------------------------------------
257!
258!             /|\
259!    --------  |  F_k+1 -------   
260!                              ----> D_k
261!             /|\              <---- E_k , A_k
262!    --------  |  F_k ---------
263!                              ----> D_k-1
264!                              <---- E_k-1 , A_k-1
265!
266!
267!
268!
269!
270!    ---------------------------
271!
272!    ----- F_lmax+1=0 ----------         \
273!            lmax     (zmax)              |
274!    ---------------------------          |
275!                                         |
276!    ---------------------------          |
277!                                         |
278!    ---------------------------          |
279!                                         |
280!    ---------------------------          |
281!                                         |
282!    ---------------------------          |
283!                                         |  E
284!    ---------------------------          |  D
285!                                         |
286!    ---------------------------          |
287!                                         |
288!    ---------------------------  \       |
289!            lalim                 |      |
290!    ---------------------------   |      |
291!                                  |      |
292!    ---------------------------   |      |
293!                                  | A    |
294!    ---------------------------   |      |
295!                                  |      |
296!    ---------------------------   |      |
297!    lmin  (=1 pour le moment)     |      |
298!    ----- F_lmin=0 ------------  /      /
299!
300!    ---------------------------
301!    //////////////////////////
302!
303!
304!=============================================================================
305!  Calculs initiaux ne faisant pas intervenir les changements de phase
306!=============================================================================
307
308!------------------------------------------------------------------
309!  1. alim_star est le profil vertical de l'alimentation a la base du
310!     panache thermique, calcule a partir de la flotabilite de l'air sec
311!  2. lmin et lalim sont les indices inferieurs et superieurs de alim_star
312!------------------------------------------------------------------
313!
314      entr_star=0. ; detr_star=0. ; alim_star=0. ; alim_star_tot=0.
315      lmin=1
316
317!-----------------------------------------------------------------------------
318!  3. wmax_sec et zmax_sec sont les vitesses et altitudes maximum d'un
319!     panache sec conservatif (e=d=0) alimente selon alim_star
320!     Il s'agit d'un calcul de type CAPE
321!     zmax_sec est utilise pour determiner la geometrie du thermique.
322!------------------------------------------------------------------------------
323!---------------------------------------------------------------------------------
324!calcul du melange et des variables dans le thermique
325!--------------------------------------------------------------------------------
326!
327      if (prt_level.ge.1) print*,'avant thermcell_plume ',lev_out
328
329!=====================================================================
330! Old version of thermcell_plume in thermcell_plume_6A.F90
331! It includes both thermcell_plume_6A and thermcell_plume_5B corresponding
332! to the 5B and 6A versions used for CMIP5 and CMIP6.
333! The latest was previously named thermcellV1_plume.
334! The new thermcell_plume is a clean version (removing obsolete
335! options) of thermcell_plume_6A.
336! The 3 versions are controled by
337! flag_thermals_ed <= 9 thermcell_plume_6A
338!                  <= 19 thermcell_plume_5B
339!                  else thermcell_plume (default 20 for convergence with 6A)
340! Fredho
341!=====================================================================
342
343      if (iflag_thermals_ed<=9) then
344!         print*,'THERM NOUVELLE/NOUVELLE Arnaud'
345         CALL thermcell_plume_6A(itap,ngrid,nlay,ptimestep,ztv,zthl,po,zl,rhobarz,&
346     &    zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot,  &
347     &    lalim,f0,detr_star,entr_star,f_star,csc,ztva,  &
348     &    ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter &
349     &    ,lev_out,lunout1,igout)
350
351      elseif (iflag_thermals_ed<=19) then
352!        print*,'THERM RIO et al 2010, version d Arnaud'
353         CALL thermcell_plume_5B(itap,ngrid,nlay,ptimestep,ztv,zthl,po,zl,rhobarz,&
354     &    zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot,  &
355     &    lalim,f0,detr_star,entr_star,f_star,csc,ztva,  &
356     &    ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter &
357     &    ,lev_out,lunout1,igout)
358      else
359         CALL thermcell_plume(itap,ngrid,nlay,ptimestep,ztv,zthl,po,zl,rhobarz,&
360     &    zlev,pplev,pphi,zpspsk,alim_star,alim_star_tot,  &
361     &    lalim,f0,detr_star,entr_star,f_star,csc,ztva,  &
362     &    ztla,zqla,zqta,zha,zw2,zw_est,ztva_est,zqsatth,lmix,lmix_bis,linter &
363     &    ,lev_out,lunout1,igout)
364      endif
365
366      if (prt_level.ge.1) print*,'apres thermcell_plume ',lev_out
367
368      call test_ltherm(ngrid,nlay,pplay,lalim,ztv,po,ztva,zqla,f_star,zw2,'thermcell_plum lalim ')
369      call test_ltherm(ngrid,nlay,pplay,lmix ,ztv,po,ztva,zqla,f_star,zw2,'thermcell_plum lmix  ')
370
371      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_plume'
372      if (prt_level.ge.10) then
373         write(lunout1,*) 'Dans thermcell_main 2'
374         write(lunout1,*) 'lmin ',lmin(igout)
375         write(lunout1,*) 'lalim ',lalim(igout)
376         write(lunout1,*) ' ig l alim_star entr_star detr_star f_star '
377         write(lunout1,'(i6,i4,4e15.5)') (igout,l,alim_star(igout,l),entr_star(igout,l),detr_star(igout,l) &
378     &    ,f_star(igout,l+1),l=1,nint(linter(igout))+5)
379      endif
380
381!-------------------------------------------------------------------------------
382! Calcul des caracteristiques du thermique:zmax,zmix,wmax
383!-------------------------------------------------------------------------------
384!
385      CALL thermcell_height(ngrid,nlay,lalim,lmin,linter,lmix,zw2,  &
386     &           zlev,lmax,zmax,zmax0,zmix,wmax)
387! Attention, w2 est transforme en sa racine carree dans cette routine
388! Le probleme vient du fait que linter et lmix sont souvent egaux a 1.
389      wmax_tmp=0.
390      do  l=1,nlay
391         wmax_tmp(:)=max(wmax_tmp(:),zw2(:,l))
392      enddo
393!     print*,"ZMAX ",lalim,lmin,linter,lmix,lmax,zmax,zmax0,zmix,wmax
394
395
396
397      call test_ltherm(ngrid,nlay,pplay,lalim,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lalim ')
398      call test_ltherm(ngrid,nlay,pplay,lmin ,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lmin  ')
399      call test_ltherm(ngrid,nlay,pplay,lmix ,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lmix  ')
400      call test_ltherm(ngrid,nlay,pplay,lmax ,ztv,po,ztva,zqla,f_star,zw2,'thermcell_heig lmax  ')
401
402      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_height'
403
404!-------------------------------------------------------------------------------
405! Fermeture,determination de f
406!-------------------------------------------------------------------------------
407!
408!
409      CALL thermcell_dry(ngrid,nlay,zlev,pphi,ztv,alim_star,  &
410    &                      lalim,lmin,zmax_sec,wmax_sec)
411
412 
413call test_ltherm(ngrid,nlay,pplay,lmin,ztv,po,ztva,zqla,f_star,zw2,'thermcell_dry  lmin  ')
414call test_ltherm(ngrid,nlay,pplay,lalim,ztv,po,ztva,zqla,f_star,zw2,'thermcell_dry  lalim ')
415
416      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_dry'
417      if (prt_level.ge.10) then
418         write(lunout1,*) 'Dans thermcell_main 1b'
419         write(lunout1,*) 'lmin ',lmin(igout)
420         write(lunout1,*) 'lalim ',lalim(igout)
421         write(lunout1,*) ' ig l alim_star entr_star detr_star f_star '
422         write(lunout1,'(i6,i4,e15.5)') (igout,l,alim_star(igout,l) &
423     &    ,l=1,lalim(igout)+4)
424      endif
425
426
427
428
429! Choix de la fonction d'alimentation utilisee pour la fermeture.
430! Apparemment sans importance
431      alim_star_clos(:,:)=alim_star(:,:)
432      alim_star_clos(:,:)=entr_star(:,:)+alim_star(:,:)
433!
434!CR Appel de la fermeture seche
435      if (iflag_thermals_closure.eq.1) then
436
437     CALL thermcell_closure(ngrid,nlay,r_aspect_thermals,ptimestep,rho,  &
438    &   zlev,lalim,alim_star_clos,zmax_sec,wmax_sec,f)
439
440!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
441! Appel avec les zmax et wmax tenant compte de la condensation
442! Semble moins bien marcher
443     else if (iflag_thermals_closure.eq.2) then
444
445     CALL thermcell_closure(ngrid,nlay,r_aspect_thermals,ptimestep,rho,  &
446    &   zlev,lalim,alim_star,zmax,wmax,f)
447
448
449     endif
450
451!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
452
453      if(prt_level.ge.1)print*,'thermcell_closure apres thermcell_closure'
454
455      if (tau_thermals>1.) then
456         lambda=exp(-ptimestep/tau_thermals)
457         f0=(1.-lambda)*f+lambda*f0
458      else
459         f0=f
460      endif
461
462! Test valable seulement en 1D mais pas genant
463      if (.not. (f0(1).ge.0.) ) then
464              abort_message = '.not. (f0(1).ge.0.)'
465              CALL abort_physic (modname,abort_message,1)
466      endif
467
468!-------------------------------------------------------------------------------
469!deduction des flux
470
471      CALL thermcell_flux2(ngrid,nlay,ptimestep,masse, &
472     &       lalim,lmax,alim_star,  &
473     &       entr_star,detr_star,f,rhobarz,zlev,zw2,fm,entr,  &
474     &       detr,zqla,lev_out,lunout1,igout)
475
476!IM 060508    &       detr,zqla,zmax,lev_out,lunout,igout)
477
478      if (prt_level.ge.1) print*,'thermcell_main apres thermcell_flux'
479      call test_ltherm(ngrid,nlay,pplay,lalim,ztv,po,ztva,zqla,f_star,zw2,'thermcell_flux lalim ')
480      call test_ltherm(ngrid,nlay,pplay,lmax ,ztv,po,ztva,zqla,f_star,zw2,'thermcell_flux lmax  ')
481
482!------------------------------------------------------------------
483!   On ne prend pas directement les profils issus des calculs precedents
484!   mais on s'autorise genereusement une relaxation vers ceci avec
485!   une constante de temps tau_thermals (typiquement 1800s).
486!------------------------------------------------------------------
487
488      if (tau_thermals>1.) then
489         lambda=exp(-ptimestep/tau_thermals)
490         fm0=(1.-lambda)*fm+lambda*fm0
491         entr0=(1.-lambda)*entr+lambda*entr0
492         detr0=(1.-lambda)*detr+lambda*detr0
493      else
494         fm0=fm
495         entr0=entr
496         detr0=detr
497      endif
498
499!------------------------------------------------------------------
500! Calcul de la fraction de l'ascendance
501!------------------------------------------------------------------
502      do ig=1,ngrid
503         fraca(ig,1)=0.
504         fraca(ig,nlay+1)=0.
505      enddo
506      do l=2,nlay
507         do ig=1,ngrid
508            if (zw2(ig,l).gt.1.e-10) then
509            fraca(ig,l)=fm(ig,l)/(rhobarz(ig,l)*zw2(ig,l))
510            else
511            fraca(ig,l)=0.
512            endif
513         enddo
514      enddo
515     
516!c------------------------------------------------------------------
517!   calcul du transport vertical
518!------------------------------------------------------------------
519      IF (iflag_thermals_down .GT. 0) THEN
520        if (debut) print*,'WARNING !!! routine thermcell_down en cours de developpement'
521        entrdn=fact_thermals_down*detr0
522        detrdn=fact_thermals_down*entr0
523        ! we want to transport potential temperature, total water and momentum
524        CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,zthl,zdthladj)
525        CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,po,pdoadj)
526        CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,zu,pduadj)
527        CALL thermcell_updown_dq(ngrid,nlay,ptimestep,lmax,entr0,detr0,entrdn,detrdn,masse,zv,pdvadj)
528      ELSE
529      !--------------------------------------------------------------
530
531        call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse,  &
532        &                    zthl,zdthladj,zta,lev_out)
533        call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse,  &
534        &                   po,pdoadj,zoa,lev_out)
535
536#ifdef ISO
537        ! C Risi: on utilise directement la meme routine
538        do ixt=1,ntiso
539          do ll=1,nlay
540            DO ig=1,ngrid
541                xtpo_tmp(ig,ll)=xtpo(ixt,ig,ll)
542                xtzo_tmp(ig,ll)=xtzo(ixt,ig,ll)
543            enddo
544          enddo
545          call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse,  &
546     &                   xtpo_tmp,xtpdoadj_tmp,xtzo_tmp,lev_out)
547          do ll=1,nlay
548            DO ig=1,ngrid
549                xtpdoadj(ixt,ig,ll)=xtpdoadj_tmp(ig,ll)
550            enddo
551          enddo
552        enddo
553#endif
554
555#ifdef ISO     
556#ifdef ISOVERIF
557      DO  ll=1,nlay
558        DO ig=1,ngrid
559          if (iso_eau.gt.0) then
560              call iso_verif_egalite(xtpo(iso_eau,ig,ll), &
561     &          po(ig,ll),'thermcell_main 594')
562              call iso_verif_egalite(xtpdoadj(iso_eau,ig,ll), &
563     &          pdoadj(ig,ll),'thermcell_main 596')
564          endif
565          if (iso_HDO.gt.0) then
566              call iso_verif_aberrant_encadre(xtpo(iso_hdo,ig,ll) &
567     &           /po(ig,ll),'thermcell_main 610')
568          endif
569        enddo
570      enddo !DO  ll=1,nlay
571      write(*,*) 'thermcell_main 600 tmp: apres thermcell_dq'
572#endif     
573#endif
574
575
576!------------------------------------------------------------------
577!  calcul du transport vertical du moment horizontal
578!------------------------------------------------------------------
579
580!IM 090508 
581      if (dvdq == 0 ) then
582
583! Calcul du transport de V tenant compte d'echange par gradient
584! de pression horizontal avec l'environnement
585
586         call thermcell_dv2(ngrid,nlay,ptimestep,fm0,entr0,masse  &
587!    &    ,fraca*dvdq,zmax &
588     &    ,fraca,zmax &
589     &    ,zu,zv,pduadj,pdvadj,zua,zva,lev_out)
590
591      else
592
593! calcul purement conservatif pour le transport de V
594         call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse  &
595     &    ,zu,pduadj,zua,lev_out)
596         call thermcell_dq(ngrid,nlay,dqimpl,ptimestep,fm0,entr0,masse  &
597     &    ,zv,pdvadj,zva,lev_out)
598
599      endif
600    ENDIF
601
602!     print*,'13 OK convect8'
603      do l=1,nlay
604         do ig=1,ngrid
605           pdtadj(ig,l)=zdthladj(ig,l)*zpspsk(ig,l) 
606         enddo
607      enddo
608
609      if (prt_level.ge.1) print*,'14 OK convect8'
610!------------------------------------------------------------------
611!   Calculs de diagnostiques pour les sorties
612!------------------------------------------------------------------
613!calcul de fraca pour les sorties
614     
615      if (sorties) then
616      if (prt_level.ge.1) print*,'14a OK convect8'
617! calcul du niveau de condensation
618! initialisation
619      do ig=1,ngrid
620         nivcon(ig)=0
621         zcon(ig)=0.
622      enddo
623!nouveau calcul
624      do ig=1,ngrid
625      CHI=zh(ig,1)/(1669.0-122.0*zo(ig,1)/zqsat(ig,1)-zh(ig,1))
626      pcon(ig)=pplay(ig,1)*(zo(ig,1)/zqsat(ig,1))**CHI
627      enddo
628!IM   do k=1,nlay
629      do k=1,nlay-1
630         do ig=1,ngrid
631         if ((pcon(ig).le.pplay(ig,k))  &
632     &      .and.(pcon(ig).gt.pplay(ig,k+1))) then
633            zcon2(ig)=zlay(ig,k)-(pcon(ig)-pplay(ig,k))/(RG*rho(ig,k))/100.
634         endif
635         enddo
636      enddo
637!IM
638      ierr=0
639      do ig=1,ngrid
640        if (pcon(ig).le.pplay(ig,nlay)) then
641           zcon2(ig)=zlay(ig,nlay)-(pcon(ig)-pplay(ig,nlay))/(RG*rho(ig,nlay))/100.
642           ierr=1
643        endif
644      enddo
645      if (ierr==1) then
646           abort_message = 'thermcellV0_main: les thermiques vont trop haut '
647           CALL abort_physic (modname,abort_message,1)
648      endif
649
650      if (prt_level.ge.1) print*,'14b OK convect8'
651      do k=nlay,1,-1
652         do ig=1,ngrid
653            if (zqla(ig,k).gt.1e-10) then
654               nivcon(ig)=k
655               zcon(ig)=zlev(ig,k)
656            endif
657         enddo
658      enddo
659      if (prt_level.ge.1) print*,'14c OK convect8'
660!calcul des moments
661!initialisation
662      do l=1,nlay
663         do ig=1,ngrid
664            q2(ig,l)=0.
665            wth2(ig,l)=0.
666            wth3(ig,l)=0.
667            ratqscth(ig,l)=0.
668            ratqsdiff(ig,l)=0.
669         enddo
670      enddo     
671      if (prt_level.ge.1) print*,'14d OK convect8'
672      if (prt_level.ge.10)write(lunout,*)                                &
673    &     'WARNING thermcell_main wth2=0. si zw2 > 1.e-10'
674      do l=1,nlay
675         do ig=1,ngrid
676            zf=fraca(ig,l)
677            zf2=zf/(1.-zf)
678!
679            thetath2(ig,l)=zf2*(ztla(ig,l)-zthl(ig,l))**2
680            if(zw2(ig,l).gt.1.e-10) then
681             wth2(ig,l)=zf2*(zw2(ig,l))**2
682            else
683             wth2(ig,l)=0.
684            endif
685            wth3(ig,l)=zf2*(1-2.*fraca(ig,l))/(1-fraca(ig,l))  &
686     &                *zw2(ig,l)*zw2(ig,l)*zw2(ig,l)
687            q2(ig,l)=zf2*(zqta(ig,l)*1000.-po(ig,l)*1000.)**2
688!test: on calcul q2/po=ratqsc
689            ratqscth(ig,l)=sqrt(max(q2(ig,l),1.e-6)/(po(ig,l)*1000.))
690         enddo
691      enddo
692!calcul des flux: q, thetal et thetav
693      do l=1,nlay
694         do ig=1,ngrid
695      wq(ig,l)=fraca(ig,l)*zw2(ig,l)*(zqta(ig,l)*1000.-po(ig,l)*1000.)
696      wthl(ig,l)=fraca(ig,l)*zw2(ig,l)*(ztla(ig,l)-zthl(ig,l))
697      wthv(ig,l)=fraca(ig,l)*zw2(ig,l)*(ztva(ig,l)-ztv(ig,l))
698         enddo
699      enddo
700
701!calcul du ratqscdiff
702      if (prt_level.ge.1) print*,'14e OK convect8'
703      var=0.
704      vardiff=0.
705      ratqsdiff(:,:)=0.
706
707      do l=1,nlay
708         do ig=1,ngrid
709            if (l<=lalim(ig)) then
710            var=var+alim_star(ig,l)*zqta(ig,l)*1000.
711            endif
712         enddo
713      enddo
714
715      if (prt_level.ge.1) print*,'14f OK convect8'
716
717      do l=1,nlay
718         do ig=1,ngrid
719            if (l<=lalim(ig)) then
720               zf=fraca(ig,l)
721               zf2=zf/(1.-zf)
722               vardiff=vardiff+alim_star(ig,l)*(zqta(ig,l)*1000.-var)**2
723            endif
724         enddo
725      enddo
726
727      if (prt_level.ge.1) print*,'14g OK convect8'
728         do l=1,nlay
729            do ig=1,ngrid
730               ratqsdiff(ig,l)=sqrt(vardiff)/(po(ig,l)*1000.)   
731            enddo
732         enddo
733      endif
734
735      if (prt_level.ge.1) print*,'thermcell_main FIN  OK'
736
737 RETURN
738      end subroutine thermcell_main
739
740!=============================================================================
741!/////////////////////////////////////////////////////////////////////////////
742!=============================================================================
743      subroutine test_ltherm(ngrid,nlay,pplay,long,ztv,po,ztva, &  ! in
744    &            zqla,f_star,zw2,comment)                          ! in
745!=============================================================================
746      USE lmdz_thermcell_ini, ONLY: prt_level
747      IMPLICIT NONE
748
749      integer i, k, ngrid,nlay
750      real, intent(in), dimension(ngrid,nlay) :: pplay,ztv,po,ztva,zqla
751      real, intent(in), dimension(ngrid,nlay) :: f_star,zw2
752      integer, intent(in), dimension(ngrid) :: long
753      real seuil
754      character*21 comment
755
756      seuil=0.25
757
758      if (prt_level.ge.1) THEN
759       print*,'WARNING !!! TEST ',comment
760      endif
761      return
762
763!  test sur la hauteur des thermiques ...
764         do i=1,ngrid
765!IMtemp           if (pplay(i,long(i)).lt.seuil*pplev(i,1)) then
766           if (prt_level.ge.10) then
767               print*,'WARNING ',comment,' au point ',i,' K= ',long(i)
768               print*,'  K  P(MB)  THV(K)     Qenv(g/kg)THVA        QLA(g/kg)   F*        W2'
769               do k=1,nlay
770                  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)
771               enddo
772           endif
773         enddo
774
775
776      return
777      end
778
779! nrlmd le 10/04/2012   Transport de la TKE par le thermique moyen pour la fermeture en ALP
780!                       On transporte pbl_tke pour donner therm_tke
781!                       Copie conforme de la subroutine DTKE dans physiq.F ecrite par Frederic Hourdin
782
783!=======================================================================
784!///////////////////////////////////////////////////////////////////////
785!=======================================================================
786
787      subroutine thermcell_tke_transport( &
788     &     ngrid,nlay,ptimestep,fm0,entr0,rg,pplev,  &   ! in
789     &     therm_tke_max)                                ! out
790      USE lmdz_thermcell_ini, ONLY: prt_level
791      implicit none
792
793!=======================================================================
794!
795!   Calcul du transport verticale dans la couche limite en presence
796!   de "thermiques" explicitement representes
797!   calcul du dq/dt une fois qu'on connait les ascendances
798!
799!=======================================================================
800
801      integer ngrid,nlay
802
803      real, intent(in) :: ptimestep
804      real, intent(in), dimension(ngrid,nlay+1) :: fm0,pplev
805      real, intent(in), dimension(ngrid,nlay) :: entr0
806      real, intent(in) :: rg
807      real, intent(out), dimension(ngrid,nlay) :: therm_tke_max
808
809      real detr0(ngrid,nlay)
810      real masse0(ngrid,nlay)
811      real masse(ngrid,nlay),fm(ngrid,nlay+1)
812      real entr(ngrid,nlay)
813      real q(ngrid,nlay)
814      integer lev_out                           ! niveau pour les print
815
816      real qa(ngrid,nlay),detr(ngrid,nlay),wqd(ngrid,nlay+1)
817      integer ig,k
818
819
820      lev_out=0
821
822
823      if (prt_level.ge.1) print*,'Q2 THERMCEL_DQ 0'
824
825!   calcul du detrainement
826      do k=1,nlay
827         detr0(:,k)=fm0(:,k)-fm0(:,k+1)+entr0(:,k)
828         masse0(:,k)=(pplev(:,k)-pplev(:,k+1))/RG
829      enddo
830
831
832! Decalage vertical des entrainements et detrainements.
833      masse(:,1)=0.5*masse0(:,1)
834      entr(:,1)=0.5*entr0(:,1)
835      detr(:,1)=0.5*detr0(:,1)
836      fm(:,1)=0.
837      do k=1,nlay-1
838         masse(:,k+1)=0.5*(masse0(:,k)+masse0(:,k+1))
839         entr(:,k+1)=0.5*(entr0(:,k)+entr0(:,k+1))
840         detr(:,k+1)=0.5*(detr0(:,k)+detr0(:,k+1))
841         fm(:,k+1)=fm(:,k)+entr(:,k)-detr(:,k)
842      enddo
843      fm(:,nlay+1)=0.
844
845
846   q(:,:)=therm_tke_max(:,:)
847!!! nrlmd le 16/09/2010
848      do ig=1,ngrid
849         qa(ig,1)=q(ig,1)
850      enddo
851!!!
852
853    if (1==1) then
854      do k=2,nlay
855         do ig=1,ngrid
856            if ((fm(ig,k+1)+detr(ig,k))*ptimestep.gt.  &
857     &         1.e-5*masse(ig,k)) then
858         qa(ig,k)=(fm(ig,k)*qa(ig,k-1)+entr(ig,k)*q(ig,k))  &
859     &         /(fm(ig,k+1)+detr(ig,k))
860            else
861               qa(ig,k)=q(ig,k)
862            endif
863            if (qa(ig,k).lt.0.) then
864!               print*,'qa<0!!!'
865            endif
866            if (q(ig,k).lt.0.) then
867!               print*,'q<0!!!'
868            endif
869         enddo
870      enddo
871
872! Calcul du flux subsident
873
874      do k=2,nlay
875         do ig=1,ngrid
876            wqd(ig,k)=fm(ig,k)*q(ig,k)
877            if (wqd(ig,k).lt.0.) then
878!               print*,'wqd<0!!!'
879            endif
880         enddo
881      enddo
882      do ig=1,ngrid
883         wqd(ig,1)=0.
884         wqd(ig,nlay+1)=0.
885      enddo
886
887! Calcul des tendances
888      do k=1,nlay
889         do ig=1,ngrid
890            q(ig,k)=q(ig,k)+(detr(ig,k)*qa(ig,k)-entr(ig,k)*q(ig,k)  &
891     &               -wqd(ig,k)+wqd(ig,k+1))  &
892     &               *ptimestep/masse(ig,k)
893         enddo
894      enddo
895
896 endif
897
898   therm_tke_max(:,:)=q(:,:)
899
900      return
901!!! fin nrlmd le 10/04/2012
902     end
903
904END MODULE lmdz_thermcell_main
Note: See TracBrowser for help on using the repository browser.