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 @ 11209

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

LB: making more efficient and more centralized use of physical/thermodynamics functions defined into "OCE/SBC/sbcblk_phy.F90".

  • Property svn:keywords set to Id
File size: 22.6 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 2016 / 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      !LB: all thermodynamics functions, rho_air, q_sat, etc... #LB
41   USE sbcblk_skin     ! cool-skin/warm layer scheme (CSWL_ECMWF)
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   !! Cool-Skin / Warm-Layer related parameters:
60   REAL(wp),    PARAMETER :: &
61      &  rdt0    = 3600.*1.5, & !: time step
62      &  rd0     = 3. ,      &  !: Depth scale [m], "d" in Eq.11 (Zeng & Beljaars 2005)
63      &  rNu0    = 0.5          !: Nu (exponent of temperature profile) Eq.11
64   !                            !: (Zeng & Beljaars 2005) !: set to 0.5 instead of
65   !                            !: 0.3 to respect a warming of +3 K in calm
66   !                            !: condition for the insolation peak of +1000W/m^2
67   INTEGER,    PARAMETER :: &
68      &  nb_itt_wl = 10         !: number of sub-itterations for solving the differential equation in warm-layer part
69   !                            !:  => use "nb_itt_wl = 1" for No itteration! => way cheaper !!!
70   !                            !:    => assumes balance between the last 2 terms of Eq.11 (lhs of eq.11 = 0)
71   !                            !:    => in that case no need for sub-itterations !
72   !                            !:    => ACCEPTABLE IN MOST CONDITIONS ! (UNLESS: sunny + very calm/low-wind conditions)
73   !                            !:  => Otherwize use "nb_itt_wl = 10"
74   !!----------------------------------------------------------------------
75CONTAINS
76
77   SUBROUTINE TURB_ECMWF( zt, zu, T_s, t_zt, q_s, q_zt, U_zu, &
78      &                   Cd, Ch, Ce, t_zu, q_zu, U_blk,      &
79      &                   Cdn, Chn, Cen,                      &
80      &                   Qsw, rad_lw, slp                   )
81      !!----------------------------------------------------------------------------------
82      !!                      ***  ROUTINE  turb_ecmwf  ***
83      !!
84      !! ** Purpose :   Computes turbulent transfert coefficients of surface
85      !!                fluxes according to IFS doc. (cycle 40)
86      !!                If relevant (zt /= zu), adjust temperature and humidity from height zt to zu
87      !!                Returns the effective bulk wind speed at 10m to be used in the bulk formulas
88      !!
89      !!                Applies the cool-skin warm-layer correction of the SST to T_s
90      !!                if the net shortwave flux at the surface (Qsw), the downwelling longwave
91      !!                radiative fluxes at the surface (rad_lw), and the sea-leve pressure (slp)
92      !!                are provided as (optional) arguments!
93      !!
94      !! INPUT :
95      !! -------
96      !!    *  zt   : height for temperature and spec. hum. of air            [m]
97      !!    *  zu   : height for wind speed (generally 10m)                   [m]
98      !!    *  U_zu : scalar wind speed at 10m                                [m/s]
99      !!    *  t_zt : potential air temperature at zt                         [K]
100      !!    *  q_zt : specific humidity of air at zt                          [kg/kg]
101      !!
102      !! INPUT/OUTPUT:
103      !! -------------
104      !!    *  T_s  : SST or skin temperature                                 [K]
105      !!    *  q_s  : SSQ aka saturation specific humidity at temp. T_s       [kg/kg]
106      !!              -> doesn't need to be given a value if skin temp computed (in case l_use_skin=True)
107      !!              -> MUST be given the correct value if not computing skint temp. (in case l_use_skin=False)
108      !!
109      !! OPTIONAL INPUT (will trigger l_use_skin=TRUE if present!):
110      !! ---------------
111      !!    *  Qsw    : net solar flux (after albedo) at the surface (>0)     [W/m^2]
112      !!    *  rad_lw : downwelling longwave radiation at the surface  (>0)   [W/m^2]
113      !!    *  slp    : sea-level pressure                                    [Pa]
114      !!
115      !! OUTPUT :
116      !! --------
117      !!    *  Cd     : drag coefficient
118      !!    *  Ch     : sensible heat coefficient
119      !!    *  Ce     : evaporation coefficient
120      !!    *  t_zu   : pot. air temperature adjusted at wind height zu       [K]
121      !!    *  q_zu   : specific humidity of air        //                    [kg/kg]
122      !!    *  U_blk  : bulk wind speed at 10m                                [m/s]
123      !!
124      !!
125      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://github.com/brodeau/aerobulk/)
126      !!----------------------------------------------------------------------------------
127      REAL(wp), INTENT(in   )                     ::   zt       ! height for t_zt and q_zt                    [m]
128      REAL(wp), INTENT(in   )                     ::   zu       ! height for U_zu                             [m]
129      REAL(wp), INTENT(inout), DIMENSION(jpi,jpj) ::   T_s      ! sea surface temperature                [Kelvin]
130      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   t_zt     ! potential air temperature              [Kelvin]
131      REAL(wp), INTENT(inout), DIMENSION(jpi,jpj) ::   q_s      ! sea surface specific humidity           [kg/kg]
132      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity at zt             [kg/kg]
133      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   U_zu     ! relative wind module at zu                [m/s]
134      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
135      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
136      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
137      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   t_zu     ! pot. air temp. adjusted at zu               [K]
138      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. humidity adjusted at zu           [kg/kg]
139      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   U_blk    ! bulk wind at 10m                          [m/s]
140      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cdn, Chn, Cen ! neutral transfer coefficients
141      !
142      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   Qsw      !             [W/m^2]
143      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   rad_lw   !             [W/m^2]
144      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   slp      !             [Pa]
145      !
146      INTEGER :: j_itt
147      LOGICAL :: l_zt_equal_zu = .FALSE.      ! if q and t are given at same height as U
148      !
149      REAL(wp), DIMENSION(jpi,jpj) ::  &
150         &  u_star, t_star, q_star, &
151         &  dt_zu, dq_zu,    &
152         &  znu_a,           & !: Nu_air, Viscosity of air
153         &  Linv,            & !: 1/L (inverse of Monin Obukhov length...
154         &  z0, z0t, z0q
155      !
156      ! Cool skin:
157      LOGICAL :: l_use_skin = .FALSE.
158      REAL(wp), DIMENSION(jpi,jpj) :: zsst    ! to back up the initial bulk SST
159
160      !
161      REAL(wp), DIMENSION(jpi,jpj) ::   func_m, func_h
162      REAL(wp), DIMENSION(jpi,jpj) ::   ztmp0, ztmp1, ztmp2
163      !!----------------------------------------------------------------------------------
164      !
165      ! Cool skin ?
166      IF( PRESENT(Qsw) .AND. PRESENT(rad_lw) .AND. PRESENT(slp) ) THEN
167         l_use_skin = .TRUE.
168      END IF
169      IF (lwp) PRINT *, ' *** LOLO(sbcblk_algo_ecmwf.F90) => l_use_skin =', l_use_skin
170
171      ! Identical first gess as in COARE, with IFS parameter values though
172      !
173      l_zt_equal_zu = .FALSE.
174      IF( ABS(zu - zt) < 0.01 )   l_zt_equal_zu = .TRUE.    ! testing "zu == zt" is risky with double precision
175
176      !! Initialization for cool skin:
177      IF( l_use_skin ) THEN
178         zsst   = T_s    ! save the bulk SST
179         T_s    = T_s - 0.25                      ! First guess of correction
180         q_s    = rdct_qsat_salt*q_sat(MAX(T_s, 200._wp), slp) ! First guess of q_s
181      END IF
182
183      !! First guess of temperature and humidity at height zu:
184      t_zu = MAX( t_zt , 0.0_wp  )   ! who knows what's given on masked-continental regions...
185      q_zu = MAX( q_zt , 1.e-6_wp)   !               "
186
187      !! Pot. temp. difference (and we don't want it to be 0!)
188      dt_zu = t_zu - T_s ;   dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
189      dq_zu = q_zu - q_s ;   dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
190
191      znu_a = visc_air(t_zt) ! Air viscosity (m^2/s) at zt given from temperature in (K)
192
193      ztmp2 = 0.5_wp*0.5_wp  ! initial guess for wind gustiness contribution
194      U_blk = SQRT(U_zu*U_zu + ztmp2)
195
196      ztmp2   = 10000._wp     ! optimization: ztmp2 == 1/z0 (with z0 first guess == 0.0001)
197      ztmp0   = LOG(zu*ztmp2)
198      ztmp1   = LOG(10.*ztmp2)
199      u_star = 0.035_wp*U_blk*ztmp1/ztmp0       ! (u* = 0.035*Un10)
200
201      z0     = charn0*u_star*u_star/grav + 0.11_wp*znu_a/u_star
202      z0     = MIN(ABS(z0), 0.001_wp)  ! (prevent FPE from stupid values from masked region later on...) !#LOLO
203      z0t    = 1._wp / ( 0.1_wp*EXP(vkarmn/(0.00115/(vkarmn/ztmp1))) )
204      z0t    = MIN(ABS(z0t), 0.001_wp)  ! (prevent FPE from stupid values from masked region later on...) !#LOLO
205
206      ztmp2  = vkarmn/ztmp0
207      Cd     = ztmp2*ztmp2    ! first guess of Cd
208
209      ztmp0 = vkarmn*vkarmn/LOG(zt/z0t)/Cd
210
211      ztmp2 = Ri_bulk( zu, T_s, t_zu, q_s, q_zu, U_blk ) ! Bulk Richardson Number (BRN)
212
213      !! First estimate of zeta_u, depending on the stability, ie sign of BRN (ztmp2):
214      ztmp1 = 0.5 + SIGN( 0.5_wp , ztmp2 )
215      func_m = ztmp0*ztmp2 ! temporary array !!
216      func_h = (1.-ztmp1) * (func_m/(1._wp+ztmp2/(-zu/(zi0*0.004_wp*Beta0**3)))) & !  BRN < 0 ! temporary array !!! func_h == zeta_u
217         &  +     ztmp1   * (func_m*(1._wp + 27._wp/9._wp*ztmp2/func_m))           !  BRN > 0
218      !#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" !
219
220      !! First guess M-O stability dependent scaling params.(u*,t*,q*) to estimate z0 and z/L
221      ztmp0  = vkarmn/(LOG(zu/z0t) - psi_h_ecmwf(func_h))
222
223      u_star = U_blk*vkarmn/(LOG(zu) - LOG(z0)  - psi_m_ecmwf(func_h))
224      t_star = dt_zu*ztmp0
225      q_star = dq_zu*ztmp0
226
227      ! What's need to be done if zt /= zu:
228      IF( .NOT. l_zt_equal_zu ) THEN
229         !! First update of values at zu (or zt for wind)
230         ztmp0 = psi_h_ecmwf(func_h) - psi_h_ecmwf(zt*func_h/zu)    ! zt*func_h/zu == zeta_t
231         ztmp1 = LOG(zt/zu) + ztmp0
232         t_zu = t_zt - t_star/vkarmn*ztmp1
233         q_zu = q_zt - q_star/vkarmn*ztmp1
234         q_zu = (0.5_wp + SIGN(0.5_wp,q_zu))*q_zu !Makes it impossible to have negative humidity :
235         !
236         dt_zu = t_zu - T_s  ; dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
237         dq_zu = q_zu - q_s  ; dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
238      END IF
239
240
241      !! => that was same first guess as in COARE...
242
243
244      !! First guess of inverse of Monin-Obukov length (1/L) :
245      Linv = One_on_L( t_zu, q_zu, u_star, t_star, q_star )
246
247      !! Functions such as  u* = U_blk*vkarmn/func_m
248      ztmp0 = zu*Linv
249      func_m = LOG(zu) - LOG(z0)  - psi_m_ecmwf(ztmp0) + psi_m_ecmwf( z0*Linv)
250      func_h = LOG(zu) - LOG(z0t) - psi_h_ecmwf(ztmp0) + psi_h_ecmwf(z0t*Linv)
251
252      !! ITERATION BLOCK
253      DO j_itt = 1, nb_itt
254
255         !! Bulk Richardson Number at z=zu (Eq. 3.25)
256         ztmp0 = Ri_bulk( zu, T_s, t_zu, q_s, q_zu, U_blk ) ! Bulk Richardson Number (BRN)
257
258         !! New estimate of the inverse of the Monin-Obukhon length (Linv == zeta/zu) :
259         Linv = ztmp0*func_m*func_m/func_h / zu     ! From Eq. 3.23, Chap.3.2.3, IFS doc - Cy40r1
260         !! Note: it is slightly different that the L we would get with the usual
261         !! expression, as in coare algorithm or in 'mod_phy.f90' (One_on_L_MO())
262         Linv = SIGN( MIN(ABS(Linv),200._wp), Linv ) ! (prevent FPE from stupid values from masked region later on...) !#LOLO
263
264         !! Update func_m with new Linv:
265         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!
266
267         !! Need to update roughness lengthes:
268         u_star = U_blk*vkarmn/func_m
269         ztmp2  = u_star*u_star
270         ztmp1  = znu_a/u_star
271         z0     = MIN( ABS( alpha_M*ztmp1 + charn0*ztmp2/grav ) , 0.001_wp)
272         z0t    = MIN( ABS( alpha_H*ztmp1                     ) , 0.001_wp)   ! eq.3.26, Chap.3, p.34, IFS doc - Cy31r1
273         z0q    = MIN( ABS( alpha_Q*ztmp1                     ) , 0.001_wp)
274
275         !! Update wind at 10m taking into acount convection-related wind gustiness:
276         !! => Chap. 3.2, IFS doc - Cy40r1, Eq.3.17 and Eq.3.18 + Eq.3.8
277         ! Only true when unstable (L<0) => when ztmp0 < 0 => - !!!
278         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)
279         !! => equivalent using Beta=1 (gustiness parameter, 1.25 for COARE, also zi0=600 in COARE..)
280         U_blk = MAX( SQRT(U_zu*U_zu + ztmp2) , 0.2_wp )              ! eq.3.17, Chap.3, p.32, IFS doc - Cy31r1
281         ! => 0.2 prevents U_blk to be 0 in stable case when U_zu=0.
282
283
284         !! Need to update "theta" and "q" at zu in case they are given at different heights
285         !! as well the air-sea differences:
286         IF( .NOT. l_zt_equal_zu ) THEN
287            !! Arrays func_m and func_h are free for a while so using them as temporary arrays...
288            func_h = psi_h_ecmwf(zu*Linv) ! temporary array !!!
289            func_m = psi_h_ecmwf(zt*Linv) ! temporary array !!!
290
291            ztmp2  = psi_h_ecmwf(z0t*Linv)
292            ztmp0  = func_h - ztmp2
293            ztmp1  = vkarmn/(LOG(zu) - LOG(z0t) - ztmp0)
294            t_star = dt_zu*ztmp1
295            ztmp2  = ztmp0 - func_m + ztmp2
296            ztmp1  = LOG(zt/zu) + ztmp2
297            t_zu   = t_zt - t_star/vkarmn*ztmp1
298
299            ztmp2  = psi_h_ecmwf(z0q*Linv)
300            ztmp0  = func_h - ztmp2
301            ztmp1  = vkarmn/(LOG(zu) - LOG(z0q) - ztmp0)
302            q_star = dq_zu*ztmp1
303            ztmp2  = ztmp0 - func_m + ztmp2
304            ztmp1  = LOG(zt/zu) + ztmp2
305            q_zu   = q_zt - q_star/vkarmn*ztmp1
306         END IF
307
308         !! Updating because of updated z0 and z0t and new Linv...
309         ztmp0 = zu*Linv
310         func_m = log(zu) - LOG(z0 ) - psi_m_ecmwf(ztmp0) + psi_m_ecmwf(z0 *Linv)
311         func_h = log(zu) - LOG(z0t) - psi_h_ecmwf(ztmp0) + psi_h_ecmwf(z0t*Linv)
312
313         !! SKIN related part
314         !! -----------------
315         IF( l_use_skin ) THEN
316            !! compute transfer coefficients at zu : lolo: verifier...
317            Cd = vkarmn*vkarmn/(func_m*func_m)
318            Ch = vkarmn*vkarmn/(func_m*func_h)
319            ztmp1 = LOG(zu) - LOG(z0q) - psi_h_ecmwf(ztmp0) + psi_h_ecmwf(z0q*Linv)   ! func_q
320            Ce = vkarmn*vkarmn/(func_m*ztmp1)
321            ! Non-Solar heat flux to the ocean:
322            ztmp1 = U_blk*MAX(rho_air(t_zu, q_zu, slp), 1._wp)     ! rho*U10
323            ztmp2 = T_s*T_s
324            ztmp1 = ztmp1 * ( Ce*rLevap*(q_zu - q_s) + Ch*rCp_dry*(t_zu - T_s) ) & ! Total turb. heat flux
325               &       +    (rad_lw - emiss_w*stefan*ztmp2*ztmp2)                  ! Net longwave flux
326            !! Updating the values of the skin temperature T_s and q_s :
327            CALL CSWL_ECMWF( Qsw, ztmp1, u_star, zsst, T_s )
328            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
329         END IF
330
331         IF( (l_use_skin).OR.(.NOT. l_zt_equal_zu) ) THEN
332            dt_zu = t_zu - T_s ;  dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
333            dq_zu = q_zu - q_s ;  dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
334         END IF
335
336      END DO !DO j_itt = 1, nb_itt
337
338      Cd = vkarmn*vkarmn/(func_m*func_m)
339      Ch = vkarmn*vkarmn/(func_m*func_h)
340      ztmp2 = log(zu/z0q) - psi_h_ecmwf(zu*Linv) + psi_h_ecmwf(z0q*Linv)   ! func_q
341      Ce = vkarmn*vkarmn/(func_m*ztmp2)
342
343      Cdn = vkarmn*vkarmn / (log(zu/z0 )*log(zu/z0 ))
344      Chn = vkarmn*vkarmn / (log(zu/z0t)*log(zu/z0t))
345      Cen = vkarmn*vkarmn / (log(zu/z0q)*log(zu/z0q))
346
347   END SUBROUTINE TURB_ECMWF
348
349
350   FUNCTION psi_m_ecmwf( pzeta )
351      !!----------------------------------------------------------------------------------
352      !! Universal profile stability function for momentum
353      !!     ECMWF / as in IFS cy31r1 documentation, available online
354      !!     at ecmwf.int
355      !!
356      !! pzeta : stability paramenter, z/L where z is altitude measurement
357      !!         and L is M-O length
358      !!
359      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://github.com/brodeau/aerobulk/)
360      !!----------------------------------------------------------------------------------
361      REAL(wp), DIMENSION(jpi,jpj) :: psi_m_ecmwf
362      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
363      !
364      INTEGER  ::   ji, jj    ! dummy loop indices
365      REAL(wp) :: zzeta, zx, ztmp, psi_unst, psi_stab, stab
366      !!----------------------------------------------------------------------------------
367      !
368      DO jj = 1, jpj
369         DO ji = 1, jpi
370            !
371            zzeta = MIN( pzeta(ji,jj) , 5._wp ) !! Very stable conditions (L positif and big!):
372            !
373            ! Unstable (Paulson 1970):
374            !   eq.3.20, Chap.3, p.33, IFS doc - Cy31r1
375            zx = SQRT(ABS(1._wp - 16._wp*zzeta))
376            ztmp = 1._wp + SQRT(zx)
377            ztmp = ztmp*ztmp
378            psi_unst = LOG( 0.125_wp*ztmp*(1._wp + zx) )   &
379               &       -2._wp*ATAN( SQRT(zx) ) + 0.5_wp*rpi
380            !
381            ! Unstable:
382            ! eq.3.22, Chap.3, p.33, IFS doc - Cy31r1
383            psi_stab = -2._wp/3._wp*(zzeta - 5._wp/0.35_wp)*EXP(-0.35_wp*zzeta) &
384               &       - zzeta - 2._wp/3._wp*5._wp/0.35_wp
385            !
386            ! Combining:
387            stab = 0.5_wp + SIGN(0.5_wp, zzeta) ! zzeta > 0 => stab = 1
388            !
389            psi_m_ecmwf(ji,jj) = (1._wp - stab) * psi_unst & ! (zzeta < 0) Unstable
390               &                +      stab  * psi_stab      ! (zzeta > 0) Stable
391            !
392         END DO
393      END DO
394      !
395   END FUNCTION psi_m_ecmwf
396
397
398   FUNCTION psi_h_ecmwf( pzeta )
399      !!----------------------------------------------------------------------------------
400      !! Universal profile stability function for temperature and humidity
401      !!     ECMWF / as in IFS cy31r1 documentation, available online
402      !!     at ecmwf.int
403      !!
404      !! pzeta : stability paramenter, z/L where z is altitude measurement
405      !!         and L is M-O length
406      !!
407      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://github.com/brodeau/aerobulk/)
408      !!----------------------------------------------------------------------------------
409      REAL(wp), DIMENSION(jpi,jpj) :: psi_h_ecmwf
410      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
411      !
412      INTEGER  ::   ji, jj     ! dummy loop indices
413      REAL(wp) ::  zzeta, zx, psi_unst, psi_stab, stab
414      !!----------------------------------------------------------------------------------
415      !
416      DO jj = 1, jpj
417         DO ji = 1, jpi
418            !
419            zzeta = MIN(pzeta(ji,jj) , 5._wp)   ! Very stable conditions (L positif and big!):
420            !
421            zx  = ABS(1._wp - 16._wp*zzeta)**.25        ! this is actually (1/phi_m)**2  !!!
422            !                                     ! eq.3.19, Chap.3, p.33, IFS doc - Cy31r1
423            ! Unstable (Paulson 1970) :
424            psi_unst = 2._wp*LOG(0.5_wp*(1._wp + zx*zx))   ! eq.3.20, Chap.3, p.33, IFS doc - Cy31r1
425            !
426            ! Stable:
427            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
428               &       - ABS(1._wp + 2._wp/3._wp*zzeta)**1.5_wp - 2._wp/3._wp*5._wp/0.35_wp + 1._wp
429            ! LB: added ABS() to avoid NaN values when unstable, which contaminates the unstable solution...
430            !
431            stab = 0.5_wp + SIGN(0.5_wp, zzeta) ! zzeta > 0 => stab = 1
432            !
433            !
434            psi_h_ecmwf(ji,jj) = (1._wp - stab) * psi_unst &   ! (zzeta < 0) Unstable
435               &                +    stab    * psi_stab        ! (zzeta > 0) Stable
436            !
437         END DO
438      END DO
439      !
440   END FUNCTION psi_h_ecmwf
441
442
443
444
445
446   !!======================================================================
447END MODULE sbcblk_algo_ecmwf
Note: See TracBrowser for help on using the repository browser.