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_coare3p6.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_coare3p6.F90 @ 11772

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

LB: solid updates+improvements of cool-skin/warm-layer capabilty of COARE and ECMWF bulk algorithms!

File size: 25.9 KB
Line 
1MODULE sbcblk_algo_coare3p6
2   !!======================================================================
3   !!                   ***  MODULE  sbcblk_algo_coare3p6  ***
4   !!
5   !!                After Fairall et al 2018 & Edson et al 2013
6   !! Computes:
7   !!   * bulk transfer coefficients C_D, C_E and C_H
8   !!   * air temp. and spec. hum. adjusted from zt (2m) to zu (10m) if needed
9   !!   * the effective bulk wind speed at 10m U_blk
10   !!   => all these are used in bulk formulas in sbcblk.F90
11   !!
12   !!       Routine turb_coare3p6 maintained and developed in AeroBulk
13   !!                     (https://github.com/brodeau/aerobulk)
14   !!
15   !! ** Author: L. Brodeau, June 2019 / AeroBulk (https://github.com/brodeau/aerobulk)
16   !!----------------------------------------------------------------------
17   !! History :  4.0  !  2016-02  (L.Brodeau)   Original code
18   !!----------------------------------------------------------------------
19
20   !!----------------------------------------------------------------------
21   !!   turb_coare3p6  : computes the bulk turbulent transfer coefficients
22   !!                   adjusts t_air and q_air from zt to zu m
23   !!                   returns the effective bulk wind speed at 10m
24   !!----------------------------------------------------------------------
25   USE oce             ! ocean dynamics and tracers
26   USE dom_oce         ! ocean space and time domain
27   USE phycst          ! physical constants
28   USE iom             ! I/O manager library
29   USE lib_mpp         ! distribued memory computing library
30   USE in_out_manager  ! I/O manager
31   USE prtctl          ! Print control
32   USE sbcwave, ONLY   :  cdn_wave ! wave module
33#if defined key_si3 || defined key_cice
34   USE sbc_ice         ! Surface boundary condition: ice fields
35#endif
36   USE lib_fortran     ! to use key_nosignedzero
37
38   USE sbc_oce         ! Surface boundary condition: ocean fields
39   USE sbcblk_phy      ! all thermodynamics functions, rho_air, q_sat, etc... !LB
40   USE sbcblk_skin_coare ! cool-skin/warm layer scheme (CSWL_ECMWF) !LB
41
42   IMPLICIT NONE
43   PRIVATE
44
45   PUBLIC :: COARE3P6_INIT, TURB_COARE3P6
46
47   !! COARE own values for given constants:
48   REAL(wp), PARAMETER :: zi0   = 600._wp     ! scale height of the atmospheric boundary layer...
49   REAL(wp), PARAMETER :: Beta0 =  1.2_wp     ! gustiness parameter
50
51   INTEGER , PARAMETER ::   nb_itt = 10             ! number of itterations
52
53   !!----------------------------------------------------------------------
54CONTAINS
55
56
57   SUBROUTINE coare3p6_init(l_use_cs, l_use_wl)
58      !!---------------------------------------------------------------------
59      !!                  ***  FUNCTION coare3p6_init  ***
60      !!
61      !! INPUT :
62      !! -------
63      !!    * l_use_cs : use the cool-skin parameterization
64      !!    * l_use_wl : use the warm-layer parameterization
65      !!---------------------------------------------------------------------
66      LOGICAL , INTENT(in) ::   l_use_cs ! use the cool-skin parameterization
67      LOGICAL , INTENT(in) ::   l_use_wl ! use the warm-layer parameterization
68      INTEGER :: ierr
69      !!---------------------------------------------------------------------
70      IF ( l_use_wl ) THEN
71         ierr = 0
72         ALLOCATE ( Tau_ac(jpi,jpj) , Qnt_ac(jpi,jpj), Hz_wl(jpi,jpj), dT_wl(jpi,jpj), STAT=ierr )
73         !IF( ierr > 0 ) STOP ' COARE3P6_INIT => allocation of Tau_ac, Qnt_ac, dT_wl & Hz_wl failed!'
74         Tau_ac(:,:) = 0._wp
75         Qnt_ac(:,:) = 0._wp
76         dT_wl(:,:)  = 0._wp
77         Hz_wl(:,:)  = Hwl_max
78      END IF
79      !!
80      IF ( l_use_cs ) THEN
81         ierr = 0
82         ALLOCATE ( dT_cs(jpi,jpj), STAT=ierr )
83         !IF( ierr > 0 ) STOP ' COARE3P6_INIT => allocation of dT_cs failed!'
84         dT_cs(:,:) = -0.25_wp  ! First guess of skin correction
85      END IF
86   END SUBROUTINE coare3p6_init
87
88
89
90   SUBROUTINE turb_coare3p6( kt, zt, zu, T_s, t_zt, q_s, q_zt, U_zu, l_use_cs, l_use_wl, &
91      &                      Cd, Ch, Ce, t_zu, q_zu, U_blk,                              &
92      &                      Cdn, Chn, Cen,                                              &
93      &                      Qsw, rad_lw, slp, pdT_cs,                                   & ! optionals for cool-skin (and warm-layer)
94      &                      pdT_wl, pHz_wl )                                                 ! optionals for warm-layer only
95      !!----------------------------------------------------------------------
96      !!                      ***  ROUTINE  turb_coare3p6  ***
97      !!
98      !! ** Purpose :   Computes turbulent transfert coefficients of surface
99      !!                fluxes according to Fairall et al. (2003)
100      !!                If relevant (zt /= zu), adjust temperature and humidity from height zt to zu
101      !!                Returns the effective bulk wind speed at zu to be used in the bulk formulas
102      !!
103      !!                Applies the cool-skin warm-layer correction of the SST to T_s
104      !!                if the net shortwave flux at the surface (Qsw), the downwelling longwave
105      !!                radiative fluxes at the surface (rad_lw), and the sea-leve pressure (slp)
106      !!                are provided as (optional) arguments!
107      !!
108      !! INPUT :
109      !! -------
110      !!    *  kt   : current time step (starts at 1)
111      !!    *  zt   : height for temperature and spec. hum. of air            [m]
112      !!    *  zu   : height for wind speed (usually 10m)                     [m]
113      !!    *  t_zt : potential air temperature at zt                         [K]
114      !!    *  q_zt : specific humidity of air at zt                          [kg/kg]
115      !!    *  U_zu : scalar wind speed at zu                                 [m/s]
116      !!    * l_use_cs : use the cool-skin parameterization
117      !!    * l_use_wl : use the warm-layer parameterization
118      !!
119      !! INPUT/OUTPUT:
120      !! -------------
121      !!    *  T_s  : always "bulk SST" as input                              [K]
122      !!              -> unchanged "bulk SST" as output if CSWL not used      [K]
123      !!              -> skin temperature as output if CSWL used              [K]
124      !!
125      !!    *  q_s  : SSQ aka saturation specific humidity at temp. T_s       [kg/kg]
126      !!              -> doesn't need to be given a value if skin temp computed (in case l_use_cs=True or l_use_wl=True)
127      !!              -> MUST be given the correct value if not computing skint temp. (in case l_use_cs=False or l_use_wl=False)
128      !!
129      !! OPTIONAL INPUT:
130      !! ---------------
131      !!    *  Qsw    : net solar flux (after albedo) at the surface (>0)     [W/m^2]
132      !!    *  rad_lw : downwelling longwave radiation at the surface  (>0)   [W/m^2]
133      !!    *  slp    : sea-level pressure                                    [Pa]
134      !!
135      !! OPTIONAL OUTPUT:
136      !! ----------------
137      !!    * pdT_cs  : SST increment "dT" for cool-skin correction           [K]
138      !!    * pdT_wl  : SST increment "dT" for warm-layer correction          [K]
139      !!    * pHz_wl  : thickness of warm-layer                               [m]
140      !!
141      !! OUTPUT :
142      !! --------
143      !!    *  Cd     : drag coefficient
144      !!    *  Ch     : sensible heat coefficient
145      !!    *  Ce     : evaporation coefficient
146      !!    *  t_zu   : pot. air temperature adjusted at wind height zu       [K]
147      !!    *  q_zu   : specific humidity of air        //                    [kg/kg]
148      !!    *  U_blk  : bulk wind speed at zu                                 [m/s]
149      !!
150      !!
151      !! ** Author: L. Brodeau, June 2019 / AeroBulk (https://github.com/brodeau/aerobulk/)
152      !!----------------------------------------------------------------------------------
153      INTEGER,  INTENT(in   )                     ::   kt       ! current time step
154      REAL(wp), INTENT(in   )                     ::   zt       ! height for t_zt and q_zt                    [m]
155      REAL(wp), INTENT(in   )                     ::   zu       ! height for U_zu                             [m]
156      REAL(wp), INTENT(inout), DIMENSION(jpi,jpj) ::   T_s      ! sea surface temperature                [Kelvin]
157      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   t_zt     ! potential air temperature              [Kelvin]
158      REAL(wp), INTENT(inout), DIMENSION(jpi,jpj) ::   q_s      ! sea surface specific humidity           [kg/kg]
159      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity at zt             [kg/kg]
160      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   U_zu     ! relative wind module at zu                [m/s]
161      LOGICAL , INTENT(in   )                     ::   l_use_cs ! use the cool-skin parameterization
162      LOGICAL , INTENT(in   )                     ::   l_use_wl ! use the warm-layer parameterization
163      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
164      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
165      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
166      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   t_zu     ! pot. air temp. adjusted at zu               [K]
167      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. humidity adjusted at zu           [kg/kg]
168      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   U_blk    ! bulk wind speed at zu                     [m/s]
169      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cdn, Chn, Cen ! neutral transfer coefficients
170      !
171      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   Qsw      !             [W/m^2]
172      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   rad_lw   !             [W/m^2]
173      REAL(wp), INTENT(in   ), OPTIONAL, DIMENSION(jpi,jpj) ::   slp      !             [Pa]
174      REAL(wp), INTENT(  out), OPTIONAL, DIMENSION(jpi,jpj) ::   pdT_cs
175      REAL(wp), INTENT(  out), OPTIONAL, DIMENSION(jpi,jpj) ::   pdT_wl   !             [K]
176      REAL(wp), INTENT(  out), OPTIONAL, DIMENSION(jpi,jpj) ::   pHz_wl   !             [m]
177      !
178      INTEGER :: j_itt
179      LOGICAL :: l_zt_equal_zu = .FALSE.      ! if q and t are given at same height as U
180      !
181      REAL(wp), DIMENSION(jpi,jpj) ::  &
182         &  u_star, t_star, q_star, &
183         &  dt_zu, dq_zu,    &
184         &  znu_a,           & !: Nu_air, Viscosity of air
185         &  z0, z0t
186      REAL(wp), DIMENSION(jpi,jpj) ::   zeta_u        ! stability parameter at height zu
187      REAL(wp), DIMENSION(jpi,jpj) ::   ztmp0, ztmp1, ztmp2
188      REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   zeta_t        ! stability parameter at height zt
189      !
190      REAL(wp), DIMENSION(:,:), ALLOCATABLE :: &
191         &                zsst,   &  ! to back up the initial bulk SST
192         &                pdTc,   &  ! SST increment "dT" for cool-skin correction           [K]
193         &                pdTw,   &  ! SST increment "dT" for warm layer correction          [K]
194         &                zHwl       ! depth of warm-layer [m]
195      CHARACTER(len=40), PARAMETER :: crtnm = 'turb_coare3p6@sbcblk_algo_coare3p6'
196      !!----------------------------------------------------------------------------------
197
198      IF ( kt == nit000 ) CALL COARE3P6_INIT(l_use_cs, l_use_wl)
199
200      l_zt_equal_zu = .FALSE.
201      IF( ABS(zu - zt) < 0.01_wp )   l_zt_equal_zu = .TRUE.    ! testing "zu == zt" is risky with double precision
202      IF( .NOT. l_zt_equal_zu )  ALLOCATE( zeta_t(jpi,jpj) )
203
204      !! Initializations for cool skin and warm layer:
205      IF ( l_use_cs ) THEN
206         IF( .NOT.(PRESENT(Qsw) .AND. PRESENT(rad_lw) .AND. PRESENT(slp)) ) THEN
207            PRINT *, ' * PROBLEM ('//TRIM(crtnm)//'): you need to provide Qsw, rad_lw & slp to use cool-skin param!'; STOP
208         END IF
209      END IF
210
211      IF ( l_use_wl ) THEN
212         IF( .NOT.(PRESENT(Qsw) .AND. PRESENT(rad_lw) .AND. PRESENT(slp)) ) THEN
213            PRINT *, ' * PROBLEM ('//TRIM(crtnm)//'): you need to provide Qsw, rad_lw & slp to use warm-layer param!'; STOP
214         END IF
215      END IF
216
217      IF ( l_use_cs .OR. l_use_wl ) THEN
218         ALLOCATE ( zsst(jpi,jpj) )
219         zsst = T_s ! backing up the bulk SST
220         IF( l_use_cs ) T_s = T_s - 0.25_wp   ! First guess of correction
221         q_s    = rdct_qsat_salt*q_sat(MAX(T_s, 200._wp), slp) ! First guess of q_s !LOLO WL too!!!
222      END IF
223
224
225      !! First guess of temperature and humidity at height zu:
226      t_zu = MAX( t_zt ,  180._wp )   ! who knows what's given on masked-continental regions...
227      q_zu = MAX( q_zt , 1.e-6_wp )   !               "
228
229      !! Pot. temp. difference (and we don't want it to be 0!)
230      dt_zu = t_zu - T_s ;   dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
231      dq_zu = q_zu - q_s ;   dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
232
233      znu_a = visc_air(t_zu) ! Air viscosity (m^2/s) at zt given from temperature in (K)
234
235      U_blk = SQRT(U_zu*U_zu + 0.5_wp*0.5_wp) ! initial guess for wind gustiness contribution
236
237      ztmp0   = LOG(    zu*10000._wp) ! optimization: 10000. == 1/z0 (with z0 first guess == 0.0001)
238      ztmp1   = LOG(10._wp*10000._wp) !       "                    "               "
239      u_star = 0.035_wp*U_blk*ztmp1/ztmp0       ! (u* = 0.035*Un10)
240
241      z0     = alfa_charn_3p6(U_zu)*u_star*u_star/grav + 0.11_wp*znu_a/u_star
242      z0     = MIN(ABS(z0), 0.001_wp)  ! (prevent FPE from stupid values from masked region later on...) !#LOLO
243      z0t    = 1._wp / ( 0.1_wp*EXP(vkarmn/(0.00115/(vkarmn/ztmp1))) )
244      z0t    = MIN(ABS(z0t), 0.001_wp)  ! (prevent FPE from stupid values from masked region later on...) !#LOLO
245
246      Cd     = (vkarmn/ztmp0)**2    ! first guess of Cd
247
248      ztmp0 = vkarmn*vkarmn/LOG(zt/z0t)/Cd
249
250      ztmp2 = Ri_bulk( zu, T_s, t_zu, q_s, q_zu, U_blk ) ! Bulk Richardson Number (BRN)
251
252      !! First estimate of zeta_u, depending on the stability, ie sign of BRN (ztmp2):
253      ztmp1 = 0.5 + SIGN( 0.5_wp , ztmp2 )
254      ztmp0 = ztmp0*ztmp2
255      zeta_u = (1._wp-ztmp1) * (ztmp0/(1._wp+ztmp2/(-zu/(zi0*0.004_wp*Beta0**3)))) & !  BRN < 0
256         &  +     ztmp1   * (ztmp0*(1._wp + 27._wp/9._wp*ztmp2/ztmp0))               !  BRN > 0
257      !#LB: should make sure that the "ztmp0" of "27./9.*ztmp2/ztmp0" is "ztmp0*ztmp2" and not "ztmp0==vkarmn*vkarmn/LOG(zt/z0t)/Cd" !
258
259      !! First guess M-O stability dependent scaling params.(u*,t*,q*) to estimate z0 and z/L
260      ztmp0  = vkarmn/(LOG(zu/z0t) - psi_h_coare(zeta_u))
261
262      u_star = U_blk*vkarmn/(LOG(zu) - LOG(z0)  - psi_m_coare(zeta_u))
263      t_star = dt_zu*ztmp0
264      q_star = dq_zu*ztmp0
265
266      ! What needs to be done if zt /= zu:
267      IF( .NOT. l_zt_equal_zu ) THEN
268         !! First update of values at zu (or zt for wind)
269         zeta_t = zt*zeta_u/zu
270         ztmp0 = psi_h_coare(zeta_u) - psi_h_coare(zeta_t)
271         ztmp1 = LOG(zt/zu) + ztmp0
272         t_zu = t_zt - t_star/vkarmn*ztmp1
273         q_zu = q_zt - q_star/vkarmn*ztmp1
274         q_zu = (0.5_wp + SIGN(0.5_wp,q_zu))*q_zu !Makes it impossible to have negative humidity :
275         !
276         dt_zu = t_zu - T_s  ; dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
277         dq_zu = q_zu - q_s  ; dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
278      END IF
279
280      !! ITERATION BLOCK
281      DO j_itt = 1, nb_itt
282
283         !!Inverse of Monin-Obukov length (1/L) :
284         ztmp0 = One_on_L(t_zu, q_zu, u_star, t_star, q_star)  ! 1/L == 1/[Monin-Obukhov length]
285         ztmp0 = SIGN( MIN(ABS(ztmp0),200._wp), ztmp0 ) ! (prevents FPE from stupid values from masked region later on...) !#LOLO
286
287         ztmp1 = u_star*u_star   ! u*^2
288
289         !! Update wind at zu with convection-related wind gustiness in unstable conditions (Fairall et al. 2003, Eq.8):
290         ztmp2 = Beta0*Beta0*ztmp1*(MAX(-zi0*ztmp0/vkarmn,0._wp))**(2._wp/3._wp) ! square of wind gustiness contribution, ztmp2 == Ug^2
291         !!   ! Only true when unstable (L<0) => when ztmp0 < 0 => explains "-" before zi0
292         U_blk = MAX(SQRT(U_zu*U_zu + ztmp2), 0.2_wp)        ! include gustiness in bulk wind speed
293         ! => 0.2 prevents U_blk to be 0 in stable case when U_zu=0.
294
295         !! Stability parameters:
296         zeta_u = zu*ztmp0
297         zeta_u = SIGN( MIN(ABS(zeta_u),50.0_wp), zeta_u )
298         IF( .NOT. l_zt_equal_zu ) THEN
299            zeta_t = zt*ztmp0
300            zeta_t = SIGN( MIN(ABS(zeta_t),50.0_wp), zeta_t )
301         END IF
302
303         !! Adjustment the wind at 10m (not needed in the current algo form):
304         !IF ( zu \= 10._wp ) U10 = U_zu + u_star/vkarmn*(LOG(10._wp/zu) - psi_m_coare(10._wp*ztmp0) + psi_m_coare(zeta_u))
305
306         !! Roughness lengthes z0, z0t (z0q = z0t) :
307         ztmp2 = u_star/vkarmn*LOG(10./z0)                                 ! Neutral wind speed at 10m
308         z0    = alfa_charn_3p6(ztmp2)*ztmp1/grav + 0.11_wp*znu_a/u_star   ! Roughness length (eq.6) [ ztmp1==u*^2 ]
309         ztmp1 = ( znu_a / (z0*u_star) )**0.72_wp     ! COARE3.6-specific! (1./Re_r)^0.72 (Re_r: roughness Reynolds number) COARE3.6-specific!
310         z0t   = MIN( 1.6E-4_wp , 5.8E-5_wp*ztmp1 )   ! COARE3.6-specific!
311
312         !! Turbulent scales at zu :
313         ztmp0   = psi_h_coare(zeta_u)
314         ztmp1   = vkarmn/(LOG(zu) - LOG(z0t) - ztmp0) ! #LOLO: in ztmp0, some use psi_h_coare(zeta_t) rather than psi_h_coare(zeta_t) ???
315
316         t_star = dt_zu*ztmp1
317         q_star = dq_zu*ztmp1
318         u_star = U_blk*vkarmn/(LOG(zu) - LOG(z0) - psi_m_coare(zeta_u))
319
320         IF( .NOT. l_zt_equal_zu ) THEN
321            !! Re-updating temperature and humidity at zu if zt /= zu :
322            ztmp1 = LOG(zt/zu) + ztmp0 - psi_h_coare(zeta_t)
323            t_zu = t_zt - t_star/vkarmn*ztmp1
324            q_zu = q_zt - q_star/vkarmn*ztmp1
325         END IF
326
327
328         IF( l_use_cs ) THEN
329            !! Cool-skin contribution
330
331            CALL UPDATE_QNSOL_TAU( zu, T_s, q_s, t_zu, q_zu, u_star, t_star, q_star, U_zu, U_blk, slp, rad_lw, &
332               &                   ztmp1, zeta_u,  Qlat=ztmp2)  ! Qnsol -> ztmp1 / Tau -> zeta_u
333
334            CALL CS_COARE( Qsw, ztmp1, u_star, zsst, ztmp2 )  ! ! Qnsol -> ztmp1 / Qlat -> ztmp2
335
336            T_s(:,:) = zsst(:,:) + dT_cs(:,:)*tmask(:,:,1)
337            IF( l_use_wl ) T_s(:,:) = T_s(:,:) + dT_wl(:,:)*tmask(:,:,1)
338            q_s(:,:) = rdct_qsat_salt*q_sat(MAX(T_s(:,:), 200._wp), slp(:,:))
339
340         END IF
341
342         IF( l_use_wl ) THEN
343            !! Warm-layer contribution
344            CALL UPDATE_QNSOL_TAU( zu, T_s, q_s, t_zu, q_zu, u_star, t_star, q_star, U_zu, U_blk, slp, rad_lw, &
345               &                   ztmp1, zeta_u)  ! Qnsol -> ztmp1 / Tau -> zeta_u
346            !! In WL_COARE or , Tau_ac and Qnt_ac must be updated at the final itteration step => add a flag to do this!
347            CALL WL_COARE( Qsw, ztmp1, zeta_u, zsst, MOD(nb_itt,j_itt) )
348
349            !! Updating T_s and q_s !!!
350            T_s(:,:) = zsst(:,:) + dT_wl(:,:)*tmask(:,:,1)
351            IF( l_use_cs ) T_s(:,:) = T_s(:,:) + dT_cs(:,:)*tmask(:,:,1)
352            q_s(:,:) = rdct_qsat_salt*q_sat(MAX(T_s(:,:), 200._wp), slp(:,:))
353         END IF
354
355         IF( l_use_cs .OR. l_use_wl .OR. (.NOT. l_zt_equal_zu) ) THEN
356            dt_zu = t_zu - T_s ;  dt_zu = SIGN( MAX(ABS(dt_zu),1.E-6_wp), dt_zu )
357            dq_zu = q_zu - q_s ;  dq_zu = SIGN( MAX(ABS(dq_zu),1.E-9_wp), dq_zu )
358         END IF
359
360      END DO !DO j_itt = 1, nb_itt
361
362      ! compute transfer coefficients at zu :
363      ztmp0 = u_star/U_blk
364      Cd   = ztmp0*ztmp0
365      Ch   = ztmp0*t_star/dt_zu
366      Ce   = ztmp0*q_star/dq_zu
367
368      ztmp1 = zu + z0
369      Cdn = vkarmn*vkarmn / (log(ztmp1/z0 )*log(ztmp1/z0 ))
370      Chn = vkarmn*vkarmn / (log(ztmp1/z0t)*log(ztmp1/z0t))
371      Cen = Chn
372
373      IF( .NOT. l_zt_equal_zu ) DEALLOCATE( zeta_t )
374
375      IF ( l_use_cs .AND. PRESENT(pdT_cs) ) pdT_cs = dT_cs
376      IF ( l_use_wl .AND. PRESENT(pdT_wl) ) pdT_wl = dT_wl
377      IF ( l_use_wl .AND. PRESENT(pHz_wl) ) pHz_wl = Hz_wl
378
379      IF ( l_use_cs .OR. l_use_wl ) DEALLOCATE ( zsst )
380
381   END SUBROUTINE turb_coare3p6
382
383
384   FUNCTION alfa_charn_3p6( pwnd )
385      !!-------------------------------------------------------------------
386      !! Computes the Charnock parameter as a function of the Neutral wind speed at 10m
387      !!  "wind speed dependent formulation"
388      !!  (Eq. 13 in Edson et al., 2013)
389      !!
390      !! Author: L. Brodeau, July 2019 / AeroBulk  (https://github.com/brodeau/aerobulk/)
391      !!-------------------------------------------------------------------
392      REAL(wp), DIMENSION(jpi,jpj) :: alfa_charn_3p6
393      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pwnd   ! neutral wind speed at 10m
394      !
395      REAL(wp), PARAMETER :: charn0_max = 0.028  !: value above which the Charnock parameter levels off for winds > 18 m/s
396      !!-------------------------------------------------------------------
397      alfa_charn_3p6 = MAX( MIN( 0.0017_wp*pwnd - 0.005_wp , charn0_max) , 0._wp )
398      !!
399   END FUNCTION alfa_charn_3p6
400
401   FUNCTION alfa_charn_3p6_wave( pus, pwsh, pwps )
402      !!-------------------------------------------------------------------
403      !! Computes the Charnock parameter as a function of wave information and u*
404      !!
405      !!  (COARE 3.6, Fairall et al., 2018)
406      !!
407      !! Author: L. Brodeau, October 2019 / AeroBulk  (https://github.com/brodeau/aerobulk/)
408      !!-------------------------------------------------------------------
409      REAL(wp), DIMENSION(jpi,jpj) :: alfa_charn_3p6_wave
410      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pus   ! friction velocity             [m/s]
411      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pwsh  ! significant wave height       [m]
412      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pwps  ! phase speed of dominant waves [m/s]
413      !!-------------------------------------------------------------------
414      alfa_charn_3p6_wave = ( pwsh*0.2_wp*(pus/pwps)**2.2_wp ) * grav/(pus*pus)
415      !!
416   END FUNCTION alfa_charn_3p6_wave
417
418
419   FUNCTION psi_m_coare( pzeta )
420      !!----------------------------------------------------------------------------------
421      !! ** Purpose: compute the universal profile stability function for momentum
422      !!             COARE 3.0, Fairall et al. 2003
423      !!             pzeta : stability paramenter, z/L where z is altitude
424      !!                     measurement and L is M-O length
425      !!       Stability function for wind speed and scalars matching Kansas and free
426      !!       convection forms with weighting f convective form, follows Fairall et
427      !!       al (1996) with profile constants from Grachev et al (2000) BLM stable
428      !!       form from Beljaars and Holtslag (1991)
429      !!
430      !! ** Author: L. Brodeau, June 2016 / AeroBulk (https://github.com/brodeau/aerobulk/)
431      !!----------------------------------------------------------------------------------
432      REAL(wp), DIMENSION(jpi,jpj) :: psi_m_coare
433      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
434      !
435      INTEGER  ::   ji, jj    ! dummy loop indices
436      REAL(wp) :: zta, zphi_m, zphi_c, zpsi_k, zpsi_c, zf, zc, zstab
437      !!----------------------------------------------------------------------------------
438      !
439      DO jj = 1, jpj
440         DO ji = 1, jpi
441            !
442            zta = pzeta(ji,jj)
443            !
444            zphi_m = ABS(1. - 15.*zta)**.25    !!Kansas unstable
445            !
446            zpsi_k = 2.*LOG((1. + zphi_m)/2.) + LOG((1. + zphi_m*zphi_m)/2.)   &
447               & - 2.*ATAN(zphi_m) + 0.5*rpi
448            !
449            zphi_c = ABS(1. - 10.15*zta)**.3333                   !!Convective
450            !
451            zpsi_c = 1.5*LOG((1. + zphi_c + zphi_c*zphi_c)/3.) &
452               &     - 1.7320508*ATAN((1. + 2.*zphi_c)/1.7320508) + 1.813799447
453            !
454            zf = zta*zta
455            zf = zf/(1. + zf)
456            zc = MIN(50._wp, 0.35_wp*zta)
457            zstab = 0.5 + SIGN(0.5_wp, zta)
458            !
459            psi_m_coare(ji,jj) = (1. - zstab) * ( (1. - zf)*zpsi_k + zf*zpsi_c ) & ! (zta < 0)
460               &                -   zstab     * ( 1. + 1.*zta     &                ! (zta > 0)
461               &                         + 0.6667*(zta - 14.28)/EXP(zc) + 8.525 )   !     "
462            !
463         END DO
464      END DO
465      !
466   END FUNCTION psi_m_coare
467
468
469   FUNCTION psi_h_coare( pzeta )
470      !!---------------------------------------------------------------------
471      !! Universal profile stability function for temperature and humidity
472      !! COARE 3.0, Fairall et al. 2003
473      !!
474      !! pzeta : stability paramenter, z/L where z is altitude measurement
475      !!         and L is M-O length
476      !!
477      !! Stability function for wind speed and scalars matching Kansas and free
478      !! convection forms with weighting f convective form, follows Fairall et
479      !! al (1996) with profile constants from Grachev et al (2000) BLM stable
480      !! form from Beljaars and Holtslag (1991)
481      !!
482      !! Author: L. Brodeau, June 2016 / AeroBulk
483      !!         (https://github.com/brodeau/aerobulk/)
484      !!----------------------------------------------------------------
485      !!
486      REAL(wp), DIMENSION(jpi,jpj) :: psi_h_coare
487      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
488      !
489      INTEGER  ::   ji, jj     ! dummy loop indices
490      REAL(wp) :: zta, zphi_h, zphi_c, zpsi_k, zpsi_c, zf, zc, zstab
491      !
492      DO jj = 1, jpj
493         DO ji = 1, jpi
494            !
495            zta = pzeta(ji,jj)
496            !
497            zphi_h = (ABS(1. - 15.*zta))**.5  !! Kansas unstable   (zphi_h = zphi_m**2 when unstable, zphi_m when stable)
498            !
499            zpsi_k = 2.*LOG((1. + zphi_h)/2.)
500            !
501            zphi_c = (ABS(1. - 34.15*zta))**.3333   !! Convective
502            !
503            zpsi_c = 1.5*LOG((1. + zphi_c + zphi_c*zphi_c)/3.) &
504               &    -1.7320508*ATAN((1. + 2.*zphi_c)/1.7320508) + 1.813799447
505            !
506            zf = zta*zta
507            zf = zf/(1. + zf)
508            zc = MIN(50._wp,0.35_wp*zta)
509            zstab = 0.5 + SIGN(0.5_wp, zta)
510            !
511            psi_h_coare(ji,jj) = (1. - zstab) * ( (1. - zf)*zpsi_k + zf*zpsi_c ) &
512               &                -   zstab     * ( (ABS(1. + 2.*zta/3.))**1.5     &
513               &                           + .6667*(zta - 14.28)/EXP(zc) + 8.525 )
514            !
515         END DO
516      END DO
517      !
518   END FUNCTION psi_h_coare
519
520   !!======================================================================
521END MODULE sbcblk_algo_coare3p6
Note: See TracBrowser for help on using the repository browser.