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_ncar.F90 in branches/2016/dev_r6711_SIMPLIF_6_aerobulk/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2016/dev_r6711_SIMPLIF_6_aerobulk/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_algo_ncar.F90 @ 6723

Last change on this file since 6723 was 6723, checked in by gm, 8 years ago

#1751 - branch SIMPLIF_6_aerobulk: add aerobulk package including NCAR, COARE and ECMWF bulk

File size: 16.4 KB
Line 
1MODULE sbcblk_algo_ncar
2   !!======================================================================
3   !!                   ***  MODULE  sbcblk_algo_ncar  ***
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 Large & Yeager 2008
11   !!
12   !!       Routine turb_ncar maintained and developed in AeroBulk
13   !!                     (http://aerobulk.sourceforge.net/)
14   !!
15   !!                         L. Brodeau, 2015
16   !!=====================================================================
17   !! History :  3.6  !  2016-02  (L.Brodeau) successor of old turb_ncar
18   !!                                         of old sbcblk_core.F90
19   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   turb_ncar  : 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 sbc_oce         ! Surface boundary condition: ocean fields
30   USE sbcwave, ONLY   :  cdn_wave ! wave module
31#if defined key_lim3 || defined key_cice
32   USE sbc_ice         ! Surface boundary condition: ice fields
33#endif
34   !
35   USE iom             ! I/O manager library
36   USE lib_mpp         ! distribued memory computing library
37   USE wrk_nemo        ! work arrays
38   USE timing          ! Timing
39   USE in_out_manager  ! I/O manager
40   USE prtctl          ! Print control
41   USE lib_fortran     ! to use key_nosignedzero
42
43
44   IMPLICIT NONE
45   PRIVATE
46
47   PUBLIC :: TURB_NCAR
48
49   !! NCAR own values for given constants:
50   REAL(wp), PARAMETER ::   rctv0 = 0.608   ! constant to obtain virtual temperature...
51   
52   !!----------------------------------------------------------------------
53CONTAINS
54
55   SUBROUTINE turb_ncar( zt, zu, sst, t_zt, ssq, q_zt, U_zu, &
56      &                  Cd, Ch, Ce, t_zu, q_zu, U_blk )
57      !!----------------------------------------------------------------------------------
58      !!                      ***  ROUTINE  turb_ncar  ***
59      !!
60      !! ** Purpose :   Computes turbulent transfert coefficients of surface
61      !!                fluxes according to Large & Yeager (2004) and Large & Yeager (2008)
62      !!                If relevant (zt /= zu), adjust temperature and humidity from height zt to zu
63      !!                Returns the effective bulk wind speed at 10m to be used in the bulk formulas
64      !!
65      !! ** Method : Monin Obukhov Similarity Theory
66      !!             + Large & Yeager (2004,2008) closure: CD_n10 = f(U_n10)
67      !!
68      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
69      !!
70      !! ** Last update: Laurent Brodeau, June 2014:
71      !!    - handles both cases zt=zu and zt/=zu
72      !!    - optimized: less 2D arrays allocated and less operations
73      !!    - better first guess of stability by checking air-sea difference of virtual temperature
74      !!       rather than temperature difference only...
75      !!    - added function "cd_neutral_10m" that uses the improved parametrization of
76      !!      Large & Yeager 2008. Drag-coefficient reduction for Cyclone conditions!
77      !!    - using code-wide physical constants defined into "phycst.mod" rather than redifining them
78      !!      => 'vkarmn' and 'grav'
79      !!
80      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
81      !!
82      !! INPUT :
83      !! -------
84      !!    *  zt   : height for temperature and spec. hum. of air            [m]
85      !!    *  zu   : height for wind speed (generally 10m)                   [m]
86      !!    *  U_zu : scalar wind speed at 10m                                [m/s]
87      !!    *  sst  : SST                                                     [K]
88      !!    *  t_zt : potential air temperature at zt                         [K]
89      !!    *  ssq  : specific humidity at saturation at SST                  [kg/kg]
90      !!    *  q_zt : specific humidity of air at zt                          [kg/kg]
91      !!
92      !!
93      !! OUTPUT :
94      !! --------
95      !!    *  Cd     : drag coefficient
96      !!    *  Ch     : sensible heat coefficient
97      !!    *  Ce     : evaporation coefficient
98      !!    *  t_zu   : pot. air temperature adjusted at wind height zu       [K]
99      !!    *  q_zu   : specific humidity of air        //                    [kg/kg]
100      !!    *  U_blk  : bulk wind at 10m                                      [m/s]
101      !!----------------------------------------------------------------------------------
102      REAL(wp), INTENT(in   )                     ::   zt       ! height for t_zt and q_zt                    [m]
103      REAL(wp), INTENT(in   )                     ::   zu       ! height for U_zu                             [m]
104      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   sst      ! sea surface temperature                [Kelvin]
105      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   t_zt     ! potential air temperature              [Kelvin]
106      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   ssq      ! sea surface specific humidity           [kg/kg]
107      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity                   [kg/kg]
108      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   U_zu     ! relative wind module at zu                [m/s]
109      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
110      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
111      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
112      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   t_zu     ! pot. air temp. adjusted at zu               [K]
113      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. humidity adjusted at zu           [kg/kg]
114      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   U_blk    ! bulk wind at 10m                          [m/s]
115      !
116      INTEGER ::   j_itt
117      LOGICAL ::   l_zt_equal_zu = .FALSE.      ! if q and t are given at same height as U
118      INTEGER , PARAMETER ::   nb_itt = 4       ! number of itterations
119      !
120      REAL(wp), DIMENSION(:,:), POINTER ::   Cx_n10        ! 10m neutral latent/sensible coefficient
121      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd_n10   ! root square of Cd_n10
122      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_u        ! stability parameter at height zu
123      REAL(wp), DIMENSION(:,:), POINTER ::   zpsi_h_u
124      REAL(wp), DIMENSION(:,:), POINTER ::   ztmp0, ztmp1, ztmp2
125      REAL(wp), DIMENSION(:,:), POINTER ::   stab          ! stability test integer
126      !!----------------------------------------------------------------------------------
127      !
128      IF( nn_timing == 1 )   CALL timing_start('turb_ncar')
129      !
130      CALL wrk_alloc( jpi,jpj,   Cx_n10, sqrt_Cd_n10, zeta_u, stab )
131      CALL wrk_alloc( jpi,jpj,   zpsi_h_u, ztmp0, ztmp1, ztmp2 )
132      !
133      l_zt_equal_zu = .FALSE.
134      IF( ABS(zu - zt) < 0.01 ) l_zt_equal_zu = .TRUE.    ! testing "zu == zt" is risky with double precision
135
136      U_blk = MAX( 0.5 , U_zu )   !  relative wind speed at zu (normally 10m), we don't want to fall under 0.5 m/s
137
138      !! First guess of stability:
139      ztmp0 = t_zt*(1. + rctv0*q_zt) - sst*(1. + rctv0*ssq) ! air-sea difference of virtual pot. temp. at zt
140      stab  = 0.5 + sign(0.5,ztmp0)                           ! stab = 1 if dTv > 0  => STABLE, 0 if unstable
141
142      !! Neutral coefficients at 10m:
143      IF( ln_cdgw ) THEN      ! wave drag case
144         cdn_wave(:,:) = cdn_wave(:,:) + rsmall * ( 1._wp - tmask(:,:,1) )
145         ztmp0   (:,:) = cdn_wave(:,:)
146      ELSE
147         ztmp0 = cd_neutral_10m( U_blk )
148      ENDIF
149
150      sqrt_Cd_n10 = SQRT( ztmp0 )
151
152      !! Initializing transf. coeff. with their first guess neutral equivalents :
153      Cd = ztmp0
154      Ce = 1.e-3*( 34.6 * sqrt_Cd_n10 )
155      Ch = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab))
156      stab = sqrt_Cd_n10   ! Temporaty array !!! stab == SQRT(Cd)
157
158      !! Initializing values at z_u with z_t values:
159      t_zu = t_zt   ;   q_zu = q_zt
160
161      !!  * Now starting iteration loop
162      DO j_itt=1, nb_itt
163         !
164         ztmp1 = t_zu - sst   ! Updating air/sea differences
165         ztmp2 = q_zu - ssq
166
167         ! Updating turbulent scales :   (L&Y 2004 eq. (7))
168         ztmp1  = Ch/stab*ztmp1    ! theta*   (stab == SQRT(Cd))
169         ztmp2  = Ce/stab*ztmp2    ! q*       (stab == SQRT(Cd))
170
171         ztmp0 = 1. + rctv0*q_zu      ! multiply this with t and you have the virtual temperature
172
173         ! Estimate the inverse of Monin-Obukov length (1/L) at height zu:
174         ztmp0 =  (grav*vkarmn/(t_zu*ztmp0)*(ztmp1*ztmp0 + rctv0*t_zu*ztmp2)) / (Cd*U_blk*U_blk)
175         !                                                      ( Cd*U_blk*U_blk is U*^2 at zu )
176
177         !! Stability parameters :
178         zeta_u   = zu*ztmp0   ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
179         zpsi_h_u = psi_h( zeta_u )
180
181         !! Shifting temperature and humidity at zu (L&Y 2004 eq. (9b-9c))
182         IF( .NOT. l_zt_equal_zu ) THEN
183            !! Array 'stab' is free for the moment so using it to store 'zeta_t'
184            stab = zt*ztmp0 ;  stab = SIGN( MIN(ABS(stab),10.0), stab )  ! Temporaty array stab == zeta_t !!!
185            stab = LOG(zt/zu) + zpsi_h_u - psi_h(stab)                   ! stab just used as temp array again!
186            t_zu = t_zt - ztmp1/vkarmn*stab    ! ztmp1 is still theta*  L&Y 2004 eq.(9b)
187            q_zu = q_zt - ztmp2/vkarmn*stab    ! ztmp2 is still q*      L&Y 2004 eq.(9c)
188            q_zu = max(0., q_zu)
189         END IF
190
191         ztmp2 = psi_m(zeta_u)
192         IF( ln_cdgw ) THEN      ! surface wave case
193            stab = vkarmn / ( vkarmn / sqrt_Cd_n10 - ztmp2 )  ! (stab == SQRT(Cd))
194            Cd      = stab * stab
195         ELSE
196            ! Update neutral wind speed at 10m and neutral Cd at 10m (L&Y 2004 eq. 9a)...
197            !   In very rare low-wind conditions, the old way of estimating the
198            !   neutral wind speed at 10m leads to a negative value that causes the code
199            !   to crash. To prevent this a threshold of 0.25m/s is imposed.
200            ztmp0 = MAX( 0.25 , U_blk/(1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - ztmp2)) ) ! U_n10 (ztmp2 == psi_m(zeta_u))
201            ztmp0 = cd_neutral_10m(ztmp0)                                               ! Cd_n10
202            sqrt_Cd_n10 = sqrt(ztmp0)
203
204            stab    = 0.5 + sign(0.5,zeta_u)                           ! update stability
205            Cx_n10  = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab))  ! L&Y 2004 eq. (6c-6d)    (Cx_n10 == Ch_n10)
206
207            !! Update of transfer coefficients:
208            ztmp1 = 1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - ztmp2)   ! L&Y 2004 eq. (10a) (ztmp2 == psi_m(zeta_u))
209            Cd      = ztmp0 / ( ztmp1*ztmp1 )
210            stab = SQRT( Cd ) ! Temporary array !!! (stab == SQRT(Cd))
211         ENDIF
212
213         ztmp0 = (LOG(zu/10.) - zpsi_h_u) / vkarmn / sqrt_Cd_n10
214         ztmp2 = stab / sqrt_Cd_n10   ! (stab == SQRT(Cd))
215         ztmp1 = 1. + Cx_n10*ztmp0    ! (Cx_n10 == Ch_n10)
216         Ch  = Cx_n10*ztmp2 / ztmp1   ! L&Y 2004 eq. (10b)
217
218         Cx_n10  = 1.e-3 * (34.6 * sqrt_Cd_n10)  ! L&Y 2004 eq. (6b)    ! Cx_n10 == Ce_n10
219         ztmp1 = 1. + Cx_n10*ztmp0
220         Ce  = Cx_n10*ztmp2 / ztmp1  ! L&Y 2004 eq. (10c)
221
222      END DO
223
224      CALL wrk_dealloc( jpi,jpj,   Cx_n10, sqrt_Cd_n10, zeta_u, stab )
225      CALL wrk_dealloc( jpi,jpj,   zpsi_h_u, ztmp0, ztmp1, ztmp2 )
226
227      IF( nn_timing == 1 )   CALL timing_stop('turb_ncar')
228
229   END SUBROUTINE turb_ncar
230
231
232   FUNCTION cd_neutral_10m( pw10 )
233      !!----------------------------------------------------------------------------------     
234      !! Estimate of the neutral drag coefficient at 10m as a function
235      !! of neutral wind  speed at 10m
236      !!
237      !! Origin: Large & Yeager 2008 eq.(11a) and eq.(11b)
238      !!
239      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
240      !!----------------------------------------------------------------------------------
241      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pw10           ! scalar wind speed at 10m (m/s)
242      REAL(wp), DIMENSION(jpi,jpj)             :: cd_neutral_10m
243      !
244      INTEGER  ::     ji, jj     ! dummy loop indices
245      REAL(wp) :: zgt33, zw, zw6 ! local scalars
246      !!----------------------------------------------------------------------------------
247      !
248      DO jj = 1, jpj
249         DO ji = 1, jpi
250            !
251            zw  = pw10(ji,jj)
252            zw6 = zw*zw*zw
253            zw6 = zw6*zw6
254            !
255            ! When wind speed > 33 m/s => Cyclone conditions => special treatment
256            zgt33 = 0.5 + SIGN( 0.5, (zw - 33.) )   ! If pw10 < 33. => 0, else => 1
257            !
258            cd_neutral_10m(ji,jj) = 1.e-3 * ( &
259               &       (1. - zgt33)*( 2.7/zw + 0.142 + zw/13.09 - 3.14807E-10*zw6) & ! wind <  33 m/s
260               &      +    zgt33   *      2.34 )                                     ! wind >= 33 m/s
261            !
262            cd_neutral_10m(ji,jj) = MAX(cd_neutral_10m(ji,jj), 1.E-6)
263            !
264         END DO
265      END DO
266      !
267   END FUNCTION cd_neutral_10m
268
269
270   FUNCTION psi_m( pzeta )
271      !!----------------------------------------------------------------------------------
272      !! Universal profile stability function for momentum
273      !!    !! Psis, L&Y 2004 eq. (8c), (8d), (8e)
274      !!     
275      !! pzet0 : stability paramenter, z/L where z is altitude measurement                                         
276      !!         and L is M-O length
277      !!
278      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
279      !!----------------------------------------------------------------------------------
280      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pzeta
281      REAL(wp), DIMENSION(jpi,jpj)             ::   psi_m
282      !
283      INTEGER  ::   ji, jj         ! dummy loop indices
284      REAL(wp) :: zx2, zx, zstab   ! local scalars
285      !!----------------------------------------------------------------------------------
286      !
287      DO jj = 1, jpj
288         DO ji = 1, jpi
289            zx2 = SQRT( ABS( 1. - 16.*pzeta(ji,jj) ) )
290            zx2 = MAX ( zx2 , 1. )
291            zx  = SQRT( zx2 )
292            zstab = 0.5 + SIGN( 0.5 , pzeta(ji,jj) )
293            !
294            psi_m(ji,jj) =        zstab  * (-5.*pzeta(ji,jj))       &          ! Stable
295               &          + (1. - zstab) * (2.*LOG((1. + zx)*0.5)   &          ! Unstable
296               &               + LOG((1. + zx2)*0.5) - 2.*ATAN(zx) + rpi*0.5)  !    "
297            !
298         END DO
299      END DO
300      !
301   END FUNCTION psi_m
302
303
304   FUNCTION psi_h( pzeta )
305      !!----------------------------------------------------------------------------------
306      !! Universal profile stability function for temperature and humidity
307      !!    !! Psis, L&Y 2004 eq. (8c), (8d), (8e)
308      !!
309      !! pzet0 : stability paramenter, z/L where z is altitude measurement                                         
310      !!         and L is M-O length
311      !!
312      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
313      !!----------------------------------------------------------------------------------
314      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pzeta
315      REAL(wp), DIMENSION(jpi,jpj)             :: psi_h
316      !
317      INTEGER  ::   ji, jj    ! dummy loop indices
318      REAL(wp) :: zx2, zstab  ! local scalars
319      !!----------------------------------------------------------------------------------
320      !
321      DO jj = 1, jpj
322         DO ji = 1, jpi
323            zx2 = SQRT( ABS( 1. - 16.*pzeta(ji,jj) ) )
324            zx2 = MAX ( zx2 , 1. )
325            zstab = 0.5 + SIGN( 0.5 , pzeta(ji,jj) )
326            !
327            psi_h(ji,jj) =         zstab  * (-5.*pzeta(ji,jj))        &  ! Stable
328               &           + (1. - zstab) * (2.*LOG( (1. + zx2)*0.5 ))   ! Unstable
329            !
330         END DO
331      END DO
332      !
333   END FUNCTION psi_h
334
335   !!======================================================================
336END MODULE sbcblk_algo_ncar
Note: See TracBrowser for help on using the repository browser.