source: LMDZ6/branches/Portage_acc/libf/phylmd/thermcell_main.F90 @ 4446

Last change on this file since 4446 was 4446, checked in by Laurent Fairhead, 15 months ago

Merged trunk revisions from 4127 to 4443 (HEAD) into branch

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