source: trunk/LMDZ.GENERIC/libf/phystd/thermcell_plume.F90 @ 2101

Last change on this file since 2101 was 2101, checked in by aboissinot, 6 years ago

Fix a bug in thermcell_alim.F90 where vertical and horizontal loops must be inverted.
In thermal plume model, arrays size declaration is standardised (no longer done thanks to dimphy module but by way of arguments).
Clean up some thermal plume model routines (remove uselesss variables...)

File size: 22.9 KB
RevLine 
[2060]1!
2!
3!
[2101]4SUBROUTINE thermcell_plume(itap,ngrid,nlay,ptimestep,ztv,                     &
5                           zthl,po,zl,rhobarz,zlev,pplev,pphi,zpspsk,         &
6                           alim_star,alim_star_tot,lalim,f0,detr_star,        &
7                           entr_star,f_star,ztva,ztla,zqla,zqta,zha,          &
8                           zw2,w_est,ztva_est,zqsatth,lmix,lmix_bis,          &
9                           lmin,lev_out,lunout1,igout)
[2060]10     
11     
12!==============================================================================
13! thermcell_plume: calcule les valeurs de qt, thetal et w dans l ascendance
[2069]14! AB : ql means "liquid water mass mixing ratio"
15!      qt means "total water mass mixing ratio"
16!      TP means "potential temperature"
17!      TRPV means "virtual potential temperature with latent heat release" 
18!      TPV means "virtual potential temperature"
19!      TR means "temperature with latent heat release"
[2060]20!==============================================================================
21     
22      USE print_control_mod, ONLY: prt_level
[2071]23      USE watercommon_h, ONLY: RLvCp, RETV, Psat_water
[2060]24      USE thermcell_mod
25     
26      IMPLICIT NONE
27     
28     
29!==============================================================================
30! Declaration
31!==============================================================================
32     
33!      inputs:
34!      -------
35     
36      INTEGER itap
37      INTEGER ngrid
[2101]38      INTEGER nlay
[2060]39      INTEGER lunout1
40      INTEGER igout
41      INTEGER lev_out                           ! niveau pour les print
42     
43      REAL ptimestep                            ! time step
[2101]44      REAL ztv(ngrid,nlay)                      ! TRPV environment
45      REAL zthl(ngrid,nlay)                     ! TP   environment
46      REAL po(ngrid,nlay)                       ! qt   environment
47      REAL zl(ngrid,nlay)                       ! ql   environment
48      REAL rhobarz(ngrid,nlay)                  ! levels density
49      REAL zlev(ngrid,nlay+1)                   ! levels altitude
50      REAL pplev(ngrid,nlay+1)                  ! levels pressure
51      REAL pphi(ngrid,nlay)                     ! geopotential
52      REAL zpspsk(ngrid,nlay)                   ! Exner function
[2060]53     
54!      outputs:
55!      --------
56     
57      INTEGER lmin(ngrid)                       ! plume base level (first unstable level)
58      INTEGER lalim(ngrid)                      ! higher alimentation level
59      INTEGER lmix(ngrid)                       ! maximum vertical speed level
60      INTEGER lmix_bis(ngrid)                   ! maximum vertical speed level (modified)
61     
[2101]62      REAL alim_star(ngrid,nlay)                ! normalized alimentation
[2060]63      REAL alim_star_tot(ngrid)                 ! integrated alimentation
64     
65      REAL f0(ngrid)                            ! previous time step mass flux norm
66     
[2101]67      REAL detr_star(ngrid,nlay)                ! normalized detrainment
68      REAL entr_star(ngrid,nlay)                ! normalized entrainment
69      REAL f_star(ngrid,nlay+1)                 ! normalized mass flux
[2060]70     
[2101]71      REAL ztva(ngrid,nlay)                     ! TRPV plume (after mixing)
72      REAL ztva_est(ngrid,nlay)                 ! TRPV plume (before mixing)
73      REAL ztla(ngrid,nlay)                     ! TP   plume
74      REAL zqla(ngrid,nlay)                     ! ql   plume (after mixing)
75      REAL zqta(ngrid,nlay)                     ! qt   plume
76      REAL zha(ngrid,nlay)                      ! TRPV plume
[2060]77     
[2101]78      REAL w_est(ngrid,nlay+1)                  ! updraft square vertical speed (before mixing)
79      REAL zw2(ngrid,nlay+1)                    ! updraft square vertical speed (after mixing)
[2060]80     
[2101]81      REAL zqsatth(ngrid,nlay)                  ! saturation vapor pressure (after mixing)
[2060]82     
83!      local:
84!      ------
85     
86      INTEGER ig, l, k
87      INTEGER lt
88      INTEGER lm
89     
[2101]90      REAL zqla_est(ngrid,nlay)                 ! ql   plume (before mixing)
91      REAL zta_est(ngrid,nlay)                  ! TR   plume (before mixing)
92      REAL zbuoy(ngrid,nlay)                    ! plume buoyancy
93      REAL zbuoyjam(ngrid,nlay)                 ! plume buoyancy (modified)
[2060]94     
95      REAL ztemp(ngrid)                         ! temperature for saturation vapor pressure computation in plume
96      REAL zqsat(ngrid)                         ! saturation vapor pressure (before mixing)
97      REAL zdz                                  ! layers height
[2101]98      REAL ztv2(ngrid,nlay)                     ! ztv + d_temp * Dirac(l=lmin)
[2060]99     
100      REAL zalpha                               !
[2101]101      REAL zdqt(ngrid,nlay)                     !
[2060]102      REAL zbetalpha                            !
103      REAL zdw2                                 !
104      REAL zdw2bis                              !
105      REAL zw2fact                              !
106      REAL zw2factbis                           !
107      REAL zw2m                                 !
108      REAL zdzbis                               !
109      REAL coefzlmel                            !
110      REAL zdz2                                 !
111      REAL zdz3                                 !
112      REAL lmel                                 !
113      REAL zlmel                                !
114      REAL zlmelup                              !
115      REAL zlmeldwn                             !
116      REAL zlt                                  !
117      REAL zltdwn                               !
118      REAL zltup                                ! useless here
119     
[2101]120      REAL psat                                 ! dummy argument for Psat_water()
[2071]121     
[2060]122      LOGICAL active(ngrid)                     ! if the plume is active at ig,l (speed and incoming mass flux > 0 or l=lmin)
123      LOGICAL activetmp(ngrid)                  ! if the plus is active at ig,l (active=true and outgoing mass flux > 0)
124      LOGICAL, SAVE :: first = .true.           ! if it is the first time step
125     
126!$OMP THREADPRIVATE(first)
127     
128!==============================================================================
129! Initialization
130!==============================================================================
131     
132      zbetalpha = betalpha / (1. + betalpha)
133     
[2069]134      ztva(:,:)        = ztv(:,:)               ! ztva     is set to the virtual potential temperature without latent heat release
135      ztva_est(:,:)    = ztva(:,:)              ! ztva_est is set to the virtual potential temperature without latent heat release
[2066]136      ztla(:,:)        = zthl(:,:)              ! ztla     is set to the potential temperature
137      zqta(:,:)        = po(:,:)                ! zqta     is set to qt
138      zqla(:,:)        = 0.                     ! zqla     is set to ql
139      zqla_est(:,:)    = 0.                     ! zqla_est is set to ql
[2069]140      zha(:,:)         = ztva(:,:)              ! zha      is set to the plume virtual potential temperature without latent heat release
[2060]141     
142      zqsat(:)         = 0.
143      zqsatth(:,:)     = 0.
144     
145      w_est(:,:)       = 0.
146      zw2(:,:)         = 0.
147     
148      zbuoy(:,:)       = 0.
149      zbuoyjam(:,:)    = 0.
150     
151      f_star(:,:)      = 0.
152      detr_star(:,:)   = 0.
153      entr_star(:,:)   = 0.
154      alim_star(:,:)   = 0.
155      alim_star_tot(:) = 0.
156     
157      lmix(:)          = 1
158      lmix_bis(:)      = 2
159      lalim(:)         = 1
160      lmin(:)          = linf
161     
[2093]162      ztv2(:,:)        = ztv(:,:)
[2101]163      ztv2(:,linf)     = ztv(:,linf) + d_temp
[2093]164     
[2060]165!==============================================================================
166! 0. Calcul de l'alimentation
167!==============================================================================
168     
169!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
170! AB : Convective plumes can go off from every layer above the linf-th and
171!      where pressure is lesser than pres_limit (cf. thermcell_plume).
172!      The second constraint is added to avoid the parametrization occurs too
173!      high when the low atmosphere is stable.
174!      However, once there is a triggered plume, it can rise as high as its
175!      velocity allows it (it can overshoot).
176!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
177      DO ig=1,ngrid
178         active(ig) = .false.
179         l = linf
[2101]180         DO WHILE ((.not.active(ig)) .and. pplev(ig,l+1).gt.pres_limit .and. l.lt.nlay)
[2093]181            IF (ztv2(ig,l).gt.ztv2(ig,l+1)) THEN
[2060]182               active(ig) = .true.
183               lmin(ig) = l
184            ENDIF
185            l = l + 1
186         ENDDO
187      ENDDO
188     
189!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2065]190! AB : On pourrait n'appeler thermcell_alim que si la plume est active
[2060]191!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2101]192      CALL thermcell_alim(ngrid,nlay,ztv2,zlev,alim_star,lalim,lmin)
[2060]193     
194!==============================================================================
195! 1. Calcul dans la premiere couche
196!==============================================================================
197     
198      DO ig=1,ngrid
199         IF (active(ig)) THEN
200           
201!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
202! AB : plume takes the environment features for every layer below lmin.
203!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
204            DO l=1,lmin(ig)
205               ztla(ig,l) = zthl(ig,l)
206               zqta(ig,l) = po(ig,l)
207               zqla(ig,l) = zl(ig,l)
208            ENDDO
209           
210            l = lmin(ig)
211            f_star(ig,l+1) = alim_star(ig,l)
212           
213            zw2(ig,l+1) = 2. * RG * (zlev(ig,l+1) - zlev(ig,l))               &
[2093]214            &                * (ztv2(ig,l) - ztv2(ig,l+1)) / ztv2(ig,l+1)
[2060]215           
216            w_est(ig,l+1) = zw2(ig,l+1)
217         ENDIF
218      ENDDO
219     
220!==============================================================================
221! 2. Boucle de calcul de la vitesse verticale dans le thermique
222!==============================================================================
223     
[2101]224      DO l=2,nlay-1
[2060]225         
226!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
227! AB : we decide here if the plume is still active or not. When the plume's
228!      first level is reached, we set active to "true". Otherwise, it is given
[2065]229!      by zw2, f_star, alim_star and entr_star.
[2060]230!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
231         DO ig=1,ngrid
232            IF (l==lmin(ig)+1) THEN
233               active(ig) = .true.
234            ENDIF
235           
236            active(ig) = active(ig)                                           &
237            &            .and. zw2(ig,l)>1.e-10                               &
238            &            .and. f_star(ig,l)+alim_star(ig,l)+entr_star(ig,l)>1.e-10
239         ENDDO
240         
241         ztemp(:) = zpspsk(:,l) * ztla(:,l-1)
[2065]242         
[2060]243         DO ig=1,ngrid
[2101]244            CALL Psat_water(ztemp(ig), pplev(ig,l), psat, zqsat(ig))
[2060]245         ENDDO
246         
247!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
248! AB : we compute thermodynamical values and speed in the plume in the layer l
249!      without mixing with environment.
250!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
251         
252         DO ig=1,ngrid
253            IF (active(ig)) THEN
254               zqla_est(ig,l) = max(0.,zqta(ig,l-1)-zqsat(ig))                   ! zqla_est set to ql   plume
255               ztva_est(ig,l) = ztla(ig,l-1)*zpspsk(ig,l)+RLvCp*zqla_est(ig,l)   ! ztva_est set to TR   plume
256               zta_est(ig,l)  = ztva_est(ig,l)                                   ! zta_est  set to TR   plume
257               ztva_est(ig,l) = ztva_est(ig,l)/zpspsk(ig,l)                      ! ztva_est set to TRP  plume
258               ztva_est(ig,l) = ztva_est(ig,l)*(1.+RETV*(zqta(ig,l-1)         &  ! ztva_est set to TRPV plume
259               &              - zqla_est(ig,l))-zqla_est(ig,l))
260               
261               zbuoy(ig,l) = RG * (ztva_est(ig,l)-ztv(ig,l)) / ztv(ig,l)
262               zdz = zlev(ig,l+1) - zlev(ig,l)
263               
264!==============================================================================
265! 3. Calcul de la flotabilite modifiee par melange avec l'air au dessus
266!==============================================================================
267               
268               lmel = fact_thermals_ed_dz * zlev(ig,l)
269               zlmel = zlev(ig,l) + lmel
270               lt = l + 1
271               zlt = zlev(ig,lt)
272               zdz2 = zlev(ig,lt) - zlev(ig,l)
273               
274               DO while (lmel.gt.zdz2)
275                  lt   = lt + 1
276                  zlt  = zlev(ig,lt)
277                  zdz2 = zlev(ig,lt) - zlev(ig,l)
278               ENDDO
279               
280!               IF (lt-l.gt.1) THEN
281!                  print *, 'WARNING: lt is greater than l+1!'
282!                  print *, 'l,lt', l, lt
283!               ENDIF
284               
285               zdz3 = zlev(ig,lt+1) - zlt
286               zltdwn = zlev(ig,lt) - zdz3 / 2
287               zlmelup = zlmel + (zdz / 2)
288               coefzlmel = Min(1.,(zlmelup - zltdwn) / zdz)
289               
290               zbuoyjam(ig,l) = 1.* RG * (coefzlmel *                         &
291               &                (ztva_est(ig,l) - ztv(ig,lt)) / ztv(ig,lt)    &
292               &              + (1. - coefzlmel) *                            &
293               &                (ztva_est(ig,l) - ztv(ig,lt-1)) / ztv(ig,lt-1)) &
294               &              + 0. * zbuoy(ig,l)
295               
296!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
297! AB : initial formulae
298!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
299!               zw2fact = fact_epsilon * 2. * zdz / (1. + betalpha)
300!               zdw2 = afact * zbuoy(ig,l) / fact_epsilon
301!               zdw2bis = afact * zbuoy(ig,l-1) / fact_epsilon
302!               w_est(ig,l+1) = Max(0.0001,exp(-zw2fact)*(w_est(ig,l)-zdw2)+zdw2)
303!               w_est(ig,l+1) = Max(0.0001,exp(-zw2fact)*(w_est(ig,l)-zdw2bis)+zdw2)
304!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2065]305! AB : own derivation for w_est (Rio 2010 formula with e=d=0)
[2060]306!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
307               zw2fact = 2. * fact_epsilon * zdz
308               zdw2 = 2. * afact * zbuoy(ig,l) * zdz
309               w_est(ig,l+1) = Max(0., exp(-zw2fact) * w_est(ig,l) + zdw2)
310               
311!               IF (w_est(ig,l+1).le.0.) THEN
312!                  print *, 'WARNING: w_est is negative!'
313!                  print *, 'l,w_est', l+1, w_est(ig,l+1)
314!               ENDIF
315            ENDIF
316         ENDDO
317         
318!==============================================================================
319! 4. Calcul de l'entrainement et du detrainement
320!==============================================================================
321         
322         DO ig=1,ngrid
323            IF (active(ig)) THEN
324               
325               zdz = zlev(ig,l+1) - zlev(ig,l)
326               
327!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2066]328! AB : The next test is added to avoid divisions by zero when w_est vanishes.
[2065]329!      Indeed, entr and detr computed here are of no importance because w_est
330!      <= 0 means it will be the last layer reached by the plume and then they
331!      will be reset in thermcell_flux.
[2060]332!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2065]333               IF (w_est(ig,l+1).eq.0.) THEN
334                  zw2m = 1.
[2066]335                  zalpha = 0.
[2060]336               ELSE
[2065]337                  zw2m = w_est(ig,l+1)
[2066]338                  zalpha = f0(ig) * f_star(ig,l) / sqrt(w_est(ig,l+1)) / rhobarz(ig,l)
[2060]339               ENDIF
340               
341!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2065]342! AB : The next test is added to avoid a division by zero if there is no water
343!      in the environment.
[2060]344!      In the case where there is no water in the env. but water in the plume
345!      (ascending from depth) we set the effect on detrainment equal to zero
346!      but at the next time step, po will be positive thanks to the mixing and
347!      then the physical effect of the water gradient will be taken on board.
348!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
349               IF (po(ig,l).lt.1.e-6) THEN
[2065]350!                  print *, 'WARNING: po=0 in layer',l,'!'
[2060]351!                  print *, 'po,zqta', po(ig,l), zqta(ig,l-1)
352                  zdqt(ig,l) = 0.0
353               ELSE
354                  zdqt(ig,l) = max(zqta(ig,l-1)-po(ig,l),0.) / po(ig,l)
355               ENDIF
356               
357!------------------------------------------------------------------------------
358! Detrainment
359!------------------------------------------------------------------------------
360               
361               detr_star(ig,l) = f_star(ig,l) * zdz * (                       &
362               &                 mix0 * 0.1 / (zalpha + 0.001)                &
363               &               + MAX(detr_min,                                &
364               &                 -afact * zbetalpha * zbuoyjam(ig,l) / zw2m   &
365               &               + detr_q_coef*(zdqt(ig,l)/zw2m)**detr_q_power) )
366               
367!               IF (detr_star(ig,l).lt.0.) THEN
368!                  print *, 'WARNING: detrainment is negative!'
369!                  print *, 'l,detr', l, detr_star(ig,l)
370!               ENDIF
371               
372!------------------------------------------------------------------------------
373! Entrainment
374!------------------------------------------------------------------------------
375               
376               entr_star(ig,l) = f_star(ig,l) * zdz * (                       &
377               &                 mix0 * 0.1 / (zalpha+0.001)                  &
378               &               + MAX(entr_min,                                &
379               &                 zbetalpha * afact * zbuoyjam(ig,l) / zw2m    &
380               &               - zbetalpha * fact_epsilon) )
381               
382!               IF (entr_star(ig,l).lt.0.) THEN
383!                  print *, 'WARNING: entrainment is negative!'
384!                  print *, 'l,entr', l, entr_star(ig,l)
385!               ENDIF
386               
387!------------------------------------------------------------------------------
388! Alimentation and entrainment
389!------------------------------------------------------------------------------
390               
391               IF (l.lt.lalim(ig)) THEN
392                  alim_star(ig,l) = max(alim_star(ig,l),entr_star(ig,l))
393                  entr_star(ig,l) = 0.
394               ENDIF
395               
396!------------------------------------------------------------------------------
397! Mass flux
398!------------------------------------------------------------------------------
399               
400               f_star(ig,l+1) = f_star(ig,l) + alim_star(ig,l)                &
401               &              + entr_star(ig,l) - detr_star(ig,l)
402               
403!               IF (f_star(ig,l+1).le.0.) THEN
404!                  print *, 'WARNING: mass flux is negative!'
405!                  print *, 'l,f_star', l+1, f_star(ig,l+1)
406!               ENDIF
407               
408            ENDIF
409         ENDDO
410         
411!==============================================================================
412! 5. Calcul de la vitesse verticale en melangeant Tl et qt du thermique
413!==============================================================================
414         
415         activetmp(:) = active(:) .and. f_star(:,l+1)>1.e-10
416         
417!------------------------------------------------------------------------------
418! Calcul du melange avec l'environnement
419!------------------------------------------------------------------------------
420         
421         DO ig=1,ngrid
422            IF (activetmp(ig)) THEN
[2065]423               ztla(ig,l) = (f_star(ig,l) * ztla(ig,l-1)                      &  ! ztla is set to TP in plume (mixed)
424               &          + (alim_star(ig,l) + entr_star(ig,l)) * zthl(ig,l)) &
425               &          / (f_star(ig,l+1) + detr_star(ig,l))
426               zqta(ig,l) = (f_star(ig,l) * zqta(ig,l-1) +                    &  ! zqta is set to qt in plume (mixed)
427               &          + (alim_star(ig,l) + entr_star(ig,l)) * po(ig,l))   &
428               &          / (f_star(ig,l+1) + detr_star(ig,l))
[2060]429            ENDIF
430         ENDDO
431         
432         ztemp(:) = zpspsk(:,l) * ztla(:,l)
433         
434         DO ig=1,ngrid
435            IF (activetmp(ig)) THEN
[2101]436               CALL Psat_water(ztemp(ig), pplev(ig,l), psat, zqsatth(ig,l))
[2060]437            ENDIF
438         ENDDO
439         
440!------------------------------------------------------------------------------
441! Calcul de la vitesse verticale zw2 apres le melange
442!------------------------------------------------------------------------------
443         
444         DO ig=1,ngrid
445            IF (activetmp(ig)) THEN
446               zqla(ig,l) = max(0.,zqta(ig,l)-zqsatth(ig,l))                     ! zqla is set to ql   plume (mixed)
447               ztva(ig,l) = ztla(ig,l) * zpspsk(ig,l)+RLvCp*zqla(ig,l)           ! ztva is set to TR   plume (mixed)
448               ztva(ig,l) = ztva(ig,l) / zpspsk(ig,l)                            ! ztva is set to TRP  plume (mixed)
449               zha(ig,l)  = ztva(ig,l)                                           ! zha  is set to TRP  plume (mixed)
450               ztva(ig,l) = ztva(ig,l) * (1. + RETV*(zqta(ig,l)-zqla(ig,l))   &  ! ztva is set to TRPV plume (mixed)
451               &                             - zqla(ig,l))
452               
453               zbuoy(ig,l) = RG * (ztva(ig,l) - ztv(ig,l)) / ztv(ig,l)
454               zdz = zlev(ig,l+1) - zlev(ig,l)
455               
456!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2065]457! AB : initial formula
[2060]458!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
459!               zw2fact = fact_epsilon * 2. * zdz / (1. + betalpha)
460!               zdw2 = afact * zbuoy(ig,l) / fact_epsilon
461!               zdw2bis = afact * zbuoy(ig,l-1) / fact_epsilon
462!               zw2(ig,l+1) = Max(0.0001,exp(-zw2fact)*(zw2(ig,l)-zdw2)+zdw2)
463!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[2065]464! AB : own derivation for zw2 (Rio 2010 formula)
[2060]465!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
466               zw2fact = 2. * (fact_epsilon * zdz  + entr_star(ig,l) / f_star(ig,l))
467               zdw2 = 2. * afact * zbuoy(ig,l) * zdz
468               zw2(ig,l+1) = Max(0., exp(-zw2fact) * zw2(ig,l) + zdw2)
469               
470!               IF (zw2(ig,l+1).le.0.) THEN
471!                  print *, 'WARNING: zw2 is negative!'
472!                  print *, 'l,zw2', l+1, zw2(ig,l+1)
473!               ENDIF
474            ENDIF
475         ENDDO
476         
477      ENDDO
478     
479!==============================================================================
480! 6. New computation of alim_star_tot
481!==============================================================================
482     
483      DO ig=1,ngrid
[2065]484         alim_star_tot(ig) = 0.
[2060]485      ENDDO
486     
487      DO ig=1,ngrid
488         DO l=1,lalim(ig)-1
489            alim_star_tot(ig) = alim_star_tot(ig) + alim_star(ig,l)
490         ENDDO
491      ENDDO
492     
493     
494RETURN
495END
Note: See TracBrowser for help on using the repository browser.