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 NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/SBC – NEMO

source: NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/SBC/sbcblk_algo_ncar.F90 @ 10345

Last change on this file since 10345 was 10345, checked in by smasson, 5 years ago

dev_r10164_HPC09_ESIWACE_PREP_MERGE: merge with trunk@10344, see #2133

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