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

Last change on this file since 4132 was 4132, checked in by Laurent Fairhead, 2 years ago

Modifications to code to start using openacc directives
LF

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