source: LMDZ6/branches/Amaury_dev/libf/phylmd/hbtm_mod.F90 @ 5143

Last change on this file since 5143 was 5143, checked in by abarral, 8 weeks ago

Put YOEGWD.h, FCTTRE.h into modules

  • 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: 29.0 KB
Line 
1module hbtm_mod
2
3  IMPLICIT NONE
4
5CONTAINS
6
7  SUBROUTINE hbtm(knon, paprs, pplay, t2m, t10m, q2m, q10m, ustar, wstar, &
8          flux_t, flux_q, u, v, t, q, pblh, cape, eauliq, ctei, pblt, therm, &
9          trmb1, trmb2, trmb3, plcl)
10    USE dimphy
11    USE lmdz_YOETHF
12    USE lmdz_fcttre, ONLY: foeew, foede, qsats, qsatl, dqsats, dqsatl, thermcep
13
14    ! ***************************************************************
15    ! *                                                             *
16    ! * HBTM2   D'apres Holstag&Boville et Troen&Mahrt              *
17    ! *                 JAS 47              BLM                     *
18    ! * Algorithme These Anne Mathieu                               *
19    ! * Critere d'Entrainement Peter Duynkerke (JAS 50)             *
20    ! * written by  : Anne MATHIEU & Alain LAHELLEC, 22/11/99       *
21    ! * features : implem. exces Mathieu                            *
22    ! ***************************************************************
23    ! * mods : decembre 99 passage th a niveau plus bas. voir fixer *
24    ! * la prise du th a z/Lambda = -.2 (max Ray)                   *
25    ! * Autre algo : entrainement ~ Theta+v =cste mais comment=>The?*
26    ! * on peut fixer q a .7qsat(cf non adiab)=>T2 et The2          *
27    ! * voir aussi //KE pblh = niveau The_e ou l = env.             *
28    ! ***************************************************************
29    ! * fin therm a la HBTM passage a forme Mathieu 12/09/2001      *
30    ! ***************************************************************
31    ! *
32
33
34    ! AM Fev 2003
35    ! Adaptation a LMDZ version couplee
36
37    ! Pour le moment on fait passer en argument les grdeurs de surface :
38    ! flux, t,q2m, t,q10m, on va utiliser systematiquement les grdeurs a 2m ms
39    ! on garde la possibilite de changer si besoin est (jusqu'a present la
40    ! forme de HB avec le 1er niveau modele etait conservee)
41
42    include "YOMCST.h"
43    REAL rlvcp, reps
44    ! Arguments:
45
46    INTEGER knon ! nombre de points a calculer
47    ! AM
48    REAL t2m(klon), t10m(klon) ! temperature a 2 et 10m
49    REAL q2m(klon), q10m(klon) ! q a 2 et 10m
50    REAL ustar(klon)
51    REAL wstar(klon) ! w*, convective velocity scale
52    REAL paprs(klon, klev + 1) ! pression a inter-couche (Pa)
53    REAL pplay(klon, klev) ! pression au milieu de couche (Pa)
54    REAL flux_t(klon, klev), flux_q(klon, klev) ! Flux
55    REAL u(klon, klev) ! vitesse U (m/s)
56    REAL v(klon, klev) ! vitesse V (m/s)
57    REAL t(klon, klev) ! temperature (K)
58    REAL q(klon, klev) ! vapeur d'eau (kg/kg)
59    ! AM      REAL cd_h(klon) ! coefficient de friction au sol pour chaleur
60    ! AM      REAL cd_m(klon) ! coefficient de friction au sol pour vitesse
61
62    INTEGER isommet
63    ! um      PARAMETER (isommet=klev) ! limite max sommet pbl
64    REAL, PARAMETER :: vk = 0.35 ! Von Karman => passer a .41 ! cf U.Olgstrom
65    REAL, PARAMETER :: ricr = 0.4
66    REAL, PARAMETER :: fak = 8.5 ! b calcul du Prandtl et de dTetas
67    REAL, PARAMETER :: fakn = 7.2 ! a
68    REAL, PARAMETER :: onet = 1.0 / 3.0
69    REAL, PARAMETER :: t_coup = 273.15
70    REAL, PARAMETER :: zkmin = 0.01
71    REAL, PARAMETER :: betam = 15.0 ! pour Phim / h dans la S.L stable
72    REAL, PARAMETER :: betah = 15.0
73
74    REAL, PARAMETER :: betas = 5.0
75    ! Phit dans la S.L. stable (mais 2 formes / z/OBL<>1
76
77    REAL, PARAMETER :: sffrac = 0.1 ! S.L. = z/h < .1
78    REAL, PARAMETER :: usmin = 1.E-12
79    REAL, PARAMETER :: binm = betam * sffrac
80    REAL, PARAMETER :: binh = betah * sffrac
81    REAL, PARAMETER :: ccon = fak * sffrac * vk
82    REAL, PARAMETER :: b1 = 70., b2 = 20.
83    REAL, PARAMETER :: zref = 2. ! Niveau de ref a 2m peut eventuellement
84    ! etre choisi a 10m
85    REAL q_star, t_star
86    REAL b212, b2sr ! Lambert correlations T' q' avec T* q*
87
88    REAL z(klon, klev)
89    ! AM      REAL pcfm(klon,klev), pcfh(klon,klev)
90    INTEGER i, k, j
91    REAL zxt
92    ! AM      REAL zxt, zxq, zxu, zxv, zxmod, taux, tauy
93    ! AM      REAL zx_alf1, zx_alf2 ! parametres pour extrapolation
94    REAL khfs(klon) ! surface kinematic heat flux [mK/s]
95    REAL kqfs(klon) ! sfc kinematic constituent flux [m/s]
96    REAL heatv(klon) ! surface virtual heat flux
97    REAL rhino(klon, klev) ! bulk Richardon no. mais en Theta_v
98    LOGICAL unstbl(klon) ! pts w/unstbl pbl (positive virtual ht flx)
99    LOGICAL stblev(klon) ! stable pbl with levels within pbl
100    LOGICAL unslev(klon) ! unstbl pbl with levels within pbl
101    LOGICAL unssrf(klon) ! unstb pbl w/lvls within srf pbl lyr
102    LOGICAL unsout(klon) ! unstb pbl w/lvls in outer pbl lyr
103    LOGICAL check(klon) ! True=>chk if Richardson no.>critcal
104    LOGICAL omegafl(klon) ! flag de prolongerment cape pour pt Omega
105    REAL pblh(klon)
106    REAL pblt(klon)
107    REAL plcl(klon)
108    ! AM      REAL cgh(klon,2:klev) ! counter-gradient term for heat [K/m]
109    ! AM      REAL cgq(klon,2:klev) ! counter-gradient term for constituents
110    ! AM      REAL cgs(klon,2:klev) ! counter-gradient star (cg/flux)
111    REAL unsobklen(klon) ! Monin-Obukhov lengh
112    ! AM      REAL ztvd, ztvu,
113    REAL zdu2
114    REAL, INTENT(OUT) :: therm(:) ! (klon) thermal virtual temperature excess
115    REAL trmb1(klon), trmb2(klon), trmb3(klon)
116    ! Algorithme thermique
117    REAL s(klon, klev) ! [P/Po]^Kappa milieux couches
118    REAL th_th(klon) ! potential temperature of thermal
119    REAL the_th(klon) ! equivalent potential temperature of thermal
120    REAL qt_th(klon) ! total water  of thermal
121    REAL tbef(klon) ! T thermique niveau precedent
122    REAL qsatbef(klon)
123    LOGICAL zsat(klon) ! le thermique est sature
124    REAL cape(klon) ! Cape du thermique
125    REAL kape(klon) ! Cape locale
126    REAL eauliq(klon) ! Eau liqu integr du thermique
127    REAL ctei(klon) ! Critere d'instab d'entrainmt des nuages de CL
128    REAL the1, the2, aa, bb, zthvd, zthvu, xintpos, qqsat
129    ! IM 091204 BEG
130    REAL a1, a2, a3
131    ! IM 091204 END
132    REAL xhis, rnum, denom, th1, th2, thv1, thv2, ql2
133    REAL dqsat_dt, qsat2, qt1, q2, t1, t2, xnull, delt_the
134    REAL delt_qt, delt_2, quadsat, spblh, reduc
135
136    REAL phiminv(klon) ! inverse phi function for momentum
137    REAL phihinv(klon) ! inverse phi function for heat
138    REAL wm(klon) ! turbulent velocity scale for momentum
139    REAL fak1(klon) ! k*ustar*pblh
140    REAL fak2(klon) ! k*wm*pblh
141    REAL fak3(klon) ! fakn*wstar/wm
142    REAL pblk(klon) ! level eddy diffusivity for momentum
143    REAL pr(klon) ! Prandtl number for eddy diffusivities
144    REAL zl(klon) ! zmzp / Obukhov length
145    REAL zh(klon) ! zmzp / pblh
146    REAL zzh(klon) ! (1-(zmzp/pblh))**2
147    REAL zm(klon) ! current level height
148    REAL zp(klon) ! current level height + one level up
149    REAL zcor, zdelta, zcvm5
150    ! AM      REAL zxqs
151    REAL fac, pblmin, zmzp, term
152
153    ! initialisations (Anne)
154    isommet = klev
155    th_th(:) = 0.
156    q_star = 0
157    t_star = 0
158    therm = 0.
159
160    b212 = sqrt(b1 * b2)
161    b2sr = sqrt(b2)
162
163    ! ============================================================
164    ! Fonctions thermo implicites
165    ! ============================================================
166    ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
167    ! Tetens : pression partielle de vap d'eau e_sat(T)
168    ! =================================================
169    ! ++ e_sat(T) = r2*exp( r3*(T-Tf)/(T-r4) ) id a r2*FOEWE
170    ! ++ avec :
171    ! ++ Tf = 273.16 K  (Temp de fusion de la glace)
172    ! ++ r2 = 611.14 Pa
173    ! ++ r3 = 17.269 (liquide) 21.875 (solide) adim
174    ! ++ r4 = 35.86             7.66           Kelvin
175    ! ++  q_sat = eps*e_sat/(p-(1-eps)*e_sat)
176    ! ++ derivé :
177    ! ++ =========
178    ! ++                   r3*(Tf-r4)*q_sat(T,p)
179    ! ++ d_qsat_dT = --------------------------------
180    ! ++             (T-r4)^2*( 1-(1-eps)*e_sat(T)/p )
181    ! ++ pour zcvm5=Lv, c'est FOEDE
182    ! ++ Rq :(1.-REPS)*esarg/Parg id a RETV*Qsat
183    ! ------------------------------------------------------------------
184
185    ! Initialisation
186    rlvcp = rlvtt / rcpd
187    reps = rd / rv
188
189
190    ! DO i = 1, klon
191    ! pcfh(i,1) = cd_h(i)
192    ! pcfm(i,1) = cd_m(i)
193    ! ENDDO
194    ! DO k = 2, klev
195    ! DO i = 1, klon
196    ! pcfh(i,k) = zkmin
197    ! pcfm(i,k) = zkmin
198    ! cgs(i,k) = 0.0
199    ! cgh(i,k) = 0.0
200    ! cgq(i,k) = 0.0
201    ! ENDDO
202    ! ENDDO
203
204    ! Calculer les hauteurs de chaque couche
205    ! (geopotentielle Int_dp/ro = Int_[Rd.T.dp/p] z = geop/g)
206    ! pourquoi ne pas utiliser Phi/RG ?
207    DO i = 1, knon
208      z(i, 1) = rd * t(i, 1) / (0.5 * (paprs(i, 1) + pplay(i, 1)))&
209              * (paprs(i, 1) - pplay(i, 1)) / rg
210      s(i, 1) = (pplay(i, 1) / paprs(i, 1))**rkappa
211    END DO
212    ! s(k) = [pplay(k)/ps]^kappa
213    ! + + + + + + + + + pplay  <-> s(k)   t  dp=pplay(k-1)-pplay(k)
214
215    ! -----------------  paprs <-> sig(k)
216
217    ! + + + + + + + + + pplay  <-> s(k-1)
218
219
220    ! + + + + + + + + + pplay  <-> s(1)   t  dp=paprs-pplay   z(1)
221
222    ! -----------------  paprs <-> sig(1)
223
224    DO k = 2, klev
225      DO i = 1, knon
226        z(i, k) = z(i, k - 1) + rd * 0.5 * (t(i, k - 1) + t(i, k)) / paprs(i, k)&
227                * (pplay(i, k - 1) - pplay(i, k)) / rg
228        s(i, k) = (pplay(i, k) / paprs(i, 1))**rkappa
229      END DO
230    END DO
231    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
232    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
233    ! +++  Determination des grandeurs de surface  +++++++++++++++++++++
234    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
235    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
236    DO i = 1, knon
237      ! AM         IF (thermcep) THEN
238      ! AM           zdelta=MAX(0.,SIGN(1.,RTT-tsol(i)))
239      ! zcvm5 = R5LES*RLVTT*(1.-zdelta) + R5IES*RLSTT*zdelta
240      ! zcvm5 = zcvm5 / RCPD / (1.0+RVTMP2*q(i,1))
241      ! AM           zxqs= r2es * FOEEW(tsol(i),zdelta)/paprs(i,1)
242      ! AM           zxqs=MIN(0.5,zxqs)
243      ! AM           zcor=1./(1.-retv*zxqs)
244      ! AM           zxqs=zxqs*zcor
245      ! AM         ELSE
246      ! AM           IF (tsol(i).LT.t_coup) THEN
247      ! AM              zxqs = qsats(tsol(i)) / paprs(i,1)
248      ! AM           ELSE
249      ! AM              zxqs = qsatl(tsol(i)) / paprs(i,1)
250      ! AM           ENDIF
251      ! AM         ENDIF
252      ! niveau de reference bulk; mais ici, c,a pourrait etre le niveau de ref
253      ! du thermique
254      ! AM        zx_alf1 = 1.0
255      ! AM        zx_alf2 = 1.0 - zx_alf1
256      ! AM        zxt = (t(i,1)+z(i,1)*RG/RCPD/(1.+RVTMP2*q(i,1)))
257      ! AM     .        *(1.+RETV*q(i,1))*zx_alf1
258      ! AM     .      + (t(i,2)+z(i,2)*RG/RCPD/(1.+RVTMP2*q(i,2)))
259      ! AM     .        *(1.+RETV*q(i,2))*zx_alf2
260      ! AM        zxu = u(i,1)*zx_alf1+u(i,2)*zx_alf2
261      ! AM        zxv = v(i,1)*zx_alf1+v(i,2)*zx_alf2
262      ! AM        zxq = q(i,1)*zx_alf1+q(i,2)*zx_alf2
263      ! AM
264      ! AMAM           zxu = u10m(i)
265      ! AMAM           zxv = v10m(i)
266      ! AMAM           zxmod = 1.0+SQRT(zxu**2+zxv**2)
267      ! AM Niveau de ref choisi a 2m
268      zxt = t2m(i)
269
270      ! ***************************************************
271      ! attention, il doit s'agir de <w'theta'>
272      ! ;Calcul de tcls virtuel et de w'theta'virtuel
273      ! ;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;;
274      ! tcls=tcls*(1+.608*qcls)
275
276      ! ;Pour avoir w'theta',
277      ! ; il faut diviser par ro.Cp
278      ! Cp=Cpd*(1+0.84*qcls)
279      ! fcs=fcs/(ro_surf*Cp)
280      ! ;On transforme w'theta' en w'thetav'
281      ! Lv=(2.501-0.00237*(tcls-273.15))*1.E6
282      ! xle=xle/(ro_surf*Lv)
283      ! fcsv=fcs+.608*xle*tcls
284      ! ***************************************************
285      ! AM        khfs(i) = (tsol(i)*(1.+RETV*q(i,1))-zxt) *zxmod*cd_h(i)
286      ! AM        kqfs(i) = (zxqs-zxq) *zxmod*cd_h(i) * beta(i)
287      ! AM
288      ! dif khfs est deja w't'_v / heatv(i) = khfs(i) + RETV*zxt*kqfs(i)
289      ! AM calcule de Ro = paprs(i,1)/Rd zxt
290      ! AM convention >0 vers le bas ds lmdz
291      khfs(i) = -flux_t(i, 1) * zxt * rd / (rcpd * paprs(i, 1))
292      kqfs(i) = -flux_q(i, 1) * zxt * rd / (paprs(i, 1))
293      ! AM   verifier que khfs et kqfs sont bien de la forme w'l'
294      heatv(i) = khfs(i) + 0.608 * zxt * kqfs(i)
295      ! a comparer aussi aux sorties de clqh : flux_T/RoCp et flux_q/RoLv
296      ! AM        heatv(i) = khfs(i)
297      ! AM ustar est en entree
298      ! AM        taux = zxu *zxmod*cd_m(i)
299      ! AM        tauy = zxv *zxmod*cd_m(i)
300      ! AM        ustar(i) = SQRT(taux**2+tauy**2)
301      ! AM        ustar(i) = MAX(SQRT(ustar(i)),0.01)
302      ! Theta et qT du thermique sans exces (interpolin vers surf)
303      ! chgt de niveau du thermique (jeudi 30/12/1999)
304      ! (interpolation lineaire avant integration phi_h)
305      ! AM        qT_th(i) = zxqs*beta(i) + 4./z(i,1)*(q(i,1)-zxqs*beta(i))
306      ! AM        qT_th(i) = max(qT_th(i),q(i,1))
307      qt_th(i) = q2m(i)
308      ! n The_th restera la Theta du thermique sans exces jusqu'a 2eme calcul
309      ! n reste a regler convention P) pour Theta
310      ! The_th(i) = tsol(i) + 4./z(i,1)*(t(i,1)-tsol(i))
311      ! -                      + RLvCp*qT_th(i)
312      ! AM        Th_th(i) = tsol(i) + 4./z(i,1)*(t(i,1)-tsol(i))
313      th_th(i) = t2m(i)
314    END DO
315
316    DO i = 1, knon
317      rhino(i, 1) = 0.0 ! Global Richardson
318      check(i) = .TRUE.
319      pblh(i) = z(i, 1) ! on initialise pblh a l'altitude du 1er niveau
320      plcl(i) = 6000.
321      ! Lambda = -u*^3 / (alpha.g.kvon.<w'Theta'v>
322      unsobklen(i) = -rg * vk * heatv(i) / (t(i, 1) * max(ustar(i), usmin)**3)
323      trmb1(i) = 0.
324      trmb2(i) = 0.
325      trmb3(i) = 0.
326    END DO
327
328
329    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
330    ! PBL height calculation:
331    ! Search for level of pbl. Scan upward until the Richardson number between
332    ! the first level and the current level exceeds the "critical" value.
333    ! (bonne idee Nu de separer le Ric et l'exces de temp du thermique)
334    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
335    fac = 100.0
336    DO k = 2, isommet
337      DO i = 1, knon
338        IF (check(i)) THEN
339          ! pourquoi / niveau 1 (au lieu du sol) et le terme en u*^2 ?
340          ! test     zdu2 =
341          ! (u(i,k)-u(i,1))**2+(v(i,k)-v(i,1))**2+fac*ustar(i)**2
342          zdu2 = u(i, k)**2 + v(i, k)**2
343          zdu2 = max(zdu2, 1.0E-20)
344          ! Theta_v environnement
345          zthvd = t(i, k) / s(i, k) * (1. + retv * q(i, k))
346
347          ! therm Theta_v sans exces (avec hypothese fausse de H&B, sinon,
348          ! passer par Theta_e et virpot)
349          ! zthvu=t(i,1)/s(i,1)*(1.+RETV*q(i,1))
350          ! AM         zthvu = Th_th(i)*(1.+RETV*q(i,1))
351          zthvu = th_th(i) * (1. + retv * qt_th(i))
352          ! Le Ri par Theta_v
353          ! AM         rhino(i,k) = (z(i,k)-z(i,1))*RG*(zthvd-zthvu)
354          ! AM     .               /(zdu2*0.5*(zthvd+zthvu))
355          ! AM On a nveau de ref a 2m ???
356          rhino(i, k) = (z(i, k) - zref) * rg * (zthvd - zthvu) / (zdu2 * 0.5&
357                  * (zthvd + zthvu))
358
359          IF (rhino(i, k)>=ricr) THEN
360            pblh(i) = z(i, k - 1) + (z(i, k - 1) - z(i, k)) * (ricr - rhino(i, k - 1))&
361                    / (rhino(i, k - 1) - rhino(i, k))
362            ! test04
363            pblh(i) = pblh(i) + 100.
364            pblt(i) = t(i, k - 1) + (t(i, k) - t(i, k - 1)) * (pblh(i) - z(i, k - 1))&
365                    / (z(i, k) - z(i, k - 1))
366            check(i) = .FALSE.
367          END IF
368        END IF
369      END DO
370    END DO
371
372
373    ! Set pbl height to maximum value where computation exceeds number of
374    ! layers allowed
375
376    DO i = 1, knon
377      IF (check(i)) pblh(i) = z(i, isommet)
378    END DO
379
380    ! Improve estimate of pbl height for the unstable points.
381    ! Find unstable points (sensible heat flux is upward):
382
383    DO i = 1, knon
384      IF (heatv(i)>0.) THEN
385        unstbl(i) = .TRUE.
386        check(i) = .TRUE.
387      ELSE
388        unstbl(i) = .FALSE.
389        check(i) = .FALSE.
390      END IF
391    END DO
392
393    ! For the unstable case, compute velocity scale and the
394    ! convective temperature excess:
395
396    DO i = 1, knon
397      IF (check(i)) THEN
398        phiminv(i) = (1. - binm * pblh(i) * unsobklen(i))**onet
399        ! ***************************************************
400        ! Wm ? et W* ? c'est la formule pour z/h < .1
401        ! ;Calcul de w* ;;
402        ! ;;;;;;;;;;;;;;;;
403        ! w_star=((g/tcls)*fcsv*z(ind))^(1/3.) [ou prendre la premiere approx
404        ! de h)
405        ! ;; CALCUL DE wm ;;
406        ! ;;;;;;;;;;;;;;;;;;
407        ! ; Ici on considerera que l'on est dans la couche de surf jusqu'a
408        ! 100 m
409        ! ; On prend svt couche de surface=0.1*h mais on ne connait pas h
410        ! ;;;;;;;;;;;Dans la couche de surface
411        ! if (z(ind) le 20) then begin
412        ! Phim=(1.-15.*(z(ind)/L))^(-1/3.)
413        ! wm=u_star/Phim
414        ! ;;;;;;;;;;;En dehors de la couche de surface
415        ! END IF ELSE IF (z(ind) gt 20) then begin
416        ! wm=(u_star^3+c1*w_star^3)^(1/3.)
417        ! END IF
418        ! ***************************************************
419        wm(i) = ustar(i) * phiminv(i)
420        ! ===================================================================
421        ! valeurs de Dominique Lambert de la campagne SEMAPHORE :
422        ! <T'^2> = 100.T*^2; <q'^2> = 20.q*^2 a 10m
423        ! <Tv'^2> = (1+1.2q).100.T* + 1.2Tv.sqrt(20*100).T*.q* +
424        ! (.608*Tv)^2*20.q*^2;
425        ! et dTetavS = sqrt(<Tv'^2>) ainsi calculee.
426        ! avec : T*=<w'T'>_s/w* et q*=<w'q'>/w*
427        ! !!! on peut donc utiliser w* pour les fluctuations <-> Lambert
428        ! (leur corellation pourrait dependre de beta par ex)
429        ! if fcsv(i,j) gt 0 then begin
430        ! dTetavs=b1*(1.+2.*.608*q_10(i,j))*(fcs(i,j)/wm(i,j))^2+$
431        ! (.608*Thetav_10(i,j))^2*b2*(xle(i,j)/wm(i,j))^2+$
432        ! 2.*.608*thetav_10(i,j)*sqrt(b1*b2)*(xle(i,j)/wm(i,j))*(fcs(i,j)
433        ! /wm(i,j))
434        ! dqs=b2*(xle(i,j)/wm(i,j))^2
435        ! theta_s(i,j)=thetav_10(i,j)+sqrt(dTetavs)
436        ! q_s(i,j)=q_10(i,j)+sqrt(dqs)
437        ! END IF else begin
438        ! Theta_s(i,j)=thetav_10(i,j)
439        ! q_s(i,j)=q_10(i,j)
440        ! endelse
441        ! ===================================================================
442
443        ! HBTM        therm(i) = heatv(i)*fak/wm(i)
444        ! forme Mathieu :
445        q_star = kqfs(i) / wm(i)
446        t_star = khfs(i) / wm(i)
447        ! IM 091204 BEG
448        IF (1==0) THEN
449          IF (t_star<0. .OR. q_star<0.) THEN
450            PRINT *, 'i t_star q_star khfs kqfs wm', i, t_star, q_star, &
451                    khfs(i), kqfs(i), wm(i)
452          END IF
453        END IF
454        ! IM 091204 END
455        ! AM Nveau cde ref 2m =>
456        ! AM        therm(i) = sqrt( b1*(1.+2.*RETV*q(i,1))*t_star**2
457        ! AM     +             + (RETV*T(i,1))**2*b2*q_star**2
458        ! AM     +             + 2.*RETV*T(i,1)*b212*q_star*t_star
459        ! AM     +                 )
460        ! IM 091204 BEG
461        a1 = b1 * (1. + 2. * retv * qt_th(i)) * t_star**2
462        a2 = (retv * th_th(i))**2 * b2 * q_star * q_star
463        a3 = 2. * retv * th_th(i) * b212 * q_star * t_star
464        aa = a1 + a2 + a3
465        IF (1==0) THEN
466          IF (aa<0.) THEN
467            PRINT *, 'i a1 a2 a3 aa', i, a1, a2, a3, aa
468            PRINT *, 'i qT_th Th_th t_star q_star RETV b1 b2 b212', i, &
469                    qt_th(i), th_th(i), t_star, q_star, retv, b1, b2, b212
470          END IF
471        END IF
472        ! IM 091204 END
473        therm(i) = sqrt(b1 * (1. + 2. * retv * qt_th(i)) * t_star**2 + (retv * th_th(&
474                i))**2 * b2 * q_star * q_star &  ! IM 101204  +             +
475                ! 2.*RETV*Th_th(i)*b212*q_star*t_star
476                + max(0., 2. * retv * th_th(i) * b212 * q_star * t_star))
477
478        ! Theta et qT du thermique (forme H&B) avec exces
479        ! (attention, on ajoute therm(i) qui est virtuelle ...)
480        ! pourquoi pas sqrt(b1)*t_star ?
481        ! dqs = b2sr*kqfs(i)/wm(i)
482        qt_th(i) = qt_th(i) + b2sr * q_star
483        ! new on differre le calcul de Theta_e
484        ! The_th(i) = The_th(i) + therm(i) + RLvCp*qT_th(i)
485        ! ou:    The_th(i) = The_th(i) + sqrt(b1)*khfs(i)/wm(i) +
486        ! RLvCp*qT_th(i)
487        rhino(i, 1) = 0.0
488      END IF
489    END DO
490
491    ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
492    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
493    ! ++ Improve pblh estimate for unstable conditions using the +++++++
494    ! ++          convective temperature excess :                +++++++
495    ! ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
496    ! +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
497
498    DO k = 2, isommet
499      DO i = 1, knon
500        IF (check(i)) THEN
501          ! test     zdu2 =
502          ! (u(i,k)-u(i,1))**2+(v(i,k)-v(i,1))**2+fac*ustar(i)**2
503          zdu2 = u(i, k)**2 + v(i, k)**2
504          zdu2 = max(zdu2, 1.0E-20)
505          ! Theta_v environnement
506          zthvd = t(i, k) / s(i, k) * (1. + retv * q(i, k))
507
508          ! et therm Theta_v (avec hypothese de constance de H&B,
509          ! zthvu=(t(i,1)+therm(i))/s(i,1)*(1.+RETV*q(i,1))
510          zthvu = th_th(i) * (1. + retv * qt_th(i)) + therm(i)
511
512
513          ! Le Ri par Theta_v
514          ! AM Niveau de ref 2m
515          ! AM         rhino(i,k) = (z(i,k)-z(i,1))*RG*(zthvd-zthvu)
516          ! AM     .               /(zdu2*0.5*(zthvd+zthvu))
517          rhino(i, k) = (z(i, k) - zref) * rg * (zthvd - zthvu) / (zdu2 * 0.5&
518                  * (zthvd + zthvu))
519
520          IF (rhino(i, k)>=ricr) THEN
521            pblh(i) = z(i, k - 1) + (z(i, k - 1) - z(i, k)) * (ricr - rhino(i, k - 1))&
522                    / (rhino(i, k - 1) - rhino(i, k))
523            ! test04
524            pblh(i) = pblh(i) + 100.
525            pblt(i) = t(i, k - 1) + (t(i, k) - t(i, k - 1)) * (pblh(i) - z(i, k - 1))&
526                    / (z(i, k) - z(i, k - 1))
527            check(i) = .FALSE.
528            ! IM 170305 BEG
529            IF (1==0) THEN
530              ! debug print -120;34       -34-        58 et    0;26 wamp
531              IF (i==950 .OR. i==192 .OR. i==624 .OR. i==118) THEN
532                PRINT *, ' i,Th_th,Therm,qT :', i, th_th(i), therm(i), &
533                        qt_th(i)
534                q_star = kqfs(i) / wm(i)
535                t_star = khfs(i) / wm(i)
536                PRINT *, 'q* t*, b1,b2,b212 ', q_star, t_star, &
537                        b1 * (1. + 2. * retv * qt_th(i)) * t_star**2, &
538                        (retv * th_th(i))**2 * b2 * q_star**2, 2. * retv * th_th(i)&
539                        * b212 * q_star * t_star
540                PRINT *, 'zdu2 ,100.*ustar(i)**2', zdu2, fac * ustar(i)**2
541              END IF
542            END IF !(1.EQ.0) THEN
543            ! IM 170305 END
544            ! q_star = kqfs(i)/wm(i)
545            ! t_star = khfs(i)/wm(i)
546            ! trmb1(i) = b1*(1.+2.*RETV*q(i,1))*t_star**2
547            ! trmb2(i) = (RETV*T(i,1))**2*b2*q_star**2
548            ! Omega now   trmb3(i) = 2.*RETV*T(i,1)*b212*q_star*t_star
549          END IF
550        END IF
551      END DO
552    END DO
553
554    ! Set pbl height to maximum value where computation exceeds number of
555    ! layers allowed
556
557    DO i = 1, knon
558      IF (check(i)) pblh(i) = z(i, isommet)
559    END DO
560
561    ! PBL height must be greater than some minimum mechanical mixing depth
562    ! Several investigators have proposed minimum mechanical mixing depth
563    ! relationships as a function of the local friction velocity, u*.  We
564    ! make use of a linear relationship of the form h = c u* where c=700.
565    ! The scaling arguments that give rise to this relationship most often
566    ! represent the coefficient c as some constant over the local coriolis
567    ! parameter.  Here we make use of the experimental results of Koracin
568    ! and Berkowicz (1988) [BLM, Vol 43] for wich they recommend 0.07/f
569    ! where f was evaluated at 39.5 N and 52 N.  Thus we use a typical mid
570    ! latitude value for f so that c = 0.07/f = 700.
571
572    DO i = 1, knon
573      pblmin = 700.0 * ustar(i)
574      pblh(i) = max(pblh(i), pblmin)
575      ! par exemple :
576      pblt(i) = t(i, 2) + (t(i, 3) - t(i, 2)) * (pblh(i) - z(i, 2)) / (z(i, 3) - z(i, 2))
577    END DO
578
579    ! ********************************************************************
580    ! pblh is now available; do preparation for diffusivity calculation :
581    ! ********************************************************************
582    DO i = 1, knon
583      check(i) = .TRUE.
584      zsat(i) = .FALSE.
585      ! omegafl utilise pour prolongement CAPE
586      omegafl(i) = .FALSE.
587      cape(i) = 0.
588      kape(i) = 0.
589      eauliq(i) = 0.
590      ctei(i) = 0.
591      pblk(i) = 0.0
592      fak1(i) = ustar(i) * pblh(i) * vk
593
594      ! Do additional preparation for unstable cases only, set temperature
595      ! and moisture perturbations depending on stability.
596      ! *** Rq: les formule sont prises dans leur forme CS ***
597      IF (unstbl(i)) THEN
598        ! AM Niveau de ref du thermique
599        ! AM          zxt=(t(i,1)-z(i,1)*0.5*RG/RCPD/(1.+RVTMP2*q(i,1)))
600        ! AM     .         *(1.+RETV*q(i,1))
601        zxt = (th_th(i) - zref * 0.5 * rg / rcpd / (1. + rvtmp2 * qt_th(i))) * &
602                (1. + retv * qt_th(i))
603        phiminv(i) = (1. - binm * pblh(i) * unsobklen(i))**onet
604        phihinv(i) = sqrt(1. - binh * pblh(i) * unsobklen(i))
605        wm(i) = ustar(i) * phiminv(i)
606        fak2(i) = wm(i) * pblh(i) * vk
607        wstar(i) = (heatv(i) * rg * pblh(i) / zxt)**onet
608        fak3(i) = fakn * wstar(i) / wm(i)
609      ELSE
610        wstar(i) = 0.
611      END IF
612      ! Computes Theta_e for thermal (all cases : to be modified)
613      ! attention ajout therm(i) = virtuelle
614      the_th(i) = th_th(i) + therm(i) + rlvcp * qt_th(i)
615      ! ou:    The_th(i) = Th_th(i) + sqrt(b1)*khfs(i)/wm(i) + RLvCp*qT_th(i)
616    END DO
617
618    ! Main level loop to compute the diffusivities and
619    ! counter-gradient terms:
620
621    DO k = 2, isommet
622
623      ! Find levels within boundary layer:
624
625      DO i = 1, knon
626        unslev(i) = .FALSE.
627        stblev(i) = .FALSE.
628        zm(i) = z(i, k - 1)
629        zp(i) = z(i, k)
630        IF (zkmin==0.0 .AND. zp(i)>pblh(i)) zp(i) = pblh(i)
631        IF (zm(i)<pblh(i)) THEN
632          zmzp = 0.5 * (zm(i) + zp(i))
633          ! debug
634          ! if (i.EQ.1864) THEN
635          ! PRINT*,'i,pblh(1864),obklen(1864)',i,pblh(i),obklen(i)
636          ! END IF
637
638          zh(i) = zmzp / pblh(i)
639          zl(i) = zmzp * unsobklen(i)
640          zzh(i) = 0.
641          IF (zh(i)<=1.0) zzh(i) = (1. - zh(i))**2
642
643          ! stblev for points zm < plbh and stable and neutral
644          ! unslev for points zm < plbh and unstable
645
646          IF (unstbl(i)) THEN
647            unslev(i) = .TRUE.
648          ELSE
649            stblev(i) = .TRUE.
650          END IF
651        END IF
652      END DO
653      ! PRINT*,'fin calcul niveaux'
654
655      ! Stable and neutral points; set diffusivities; counter-gradient
656      ! terms zero for stable case:
657
658      DO i = 1, knon
659        IF (stblev(i)) THEN
660          IF (zl(i)<=1.) THEN
661            pblk(i) = fak1(i) * zh(i) * zzh(i) / (1. + betas * zl(i))
662          ELSE
663            pblk(i) = fak1(i) * zh(i) * zzh(i) / (betas + zl(i))
664          END IF
665          ! pcfm(i,k) = pblk(i)
666          ! pcfh(i,k) = pcfm(i,k)
667        END IF
668      END DO
669
670      ! unssrf, unstable within surface layer of pbl
671      ! unsout, unstable within outer   layer of pbl
672
673      DO i = 1, knon
674        unssrf(i) = .FALSE.
675        unsout(i) = .FALSE.
676        IF (unslev(i)) THEN
677          IF (zh(i)<sffrac) THEN
678            unssrf(i) = .TRUE.
679          ELSE
680            unsout(i) = .TRUE.
681          END IF
682        END IF
683      END DO
684
685      ! Unstable for surface layer; counter-gradient terms zero
686
687      DO i = 1, knon
688        IF (unssrf(i)) THEN
689          term = (1. - betam * zl(i))**onet
690          pblk(i) = fak1(i) * zh(i) * zzh(i) * term
691          pr(i) = term / sqrt(1. - betah * zl(i))
692        END IF
693      END DO
694      ! PRINT*,'fin counter-gradient terms zero'
695
696      ! Unstable for outer layer; counter-gradient terms non-zero:
697
698      DO i = 1, knon
699        IF (unsout(i)) THEN
700          pblk(i) = fak2(i) * zh(i) * zzh(i)
701          ! cgs(i,k) = fak3(i)/(pblh(i)*wm(i))
702          ! cgh(i,k) = khfs(i)*cgs(i,k)
703          pr(i) = phiminv(i) / phihinv(i) + ccon * fak3(i) / fak
704          ! cgq(i,k) = kqfs(i)*cgs(i,k)
705        END IF
706      END DO
707      ! PRINT*,'fin counter-gradient terms non zero'
708
709      ! For all unstable layers, compute diffusivities and ctrgrad ter m
710
711      ! DO i = 1, knon
712      ! IF (unslev(i)) THEN
713      ! pcfm(i,k) = pblk(i)
714      ! pcfh(i,k) = pblk(i)/pr(i)
715      ! etc cf original
716      ! ENDIF
717      ! ENDDO
718
719      ! For all layers, compute integral info and CTEI
720
721      DO i = 1, knon
722        IF (check(i) .OR. omegafl(i)) THEN
723          IF (.NOT. zsat(i)) THEN
724            ! Th2 = The_th(i) - RLvCp*qT_th(i)
725            th2 = th_th(i)
726            t2 = th2 * s(i, k)
727            ! thermodyn functions
728            zdelta = max(0., sign(1., rtt - t2))
729            qqsat = r2es * foeew(t2, zdelta) / pplay(i, k)
730            qqsat = min(0.5, qqsat)
731            zcor = 1. / (1. - retv * qqsat)
732            qqsat = qqsat * zcor
733
734            IF (qqsat<qt_th(i)) THEN
735              ! on calcule lcl
736              IF (k==2) THEN
737                plcl(i) = z(i, k)
738              ELSE
739                plcl(i) = z(i, k - 1) + (z(i, k - 1) - z(i, k)) * (qt_th(i)&
740                        - qsatbef(i)) / (qsatbef(i) - qqsat)
741              END IF
742              zsat(i) = .TRUE.
743              tbef(i) = t2
744            END IF
745
746            qsatbef(i) = qqsat ! bug dans la version orig ???
747          END IF
748          ! amn ???? cette ligne a deja ete faite normalement ?
749        END IF
750        ! PRINT*,'hbtm2 i,k=',i,k
751      END DO
752    END DO ! end of level loop
753    ! IM 170305 BEG
754    IF (1==0) THEN
755      PRINT *, 'hbtm2  ok'
756    END IF !(1.EQ.0) THEN
757    ! IM 170305 END
758
759  END SUBROUTINE hbtm
760
761END MODULE hbtm_mod
Note: See TracBrowser for help on using the repository browser.