source: LMDZ6/trunk/libf/phylmd/calltherm.F90 @ 5185

Last change on this file since 5185 was 4843, checked in by crio, 10 months ago

Nouvelle formulation du strig et correction thermiques montent trop haut

  • 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: 18.5 KB
RevLine 
[878]1!
[1403]2! $Id: calltherm.F90 4843 2024-03-04 17:58:03Z abarral $
[878]3!
4      subroutine calltherm(dtime  &
5     &      ,pplay,paprs,pphi,weak_inversion  &
[4690]6     &      ,u_seri_,v_seri_,t_seri_,q_seri_,t_env,q_env,zqsat,debut  &
[878]7     &      ,d_u_ajs,d_v_ajs,d_t_ajs,d_q_ajs  &
[4837]8     &    ,fm_therm,entr_therm,detr_therm,zqasc,clwcon0,lmax,ratqscth,&
[4089]9     &       ratqsdiff,zqsatth,ale_bl,alp_bl,lalim_conv,wght_th, &
[1638]10     &       zmax0,f0,zw2,fraca,ztv,zpspsk,ztla,zthl &
11!!! nrlmd le 10/04/2012
12     &      ,pbl_tke,pctsrf,omega,airephy &
13     &      ,zlcl_th,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &
[4843]14     &      ,n2,s2,strig,zcong,ale_bl_stat &
[1638]15     &      ,therm_tke_max,env_tke_max &
16     &      ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &
17     &      ,alp_bl_conv,alp_bl_stat &
18!!! fin nrlmd le 10/04/2012
[4089]19     &      ,zqla,ztva &
20#ifdef ISO         
21     &      ,xt_seri,d_xt_ajs &
22#ifdef DIAGISO         
23     &      ,q_the,xt_the &
24#endif
25#endif         
26     &   )
[878]27
[940]28      USE dimphy
[1785]29      USE indice_sol_mod
[2311]30      USE print_control_mod, ONLY: prt_level,lunout
[4590]31      USE lmdz_thermcell_alp, ONLY: thermcell_alp
32      USE lmdz_thermcell_main, ONLY: thermcell_main
33      USE lmdz_thermcell_old, ONLY: thermcell, thermcell_2002, thermcell_eau, calcul_sec, thermcell_sec
[4089]34#ifdef ISO
[4143]35      use infotrac_phy, ONLY: ntiso
[4089]36#ifdef ISOVERIF
37      USE isotopes_mod, ONLY: iso_eau,iso_HDO
38      USE isotopes_verif_mod, ONLY: iso_verif_aberrant_enc_vect2D, &
39        iso_verif_egalite_vect2D
40#endif   
41#endif
[1790]42
[878]43      implicit none
[4089]44      include "clesphys.h"
45      include "thermcell_old.h"
[878]46
[1638]47
[973]48!IM 140508
[1776]49      INTEGER, SAVE ::  itap
50!$OMP THREADPRIVATE(itap)
[878]51      REAL dtime
52      LOGICAL debut
[973]53      LOGICAL logexpr0, logexpr2(klon,klev), logexpr1(klon)
54      REAL fact(klon)
55      INTEGER nbptspb
56
[4678]57      REAL, DIMENSION(klon,klev), INTENT(IN) :: u_seri_,v_seri_
[4837]58      REAL, DIMENSION(klon,klev), INTENT(IN) :: t_seri_,q_seri_
59      REAL, DIMENSION(klon,klev), INTENT(IN) :: t_env,q_env
[4678]60      REAL, DIMENSION(klon,klev) :: u_seri,v_seri
61      REAL, DIMENSION(klon,klev) :: t_seri,q_seri
62      REAL, DIMENSION(klon,klev) :: qmemoire
[878]63      REAL weak_inversion(klon)
64      REAL paprs(klon,klev+1)
65      REAL pplay(klon,klev)
66      REAL pphi(klon,klev)
67      real zlev(klon,klev+1)
[879]68!test: on sort lentr et a* pour alimenter KE
69      REAL wght_th(klon,klev)
70      INTEGER lalim_conv(klon)
[1026]71      REAL zw2(klon,klev+1),fraca(klon,klev+1)
[878]72
73!FH Update Thermiques
74      REAL d_t_ajs(klon,klev), d_q_ajs(klon,klev)
75      REAL d_u_ajs(klon,klev),d_v_ajs(klon,klev)
[973]76      real fm_therm(klon,klev+1)
77      real entr_therm(klon,klev),detr_therm(klon,klev)
[878]78
79!********************************************************
80!     declarations
[1403]81      LOGICAL flag_bidouille_stratocu
[973]82      real fmc_therm(klon,klev+1),zqasc(klon,klev)
[878]83      real zqla(klon,klev)
[1790]84      real ztv(klon,klev),ztva(klon,klev)
[1403]85      real zpspsk(klon,klev)
86      real ztla(klon,klev)
87      real zthl(klon,klev)
[878]88      real wmax_sec(klon)
[4843]89      real zcong(klon)
[878]90      real zmax_sec(klon)
91      real f_sec(klon)
[973]92      real detrc_therm(klon,klev)
93! FH WARNING : il semble que ces save ne servent a rien
94!     save fmc_therm, detrc_therm
[878]95      real clwcon0(klon,klev)
96      real zqsat(klon,klev)
97      real zw_sec(klon,klev+1)
98      integer lmix_sec(klon)
99      integer lmax(klon)
100      real ratqscth(klon,klev)
101      real ratqsdiff(klon,klev)
102      real zqsatth(klon,klev) 
[879]103!nouvelles variables pour la convection
[4089]104      real ale_bl(klon)
105      real alp_bl(klon)
106      real ale(klon)
107      real alp(klon)
[879]108!RC
[927]109      !on garde le zmax du pas de temps precedent
110      real zmax0(klon), f0(klon)
[1638]111
112!!! nrlmd le 10/04/2012
113      real pbl_tke(klon,klev+1,nbsrf)
114      real pctsrf(klon,nbsrf)
115      real omega(klon,klev)
116      real airephy(klon)
117      real zlcl_th(klon),fraca0(klon),w0(klon),w_conv(klon)
118      real therm_tke_max0(klon),env_tke_max0(klon)
[4843]119      real n2(klon),s2(klon),strig(klon)
[1638]120      real ale_bl_stat(klon)
121      real therm_tke_max(klon,klev),env_tke_max(klon,klev)
122      real alp_bl_det(klon),alp_bl_fluct_m(klon),alp_bl_fluct_tke(klon),alp_bl_conv(klon),alp_bl_stat(klon)
123!!! fin nrlmd le 10/04/2012
124
[878]125!********************************************************
126
[4089]127      real, dimension(klon) :: pcon
128      real, dimension(klon,klev) :: rhobarz,wth3
129      integer,dimension(klon) :: lalim
130      real, dimension(klon,klev+1) :: fm
131      real, dimension(klon,klev) :: alim_star
132      real, dimension(klon) :: zmax
[878]133
[4089]134
135
136
[878]137! variables locales
138      REAL d_t_the(klon,klev), d_q_the(klon,klev)
139      REAL d_u_the(klon,klev),d_v_the(klon,klev)
140!
[973]141      real zfm_therm(klon,klev+1),zdt
142      real zentr_therm(klon,klev),zdetr_therm(klon,klev)
143! FH A VERIFIER : SAVE INUTILES
[940]144!      save zentr_therm,zfm_therm
[973]145
[1403]146      character (len=20) :: modname='calltherm'
147      character (len=80) :: abort_message
148
[4089]149      integer i,k,isplit
[973]150      logical, save :: first=.true.
[4089]151      logical :: new_thermcell
152
153#ifdef ISO
[4143]154      REAL xt_seri(ntiso,klon,klev),xtmemoire(ntiso,klon,klev)
155      REAL d_xt_ajs(ntiso,klon,klev)
156      real d_xt_the(ntiso,klon,klev)
[4089]157#ifdef DIAGISO
158      real q_the(klon,klev)
[4143]159      real xt_the(ntiso,klon,klev)
[4089]160#endif
161      real qprec(klon,klev)
162      integer ixt
163#endif
164
165
[987]166!$OMP THREADPRIVATE(first)
[878]167!********************************************************
[973]168      if (first) then
169        itap=0
170        first=.false.
171      endif
[878]172
[4678]173      u_seri(:,:)=u_seri_(:,:)
174      v_seri(:,:)=v_seri_(:,:)
175      t_seri(:,:)=t_seri_(:,:)
176      q_seri(:,:)=q_seri_(:,:)
177
[973]178! Incrementer le compteur de la physique
179     itap   = itap + 1
180
[878]181!  Modele du thermique
182!  ===================
183!         print*,'thermiques: WARNING on passe t au lieu de t_seri'
184
185
[973]186! On prend comme valeur initiale des thermiques la valeur du pas
187! de temps precedent
188         zfm_therm(:,:)=fm_therm(:,:)
189         zdetr_therm(:,:)=detr_therm(:,:)
190         zentr_therm(:,:)=entr_therm(:,:)
191
192! On reinitialise les flux de masse a zero pour le cumul en
193! cas de splitting
[878]194         fm_therm(:,:)=0.
195         entr_therm(:,:)=0.
[973]196         detr_therm(:,:)=0.
197
[4089]198         ale_bl(:)=0.
199         alp_bl(:)=0.
[938]200         if (prt_level.ge.10) then
201          print*,'thermV4 nsplit: ',nsplit_thermals,' weak_inversion'
202         endif
[878]203
204!   tests sur les valeurs negatives de l'eau
[973]205         logexpr0=prt_level.ge.10
206         nbptspb=0
[878]207         do k=1,klev
208            do i=1,klon
[1146]209! Attention teste abderr 19-03-09
210!               logexpr2(i,k)=.not.q_seri(i,k).ge.0.
211                logexpr2(i,k)=.not.q_seri(i,k).ge.1.e-15
[973]212               if (logexpr2(i,k)) then
[4089]213#ifdef ISO
214                qprec(i,k)=q_seri(i,k)
215#endif
[938]216                q_seri(i,k)=1.e-15
[973]217                nbptspb=nbptspb+1
[4089]218#ifdef ISO
[4143]219                do ixt=1,ntiso
[4089]220                  xt_seri(ixt,i,k)=1.e-15*(xt_seri(ixt,i,k)/qprec(i,k))
221                  ! xt_seri(ixt,i,k)=1.e-15*(Rdefault(index_iso(ixt)))
222                enddo
223#endif
[878]224               endif
[973]225!               if (logexpr0) &
226!    &             print*,'WARN eau<0 avant therm i=',i,'  k=',k  &
227!    &         ,' dq,q',d_q_the(i,k),q_seri(i,k)
[878]228            enddo
229         enddo
[973]230         if(nbptspb.GT.0) print*,'Number of points with q_seri(i,k)<=0 ',nbptspb   
[878]231
[4089]232
233         new_thermcell=iflag_thermals>=15.and.iflag_thermals<=18
234#ifdef ISO
235      if (.not.new_thermcell) then
236           CALL abort_gcm('calltherm 234','isos pas prevus ici',1)
237      endif
238#ifdef ISOVERIF
239      if (iso_eau.gt.0) then
240       call iso_verif_egalite_vect2D( &
241     &           xt_seri,q_seri, &
[4143]242     &           'calltherm 174',ntiso,klon,klev)
[4089]243      endif !if (iso_eau.gt.0) then
244#endif   
245#endif
[1403]246         zdt=dtime/REAL(nsplit_thermals)
[4089]247
248
[878]249         do isplit=1,nsplit_thermals
250
[1943]251          if (iflag_thermals>=1000) then
252            CALL thermcell_2002(klon,klev,zdt,iflag_thermals   &
[878]253     &      ,pplay,paprs,pphi  &
254     &      ,u_seri,v_seri,t_seri,q_seri  &
255     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
[1943]256     &      ,zfm_therm,zentr_therm,fraca,zw2  &
[878]257     &      ,r_aspect_thermals,30.,w2di_thermals  &
[1517]258     &      ,tau_thermals)
[878]259          else if (iflag_thermals.eq.2) then
260            CALL thermcell_sec(klon,klev,zdt  &
261     &      ,pplay,paprs,pphi,zlev  &
262     &      ,u_seri,v_seri,t_seri,q_seri  &
263     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
264     &      ,zfm_therm,zentr_therm  &
265     &      ,r_aspect_thermals,30.,w2di_thermals  &
[1517]266     &      ,tau_thermals)
[878]267          else if (iflag_thermals.eq.3) then
268            CALL thermcell(klon,klev,zdt  &
269     &      ,pplay,paprs,pphi  &
270     &      ,u_seri,v_seri,t_seri,q_seri  &
271     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
272     &      ,zfm_therm,zentr_therm  &
273     &      ,r_aspect_thermals,l_mix_thermals,w2di_thermals  &
[1517]274     &      ,tau_thermals)
[878]275          else if (iflag_thermals.eq.10) then
276            CALL thermcell_eau(klon,klev,zdt  &
277     &      ,pplay,paprs,pphi  &
278     &      ,u_seri,v_seri,t_seri,q_seri  &
279     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
280     &      ,zfm_therm,zentr_therm  &
281     &      ,r_aspect_thermals,l_mix_thermals,w2di_thermals  &
[1517]282     &      ,tau_thermals)
[878]283          else if (iflag_thermals.eq.11) then
[1403]284              abort_message = 'cas non prevu dans calltherm'
[2311]285              CALL abort_physic (modname,abort_message,1)
[878]286          else if (iflag_thermals.eq.12) then
287            CALL calcul_sec(klon,klev,zdt  &
288     &      ,pplay,paprs,pphi,zlev  &
289     &      ,u_seri,v_seri,t_seri,q_seri  &
290     &      ,zmax_sec,wmax_sec,zw_sec,lmix_sec  &
291     &      ,r_aspect_thermals,l_mix_thermals,w2di_thermals  &
[973]292     &      ,tau_thermals)
[1403]293          else if (iflag_thermals==13.or.iflag_thermals==14) then
[4089]294              abort_message = 'thermcellV0_main enleve svn>2084'
295              CALL abort_physic (modname,abort_message,1)
296          else if (new_thermcell) then
[1403]297            CALL thermcell_main(itap,klon,klev,zdt  &
298     &      ,pplay,paprs,pphi,debut  &
[4690]299     &      ,u_seri,v_seri,t_seri,q_seri,t_env,q_env  &
[1403]300     &      ,d_u_the,d_v_the,d_t_the,d_q_the  &
301     &      ,zfm_therm,zentr_therm,zdetr_therm,zqasc,zqla,lmax  &
302     &      ,ratqscth,ratqsdiff,zqsatth  &
303     &      ,zmax0,f0,zw2,fraca,ztv,zpspsk &
[4089]304     &      ,ztla,zthl,ztva &
[4843]305     &      ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax,zcong &
[4089]306#ifdef ISO         
307     &      ,xt_seri,d_xt_the &
308#endif         
309     &   )
310
311            CALL thermcell_alp(klon,klev,zdt  &                      ! in
312     &        ,pplay,paprs  &                                        ! in
313     &        ,zfm_therm,zentr_therm,lmax  &                         ! in
314     &        ,pbl_tke,pctsrf,omega,airephy &                        ! in
315     &        ,zw2,fraca &                                           ! in
316     &        ,pcon,rhobarz,wth3,wmax_sec,lalim,fm,alim_star,zmax &  ! in
[4843]317     &        ,zcong,ale,alp,lalim_conv,wght_th &                          ! out
[4089]318     &        ,zlcl_th,fraca0,w0,w_conv,therm_tke_max0,env_tke_max0 &! out
[4843]319     &        ,n2,s2,strig,ale_bl_stat &                                   ! out
[4089]320     &        ,therm_tke_max,env_tke_max &                           ! out
321     &        ,alp_bl_det,alp_bl_fluct_m,alp_bl_fluct_tke &          ! out
322     &        ,alp_bl_conv,alp_bl_stat &                             ! out
323     &        )
324
[1403]325           if (prt_level.gt.10) write(lunout,*)'Apres thermcell_main OK'
326         else
327           abort_message = 'Cas des thermiques non prevu'
[2311]328           CALL abort_physic (modname,abort_message,1)
[878]329         endif
330
[1638]331! Attention : les noms sont contre intuitif.
332! flag_bidouille_stratocu est .true. si on ne fait pas de bidouille.
333! Il aurait mieux valu avoir un nobidouille_stratocu
334! Et pour simplifier :
335! nobidouille_stratocu=.not.(iflag_thermals==13.or.iflag_thermals=15)
336! Ce serait bien de changer, mai en prenant le temps de vérifier que ca
337! fait bien ce qu'on croit.
[878]338
[1738]339       flag_bidouille_stratocu=iflag_thermals<=12.or.iflag_thermals==14.or.iflag_thermals==16.or.iflag_thermals==18
[1638]340
[1943]341! Calcul a posteriori du niveau max des thermiques pour les schémas qui
342! ne la sortent pas.
343      if (iflag_thermals<=12.or.iflag_thermals>=1000) then
344         lmax(:)=1
[1638]345         do k=1,klev-1
346            zdetr_therm(:,k)=zentr_therm(:,k)+zfm_therm(:,k)-zfm_therm(:,k+1)
347         enddo
[1943]348         do k=1,klev-1
349            do i=1,klon
350               if (zfm_therm(i,k+1)>0.) lmax(i)=k
351            enddo
352         enddo
[1638]353      endif
354
[973]355      fact(:)=0.
[878]356      DO i=1,klon
[1403]357       logexpr1(i)=flag_bidouille_stratocu.or.weak_inversion(i).gt.0.5
358       IF(logexpr1(i)) fact(i)=1./REAL(nsplit_thermals)
[973]359      ENDDO
[878]360
[973]361     DO k=1,klev
[878]362!  transformation de la derivee en tendance
[973]363            d_t_the(:,k)=d_t_the(:,k)*dtime*fact(:)
364            d_u_the(:,k)=d_u_the(:,k)*dtime*fact(:)
365            d_v_the(:,k)=d_v_the(:,k)*dtime*fact(:)
366            d_q_the(:,k)=d_q_the(:,k)*dtime*fact(:)
367            fm_therm(:,k)=fm_therm(:,k)  &
368     &      +zfm_therm(:,k)*fact(:)
369            entr_therm(:,k)=entr_therm(:,k)  &
370     &       +zentr_therm(:,k)*fact(:)
[1026]371            detr_therm(:,k)=detr_therm(:,k)  &
372     &       +zdetr_therm(:,k)*fact(:)
[4089]373#ifdef ISO
[4143]374            do ixt=1,ntiso
[4089]375              d_xt_the(ixt,:,k)=d_xt_the(ixt,:,k)*dtime*fact(:)
376            enddo
377#endif
[973]378      ENDDO
379       fm_therm(:,klev+1)=0.
[878]380
381
382
383!  accumulation de la tendance
[973]384            d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:)
385            d_u_ajs(:,:)=d_u_ajs(:,:)+d_u_the(:,:)
386            d_v_ajs(:,:)=d_v_ajs(:,:)+d_v_the(:,:)
387            d_q_ajs(:,:)=d_q_ajs(:,:)+d_q_the(:,:)
[4089]388#ifdef ISO
389            d_xt_ajs(:,:,:)=d_xt_ajs(:,:,:)+d_xt_the(:,:,:)
390#endif
[878]391
392!  incrementation des variables meteo
[973]393            t_seri(:,:) = t_seri(:,:) + d_t_the(:,:)
394            u_seri(:,:) = u_seri(:,:) + d_u_the(:,:)
395            v_seri(:,:) = v_seri(:,:) + d_v_the(:,:)
396            qmemoire(:,:)=q_seri(:,:)
397            q_seri(:,:) = q_seri(:,:) + d_q_the(:,:)
[4089]398#ifdef ISO
399            xtmemoire(:,:,:)=xt_seri(:,:,:)
400            xt_seri(:,:,:) = xt_seri(:,:,:) + d_xt_the(:,:,:)
401#ifdef ISOVERIF
402!      write(*,*) 'calltherm 350 tmp: ajout d_xt_the'
403      if (iso_HDO.gt.0) then
404!      i=479
405!      k=4
406!      write(*,*) 'xt_seri(iso_hdo,i,k),q_seri(i,k)=', &
407!     &   xt_seri(iso_hdo,i,k),q_seri(i,k)
408!      write(*,*) 'd_xt_the(iso_hdo,i,k),d_q_the(i,k)=', &
409!     &   d_xt_the(iso_hdo,i,k),d_q_the(i,k)
410      call iso_verif_aberrant_enc_vect2D( &
411     &        xt_seri,q_seri, &
[4143]412     &        'calltherm 353, apres ajout d_xt_the',ntiso,klon,klev)
[4089]413      endif     
414#endif
415#endif
[1403]416           if (prt_level.gt.10) write(lunout,*)'Apres apres thermcell_main OK'
[878]417
[879]418       DO i=1,klon
419            fm_therm(i,klev+1)=0.
[4089]420            ale_bl(i)=ale_bl(i)+ale(i)/REAL(nsplit_thermals)
421!            write(22,*)'ALE CALLTHERM',ale_bl(i),ale(i)
422            alp_bl(i)=alp_bl(i)+alp(i)/REAL(nsplit_thermals)
423!            write(23,*)'ALP CALLTHERM',alp_bl(i),alp(i)
424        if(prt_level.GE.10) print*,'calltherm i alp_bl alp ale_bl ale',i,alp_bl(i),alp(i),ale_bl(i),ale(i)
[879]425       ENDDO
426
[973]427!IM 060508 marche pas comme cela !!!        enddo ! isplit
428
[878]429!   tests sur les valeurs negatives de l'eau
[973]430         nbptspb=0
[878]431            DO k = 1, klev
432            DO i = 1, klon
[973]433               logexpr2(i,k)=.not.q_seri(i,k).ge.0.
434               if (logexpr2(i,k)) then
435                q_seri(i,k)=1.e-15
436                nbptspb=nbptspb+1
[4089]437#ifdef ISO
[4143]438                do ixt=1,ntiso
[4089]439                  xt_seri(ixt,i,k)=1.e-15*(xtmemoire(ixt,i,k)/qmemoire(i,k))
440                enddo
441#endif
[973]442!                if (prt_level.ge.10) then
443!                  print*,'WARN eau<0 apres therm i=',i,'  k=',k  &
444!    &         ,' dq,q',d_q_the(i,k),q_seri(i,k),  &
445!    &         'fm=',zfm_therm(i,k),'entr=',entr_therm(i,k)
[938]446                 endif
[878]447            ENDDO
448            ENDDO
[4089]449#ifdef ISO
450#ifdef ISOVERIF
451      if (iso_HDO.gt.0) then
452      call iso_verif_aberrant_enc_vect2D( &
453     &        xt_seri,q_seri, &
[4143]454     &        'calltherm 393, apres bidouille q<0',ntiso,klon,klev)
[4089]455      endif     
456#endif
457#endif
458
[973]459        IF(nbptspb.GT.0) print*,'Number of points with q_seri(i,k)<=0 ',nbptspb   
[878]460! tests sur les valeurs de la temperature
[973]461        nbptspb=0
[878]462            DO k = 1, klev
463            DO i = 1, klon
[973]464               logexpr2(i,k)=t_seri(i,k).lt.50..or.t_seri(i,k).gt.370.
465               if (logexpr2(i,k)) nbptspb=nbptspb+1
466!              if ((t_seri(i,k).lt.50.) .or.  &
467!    &              (t_seri(i,k).gt.370.)) then
468!                 print*,'WARN temp apres therm i=',i,'  k=',k  &
469!    &         ,' t_seri',t_seri(i,k)
[878]470!              CALL abort
[973]471!              endif
[878]472            ENDDO
473            ENDDO
[973]474        IF(nbptspb.GT.0) print*,'Number of points with q_seri(i,k)<=0 ',nbptspb
[878]475         enddo ! isplit
476
477!
478!***************************************************************
479!     calcul du flux ascencant conservatif
480!            print*,'<<<<calcul flux ascendant conservatif'
481
482      fmc_therm=0.
483               do k=1,klev
484            do i=1,klon
485                  if (entr_therm(i,k).gt.0.) then
486                     fmc_therm(i,k+1)=fmc_therm(i,k)+entr_therm(i,k)
487                  else
488                     fmc_therm(i,k+1)=fmc_therm(i,k)
489                  endif
490                  detrc_therm(i,k)=(fmc_therm(i,k+1)-fm_therm(i,k+1))  &
491     &                 -(fmc_therm(i,k)-fm_therm(i,k))
492               enddo
493            enddo
494     
495     
496!****************************************************************
497!     calcul de l'humidite dans l'ascendance
498!      print*,'<<<<calcul de lhumidite dans thermique'
499!CR:on ne le calcule que pour le cas sec
500      if (iflag_thermals.le.11) then     
501      do i=1,klon
502         zqasc(i,1)=q_seri(i,1)
503         do k=2,klev
504            if (fmc_therm(i,k+1).gt.1.e-6) then
505               zqasc(i,k)=(fmc_therm(i,k)*zqasc(i,k-1)  &
506     &              +entr_therm(i,k)*q_seri(i,k))/fmc_therm(i,k+1)
507!CR:test on asseche le thermique
508!               zqasc(i,k)=zqasc(i,k)/2.
509!            else
510!               zqasc(i,k)=q_seri(i,k)
511            endif
512         enddo
513       enddo
514     
515
516!     calcul de l'eau condensee dans l'ascendance
517!             print*,'<<<<calcul de leau condensee dans thermique'
518             do i=1,klon
519                do k=1,klev
520                   clwcon0(i,k)=zqasc(i,k)-zqsat(i,k)
521                   if (clwcon0(i,k).lt.0. .or.   &
522     &             (fm_therm(i,k+1)+detrc_therm(i,k)).lt.1.e-6) then
523                      clwcon0(i,k)=0.
524                   endif
525                enddo
526             enddo
527       else
528              do i=1,klon
529                do k=1,klev
530                   clwcon0(i,k)=zqla(i,k) 
531                   if (clwcon0(i,k).lt.0. .or.   &
532     &             (fm_therm(i,k+1)+detrc_therm(i,k)).lt.1.e-6) then
533                   clwcon0(i,k)=0.
534                   endif
535                enddo
536             enddo
537       endif
538!*******************************************************************   
539
540
[1428]541!jyg  Protection contre les temperatures nulles
542          do i=1,klon
543             do k=1,klev
544                if (ztla(i,k) .lt. 1.e-10) fraca(i,k) =0.
545             enddo
546          enddo
547
548
[878]549      return
550
551      end
Note: See TracBrowser for help on using the repository browser.