New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
sbcblk_algo_ecmwf.F90 in NEMO/branches/2019/dev_r11085_ASINTER-05_Brodeau_Advanced_Bulk/src/OCE/SBC – NEMO

source: NEMO/branches/2019/dev_r11085_ASINTER-05_Brodeau_Advanced_Bulk/src/OCE/SBC/sbcblk_algo_ecmwf.F90 @ 11217

Last change on this file since 11217 was 11215, checked in by laurent, 5 years ago

LB: inclusion of "cool-skin/warm-layer" support (a la ECMWF) into COARE 3.0 algorithm.

  • Property svn:keywords set to Id
File size: 21.4 KB
Line 
1MODULE sbcblk_algo_ecmwf
2   !!======================================================================
3   !!                   ***  MODULE  sbcblk_algo_ecmwf  ***
4   !! Computes:
5   !!   * bulk transfer coefficients C_D, C_E and C_H
6   !!   * air temp. and spec. hum. adjusted from zt (2m) to zu (10m) if needed
7   !!   * the effective bulk wind speed at 10m U_blk
8   !!   => all these are used in bulk formulas in sbcblk.F90
9   !!
10   !!    Using the bulk formulation/param. of IFS of ECMWF (cycle 40r1)
11   !!         based on IFS doc (avaible online on the ECMWF's website)
12   !!
13   !!       Routine turb_ecmwf maintained and developed in AeroBulk
14   !!                     (https://github.com/brodeau/aerobulk)
15   !!
16   !! ** Author: L. Brodeau, June 2019 / AeroBulk (https://github.com/brodeau/aerobulk)
17   !!----------------------------------------------------------------------
18   !! History :  4.0  !  2016-02  (L.Brodeau)   Original code
19   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   turb_ecmwf  : computes the bulk turbulent transfer coefficients
23   !!                   adjusts t_air and q_air from zt to zu m
24   !!                   returns the effective bulk wind speed at 10m
25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and tracers
27   USE dom_oce         ! ocean space and time domain
28   USE phycst          ! physical constants
29   USE iom             ! I/O manager library
30   USE lib_mpp         ! distribued memory computing library
31   USE in_out_manager  ! I/O manager
32   USE prtctl          ! Print control
33   USE sbcwave, ONLY   :  cdn_wave ! wave module
34#if defined key_si3 || defined key_cice
35   USE sbc_ice         ! Surface boundary condition: ice fields
36#endif
37   USE lib_fortran     ! to use key_nosignedzero
38
39   USE sbc_oce         ! Surface boundary condition: ocean fields
40   USE sbcblk_phy      ! all thermodynamics functions, rho_air, q_sat, etc... !LB
41   USE sbcblk_skin     ! cool-skin/warm layer scheme (CSWL_ECMWF) !LB
42
43   IMPLICIT NONE
44   PRIVATE
45
46   PUBLIC ::   TURB_ECMWF   ! called by sbcblk.F90
47
48   !                   !! ECMWF own values for given constants, taken form IFS documentation...
49   REAL(wp), PARAMETER ::   charn0 = 0.018    ! Charnock constant (pretty high value here !!!
50   !                                          !    =>  Usually 0.011 for moderate winds)
51   REAL(wp), PARAMETER ::   zi0     = 1000.   ! scale height of the atmospheric boundary layer...1
52   REAL(wp), PARAMETER ::   Beta0    = 1.     ! gustiness parameter ( = 1.25 in COAREv3)
53   REAL(wp), PARAMETER ::   alpha_M = 0.11    ! For roughness length (smooth surface term)
54   REAL(wp), PARAMETER ::   alpha_H = 0.40    ! (Chapter 3, p.34, IFS doc Cy31r1)
55   REAL(wp), PARAMETER ::   alpha_Q = 0.62    !
56
57   INTEGER , PARAMETER ::   nb_itt = 5        ! number of itterations
58
59   !!----------------------------------------------------------------------
60CONTAINS
61
62   SUBROUTINE TURB_ECMWF( zt, zu, T_s, t_zt, q_s, q_zt, U_zu, &
63      &                   Cd, Ch, Ce, t_zu, q_zu, U_blk,      &
64      &                   Cdn, Chn, Cen,                      &
65      &                   Qsw, rad_lw, slp                   )
66      !!----------------------------------------------------------------------------------
67      !!                      ***  ROUTINE  turb_ecmwf  ***
68      !!
69      !! ** Purpose :   Computes turbulent transfert coefficients of surface
70      !!                fluxes according to IFS doc. (cycle 40)
71      !!                If relevant (zt /= zu), adjust temperature and humidity from height zt to zu
72      !!                Returns the effective bulk wind speed at 10m to be used in the bulk formulas
73      !!
74      !!                Applies the cool-skin warm-layer correction of the SST to T_s
75      !!                if the net shortwave flux at the surface (Qsw), the downwelling longwave
76      !!                radiative fluxes at the surface (rad_lw), and the sea-leve pressure (slp)
77      !!                are provided as (optional) arguments!
78      !!
79      !! INPUT :
80      !! -------
81      !!    *  zt   : height for temperature and spec. hum. of air            [m]
82      !!    *  zu   : height for wind speed (generally 10m)                   [m]
83      !!    *  U_zu : scalar wind speed at 10m                                [m/s]
84      !!    *  t_zt : potential air temperature at zt                         [K]
85      !!    *  q_zt : specific humidity of air at zt                          [kg/kg]
86      !!
87      !! INPUT/OUTPUT:
88      !! -------------
89      !!    *  T_s  : SST or skin temperature                                 [K]
90      !!    *  q_s  : SSQ aka saturation specific humidity at temp. T_s       [kg/kg]
91      !!              -> doesn't need to be given a value if skin temp computed (in case l_use_skin=True)
92      !!              -> MUST be given the correct value if not computing skint temp. (in case l_use_skin=False)
93      !!
94      !! OPTIONAL INPUT (will trigger l_use_skin=TRUE if present!):
95      !! ---------------
96      !!    *  Qsw    : net solar flux (after albedo) at the surface (>0)     [W/m^2]
97      !!    *  rad_lw : downwelling longwave radiation at the surface  (>0)   [W/m^2]
98      !!    *  slp    : sea-level pressure                                    [Pa]
99      !!
100      !! OUTPUT :
101      !! --------
102      !!    *  Cd     : drag coefficient
103      !!    *  Ch     : sensible heat coefficient
104      !!    *  Ce     : evaporation coefficient
105      !!    *  t_zu   : pot. air temperature adjusted at wind height zu       [K]
106      !!    *  q_zu   : specific humidity of air        //                    [kg/kg]
107      !!    *  U_blk  : bulk wind speed at 10m                                [m/s]
108      !!
109      !!
110      !! ** Author: L. Brodeau, June 2019 / AeroBulk (https://github.com/brodeau/aerobulk/)
111      !!----------------------------------------------------------------------------------
112      REAL(wp), INTENT(in   )                     ::   zt       ! height for t_zt and q_zt                    [m]
113      REAL(wp), INTENT(in   )                     ::   zu       ! height for U_zu                             [m]
114      REAL(wp), INTENT(inout), DIMENSION(jpi,jpj) ::   T_s      ! sea surface temperature                [Kelvin]
115      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   t_zt     ! potential air temperature              [Kelvin]
116      REAL(wp), INTENT(inout), DIMENSION(jpi,jpj) ::   q_s      ! sea surface specific humidity           [kg/kg]
117      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity at zt             [kg/kg]
118      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   U_zu     ! relative wind module at zu                [m/s]
119      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
120      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
121      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
122      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   t_zu     ! pot. air temp. adjusted at zu               [K]
123      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. humidity adjusted at zu           [kg/kg]
124      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   U_blk    ! bulk wind at 10m                          [m/s]
125      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cdn, Chn, Cen ! neutral transfer coefficients
126      !
127      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   Qsw      !             [W/m^2]
128      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   rad_lw   !             [W/m^2]
129      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   slp      !             [Pa]
130      !
131      INTEGER :: j_itt
132      LOGICAL :: l_zt_equal_zu = .FALSE.      ! if q and t are given at same height as U
133      !
134      REAL(wp), DIMENSION(jpi,jpj) ::  &
135         &  u_star, t_star, q_star, &
136         &  dt_zu, dq_zu,    &
137         &  znu_a,           & !: Nu_air, Viscosity of air
138         &  Linv,            & !: 1/L (inverse of Monin Obukhov length...
139         &  z0, z0t, z0q
140      !
141      ! Cool skin:
142      LOGICAL :: l_use_skin = .FALSE.
143      REAL(wp), DIMENSION(jpi,jpj) :: zsst    ! to back up the initial bulk SST
144      !
145      REAL(wp), DIMENSION(jpi,jpj) ::   func_m, func_h
146      REAL(wp), DIMENSION(jpi,jpj) ::   ztmp0, ztmp1, ztmp2
147      !!----------------------------------------------------------------------------------
148
149      ! Cool skin ?
150      IF( PRESENT(Qsw) .AND. PRESENT(rad_lw) .AND. PRESENT(slp) ) THEN
151         l_use_skin = .TRUE.
152      END IF
153      IF (lwp) PRINT *, ' *** LOLO(sbcblk_algo_ecmwf.F90) => l_use_skin =', l_use_skin
154
155      ! Identical first gess as in COARE, with IFS parameter values though...
156      !
157      l_zt_equal_zu = .FALSE.
158      IF( ABS(zu - zt) < 0.01 )   l_zt_equal_zu = .TRUE.    ! testing "zu == zt" is risky with double precision
159
160      !! Initialization for cool skin:
161      IF( l_use_skin ) THEN
162         zsst   = T_s    ! save the bulk SST
163         T_s    = T_s - 0.25                      ! First guess of correction
164         q_s    = rdct_qsat_salt*q_sat(MAX(T_s, 200._wp), slp) ! First guess of q_s
165      END IF
166
167      !! First guess of temperature and humidity at height zu:
168      t_zu = MAX( t_zt ,  180._wp )   ! who knows what's given on masked-continental regions...
169      q_zu = MAX( q_zt , 1.e-6_wp )   !               "
170
171      !! Pot. temp. difference (and we don't want it to be 0!)
172      dt_zu = t_zu - T_s ;   dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
173      dq_zu = q_zu - q_s ;   dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
174
175      znu_a = visc_air(t_zu) ! Air viscosity (m^2/s) at zt given from temperature in (K)
176
177      ztmp2 = 0.5_wp*0.5_wp  ! initial guess for wind gustiness contribution
178      U_blk = SQRT(U_zu*U_zu + ztmp2)
179
180      ztmp2   = 10000._wp     ! optimization: ztmp2 == 1/z0 (with z0 first guess == 0.0001)
181      ztmp0   = LOG(zu*ztmp2)
182      ztmp1   = LOG(10.*ztmp2)
183      u_star = 0.035_wp*U_blk*ztmp1/ztmp0       ! (u* = 0.035*Un10)
184
185      z0     = charn0*u_star*u_star/grav + 0.11_wp*znu_a/u_star
186      z0     = MIN(ABS(z0), 0.001_wp)  ! (prevent FPE from stupid values from masked region later on...) !#LOLO
187      z0t    = 1._wp / ( 0.1_wp*EXP(vkarmn/(0.00115/(vkarmn/ztmp1))) )
188      z0t    = MIN(ABS(z0t), 0.001_wp)  ! (prevent FPE from stupid values from masked region later on...) !#LOLO
189
190      ztmp2  = vkarmn/ztmp0
191      Cd     = ztmp2*ztmp2    ! first guess of Cd
192
193      ztmp0 = vkarmn*vkarmn/LOG(zt/z0t)/Cd
194
195      ztmp2 = Ri_bulk( zu, T_s, t_zu, q_s, q_zu, U_blk ) ! Bulk Richardson Number (BRN)
196
197      !! First estimate of zeta_u, depending on the stability, ie sign of BRN (ztmp2):
198      ztmp1 = 0.5 + SIGN( 0.5_wp , ztmp2 )
199      func_m = ztmp0*ztmp2 ! temporary array !!
200      func_h = (1._wp-ztmp1) * (func_m/(1._wp+ztmp2/(-zu/(zi0*0.004_wp*Beta0**3)))) & !  BRN < 0 ! temporary array !!! func_h == zeta_u
201         &  +     ztmp1   * (func_m*(1._wp + 27._wp/9._wp*ztmp2/func_m))              !  BRN > 0
202      !#LB: should make sure that the "func_m" of "27./9.*ztmp2/func_m" is "ztmp0*ztmp2" and not "ztmp0==vkarmn*vkarmn/LOG(zt/z0t)/Cd" !
203
204      !! First guess M-O stability dependent scaling params.(u*,t*,q*) to estimate z0 and z/L
205      ztmp0  = vkarmn/(LOG(zu/z0t) - psi_h_ecmwf(func_h))
206
207      u_star = U_blk*vkarmn/(LOG(zu) - LOG(z0)  - psi_m_ecmwf(func_h))
208      t_star = dt_zu*ztmp0
209      q_star = dq_zu*ztmp0
210
211      ! What's need to be done if zt /= zu:
212      IF( .NOT. l_zt_equal_zu ) THEN
213         !! First update of values at zu (or zt for wind)
214         ztmp0 = psi_h_ecmwf(func_h) - psi_h_ecmwf(zt*func_h/zu)    ! zt*func_h/zu == zeta_t
215         ztmp1 = LOG(zt/zu) + ztmp0
216         t_zu = t_zt - t_star/vkarmn*ztmp1
217         q_zu = q_zt - q_star/vkarmn*ztmp1
218         q_zu = (0.5_wp + SIGN(0.5_wp,q_zu))*q_zu !Makes it impossible to have negative humidity :
219         !
220         dt_zu = t_zu - T_s  ; dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
221         dq_zu = q_zu - q_s  ; dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
222      END IF
223
224
225      !! => that was same first guess as in COARE...
226
227
228      !! First guess of inverse of Monin-Obukov length (1/L) :
229      Linv = One_on_L( t_zu, q_zu, u_star, t_star, q_star )
230
231      !! Functions such as  u* = U_blk*vkarmn/func_m
232      ztmp0 = zu*Linv
233      func_m = LOG(zu) - LOG(z0)  - psi_m_ecmwf(ztmp0) + psi_m_ecmwf( z0*Linv)
234      func_h = LOG(zu) - LOG(z0t) - psi_h_ecmwf(ztmp0) + psi_h_ecmwf(z0t*Linv)
235
236      !! ITERATION BLOCK
237      DO j_itt = 1, nb_itt
238
239         !! Bulk Richardson Number at z=zu (Eq. 3.25)
240         ztmp0 = Ri_bulk( zu, T_s, t_zu, q_s, q_zu, U_blk ) ! Bulk Richardson Number (BRN)
241
242         !! New estimate of the inverse of the Monin-Obukhon length (Linv == zeta/zu) :
243         Linv = ztmp0*func_m*func_m/func_h / zu     ! From Eq. 3.23, Chap.3.2.3, IFS doc - Cy40r1
244         !! Note: it is slightly different that the L we would get with the usual
245         Linv = SIGN( MIN(ABS(Linv),200._wp), Linv ) ! (prevent FPE from stupid values from masked region later on...) !#LOLO
246
247         !! Update func_m with new Linv:
248         func_m = LOG(zu) -LOG(z0) - psi_m_ecmwf(zu*Linv) + psi_m_ecmwf(z0*Linv) ! LB: should be "zu+z0" rather than "zu" alone, but z0 is tiny wrt zu!
249
250         !! Need to update roughness lengthes:
251         u_star = U_blk*vkarmn/func_m
252         ztmp2  = u_star*u_star
253         ztmp1  = znu_a/u_star
254         z0     = MIN( ABS( alpha_M*ztmp1 + charn0*ztmp2/grav ) , 0.001_wp)
255         z0t    = MIN( ABS( alpha_H*ztmp1                     ) , 0.001_wp)   ! eq.3.26, Chap.3, p.34, IFS doc - Cy31r1
256         z0q    = MIN( ABS( alpha_Q*ztmp1                     ) , 0.001_wp)
257
258         !! Update wind at 10m taking into acount convection-related wind gustiness:
259         !! => Chap. 3.2, IFS doc - Cy40r1, Eq.3.17 and Eq.3.18 + Eq.3.8
260         ! Only true when unstable (L<0) => when ztmp0 < 0 => - !!!
261         ztmp2 = ztmp2 * ( MAX(-zi0*Linv/vkarmn , 0._wp))**(2._wp/3._wp) ! => w*^2  (combining Eq. 3.8 and 3.18, hap.3, IFS doc - Cy31r1)
262         !! => equivalent using Beta=1 (gustiness parameter, 1.25 for COARE, also zi0=600 in COARE..)
263         U_blk = MAX( SQRT(U_zu*U_zu + ztmp2) , 0.2_wp )              ! eq.3.17, Chap.3, p.32, IFS doc - Cy31r1
264         ! => 0.2 prevents U_blk to be 0 in stable case when U_zu=0.
265
266
267         !! Need to update "theta" and "q" at zu in case they are given at different heights
268         !! as well the air-sea differences:
269         IF( .NOT. l_zt_equal_zu ) THEN
270            !! Arrays func_m and func_h are free for a while so using them as temporary arrays...
271            func_h = psi_h_ecmwf(zu*Linv) ! temporary array !!!
272            func_m = psi_h_ecmwf(zt*Linv) ! temporary array !!!
273
274            ztmp2  = psi_h_ecmwf(z0t*Linv)
275            ztmp0  = func_h - ztmp2
276            ztmp1  = vkarmn/(LOG(zu) - LOG(z0t) - ztmp0)
277            t_star = dt_zu*ztmp1
278            ztmp2  = ztmp0 - func_m + ztmp2
279            ztmp1  = LOG(zt/zu) + ztmp2
280            t_zu   = t_zt - t_star/vkarmn*ztmp1
281
282            ztmp2  = psi_h_ecmwf(z0q*Linv)
283            ztmp0  = func_h - ztmp2
284            ztmp1  = vkarmn/(LOG(zu) - LOG(z0q) - ztmp0)
285            q_star = dq_zu*ztmp1
286            ztmp2  = ztmp0 - func_m + ztmp2
287            ztmp1  = LOG(zt/zu) + ztmp2
288            q_zu   = q_zt - q_star/vkarmn*ztmp1
289         END IF
290
291         !! Updating because of updated z0 and z0t and new Linv...
292         ztmp0 = zu*Linv
293         func_m = log(zu) - LOG(z0 ) - psi_m_ecmwf(ztmp0) + psi_m_ecmwf(z0 *Linv)
294         func_h = log(zu) - LOG(z0t) - psi_h_ecmwf(ztmp0) + psi_h_ecmwf(z0t*Linv)
295
296         !! SKIN related part
297         !! -----------------
298         IF( l_use_skin ) THEN
299            !! compute transfer coefficients at zu : lolo: verifier...
300            Ch = vkarmn*vkarmn/(func_m*func_h)
301            ztmp1 = LOG(zu) - LOG(z0q) - psi_h_ecmwf(ztmp0) + psi_h_ecmwf(z0q*Linv)   ! func_q
302            Ce = vkarmn*vkarmn/(func_m*ztmp1)
303            ! Non-Solar heat flux to the ocean:
304            ztmp1 = U_blk*MAX(rho_air(t_zu, q_zu, slp), 1._wp)     ! rho*U10
305            ztmp2 = T_s*T_s
306            ztmp1 = ztmp1 * ( Ce*rLevap*(q_zu - q_s) + Ch*rCp_dry*(t_zu - T_s) ) & ! Total turb. heat flux
307               &       +    (rad_lw - emiss_w*stefan*ztmp2*ztmp2)                  ! Net longwave flux
308            !! Updating the values of the skin temperature T_s and q_s :
309            CALL CSWL_ECMWF( Qsw, ztmp1, u_star, zsst, T_s )
310            q_s = rdct_qsat_salt*q_sat(MAX(T_s, 200._wp), slp)  ! 200 -> just to avoid numerics problem on masked regions if silly values are given
311         END IF
312
313         IF( (l_use_skin).OR.(.NOT. l_zt_equal_zu) ) THEN
314            dt_zu = t_zu - T_s ;  dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
315            dq_zu = q_zu - q_s ;  dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
316         END IF
317
318      END DO !DO j_itt = 1, nb_itt
319
320      Cd = vkarmn*vkarmn/(func_m*func_m)
321      Ch = vkarmn*vkarmn/(func_m*func_h)
322      ztmp2 = log(zu/z0q) - psi_h_ecmwf(zu*Linv) + psi_h_ecmwf(z0q*Linv)   ! func_q
323      Ce = vkarmn*vkarmn/(func_m*ztmp2)
324
325      Cdn = vkarmn*vkarmn / (log(zu/z0 )*log(zu/z0 ))
326      Chn = vkarmn*vkarmn / (log(zu/z0t)*log(zu/z0t))
327      Cen = vkarmn*vkarmn / (log(zu/z0q)*log(zu/z0q))
328
329   END SUBROUTINE TURB_ECMWF
330
331
332   FUNCTION psi_m_ecmwf( pzeta )
333      !!----------------------------------------------------------------------------------
334      !! Universal profile stability function for momentum
335      !!     ECMWF / as in IFS cy31r1 documentation, available online
336      !!     at ecmwf.int
337      !!
338      !! pzeta : stability paramenter, z/L where z is altitude measurement
339      !!         and L is M-O length
340      !!
341      !! ** Author: L. Brodeau, June 2016 / AeroBulk (https://github.com/brodeau/aerobulk/)
342      !!----------------------------------------------------------------------------------
343      REAL(wp), DIMENSION(jpi,jpj) :: psi_m_ecmwf
344      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
345      !
346      INTEGER  ::   ji, jj    ! dummy loop indices
347      REAL(wp) :: zzeta, zx, ztmp, psi_unst, psi_stab, stab
348      !!----------------------------------------------------------------------------------
349      !
350      DO jj = 1, jpj
351         DO ji = 1, jpi
352            !
353            zzeta = MIN( pzeta(ji,jj) , 5._wp ) !! Very stable conditions (L positif and big!):
354            !
355            ! Unstable (Paulson 1970):
356            !   eq.3.20, Chap.3, p.33, IFS doc - Cy31r1
357            zx = SQRT(ABS(1._wp - 16._wp*zzeta))
358            ztmp = 1._wp + SQRT(zx)
359            ztmp = ztmp*ztmp
360            psi_unst = LOG( 0.125_wp*ztmp*(1._wp + zx) )   &
361               &       -2._wp*ATAN( SQRT(zx) ) + 0.5_wp*rpi
362            !
363            ! Unstable:
364            ! eq.3.22, Chap.3, p.33, IFS doc - Cy31r1
365            psi_stab = -2._wp/3._wp*(zzeta - 5._wp/0.35_wp)*EXP(-0.35_wp*zzeta) &
366               &       - zzeta - 2._wp/3._wp*5._wp/0.35_wp
367            !
368            ! Combining:
369            stab = 0.5_wp + SIGN(0.5_wp, zzeta) ! zzeta > 0 => stab = 1
370            !
371            psi_m_ecmwf(ji,jj) = (1._wp - stab) * psi_unst & ! (zzeta < 0) Unstable
372               &                +      stab  * psi_stab      ! (zzeta > 0) Stable
373            !
374         END DO
375      END DO
376      !
377   END FUNCTION psi_m_ecmwf
378
379
380   FUNCTION psi_h_ecmwf( pzeta )
381      !!----------------------------------------------------------------------------------
382      !! Universal profile stability function for temperature and humidity
383      !!     ECMWF / as in IFS cy31r1 documentation, available online
384      !!     at ecmwf.int
385      !!
386      !! pzeta : stability paramenter, z/L where z is altitude measurement
387      !!         and L is M-O length
388      !!
389      !! ** Author: L. Brodeau, June 2016 / AeroBulk (https://github.com/brodeau/aerobulk/)
390      !!----------------------------------------------------------------------------------
391      REAL(wp), DIMENSION(jpi,jpj) :: psi_h_ecmwf
392      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
393      !
394      INTEGER  ::   ji, jj     ! dummy loop indices
395      REAL(wp) ::  zzeta, zx, psi_unst, psi_stab, stab
396      !!----------------------------------------------------------------------------------
397      !
398      DO jj = 1, jpj
399         DO ji = 1, jpi
400            !
401            zzeta = MIN(pzeta(ji,jj) , 5._wp)   ! Very stable conditions (L positif and big!):
402            !
403            zx  = ABS(1._wp - 16._wp*zzeta)**.25        ! this is actually (1/phi_m)**2  !!!
404            !                                     ! eq.3.19, Chap.3, p.33, IFS doc - Cy31r1
405            ! Unstable (Paulson 1970) :
406            psi_unst = 2._wp*LOG(0.5_wp*(1._wp + zx*zx))   ! eq.3.20, Chap.3, p.33, IFS doc - Cy31r1
407            !
408            ! Stable:
409            psi_stab = -2._wp/3._wp*(zzeta - 5._wp/0.35_wp)*EXP(-0.35_wp*zzeta) & ! eq.3.22, Chap.3, p.33, IFS doc - Cy31r1
410               &       - ABS(1._wp + 2._wp/3._wp*zzeta)**1.5_wp - 2._wp/3._wp*5._wp/0.35_wp + 1._wp
411            ! LB: added ABS() to avoid NaN values when unstable, which contaminates the unstable solution...
412            !
413            stab = 0.5_wp + SIGN(0.5_wp, zzeta) ! zzeta > 0 => stab = 1
414            !
415            !
416            psi_h_ecmwf(ji,jj) = (1._wp - stab) * psi_unst &   ! (zzeta < 0) Unstable
417               &                +    stab    * psi_stab        ! (zzeta > 0) Stable
418            !
419         END DO
420      END DO
421      !
422   END FUNCTION psi_h_ecmwf
423
424
425
426
427
428   !!======================================================================
429END MODULE sbcblk_algo_ecmwf
Note: See TracBrowser for help on using the repository browser.