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.F90 in NEMO/branches/2020/tickets_icb_1900/src/OCE/SBC – NEMO

source: NEMO/branches/2020/tickets_icb_1900/src/OCE/SBC/sbcblk.F90 @ 14016

Last change on this file since 14016 was 14016, checked in by mathiot, 3 years ago

ticket 1900: upgrade to trunk@r14015 (head trunk at 16h27)

  • Property svn:keywords set to Id
File size: 83.6 KB
Line 
1MODULE sbcblk
2   !!======================================================================
3   !!                       ***  MODULE  sbcblk  ***
4   !! Ocean forcing:  momentum, heat and freshwater flux formulation
5   !!                         Aerodynamic Bulk Formulas
6   !!                        SUCCESSOR OF "sbcblk_core"
7   !!=====================================================================
8   !! History :  1.0  !  2004-08  (U. Schweckendiek)  Original CORE code
9   !!            2.0  !  2005-04  (L. Brodeau, A.M. Treguier)  improved CORE bulk and its user interface
10   !!            3.0  !  2006-06  (G. Madec)  sbc rewritting
11   !!             -   !  2006-12  (L. Brodeau)  Original code for turb_core
12   !!            3.2  !  2009-04  (B. Lemaire)  Introduce iom_put
13   !!            3.3  !  2010-10  (S. Masson)  add diurnal cycle
14   !!            3.4  !  2011-11  (C. Harris)  Fill arrays required by CICE
15   !!            3.7  !  2014-06  (L. Brodeau)  simplification and optimization of CORE bulk
16   !!            4.0  !  2016-06  (L. Brodeau)  sbcblk_core becomes sbcblk and is not restricted to the CORE algorithm anymore
17   !!                 !                        ==> based on AeroBulk (https://github.com/brodeau/aerobulk/)
18   !!            4.0  !  2016-10  (G. Madec)  introduce a sbc_blk_init routine
19   !!            4.0  !  2016-10  (M. Vancoppenolle)  Introduce conduction flux emulator (M. Vancoppenolle)
20   !!            4.0  !  2019-03  (F. Lemarié & G. Samson)  add ABL compatibility (ln_abl=TRUE)
21   !!----------------------------------------------------------------------
22
23   !!----------------------------------------------------------------------
24   !!   sbc_blk_init  : initialisation of the chosen bulk formulation as ocean surface boundary condition
25   !!   sbc_blk       : bulk formulation as ocean surface boundary condition
26   !!   blk_oce_1     : computes pieces of momentum, heat and freshwater fluxes over ocean for ABL model  (ln_abl=TRUE)
27   !!   blk_oce_2     : finalizes momentum, heat and freshwater fluxes computation over ocean after the ABL step  (ln_abl=TRUE)
28   !!             sea-ice case only :
29   !!   blk_ice_1   : provide the air-ice stress
30   !!   blk_ice_2   : provide the heat and mass fluxes at air-ice interface
31   !!   blk_ice_qcn   : provide ice surface temperature and snow/ice conduction flux (emulating conduction flux)
32   !!   Cdn10_Lupkes2012 : Lupkes et al. (2012) air-ice drag
33   !!   Cdn10_Lupkes2015 : Lupkes et al. (2015) air-ice drag
34   !!----------------------------------------------------------------------
35   USE oce            ! ocean dynamics and tracers
36   USE dom_oce        ! ocean space and time domain
37   USE phycst         ! physical constants
38   USE fldread        ! read input fields
39   USE sbc_oce        ! Surface boundary condition: ocean fields
40   USE cyclone        ! Cyclone 10m wind form trac of cyclone centres
41   USE sbcdcy         ! surface boundary condition: diurnal cycle
42   USE sbcwave , ONLY :   cdn_wave ! wave module
43   USE sbc_ice        ! Surface boundary condition: ice fields
44   USE lib_fortran    ! to use key_nosignedzero
45#if defined key_si3
46   USE ice     , ONLY :   u_ice, v_ice, jpl, a_i_b, at_i_b, t_su, rn_cnd_s, hfx_err_dif, nn_qtrice
47   USE icevar         ! for CALL ice_var_snwblow
48#endif
49   USE sbcblk_algo_ncar     ! => turb_ncar     : NCAR - CORE (Large & Yeager, 2009)
50   USE sbcblk_algo_coare3p0 ! => turb_coare3p0 : COAREv3.0 (Fairall et al. 2003)
51   USE sbcblk_algo_coare3p6 ! => turb_coare3p6 : COAREv3.6 (Fairall et al. 2018 + Edson et al. 2013)
52   USE sbcblk_algo_ecmwf    ! => turb_ecmwf    : ECMWF (IFS cycle 45r1)
53   !
54   USE iom            ! I/O manager library
55   USE in_out_manager ! I/O manager
56   USE lib_mpp        ! distribued memory computing library
57   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
58   USE prtctl         ! Print control
59
60   USE sbcblk_phy     ! a catalog of functions for physical/meteorological parameters in the marine boundary layer, rho_air, q_sat, etc...
61
62
63   IMPLICIT NONE
64   PRIVATE
65
66   PUBLIC   sbc_blk_init  ! called in sbcmod
67   PUBLIC   sbc_blk       ! called in sbcmod
68   PUBLIC   blk_oce_1     ! called in sbcabl
69   PUBLIC   blk_oce_2     ! called in sbcabl
70#if defined key_si3
71   PUBLIC   blk_ice_1     ! routine called in icesbc
72   PUBLIC   blk_ice_2     ! routine called in icesbc
73   PUBLIC   blk_ice_qcn   ! routine called in icesbc
74#endif
75
76   INTEGER , PUBLIC, PARAMETER ::   jp_wndi  =  1   ! index of 10m wind velocity (i-component) (m/s)    at T-point
77   INTEGER , PUBLIC, PARAMETER ::   jp_wndj  =  2   ! index of 10m wind velocity (j-component) (m/s)    at T-point
78   INTEGER , PUBLIC, PARAMETER ::   jp_tair  =  3   ! index of 10m air temperature             (Kelvin)
79   INTEGER , PUBLIC, PARAMETER ::   jp_humi  =  4   ! index of specific humidity               ( % )
80   INTEGER , PUBLIC, PARAMETER ::   jp_qsr   =  5   ! index of solar heat                      (W/m2)
81   INTEGER , PUBLIC, PARAMETER ::   jp_qlw   =  6   ! index of Long wave                       (W/m2)
82   INTEGER , PUBLIC, PARAMETER ::   jp_prec  =  7   ! index of total precipitation (rain+snow) (Kg/m2/s)
83   INTEGER , PUBLIC, PARAMETER ::   jp_snow  =  8   ! index of snow (solid prcipitation)       (kg/m2/s)
84   INTEGER , PUBLIC, PARAMETER ::   jp_slp   =  9   ! index of sea level pressure              (Pa)
85   INTEGER , PUBLIC, PARAMETER ::   jp_uoatm = 10   ! index of surface current (i-component)
86   !                                                !          seen by the atmospheric forcing (m/s) at T-point
87   INTEGER , PUBLIC, PARAMETER ::   jp_voatm = 11   ! index of surface current (j-component)
88   !                                                !          seen by the atmospheric forcing (m/s) at T-point
89   INTEGER , PUBLIC, PARAMETER ::   jp_cc    = 12   ! index of cloud cover                     (-)      range:0-1
90   INTEGER , PUBLIC, PARAMETER ::   jp_hpgi  = 13   ! index of ABL geostrophic wind or hpg (i-component) (m/s) at T-point
91   INTEGER , PUBLIC, PARAMETER ::   jp_hpgj  = 14   ! index of ABL geostrophic wind or hpg (j-component) (m/s) at T-point
92   INTEGER , PUBLIC, PARAMETER ::   jpfld    = 14   ! maximum number of files to read
93
94   ! Warning: keep this structure allocatable for Agrif...
95   TYPE(FLD), PUBLIC, ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input atmospheric fields (file informations, fields read)
96
97   !                           !!* Namelist namsbc_blk : bulk parameters
98   LOGICAL  ::   ln_NCAR        ! "NCAR"      algorithm   (Large and Yeager 2008)
99   LOGICAL  ::   ln_COARE_3p0   ! "COARE 3.0" algorithm   (Fairall et al. 2003)
100   LOGICAL  ::   ln_COARE_3p6   ! "COARE 3.6" algorithm   (Edson et al. 2013)
101   LOGICAL  ::   ln_ECMWF       ! "ECMWF"     algorithm   (IFS cycle 45r1)
102   !
103   LOGICAL  ::   ln_Cd_L12      ! ice-atm drag = F( ice concentration )                        (Lupkes et al. JGR2012)
104   LOGICAL  ::   ln_Cd_L15      ! ice-atm drag = F( ice concentration, atmospheric stability ) (Lupkes et al. JGR2015)
105   !
106   LOGICAL  ::   ln_crt_fbk     ! Add surface current feedback to the wind stress computation  (Renault et al. 2020)
107   REAL(wp) ::   rn_stau_a      ! Alpha and Beta coefficients of Renault et al. 2020, eq. 10: Stau = Alpha * Wnd + Beta
108   REAL(wp) ::   rn_stau_b      !
109   !
110   REAL(wp)         ::   rn_pfac   ! multiplication factor for precipitation
111   REAL(wp), PUBLIC ::   rn_efac   ! multiplication factor for evaporation
112   REAL(wp)         ::   rn_zqt    ! z(q,t) : height of humidity and temperature measurements
113   REAL(wp)         ::   rn_zu     ! z(u)   : height of wind measurements
114   !
115   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) ::   Cdn_oce, Chn_oce, Cen_oce  ! neutral coeffs over ocean (L15 bulk scheme and ABL)
116   REAL(wp),         ALLOCATABLE, DIMENSION(:,:) ::   Cd_ice , Ch_ice , Ce_ice   ! transfert coefficients over ice
117   REAL(wp),         ALLOCATABLE, DIMENSION(:,:) ::   t_zu, q_zu                 ! air temp. and spec. hum. at wind speed height (L15 bulk scheme)
118
119   LOGICAL  ::   ln_skin_cs     ! use the cool-skin (only available in ECMWF and COARE algorithms) !LB
120   LOGICAL  ::   ln_skin_wl     ! use the warm-layer parameterization (only available in ECMWF and COARE algorithms) !LB
121   LOGICAL  ::   ln_humi_sph    ! humidity read in files ("sn_humi") is specific humidity [kg/kg] if .true. !LB
122   LOGICAL  ::   ln_humi_dpt    ! humidity read in files ("sn_humi") is dew-point temperature [K] if .true. !LB
123   LOGICAL  ::   ln_humi_rlh    ! humidity read in files ("sn_humi") is relative humidity     [%] if .true. !LB
124   LOGICAL  ::   ln_tpot        !!GS: flag to compute or not potential temperature
125   !
126   INTEGER  ::   nhumi          ! choice of the bulk algorithm
127   !                            ! associated indices:
128   INTEGER, PARAMETER :: np_humi_sph = 1
129   INTEGER, PARAMETER :: np_humi_dpt = 2
130   INTEGER, PARAMETER :: np_humi_rlh = 3
131
132   INTEGER  ::   nblk           ! choice of the bulk algorithm
133   !                            ! associated indices:
134   INTEGER, PARAMETER ::   np_NCAR      = 1   ! "NCAR" algorithm        (Large and Yeager 2008)
135   INTEGER, PARAMETER ::   np_COARE_3p0 = 2   ! "COARE 3.0" algorithm   (Fairall et al. 2003)
136   INTEGER, PARAMETER ::   np_COARE_3p6 = 3   ! "COARE 3.6" algorithm   (Edson et al. 2013)
137   INTEGER, PARAMETER ::   np_ECMWF     = 4   ! "ECMWF" algorithm       (IFS cycle 45r1)
138
139   !! * Substitutions
140#  include "do_loop_substitute.h90"
141   !!----------------------------------------------------------------------
142   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
143   !! $Id$
144   !! Software governed by the CeCILL license (see ./LICENSE)
145   !!----------------------------------------------------------------------
146CONTAINS
147
148   INTEGER FUNCTION sbc_blk_alloc()
149      !!-------------------------------------------------------------------
150      !!             ***  ROUTINE sbc_blk_alloc ***
151      !!-------------------------------------------------------------------
152      ALLOCATE( t_zu(jpi,jpj)   , q_zu(jpi,jpj)   ,                                      &
153         &      Cdn_oce(jpi,jpj), Chn_oce(jpi,jpj), Cen_oce(jpi,jpj),                    &
154         &      Cd_ice (jpi,jpj), Ch_ice (jpi,jpj), Ce_ice (jpi,jpj), STAT=sbc_blk_alloc )
155      !
156      CALL mpp_sum ( 'sbcblk', sbc_blk_alloc )
157      IF( sbc_blk_alloc /= 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_alloc: failed to allocate arrays' )
158   END FUNCTION sbc_blk_alloc
159
160
161   SUBROUTINE sbc_blk_init
162      !!---------------------------------------------------------------------
163      !!                    ***  ROUTINE sbc_blk_init  ***
164      !!
165      !! ** Purpose :   choose and initialize a bulk formulae formulation
166      !!
167      !! ** Method  :
168      !!
169      !!----------------------------------------------------------------------
170      INTEGER  ::   jfpr                  ! dummy loop indice and argument
171      INTEGER  ::   ios, ierror, ioptio   ! Local integer
172      !!
173      CHARACTER(len=100)            ::   cn_dir                ! Root directory for location of atmospheric forcing files
174      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i                 ! array of namelist informations on the fields to read
175      TYPE(FLD_N) ::   sn_wndi, sn_wndj , sn_humi, sn_qsr      ! informations about the fields to be read
176      TYPE(FLD_N) ::   sn_qlw , sn_tair , sn_prec, sn_snow     !       "                        "
177      TYPE(FLD_N) ::   sn_slp , sn_uoatm, sn_voatm             !       "                        "
178      TYPE(FLD_N) ::   sn_cc, sn_hpgi, sn_hpgj                 !       "                        "
179      INTEGER     ::   ipka                                    ! number of levels in the atmospheric variable
180      NAMELIST/namsbc_blk/ sn_wndi, sn_wndj, sn_humi, sn_qsr, sn_qlw ,                &   ! input fields
181         &                 sn_tair, sn_prec, sn_snow, sn_slp, sn_uoatm, sn_voatm,     &
182         &                 sn_cc, sn_hpgi, sn_hpgj,                                   &
183         &                 ln_NCAR, ln_COARE_3p0, ln_COARE_3p6, ln_ECMWF,             &   ! bulk algorithm
184         &                 cn_dir , rn_zqt, rn_zu,                                    &
185         &                 rn_pfac, rn_efac, ln_Cd_L12, ln_Cd_L15, ln_tpot,           &
186         &                 ln_crt_fbk, rn_stau_a, rn_stau_b,                          &   ! current feedback
187         &                 ln_skin_cs, ln_skin_wl, ln_humi_sph, ln_humi_dpt, ln_humi_rlh  ! cool-skin / warm-layer !LB
188      !!---------------------------------------------------------------------
189      !
190      !                                      ! allocate sbc_blk_core array
191      IF( sbc_blk_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'sbc_blk : unable to allocate standard arrays' )
192      !
193      !                             !** read bulk namelist
194      READ  ( numnam_ref, namsbc_blk, IOSTAT = ios, ERR = 901)
195901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namsbc_blk in reference namelist' )
196      !
197      READ  ( numnam_cfg, namsbc_blk, IOSTAT = ios, ERR = 902 )
198902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namsbc_blk in configuration namelist' )
199      !
200      IF(lwm) WRITE( numond, namsbc_blk )
201      !
202      !                             !** initialization of the chosen bulk formulae (+ check)
203      !                                   !* select the bulk chosen in the namelist and check the choice
204      ioptio = 0
205      IF( ln_NCAR      ) THEN
206         nblk =  np_NCAR        ;   ioptio = ioptio + 1
207      ENDIF
208      IF( ln_COARE_3p0 ) THEN
209         nblk =  np_COARE_3p0   ;   ioptio = ioptio + 1
210      ENDIF
211      IF( ln_COARE_3p6 ) THEN
212         nblk =  np_COARE_3p6   ;   ioptio = ioptio + 1
213      ENDIF
214      IF( ln_ECMWF     ) THEN
215         nblk =  np_ECMWF       ;   ioptio = ioptio + 1
216      ENDIF
217      IF( ioptio /= 1 )   CALL ctl_stop( 'sbc_blk_init: Choose one and only one bulk algorithm' )
218
219      !                             !** initialization of the cool-skin / warm-layer parametrization
220      IF( ln_skin_cs .OR. ln_skin_wl ) THEN
221         !! Some namelist sanity tests:
222         IF( ln_NCAR )      &
223            & CALL ctl_stop( 'sbc_blk_init: Cool-skin/warm-layer param. not compatible with NCAR algorithm' )
224         IF( nn_fsbc /= 1 ) &
225            & CALL ctl_stop( 'sbc_blk_init: Please set "nn_fsbc" to 1 when using cool-skin/warm-layer param.')
226      END IF
227
228      IF( ln_skin_wl ) THEN
229         !! Check if the frequency of downwelling solar flux input makes sense and if ln_dm2dc=T if it is daily!
230         IF( (sn_qsr%freqh  < 0.).OR.(sn_qsr%freqh  > 24.) ) &
231            & CALL ctl_stop( 'sbc_blk_init: Warm-layer param. (ln_skin_wl) not compatible with freq. of solar flux > daily' )
232         IF( (sn_qsr%freqh == 24.).AND.(.NOT. ln_dm2dc) ) &
233            & CALL ctl_stop( 'sbc_blk_init: Please set ln_dm2dc=T for warm-layer param. (ln_skin_wl) to work properly' )
234      END IF
235
236      ioptio = 0
237      IF( ln_humi_sph ) THEN
238         nhumi =  np_humi_sph    ;   ioptio = ioptio + 1
239      ENDIF
240      IF( ln_humi_dpt ) THEN
241         nhumi =  np_humi_dpt    ;   ioptio = ioptio + 1
242      ENDIF
243      IF( ln_humi_rlh ) THEN
244         nhumi =  np_humi_rlh    ;   ioptio = ioptio + 1
245      ENDIF
246      IF( ioptio /= 1 )   CALL ctl_stop( 'sbc_blk_init: Choose one and only one type of air humidity' )
247      !
248      IF( ln_dm2dc ) THEN                 !* check: diurnal cycle on Qsr
249         IF( sn_qsr%freqh /= 24. )   CALL ctl_stop( 'sbc_blk_init: ln_dm2dc=T only with daily short-wave input' )
250         IF( sn_qsr%ln_tint ) THEN
251            CALL ctl_warn( 'sbc_blk_init: ln_dm2dc=T daily qsr time interpolation done by sbcdcy module',   &
252               &           '              ==> We force time interpolation = .false. for qsr' )
253            sn_qsr%ln_tint = .false.
254         ENDIF
255      ENDIF
256      !                                   !* set the bulk structure
257      !                                      !- store namelist information in an array
258      !
259      slf_i(jp_wndi ) = sn_wndi    ;   slf_i(jp_wndj ) = sn_wndj
260      slf_i(jp_qsr  ) = sn_qsr     ;   slf_i(jp_qlw  ) = sn_qlw
261      slf_i(jp_tair ) = sn_tair    ;   slf_i(jp_humi ) = sn_humi
262      slf_i(jp_prec ) = sn_prec    ;   slf_i(jp_snow ) = sn_snow
263      slf_i(jp_slp  ) = sn_slp     ;   slf_i(jp_cc   ) = sn_cc
264      slf_i(jp_uoatm) = sn_uoatm   ;   slf_i(jp_voatm) = sn_voatm
265      slf_i(jp_hpgi ) = sn_hpgi    ;   slf_i(jp_hpgj ) = sn_hpgj
266      !
267      IF( .NOT. ln_abl ) THEN   ! force to not use jp_hpgi and jp_hpgj, should already be done in namelist_* but we never know...
268         slf_i(jp_hpgi)%clname = 'NOT USED'
269         slf_i(jp_hpgj)%clname = 'NOT USED'
270      ENDIF
271      !
272      !                                      !- allocate the bulk structure
273      ALLOCATE( sf(jpfld), STAT=ierror )
274      IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_init: unable to allocate sf structure' )
275      !
276      !                                      !- fill the bulk structure with namelist informations
277      CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_init', 'surface boundary condition -- bulk formulae', 'namsbc_blk' )
278      !
279      DO jfpr= 1, jpfld
280         !
281         IF(   ln_abl    .AND.                                                      &
282            &    ( jfpr == jp_wndi .OR. jfpr == jp_wndj .OR. jfpr == jp_humi .OR.   &
283            &      jfpr == jp_hpgi .OR. jfpr == jp_hpgj .OR. jfpr == jp_tair     )  ) THEN
284            ipka = jpka   ! ABL: some fields are 3D input
285         ELSE
286            ipka = 1
287         ENDIF
288         !
289         ALLOCATE( sf(jfpr)%fnow(jpi,jpj,ipka) )
290         !
291         IF( TRIM(sf(jfpr)%clrootname) == 'NOT USED' ) THEN    !--  not used field  --!   (only now allocated and set to default)
292            IF(     jfpr == jp_slp ) THEN
293               sf(jfpr)%fnow(:,:,1:ipka) = 101325._wp   ! use standard pressure in Pa
294            ELSEIF( jfpr == jp_prec .OR. jfpr == jp_snow .OR. jfpr == jp_uoatm .OR. jfpr == jp_voatm ) THEN
295               sf(jfpr)%fnow(:,:,1:ipka) = 0._wp        ! no precip or no snow or no surface currents
296            ELSEIF( jfpr == jp_hpgi .OR. jfpr == jp_hpgj ) THEN
297               IF( .NOT. ln_abl ) THEN
298                  DEALLOCATE( sf(jfpr)%fnow )   ! deallocate as not used in this case
299               ELSE
300                  sf(jfpr)%fnow(:,:,1:ipka) = 0._wp
301               ENDIF
302            ELSEIF( jfpr == jp_cc  ) THEN
303               sf(jp_cc)%fnow(:,:,1:ipka) = pp_cldf
304            ELSE
305               WRITE(ctmp1,*) 'sbc_blk_init: no default value defined for field number', jfpr
306               CALL ctl_stop( ctmp1 )
307            ENDIF
308         ELSE                                                  !-- used field  --!
309            IF( sf(jfpr)%ln_tint )   ALLOCATE( sf(jfpr)%fdta(jpi,jpj,ipka,2) )   ! allocate array for temporal interpolation
310            !
311            IF( sf(jfpr)%freqh > 0. .AND. MOD( NINT(3600. * sf(jfpr)%freqh), nn_fsbc * NINT(rn_Dt) ) /= 0 )   &
312         &  CALL ctl_warn( 'sbc_blk_init: sbcmod timestep rn_Dt*nn_fsbc is NOT a submultiple of atmospheric forcing frequency.',   &
313         &                 '               This is not ideal. You should consider changing either rn_Dt or nn_fsbc value...' )
314         ENDIF
315      END DO
316      !
317      IF( ln_abl ) THEN       ! ABL: read 3D fields for wind, temperature, humidity and pressure gradient
318         rn_zqt = ght_abl(2)          ! set the bulk altitude to ABL first level
319         rn_zu  = ght_abl(2)
320         IF(lwp) WRITE(numout,*)
321         IF(lwp) WRITE(numout,*) '   ABL formulation: overwrite rn_zqt & rn_zu with ABL first level altitude'
322      ENDIF
323      !
324      ! set transfer coefficients to default sea-ice values
325      Cd_ice(:,:) = rCd_ice
326      Ch_ice(:,:) = rCd_ice
327      Ce_ice(:,:) = rCd_ice
328      !
329      IF(lwp) THEN                     !** Control print
330         !
331         WRITE(numout,*)                  !* namelist
332         WRITE(numout,*) '   Namelist namsbc_blk (other than data information):'
333         WRITE(numout,*) '      "NCAR"      algorithm   (Large and Yeager 2008)     ln_NCAR      = ', ln_NCAR
334         WRITE(numout,*) '      "COARE 3.0" algorithm   (Fairall et al. 2003)       ln_COARE_3p0 = ', ln_COARE_3p0
335         WRITE(numout,*) '      "COARE 3.6" algorithm (Fairall 2018 + Edson al 2013)ln_COARE_3p6 = ', ln_COARE_3p6
336         WRITE(numout,*) '      "ECMWF"     algorithm   (IFS cycle 45r1)            ln_ECMWF     = ', ln_ECMWF
337         WRITE(numout,*) '      Air temperature and humidity reference height (m)   rn_zqt       = ', rn_zqt
338         WRITE(numout,*) '      Wind vector reference height (m)                    rn_zu        = ', rn_zu
339         WRITE(numout,*) '      factor applied on precipitation (total & snow)      rn_pfac      = ', rn_pfac
340         WRITE(numout,*) '      factor applied on evaporation                       rn_efac      = ', rn_efac
341         WRITE(numout,*) '         (form absolute (=0) to relative winds(=1))'
342         WRITE(numout,*) '      use ice-atm drag from Lupkes2012                    ln_Cd_L12    = ', ln_Cd_L12
343         WRITE(numout,*) '      use ice-atm drag from Lupkes2015                    ln_Cd_L15    = ', ln_Cd_L15
344         WRITE(numout,*) '      use surface current feedback on wind stress         ln_crt_fbk   = ', ln_crt_fbk
345         IF(ln_crt_fbk) THEN
346         WRITE(numout,*) '         Renault et al. 2020, eq. 10: Stau = Alpha * Wnd + Beta'
347         WRITE(numout,*) '            Alpha                                         rn_stau_a    = ', rn_stau_a
348         WRITE(numout,*) '            Beta                                          rn_stau_b    = ', rn_stau_b
349         ENDIF
350         !
351         WRITE(numout,*)
352         SELECT CASE( nblk )              !* Print the choice of bulk algorithm
353         CASE( np_NCAR      )   ;   WRITE(numout,*) '   ==>>>   "NCAR" algorithm        (Large and Yeager 2008)'
354         CASE( np_COARE_3p0 )   ;   WRITE(numout,*) '   ==>>>   "COARE 3.0" algorithm   (Fairall et al. 2003)'
355         CASE( np_COARE_3p6 )   ;   WRITE(numout,*) '   ==>>>   "COARE 3.6" algorithm (Fairall 2018+Edson et al. 2013)'
356         CASE( np_ECMWF     )   ;   WRITE(numout,*) '   ==>>>   "ECMWF" algorithm       (IFS cycle 45r1)'
357         END SELECT
358         !
359         WRITE(numout,*)
360         WRITE(numout,*) '      use cool-skin  parameterization (SSST)  ln_skin_cs  = ', ln_skin_cs
361         WRITE(numout,*) '      use warm-layer parameterization (SSST)  ln_skin_wl  = ', ln_skin_wl
362         !
363         WRITE(numout,*)
364         SELECT CASE( nhumi )              !* Print the choice of air humidity
365         CASE( np_humi_sph )   ;   WRITE(numout,*) '   ==>>>   air humidity is SPECIFIC HUMIDITY     [kg/kg]'
366         CASE( np_humi_dpt )   ;   WRITE(numout,*) '   ==>>>   air humidity is DEW-POINT TEMPERATURE [K]'
367         CASE( np_humi_rlh )   ;   WRITE(numout,*) '   ==>>>   air humidity is RELATIVE HUMIDITY     [%]'
368         END SELECT
369         !
370      ENDIF
371      !
372   END SUBROUTINE sbc_blk_init
373
374
375   SUBROUTINE sbc_blk( kt )
376      !!---------------------------------------------------------------------
377      !!                    ***  ROUTINE sbc_blk  ***
378      !!
379      !! ** Purpose :   provide at each time step the surface ocean fluxes
380      !!              (momentum, heat, freshwater and runoff)
381      !!
382      !! ** Method  :
383      !!              (1) READ each fluxes in NetCDF files:
384      !!      the wind velocity (i-component) at z=rn_zu  (m/s) at T-point
385      !!      the wind velocity (j-component) at z=rn_zu  (m/s) at T-point
386      !!      the specific humidity           at z=rn_zqt (kg/kg)
387      !!      the air temperature             at z=rn_zqt (Kelvin)
388      !!      the solar heat                              (W/m2)
389      !!      the Long wave                               (W/m2)
390      !!      the total precipitation (rain+snow)         (Kg/m2/s)
391      !!      the snow (solid precipitation)              (kg/m2/s)
392      !!      ABL dynamical forcing (i/j-components of either hpg or geostrophic winds)
393      !!              (2) CALL blk_oce_1 and blk_oce_2
394      !!
395      !!      C A U T I O N : never mask the surface stress fields
396      !!                      the stress is assumed to be in the (i,j) mesh referential
397      !!
398      !! ** Action  :   defined at each time-step at the air-sea interface
399      !!              - utau, vtau  i- and j-component of the wind stress
400      !!              - taum        wind stress module at T-point
401      !!              - wndm        wind speed  module at T-point over free ocean or leads in presence of sea-ice
402      !!              - qns, qsr    non-solar and solar heat fluxes
403      !!              - emp         upward mass flux (evapo. - precip.)
404      !!              - sfx         salt flux due to freezing/melting (non-zero only if ice is present)
405      !!
406      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
407      !!                   Brodeau et al. Ocean Modelling 2010
408      !!----------------------------------------------------------------------
409      INTEGER, INTENT(in) ::   kt   ! ocean time step
410      !!----------------------------------------------------------------------
411      REAL(wp), DIMENSION(jpi,jpj) ::   zssq, zcd_du, zsen, zevp
412      REAL(wp) :: ztmp
413      !!----------------------------------------------------------------------
414      !
415      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
416
417      ! Sanity/consistence test on humidity at first time step to detect potential screw-up:
418      IF( kt == nit000 ) THEN
419         IF(lwp) WRITE(numout,*) ''
420#if defined key_agrif
421         IF(lwp) WRITE(numout,*) ' === AGRIF => Sanity/consistence test on air humidity SKIPPED! :( ==='
422#else
423         ztmp = SUM(tmask(:,:,1)) ! number of ocean points on local proc domain
424         IF( ztmp > 8._wp ) THEN ! test only on proc domains with at least 8 ocean points!
425            ztmp = SUM(sf(jp_humi)%fnow(:,:,1)*tmask(:,:,1))/ztmp ! mean humidity over ocean on proc
426            SELECT CASE( nhumi )
427            CASE( np_humi_sph ) ! specific humidity => expect: 0. <= something < 0.065 [kg/kg] (0.061 is saturation at 45degC !!!)
428               IF(  (ztmp < 0._wp) .OR. (ztmp > 0.065)  ) ztmp = -1._wp
429            CASE( np_humi_dpt ) ! dew-point temperature => expect: 110. <= something < 320. [K]
430               IF( (ztmp < 110._wp).OR.(ztmp > 320._wp) ) ztmp = -1._wp
431            CASE( np_humi_rlh ) ! relative humidity => expect: 0. <= something < 100. [%]
432               IF(  (ztmp < 0._wp) .OR.(ztmp > 100._wp) ) ztmp = -1._wp
433            END SELECT
434            IF(ztmp < 0._wp) THEN
435               IF (lwp) WRITE(numout,'("   Mean humidity value found on proc #",i6.6," is: ",f10.5)') narea, ztmp
436               CALL ctl_stop( 'STOP', 'Something is wrong with air humidity!!!', &
437                  &   ' ==> check the unit in your input files'       , &
438                  &   ' ==> check consistence of namelist choice: specific? relative? dew-point?', &
439                  &   ' ==> ln_humi_sph -> [kg/kg] | ln_humi_rlh -> [%] | ln_humi_dpt -> [K] !!!' )
440            END IF
441         END IF
442         IF(lwp) WRITE(numout,*) ' === Sanity/consistence test on air humidity sucessfuly passed! ==='
443#endif
444         IF(lwp) WRITE(numout,*) ''
445      END IF !IF( kt == nit000 )
446      !                                            ! compute the surface ocean fluxes using bulk formulea
447      IF( MOD( kt - 1, nn_fsbc ) == 0 ) THEN
448         CALL blk_oce_1( kt, sf(jp_wndi )%fnow(:,:,1), sf(jp_wndj )%fnow(:,:,1),   &   !   <<= in
449            &                sf(jp_tair )%fnow(:,:,1), sf(jp_humi )%fnow(:,:,1),   &   !   <<= in
450            &                sf(jp_slp  )%fnow(:,:,1), sst_m, ssu_m, ssv_m,        &   !   <<= in
451            &                sf(jp_uoatm)%fnow(:,:,1), sf(jp_voatm)%fnow(:,:,1),   &   !   <<= in
452            &                sf(jp_qsr  )%fnow(:,:,1), sf(jp_qlw  )%fnow(:,:,1),   &   !   <<= in (wl/cs)
453            &                tsk_m, zssq, zcd_du, zsen, zevp )                         !   =>> out
454
455         CALL blk_oce_2(     sf(jp_tair )%fnow(:,:,1), sf(jp_qsr  )%fnow(:,:,1),   &   !   <<= in
456            &                sf(jp_qlw  )%fnow(:,:,1), sf(jp_prec )%fnow(:,:,1),   &   !   <<= in
457            &                sf(jp_snow )%fnow(:,:,1), tsk_m,                      &   !   <<= in
458            &                zsen, zevp )                                              !   <=> in out
459      ENDIF
460      !
461#if defined key_cice
462      IF( MOD( kt - 1, nn_fsbc ) == 0 )   THEN
463         qlw_ice(:,:,1)   = sf(jp_qlw )%fnow(:,:,1)
464         IF( ln_dm2dc ) THEN
465            qsr_ice(:,:,1) = sbc_dcy( sf(jp_qsr)%fnow(:,:,1) )
466         ELSE
467            qsr_ice(:,:,1) =          sf(jp_qsr)%fnow(:,:,1)
468         ENDIF
469         tatm_ice(:,:)    = sf(jp_tair)%fnow(:,:,1)
470
471         SELECT CASE( nhumi )
472         CASE( np_humi_sph )
473            qatm_ice(:,:) =           sf(jp_humi)%fnow(:,:,1)
474         CASE( np_humi_dpt )
475            qatm_ice(:,:) = q_sat(    sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1) )
476         CASE( np_humi_rlh )
477            qatm_ice(:,:) = q_air_rh( 0.01_wp*sf(jp_humi)%fnow(:,:,1), sf(jp_tair)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1)) !LB: 0.01 => RH is % percent in file
478         END SELECT
479
480         tprecip(:,:)     = sf(jp_prec)%fnow(:,:,1) * rn_pfac
481         sprecip(:,:)     = sf(jp_snow)%fnow(:,:,1) * rn_pfac
482         wndi_ice(:,:)    = sf(jp_wndi)%fnow(:,:,1)
483         wndj_ice(:,:)    = sf(jp_wndj)%fnow(:,:,1)
484      ENDIF
485#endif
486      !
487   END SUBROUTINE sbc_blk
488
489
490   SUBROUTINE blk_oce_1( kt, pwndi, pwndj, ptair, phumi,         &  ! inp
491      &                      pslp , pst  , pu   , pv,            &  ! inp
492      &                      puatm, pvatm, pqsr , pqlw ,         &  ! inp
493      &                      ptsk , pssq , pcd_du, psen, pevp   )   ! out
494      !!---------------------------------------------------------------------
495      !!                     ***  ROUTINE blk_oce_1  ***
496      !!
497      !! ** Purpose :   if ln_blk=T, computes surface momentum, heat and freshwater fluxes
498      !!                if ln_abl=T, computes Cd x |U|, Ch x |U|, Ce x |U| for ABL integration
499      !!
500      !! ** Method  :   bulk formulae using atmospheric fields from :
501      !!                if ln_blk=T, atmospheric fields read in sbc_read
502      !!                if ln_abl=T, the ABL model at previous time-step
503      !!
504      !! ** Outputs : - pssq    : surface humidity used to compute latent heat flux (kg/kg)
505      !!              - pcd_du  : Cd x |dU| at T-points  (m/s)
506      !!              - psen    : Ch x |dU| at T-points  (m/s)
507      !!              - pevp    : Ce x |dU| at T-points  (m/s)
508      !!---------------------------------------------------------------------
509      INTEGER , INTENT(in   )                 ::   kt     ! time step index
510      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pwndi  ! atmospheric wind at U-point              [m/s]
511      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pwndj  ! atmospheric wind at V-point              [m/s]
512      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   phumi  ! specific humidity at T-points            [kg/kg]
513      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   ptair  ! potential temperature at T-points        [Kelvin]
514      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pslp   ! sea-level pressure                       [Pa]
515      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pst    ! surface temperature                      [Celsius]
516      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pu     ! surface current at U-point (i-component) [m/s]
517      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pv     ! surface current at V-point (j-component) [m/s]
518      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   puatm  ! surface current seen by the atm at T-point (i-component) [m/s]
519      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pvatm  ! surface current seen by the atm at T-point (j-component) [m/s]
520      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pqsr   !
521      REAL(wp), INTENT(in   ), DIMENSION(:,:) ::   pqlw   !
522      REAL(wp), INTENT(  out), DIMENSION(:,:) ::   ptsk   ! skin temp. (or SST if CS & WL not used)  [Celsius]
523      REAL(wp), INTENT(  out), DIMENSION(:,:) ::   pssq   ! specific humidity at pst                 [kg/kg]
524      REAL(wp), INTENT(  out), DIMENSION(:,:) ::   pcd_du ! Cd x |dU| at T-points                    [m/s]
525      REAL(wp), INTENT(  out), DIMENSION(:,:) ::   psen   ! Ch x |dU| at T-points                    [m/s]
526      REAL(wp), INTENT(  out), DIMENSION(:,:) ::   pevp   ! Ce x |dU| at T-points                    [m/s]
527      !
528      INTEGER  ::   ji, jj               ! dummy loop indices
529      REAL(wp) ::   zztmp                ! local variable
530      REAL(wp) ::   zstmax, zstau
531#if defined key_cyclone
532      REAL(wp), DIMENSION(jpi,jpj) ::   zwnd_i, zwnd_j    ! wind speed components at T-point
533#endif
534      REAL(wp), DIMENSION(jpi,jpj) ::   ztau_i, ztau_j    ! wind stress components at T-point
535      REAL(wp), DIMENSION(jpi,jpj) ::   zU_zu             ! bulk wind speed at height zu  [m/s]
536      REAL(wp), DIMENSION(jpi,jpj) ::   ztpot             ! potential temperature of air at z=rn_zqt [K]
537      REAL(wp), DIMENSION(jpi,jpj) ::   zqair             ! specific humidity     of air at z=rn_zqt [kg/kg]
538      REAL(wp), DIMENSION(jpi,jpj) ::   zcd_oce           ! momentum transfert coefficient over ocean
539      REAL(wp), DIMENSION(jpi,jpj) ::   zch_oce           ! sensible heat transfert coefficient over ocean
540      REAL(wp), DIMENSION(jpi,jpj) ::   zce_oce           ! latent   heat transfert coefficient over ocean
541      REAL(wp), DIMENSION(jpi,jpj) ::   zqla              ! latent heat flux
542      REAL(wp), DIMENSION(jpi,jpj) ::   zztmp1, zztmp2
543      !!---------------------------------------------------------------------
544      !
545      ! local scalars ( place there for vector optimisation purposes)
546      !                           ! Temporary conversion from Celcius to Kelvin (and set minimum value far above 0 K)
547      ptsk(:,:) = pst(:,:) + rt0  ! by default: skin temperature = "bulk SST" (will remain this way if NCAR algorithm used!)
548
549      ! --- cloud cover --- !
550      cloud_fra(:,:) = sf(jp_cc)%fnow(:,:,1)
551
552      ! ----------------------------------------------------------------------------- !
553      !      0   Wind components and module at T-point relative to the moving ocean   !
554      ! ----------------------------------------------------------------------------- !
555
556      ! ... components ( U10m - U_oce ) at T-point (unmasked)
557#if defined key_cyclone
558      zwnd_i(:,:) = 0._wp
559      zwnd_j(:,:) = 0._wp
560      CALL wnd_cyc( kt, zwnd_i, zwnd_j )    ! add analytical tropical cyclone (Vincent et al. JGR 2012)
561      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
562         zwnd_i(ji,jj) = pwndi(ji,jj) + zwnd_i(ji,jj)
563         zwnd_j(ji,jj) = pwndj(ji,jj) + zwnd_j(ji,jj)
564         ! ... scalar wind at T-point (not masked)
565         wndm(ji,jj) = SQRT( zwnd_i(ji,jj) * zwnd_i(ji,jj) + zwnd_j(ji,jj) * zwnd_j(ji,jj) )
566      END_2D
567#else
568      ! ... scalar wind module at T-point (not masked)
569      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
570         wndm(ji,jj) = SQRT(  pwndi(ji,jj) * pwndi(ji,jj) + pwndj(ji,jj) * pwndj(ji,jj)  )
571      END_2D
572#endif
573      ! ----------------------------------------------------------------------------- !
574      !      I   Solar FLUX                                                           !
575      ! ----------------------------------------------------------------------------- !
576
577      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
578      zztmp = 1. - albo
579      IF( ln_dm2dc ) THEN
580         qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
581      ELSE
582         qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
583      ENDIF
584
585
586      ! ----------------------------------------------------------------------------- !
587      !     II   Turbulent FLUXES                                                     !
588      ! ----------------------------------------------------------------------------- !
589
590      ! specific humidity at SST
591      pssq(:,:) = rdct_qsat_salt * q_sat( ptsk(:,:), pslp(:,:) )
592
593      IF( ln_skin_cs .OR. ln_skin_wl ) THEN
594         !! Backup "bulk SST" and associated spec. hum.
595         zztmp1(:,:) = ptsk(:,:)
596         zztmp2(:,:) = pssq(:,:)
597      ENDIF
598
599      ! specific humidity of air at "rn_zqt" m above the sea
600      SELECT CASE( nhumi )
601      CASE( np_humi_sph )
602         zqair(:,:) = phumi(:,:)      ! what we read in file is already a spec. humidity!
603      CASE( np_humi_dpt )
604         !IF(lwp) WRITE(numout,*) ' *** blk_oce => computing q_air out of d_air and slp !' !LBrm
605         zqair(:,:) = q_sat( phumi(:,:), pslp(:,:) )
606      CASE( np_humi_rlh )
607         !IF(lwp) WRITE(numout,*) ' *** blk_oce => computing q_air out of RH, t_air and slp !' !LBrm
608         zqair(:,:) = q_air_rh( 0.01_wp*phumi(:,:), ptair(:,:), pslp(:,:) ) !LB: 0.01 => RH is % percent in file
609      END SELECT
610      !
611      ! potential temperature of air at "rn_zqt" m above the sea
612      IF( ln_abl ) THEN
613         ztpot = ptair(:,:)
614      ELSE
615         ! Estimate of potential temperature at z=rn_zqt, based on adiabatic lapse-rate
616         !    (see Josey, Gulev & Yu, 2013) / doi=10.1016/B978-0-12-391851-2.00005-2
617         !    (since reanalysis products provide T at z, not theta !)
618         !#LB: because AGRIF hates functions that return something else than a scalar, need to
619         !     use scalar version of gamma_moist() ...
620         IF( ln_tpot ) THEN
621            DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
622               ztpot(ji,jj) = ptair(ji,jj) + gamma_moist( ptair(ji,jj), zqair(ji,jj) ) * rn_zqt
623            END_2D
624         ELSE
625            ztpot = ptair(:,:)
626         ENDIF
627      ENDIF
628
629      !! Time to call the user-selected bulk parameterization for
630      !!  ==  transfer coefficients  ==!   Cd, Ch, Ce at T-point, and more...
631      SELECT CASE( nblk )
632
633      CASE( np_NCAR      )
634         CALL turb_ncar    ( rn_zqt, rn_zu, ptsk, ztpot, pssq, zqair, wndm,                              &
635            &                zcd_oce, zch_oce, zce_oce, t_zu, q_zu, zU_zu, cdn_oce, chn_oce, cen_oce )
636
637      CASE( np_COARE_3p0 )
638         CALL turb_coare3p0 ( kt, rn_zqt, rn_zu, ptsk, ztpot, pssq, zqair, wndm, ln_skin_cs, ln_skin_wl, &
639            &                zcd_oce, zch_oce, zce_oce, t_zu, q_zu, zU_zu, cdn_oce, chn_oce, cen_oce,   &
640            &                Qsw=qsr(:,:), rad_lw=pqlw(:,:), slp=pslp(:,:) )
641
642      CASE( np_COARE_3p6 )
643         CALL turb_coare3p6 ( kt, rn_zqt, rn_zu, ptsk, ztpot, pssq, zqair, wndm, ln_skin_cs, ln_skin_wl, &
644            &                zcd_oce, zch_oce, zce_oce, t_zu, q_zu, zU_zu, cdn_oce, chn_oce, cen_oce,   &
645            &                Qsw=qsr(:,:), rad_lw=pqlw(:,:), slp=pslp(:,:) )
646
647      CASE( np_ECMWF     )
648         CALL turb_ecmwf   ( kt, rn_zqt, rn_zu, ptsk, ztpot, pssq, zqair, wndm, ln_skin_cs, ln_skin_wl,  &
649            &                zcd_oce, zch_oce, zce_oce, t_zu, q_zu, zU_zu, cdn_oce, chn_oce, cen_oce,   &
650            &                Qsw=qsr(:,:), rad_lw=pqlw(:,:), slp=pslp(:,:) )
651
652      CASE DEFAULT
653         CALL ctl_stop( 'STOP', 'sbc_oce: non-existing bulk formula selected' )
654
655      END SELECT
656     
657      IF( iom_use('Cd_oce') )   CALL iom_put("Cd_oce",   zcd_oce * tmask(:,:,1))
658      IF( iom_use('Ce_oce') )   CALL iom_put("Ce_oce",   zce_oce * tmask(:,:,1))
659      IF( iom_use('Ch_oce') )   CALL iom_put("Ch_oce",   zch_oce * tmask(:,:,1))
660      !! LB: mainly here for debugging purpose:
661      IF( iom_use('theta_zt') ) CALL iom_put("theta_zt", (ztpot-rt0) * tmask(:,:,1)) ! potential temperature at z=zt
662      IF( iom_use('q_zt') )     CALL iom_put("q_zt",     zqair       * tmask(:,:,1)) ! specific humidity       "
663      IF( iom_use('theta_zu') ) CALL iom_put("theta_zu", (t_zu -rt0) * tmask(:,:,1)) ! potential temperature at z=zu
664      IF( iom_use('q_zu') )     CALL iom_put("q_zu",     q_zu        * tmask(:,:,1)) ! specific humidity       "
665      IF( iom_use('ssq') )      CALL iom_put("ssq",      pssq        * tmask(:,:,1)) ! saturation specific humidity at z=0
666      IF( iom_use('wspd_blk') ) CALL iom_put("wspd_blk", zU_zu       * tmask(:,:,1)) ! bulk wind speed at z=zu
667     
668      IF( ln_skin_cs .OR. ln_skin_wl ) THEN
669         !! ptsk and pssq have been updated!!!
670         !!
671         !! In the presence of sea-ice we forget about the cool-skin/warm-layer update of ptsk and pssq:
672         WHERE ( fr_i(:,:) > 0.001_wp )
673            ! sea-ice present, we forget about the update, using what we backed up before call to turb_*()
674            ptsk(:,:) = zztmp1(:,:)
675            pssq(:,:) = zztmp2(:,:)
676         END WHERE
677      END IF
678
679      !  Turbulent fluxes over ocean  => BULK_FORMULA @ sbcblk_phy.F90
680      ! -------------------------------------------------------------
681
682      IF( ln_abl ) THEN         !==  ABL formulation  ==!   multiplication by rho_air and turbulent fluxes computation done in ablstp
683         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
684            zztmp = zU_zu(ji,jj)
685            wndm(ji,jj)   = zztmp                   ! Store zU_zu in wndm to compute ustar2 in ablmod
686            pcd_du(ji,jj) = zztmp * zcd_oce(ji,jj)
687            psen(ji,jj)   = zztmp * zch_oce(ji,jj)
688            pevp(ji,jj)   = zztmp * zce_oce(ji,jj)
689            rhoa(ji,jj)   = rho_air( ptair(ji,jj), phumi(ji,jj), pslp(ji,jj) )
690         END_2D
691      ELSE                      !==  BLK formulation  ==!   turbulent fluxes computation
692         CALL BULK_FORMULA( rn_zu, ptsk(:,:), pssq(:,:), t_zu(:,:), q_zu(:,:), &
693            &               zcd_oce(:,:), zch_oce(:,:), zce_oce(:,:),          &
694            &               wndm(:,:), zU_zu(:,:), pslp(:,:),                  &
695            &               taum(:,:), psen(:,:), zqla(:,:),                   &
696            &               pEvap=pevp(:,:), prhoa=rhoa(:,:), pfact_evap=rn_efac )
697
698         zqla(:,:) = zqla(:,:) * tmask(:,:,1)
699         psen(:,:) = psen(:,:) * tmask(:,:,1)
700         taum(:,:) = taum(:,:) * tmask(:,:,1)
701         pevp(:,:) = pevp(:,:) * tmask(:,:,1)
702
703         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
704            IF( wndm(ji,jj) > 0._wp ) THEN
705               zztmp = taum(ji,jj) / wndm(ji,jj)
706#if defined key_cyclone
707               ztau_i(ji,jj) = zztmp * zwnd_i(ji,jj)
708               ztau_j(ji,jj) = zztmp * zwnd_j(ji,jj)
709#else
710               ztau_i(ji,jj) = zztmp * pwndi(ji,jj)
711               ztau_j(ji,jj) = zztmp * pwndj(ji,jj)
712#endif
713            ELSE
714               ztau_i(ji,jj) = 0._wp
715               ztau_j(ji,jj) = 0._wp                 
716            ENDIF
717         END_2D
718
719         IF( ln_crt_fbk ) THEN   ! aply eq. 10 and 11 of Renault et al. 2020 (doi: 10.1029/2019MS001715)
720            zstmax = MIN( rn_stau_a * 3._wp + rn_stau_b, 0._wp )   ! set the max value of Stau corresponding to a wind of 3 m/s (<0)
721            DO_2D( 0, 1, 0, 1 )   ! end at jpj and jpi, as ztau_j(ji,jj+1) ztau_i(ji+1,jj) used in the next loop
722               zstau = MIN( rn_stau_a * wndm(ji,jj) + rn_stau_b, zstmax )   ! stau (<0) must be smaller than zstmax
723               ztau_i(ji,jj) = ztau_i(ji,jj) + zstau * ( 0.5_wp * ( pu(ji-1,jj  ) + pu(ji,jj) ) - puatm(ji,jj) )
724               ztau_j(ji,jj) = ztau_j(ji,jj) + zstau * ( 0.5_wp * ( pv(ji  ,jj-1) + pv(ji,jj) ) - pvatm(ji,jj) )
725               taum(ji,jj) = SQRT( ztau_i(ji,jj) * ztau_i(ji,jj) + ztau_j(ji,jj) * ztau_j(ji,jj) )
726            END_2D
727         ENDIF
728
729         ! ... utau, vtau at U- and V_points, resp.
730         !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
731         !     Note that coastal wind stress is not used in the code... so this extra care has no effect
732         DO_2D( 0, 0, 0, 0 )              ! start loop at 2, in case ln_crt_fbk = T
733            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( ztau_i(ji,jj) + ztau_i(ji+1,jj  ) ) &
734               &              * MAX(tmask(ji,jj,1),tmask(ji+1,jj,1))
735            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( ztau_j(ji,jj) + ztau_j(ji  ,jj+1) ) &
736               &              * MAX(tmask(ji,jj,1),tmask(ji,jj+1,1))
737         END_2D
738
739         IF( ln_crt_fbk ) THEN
740            CALL lbc_lnk_multi( 'sbcblk', utau, 'U', -1., vtau, 'V', -1., taum, 'T', -1. )
741         ELSE
742            CALL lbc_lnk_multi( 'sbcblk', utau, 'U', -1., vtau, 'V', -1. )
743         ENDIF
744
745         CALL iom_put( "taum_oce", taum )   ! output wind stress module
746
747         IF(sn_cfctl%l_prtctl) THEN
748            CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce_1: wndm   : ')
749            CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce_1: utau   : ', mask1=umask,   &
750               &          tab2d_2=vtau  , clinfo2='            vtau   : ', mask2=vmask )
751         ENDIF
752         !
753      ENDIF !IF( ln_abl )
754     
755      ptsk(:,:) = ( ptsk(:,:) - rt0 ) * tmask(:,:,1)  ! Back to Celsius
756           
757      IF( ln_skin_cs .OR. ln_skin_wl ) THEN
758         CALL iom_put( "t_skin" ,  ptsk        )  ! T_skin in Celsius
759         CALL iom_put( "dt_skin" , ptsk - pst  )  ! T_skin - SST temperature difference...
760      ENDIF
761
762      IF(sn_cfctl%l_prtctl) THEN
763         CALL prt_ctl( tab2d_1=pevp  , clinfo1=' blk_oce_1: pevp   : ' )
764         CALL prt_ctl( tab2d_1=psen  , clinfo1=' blk_oce_1: psen   : ' )
765         CALL prt_ctl( tab2d_1=pssq  , clinfo1=' blk_oce_1: pssq   : ' )
766      ENDIF
767      !
768   END SUBROUTINE blk_oce_1
769
770
771   SUBROUTINE blk_oce_2( ptair, pqsr, pqlw, pprec,   &   ! <<= in
772      &                  psnow, ptsk, psen, pevp     )   ! <<= in
773      !!---------------------------------------------------------------------
774      !!                     ***  ROUTINE blk_oce_2  ***
775      !!
776      !! ** Purpose :   finalize the momentum, heat and freshwater fluxes computation
777      !!                at the ocean surface at each time step knowing Cd, Ch, Ce and
778      !!                atmospheric variables (from ABL or external data)
779      !!
780      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
781      !!              - vtau    : j-component of the stress at V-point  (N/m2)
782      !!              - taum    : Wind stress module at T-point         (N/m2)
783      !!              - wndm    : Wind speed module at T-point          (m/s)
784      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
785      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
786      !!              - emp     : evaporation minus precipitation       (kg/m2/s)
787      !!---------------------------------------------------------------------
788      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptair
789      REAL(wp), INTENT(in), DIMENSION(:,:) ::   pqsr
790      REAL(wp), INTENT(in), DIMENSION(:,:) ::   pqlw
791      REAL(wp), INTENT(in), DIMENSION(:,:) ::   pprec
792      REAL(wp), INTENT(in), DIMENSION(:,:) ::   psnow
793      REAL(wp), INTENT(in), DIMENSION(:,:) ::   ptsk   ! SKIN surface temperature   [Celsius]
794      REAL(wp), INTENT(in), DIMENSION(:,:) ::   psen
795      REAL(wp), INTENT(in), DIMENSION(:,:) ::   pevp
796      !
797      INTEGER  ::   ji, jj               ! dummy loop indices
798      REAL(wp) ::   zztmp,zz1,zz2,zz3    ! local variable
799      REAL(wp), DIMENSION(jpi,jpj) ::   ztskk             ! skin temp. in Kelvin
800      REAL(wp), DIMENSION(jpi,jpj) ::   zqlw              ! long wave and sensible heat fluxes     
801      REAL(wp), DIMENSION(jpi,jpj) ::   zqla              ! latent heat fluxes and evaporation
802      !!---------------------------------------------------------------------
803      !
804      ! local scalars ( place there for vector optimisation purposes)
805
806
807      ztskk(:,:) = ptsk(:,:) + rt0  ! => ptsk in Kelvin rather than Celsius
808     
809      ! ----------------------------------------------------------------------------- !
810      !     III    Net longwave radiative FLUX                                        !
811      ! ----------------------------------------------------------------------------- !
812
813      !! LB: now moved after Turbulent fluxes because must use the skin temperature rather that the SST
814      !! (ztskk is skin temperature if ln_skin_cs==.TRUE. .OR. ln_skin_wl==.TRUE.)
815      zqlw(:,:) = emiss_w * ( pqlw(:,:) - stefan*ztskk(:,:)*ztskk(:,:)*ztskk(:,:)*ztskk(:,:) ) * tmask(:,:,1)   ! Net radiative longwave flux
816
817      !  Latent flux over ocean
818      ! -----------------------
819
820      ! use scalar version of L_vap() for AGRIF compatibility
821      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
822         zqla(ji,jj) = - L_vap( ztskk(ji,jj) ) * pevp(ji,jj)    ! Latent Heat flux !!GS: possibility to add a global qla to avoid recomputation after abl update
823      END_2D
824
825      IF(sn_cfctl%l_prtctl) THEN
826         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce_2: zqla   : ' )
827         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce_2: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
828
829      ENDIF
830
831      ! ----------------------------------------------------------------------------- !
832      !     IV    Total FLUXES                                                       !
833      ! ----------------------------------------------------------------------------- !
834      !
835      emp (:,:) = (  pevp(:,:)                                       &   ! mass flux (evap. - precip.)
836         &         - pprec(:,:) * rn_pfac  ) * tmask(:,:,1)
837      !
838      qns(:,:) = zqlw(:,:) + psen(:,:) + zqla(:,:)                   &   ! Downward Non Solar
839         &     - psnow(:,:) * rn_pfac * rLfus                        &   ! remove latent melting heat for solid precip
840         &     - pevp(:,:) * ptsk(:,:) * rcp                         &   ! remove evap heat content at SST
841         &     + ( pprec(:,:) - psnow(:,:) ) * rn_pfac               &   ! add liquid precip heat content at Tair
842         &     * ( ptair(:,:) - rt0 ) * rcp                          &
843         &     + psnow(:,:) * rn_pfac                                &   ! add solid  precip heat content at min(Tair,Tsnow)
844         &     * ( MIN( ptair(:,:), rt0 ) - rt0 ) * rcpi
845      qns(:,:) = qns(:,:) * tmask(:,:,1)
846      !
847#if defined key_si3
848      qns_oce(:,:) = zqlw(:,:) + psen(:,:) + zqla(:,:)                             ! non solar without emp (only needed by SI3)
849      qsr_oce(:,:) = qsr(:,:)
850#endif
851      !
852      CALL iom_put( "rho_air"  , rhoa*tmask(:,:,1) )       ! output air density [kg/m^3]
853      CALL iom_put( "evap_oce" , pevp )                    ! evaporation
854      CALL iom_put( "qlw_oce"  , zqlw )                    ! output downward longwave heat over the ocean
855      CALL iom_put( "qsb_oce"  , psen )                    ! output downward sensible heat over the ocean
856      CALL iom_put( "qla_oce"  , zqla )                    ! output downward latent   heat over the ocean
857      tprecip(:,:) = pprec(:,:) * rn_pfac * tmask(:,:,1)   ! output total precipitation [kg/m2/s]
858      sprecip(:,:) = psnow(:,:) * rn_pfac * tmask(:,:,1)   ! output solid precipitation [kg/m2/s]
859      CALL iom_put( 'snowpre', sprecip )                   ! Snow
860      CALL iom_put( 'precip' , tprecip )                   ! Total precipitation
861      !
862      IF ( nn_ice == 0 ) THEN
863         CALL iom_put( "qemp_oce" , qns-zqlw-psen-zqla )   ! output downward heat content of E-P over the ocean
864         CALL iom_put( "qns_oce"  ,   qns  )               ! output downward non solar heat over the ocean
865         CALL iom_put( "qsr_oce"  ,   qsr  )               ! output downward solar heat over the ocean
866         CALL iom_put( "qt_oce"   ,   qns+qsr )            ! output total downward heat over the ocean
867      ENDIF
868      !
869      IF(sn_cfctl%l_prtctl) THEN
870         CALL prt_ctl(tab2d_1=zqlw , clinfo1=' blk_oce_2: zqlw  : ')
871         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_2: zqla  : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
872         CALL prt_ctl(tab2d_1=emp  , clinfo1=' blk_oce_2: emp   : ')
873      ENDIF
874      !
875   END SUBROUTINE blk_oce_2
876
877
878#if defined key_si3
879   !!----------------------------------------------------------------------
880   !!   'key_si3'                                       SI3 sea-ice model
881   !!----------------------------------------------------------------------
882   !!   blk_ice_1   : provide the air-ice stress
883   !!   blk_ice_2   : provide the heat and mass fluxes at air-ice interface
884   !!   blk_ice_qcn : provide ice surface temperature and snow/ice conduction flux (emulating conduction flux)
885   !!   Cdn10_Lupkes2012 : Lupkes et al. (2012) air-ice drag
886   !!   Cdn10_Lupkes2015 : Lupkes et al. (2015) air-ice drag
887   !!----------------------------------------------------------------------
888
889   SUBROUTINE blk_ice_1( pwndi, pwndj, ptair, phumi, pslp , puice, pvice, ptsui,  &   ! inputs
890      &                  putaui, pvtaui, pseni, pevpi, pssqi, pcd_dui             )   ! optional outputs
891      !!---------------------------------------------------------------------
892      !!                     ***  ROUTINE blk_ice_1  ***
893      !!
894      !! ** Purpose :   provide the surface boundary condition over sea-ice
895      !!
896      !! ** Method  :   compute momentum using bulk formulation
897      !!                formulea, ice variables and read atmospheric fields.
898      !!                NB: ice drag coefficient is assumed to be a constant
899      !!---------------------------------------------------------------------
900      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   pslp    ! sea-level pressure [Pa]
901      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   pwndi   ! atmospheric wind at T-point [m/s]
902      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   pwndj   ! atmospheric wind at T-point [m/s]
903      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   ptair   ! atmospheric wind at T-point [m/s]
904      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   phumi   ! atmospheric wind at T-point [m/s]
905      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   puice   ! sea-ice velocity on I or C grid [m/s]
906      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   pvice   ! "
907      REAL(wp) , INTENT(in   ), DIMENSION(:,:  ) ::   ptsui   ! sea-ice surface temperature [K]
908      REAL(wp) , INTENT(  out), DIMENSION(:,:  ), OPTIONAL ::   putaui  ! if ln_blk
909      REAL(wp) , INTENT(  out), DIMENSION(:,:  ), OPTIONAL ::   pvtaui  ! if ln_blk
910      REAL(wp) , INTENT(  out), DIMENSION(:,:  ), OPTIONAL ::   pseni   ! if ln_abl
911      REAL(wp) , INTENT(  out), DIMENSION(:,:  ), OPTIONAL ::   pevpi   ! if ln_abl
912      REAL(wp) , INTENT(  out), DIMENSION(:,:  ), OPTIONAL ::   pssqi   ! if ln_abl
913      REAL(wp) , INTENT(  out), DIMENSION(:,:  ), OPTIONAL ::   pcd_dui ! if ln_abl
914      !
915      INTEGER  ::   ji, jj    ! dummy loop indices
916      REAL(wp) ::   zootm_su                      ! sea-ice surface mean temperature
917      REAL(wp) ::   zztmp1, zztmp2                ! temporary arrays
918      REAL(wp), DIMENSION(jpi,jpj) ::   zcd_dui   ! transfer coefficient for momentum      (tau)
919      !!---------------------------------------------------------------------
920      !
921
922      ! ------------------------------------------------------------ !
923      !    Wind module relative to the moving ice ( U10m - U_ice )   !
924      ! ------------------------------------------------------------ !
925      ! C-grid ice dynamics :   U & V-points (same as ocean)
926      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
927         wndm_ice(ji,jj) = SQRT( pwndi(ji,jj) * pwndi(ji,jj) + pwndj(ji,jj) * pwndj(ji,jj) )
928      END_2D
929      !
930      ! Make ice-atm. drag dependent on ice concentration
931      IF    ( ln_Cd_L12 ) THEN   ! calculate new drag from Lupkes(2012) equations
932         CALL Cdn10_Lupkes2012( Cd_ice )
933         Ch_ice(:,:) = Cd_ice(:,:)       ! momentum and heat transfer coef. are considered identical
934         Ce_ice(:,:) = Cd_ice(:,:)
935      ELSEIF( ln_Cd_L15 ) THEN   ! calculate new drag from Lupkes(2015) equations
936         CALL Cdn10_Lupkes2015( ptsui, pslp, Cd_ice, Ch_ice )
937         Ce_ice(:,:) = Ch_ice(:,:)       ! sensible and latent heat transfer coef. are considered identical
938      ENDIF
939     
940      IF( iom_use('Cd_ice') ) CALL iom_put("Cd_ice", Cd_ice)
941      IF( iom_use('Ce_ice') ) CALL iom_put("Ce_ice", Ce_ice)
942      IF( iom_use('Ch_ice') ) CALL iom_put("Ch_ice", Ch_ice)
943     
944      ! local scalars ( place there for vector optimisation purposes)
945      zcd_dui(:,:) = wndm_ice(:,:) * Cd_ice(:,:)
946
947      IF( ln_blk ) THEN
948         ! ---------------------------------------------------- !
949         !    Wind stress relative to nonmoving ice ( U10m )    !
950         ! ---------------------------------------------------- !
951         ! supress moving ice in wind stress computation as we don't know how to do it properly...
952         DO_2D( 0, 1, 0, 1 )    ! at T point
953            putaui(ji,jj) = rhoa(ji,jj) * zcd_dui(ji,jj) * pwndi(ji,jj)
954            pvtaui(ji,jj) = rhoa(ji,jj) * zcd_dui(ji,jj) * pwndj(ji,jj)
955         END_2D
956         !
957         DO_2D( 0, 0, 0, 0 )    ! U & V-points (same as ocean).
958            ! take care of the land-sea mask to avoid "pollution" of coastal stress. p[uv]taui used in frazil and  rheology
959            zztmp1 = 0.5_wp * ( 2. - umask(ji,jj,1) ) * MAX( tmask(ji,jj,1),tmask(ji+1,jj  ,1) )
960            zztmp2 = 0.5_wp * ( 2. - vmask(ji,jj,1) ) * MAX( tmask(ji,jj,1),tmask(ji  ,jj+1,1) )
961            putaui(ji,jj) = zztmp1 * ( putaui(ji,jj) + putaui(ji+1,jj  ) )
962            pvtaui(ji,jj) = zztmp2 * ( pvtaui(ji,jj) + pvtaui(ji  ,jj+1) )
963         END_2D
964         CALL lbc_lnk_multi( 'sbcblk', putaui, 'U', -1._wp, pvtaui, 'V', -1._wp )
965         !
966         IF(sn_cfctl%l_prtctl)  CALL prt_ctl( tab2d_1=putaui  , clinfo1=' blk_ice: putaui : '   &
967            &                               , tab2d_2=pvtaui  , clinfo2='          pvtaui : ' )
968      ELSE ! ln_abl
969         zztmp1 = 11637800.0_wp
970         zztmp2 =    -5897.8_wp
971         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
972            pcd_dui(ji,jj) = zcd_dui (ji,jj)
973            pseni  (ji,jj) = wndm_ice(ji,jj) * Ch_ice(ji,jj)
974            pevpi  (ji,jj) = wndm_ice(ji,jj) * Ce_ice(ji,jj)
975            zootm_su       = zztmp2 / ptsui(ji,jj)   ! ptsui is in K (it can't be zero ??)
976            pssqi  (ji,jj) = zztmp1 * EXP( zootm_su ) / rhoa(ji,jj)
977         END_2D
978      ENDIF
979      !
980      IF(sn_cfctl%l_prtctl)  CALL prt_ctl(tab2d_1=wndm_ice  , clinfo1=' blk_ice: wndm_ice : ')
981      !
982   END SUBROUTINE blk_ice_1
983
984
985   SUBROUTINE blk_ice_2( ptsu, phs, phi, palb, ptair, phumi, pslp, pqlw, pprec, psnow  )
986      !!---------------------------------------------------------------------
987      !!                     ***  ROUTINE blk_ice_2  ***
988      !!
989      !! ** Purpose :   provide the heat and mass fluxes at air-ice interface
990      !!
991      !! ** Method  :   compute heat and freshwater exchanged
992      !!                between atmosphere and sea-ice using bulk formulation
993      !!                formulea, ice variables and read atmmospheric fields.
994      !!
995      !! caution : the net upward water flux has with mm/day unit
996      !!---------------------------------------------------------------------
997      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   ptsu   ! sea ice surface temperature [K]
998      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   phs    ! snow thickness
999      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   phi    ! ice thickness
1000      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   palb   ! ice albedo (all skies)
1001      REAL(wp), DIMENSION(:,:  ), INTENT(in)  ::   ptair
1002      REAL(wp), DIMENSION(:,:  ), INTENT(in)  ::   phumi
1003      REAL(wp), DIMENSION(:,:  ), INTENT(in)  ::   pslp
1004      REAL(wp), DIMENSION(:,:  ), INTENT(in)  ::   pqlw
1005      REAL(wp), DIMENSION(:,:  ), INTENT(in)  ::   pprec
1006      REAL(wp), DIMENSION(:,:  ), INTENT(in)  ::   psnow
1007      !!
1008      INTEGER  ::   ji, jj, jl               ! dummy loop indices
1009      REAL(wp) ::   zst3                     ! local variable
1010      REAL(wp) ::   zcoef_dqlw, zcoef_dqla   !   -      -
1011      REAL(wp) ::   zztmp, zztmp2, z1_rLsub  !   -      -
1012      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   z1_st         ! inverse of surface temperature
1013      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   z_qlw         ! long wave heat flux over ice
1014      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   z_qsb         ! sensible  heat flux over ice
1015      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   z_dqlw        ! long wave heat sensitivity over ice
1016      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   z_dqsb        ! sensible  heat sensitivity over ice
1017      REAL(wp), DIMENSION(jpi,jpj)     ::   zevap, zsnw   ! evaporation and snw distribution after wind blowing (SI3)
1018      REAL(wp), DIMENSION(jpi,jpj)     ::   zqair         ! specific humidity of air at z=rn_zqt [kg/kg] !LB
1019      REAL(wp), DIMENSION(jpi,jpj)     ::   ztmp, ztmp2
1020      REAL(wp), DIMENSION(jpi,jpj)     ::   ztri
1021      !!---------------------------------------------------------------------
1022      !
1023      zcoef_dqlw = 4._wp * 0.95_wp * stefan             ! local scalars
1024      zcoef_dqla = -rLsub * 11637800._wp * (-5897.8_wp) !LB: BAD!
1025      !
1026      SELECT CASE( nhumi )
1027      CASE( np_humi_sph )
1028         zqair(:,:) =  phumi(:,:)      ! what we read in file is already a spec. humidity!
1029      CASE( np_humi_dpt )
1030         zqair(:,:) = q_sat( phumi(:,:), pslp )
1031      CASE( np_humi_rlh )
1032         zqair(:,:) = q_air_rh( 0.01_wp*phumi(:,:), ptair(:,:), pslp(:,:) ) !LB: 0.01 => RH is % percent in file
1033      END SELECT
1034      !
1035      zztmp = 1. / ( 1. - albo )
1036      WHERE( ptsu(:,:,:) /= 0._wp )
1037         z1_st(:,:,:) = 1._wp / ptsu(:,:,:)
1038      ELSEWHERE
1039         z1_st(:,:,:) = 0._wp
1040      END WHERE
1041      !                                     ! ========================== !
1042      DO jl = 1, jpl                        !  Loop over ice categories  !
1043         !                                  ! ========================== !
1044         DO jj = 1 , jpj
1045            DO ji = 1, jpi
1046               ! ----------------------------!
1047               !      I   Radiative FLUXES   !
1048               ! ----------------------------!
1049               zst3 = ptsu(ji,jj,jl) * ptsu(ji,jj,jl) * ptsu(ji,jj,jl)
1050               ! Short Wave (sw)
1051               qsr_ice(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
1052               ! Long  Wave (lw)
1053               z_qlw(ji,jj,jl) = 0.95 * ( pqlw(ji,jj) - stefan * ptsu(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
1054               ! lw sensitivity
1055               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3
1056
1057               ! ----------------------------!
1058               !     II    Turbulent FLUXES  !
1059               ! ----------------------------!
1060
1061               ! ... turbulent heat fluxes with Ch_ice recalculated in blk_ice_1
1062               ! Sensible Heat
1063               z_qsb(ji,jj,jl) = rhoa(ji,jj) * rCp_air * Ch_ice(ji,jj) * wndm_ice(ji,jj) * (ptsu(ji,jj,jl) - ptair(ji,jj))
1064               ! Latent Heat
1065               zztmp2 = EXP( -5897.8 * z1_st(ji,jj,jl) )
1066               qla_ice(ji,jj,jl) = rn_efac * MAX( 0.e0, rhoa(ji,jj) * rLsub  * Ce_ice(ji,jj) * wndm_ice(ji,jj) *  &
1067                  &                ( 11637800. * zztmp2 / rhoa(ji,jj) - zqair(ji,jj) ) )
1068               ! Latent heat sensitivity for ice (Dqla/Dt)
1069               IF( qla_ice(ji,jj,jl) > 0._wp ) THEN
1070                  dqla_ice(ji,jj,jl) = rn_efac * zcoef_dqla * Ce_ice(ji,jj) * wndm_ice(ji,jj) *  &
1071                     &                 z1_st(ji,jj,jl) * z1_st(ji,jj,jl) * zztmp2
1072               ELSE
1073                  dqla_ice(ji,jj,jl) = 0._wp
1074               ENDIF
1075
1076               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
1077               z_dqsb(ji,jj,jl) = rhoa(ji,jj) * rCp_air * Ch_ice(ji,jj) * wndm_ice(ji,jj)
1078
1079               ! ----------------------------!
1080               !     III    Total FLUXES     !
1081               ! ----------------------------!
1082               ! Downward Non Solar flux
1083               qns_ice (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - qla_ice (ji,jj,jl)
1084               ! Total non solar heat flux sensitivity for ice
1085               dqns_ice(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + dqla_ice(ji,jj,jl) )
1086            END DO
1087            !
1088         END DO
1089         !
1090      END DO
1091      !
1092      tprecip(:,:) = pprec(:,:) * rn_pfac * tmask(:,:,1)  ! total precipitation [kg/m2/s]
1093      sprecip(:,:) = psnow(:,:) * rn_pfac * tmask(:,:,1)  ! solid precipitation [kg/m2/s]
1094      CALL iom_put( 'snowpre', sprecip )                  ! Snow precipitation
1095      CALL iom_put( 'precip' , tprecip )                  ! Total precipitation
1096
1097      ! --- evaporation --- !
1098      z1_rLsub = 1._wp / rLsub
1099      evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * z1_rLsub    ! sublimation
1100      devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * z1_rLsub    ! d(sublimation)/dT
1101      zevap    (:,:)   = emp(:,:) + tprecip(:,:)   ! evaporation over ocean  !LB: removed rn_efac here, correct???
1102
1103      ! --- evaporation minus precipitation --- !
1104      zsnw(:,:) = 0._wp
1105      CALL ice_var_snwblow( (1.-at_i_b(:,:)), zsnw )  ! snow distribution over ice after wind blowing
1106      emp_oce(:,:) = ( 1._wp - at_i_b(:,:) ) * zevap(:,:) - ( tprecip(:,:) - sprecip(:,:) ) - sprecip(:,:) * (1._wp - zsnw )
1107      emp_ice(:,:) = SUM( a_i_b(:,:,:) * evap_ice(:,:,:), dim=3 ) - sprecip(:,:) * zsnw
1108      emp_tot(:,:) = emp_oce(:,:) + emp_ice(:,:)
1109
1110      ! --- heat flux associated with emp --- !
1111      qemp_oce(:,:) = - ( 1._wp - at_i_b(:,:) ) * zevap(:,:) * sst_m(:,:) * rcp                  & ! evap at sst
1112         &          + ( tprecip(:,:) - sprecip(:,:) ) * ( ptair(:,:) - rt0 ) * rcp               & ! liquid precip at Tair
1113         &          +   sprecip(:,:) * ( 1._wp - zsnw ) *                                        & ! solid precip at min(Tair,Tsnow)
1114         &              ( ( MIN( ptair(:,:), rt0 ) - rt0 ) * rcpi * tmask(:,:,1) - rLfus )
1115      qemp_ice(:,:) =   sprecip(:,:) * zsnw *                                                    & ! solid precip (only)
1116         &              ( ( MIN( ptair(:,:), rt0 ) - rt0 ) * rcpi * tmask(:,:,1) - rLfus )
1117
1118      ! --- total solar and non solar fluxes --- !
1119      qns_tot(:,:) = ( 1._wp - at_i_b(:,:) ) * qns_oce(:,:) + SUM( a_i_b(:,:,:) * qns_ice(:,:,:), dim=3 )  &
1120         &           + qemp_ice(:,:) + qemp_oce(:,:)
1121      qsr_tot(:,:) = ( 1._wp - at_i_b(:,:) ) * qsr_oce(:,:) + SUM( a_i_b(:,:,:) * qsr_ice(:,:,:), dim=3 )
1122
1123      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
1124      qprec_ice(:,:) = rhos * ( ( MIN( ptair(:,:), rt0 ) - rt0 ) * rcpi * tmask(:,:,1) - rLfus )
1125
1126      ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) ---
1127      DO jl = 1, jpl
1128         qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * rcpi * tmask(:,:,1) )
1129         !                         ! But we do not have Tice => consider it at 0degC => evap=0
1130      END DO
1131
1132      ! --- shortwave radiation transmitted thru the surface scattering layer (W/m2) --- !
1133      IF( nn_qtrice == 0 ) THEN
1134         ! formulation derived from Grenfell and Maykut (1977), where transmission rate
1135         !    1) depends on cloudiness
1136         !    2) is 0 when there is any snow
1137         !    3) tends to 1 for thin ice
1138         ztri(:,:) = 0.18 * ( 1.0 - cloud_fra(:,:) ) + 0.35 * cloud_fra(:,:)  ! surface transmission when hi>10cm
1139         DO jl = 1, jpl
1140            WHERE    ( phs(:,:,jl) <= 0._wp .AND. phi(:,:,jl) <  0.1_wp )     ! linear decrease from hi=0 to 10cm 
1141               qtr_ice_top(:,:,jl) = qsr_ice(:,:,jl) * ( ztri(:,:) + ( 1._wp - ztri(:,:) ) * ( 1._wp - phi(:,:,jl) * 10._wp ) )
1142            ELSEWHERE( phs(:,:,jl) <= 0._wp .AND. phi(:,:,jl) >= 0.1_wp )     ! constant (ztri) when hi>10cm
1143               qtr_ice_top(:,:,jl) = qsr_ice(:,:,jl) * ztri(:,:)
1144            ELSEWHERE                                                         ! zero when hs>0
1145               qtr_ice_top(:,:,jl) = 0._wp 
1146            END WHERE
1147         ENDDO
1148      ELSEIF( nn_qtrice == 1 ) THEN
1149         ! formulation is derived from the thesis of M. Lebrun (2019).
1150         !    It represents the best fit using several sets of observations
1151         !    It comes with snow conductivities adapted to freezing/melting conditions (see icethd_zdf_bl99.F90)
1152         qtr_ice_top(:,:,:) = 0.3_wp * qsr_ice(:,:,:)
1153      ENDIF
1154      !
1155      IF( iom_use('evap_ao_cea') .OR. iom_use('hflx_evap_cea') ) THEN
1156         ztmp(:,:) = zevap(:,:) * ( 1._wp - at_i_b(:,:) )
1157         IF( iom_use('evap_ao_cea'  ) )  CALL iom_put( 'evap_ao_cea'  , ztmp(:,:) * tmask(:,:,1) )   ! ice-free oce evap (cell average)
1158         IF( iom_use('hflx_evap_cea') )  CALL iom_put( 'hflx_evap_cea', ztmp(:,:) * sst_m(:,:) * rcp * tmask(:,:,1) )   ! heat flux from evap (cell average)
1159      ENDIF
1160      IF( iom_use('hflx_rain_cea') ) THEN
1161         ztmp(:,:) = rcp * ( SUM( (ptsu-rt0) * a_i_b, dim=3 ) + sst_m(:,:) * ( 1._wp - at_i_b(:,:) ) )
1162         IF( iom_use('hflx_rain_cea') )  CALL iom_put( 'hflx_rain_cea', ( tprecip(:,:) - sprecip(:,:) ) * ztmp(:,:) )   ! heat flux from rain (cell average)
1163      ENDIF
1164      IF( iom_use('hflx_snow_cea') .OR. iom_use('hflx_snow_ao_cea') .OR. iom_use('hflx_snow_ai_cea')  )  THEN
1165         WHERE( SUM( a_i_b, dim=3 ) > 1.e-10 )
1166            ztmp(:,:) = rcpi * SUM( (ptsu-rt0) * a_i_b, dim=3 ) / SUM( a_i_b, dim=3 )
1167         ELSEWHERE
1168            ztmp(:,:) = rcp * sst_m(:,:)
1169         ENDWHERE
1170         ztmp2(:,:) = sprecip(:,:) * ( ztmp(:,:) - rLfus )
1171         IF( iom_use('hflx_snow_cea')    ) CALL iom_put('hflx_snow_cea'   , ztmp2(:,:) ) ! heat flux from snow (cell average)
1172         IF( iom_use('hflx_snow_ao_cea') ) CALL iom_put('hflx_snow_ao_cea', ztmp2(:,:) * ( 1._wp - zsnw(:,:) ) ) ! heat flux from snow (over ocean)
1173         IF( iom_use('hflx_snow_ai_cea') ) CALL iom_put('hflx_snow_ai_cea', ztmp2(:,:) *           zsnw(:,:)   ) ! heat flux from snow (over ice)
1174      ENDIF
1175      !
1176      IF(sn_cfctl%l_prtctl) THEN
1177         CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' blk_ice: qla_ice  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb    : ', kdim=jpl)
1178         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice: z_qlw    : ', tab3d_2=dqla_ice, clinfo2=' dqla_ice : ', kdim=jpl)
1179         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice: z_dqsb   : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw   : ', kdim=jpl)
1180         CALL prt_ctl(tab3d_1=dqns_ice, clinfo1=' blk_ice: dqns_ice : ', tab3d_2=qsr_ice , clinfo2=' qsr_ice  : ', kdim=jpl)
1181         CALL prt_ctl(tab3d_1=ptsu    , clinfo1=' blk_ice: ptsu     : ', tab3d_2=qns_ice , clinfo2=' qns_ice  : ', kdim=jpl)
1182         CALL prt_ctl(tab2d_1=tprecip , clinfo1=' blk_ice: tprecip  : ', tab2d_2=sprecip , clinfo2=' sprecip  : ')
1183      ENDIF
1184      !
1185   END SUBROUTINE blk_ice_2
1186
1187
1188   SUBROUTINE blk_ice_qcn( ld_virtual_itd, ptsu, ptb, phs, phi )
1189      !!---------------------------------------------------------------------
1190      !!                     ***  ROUTINE blk_ice_qcn  ***
1191      !!
1192      !! ** Purpose :   Compute surface temperature and snow/ice conduction flux
1193      !!                to force sea ice / snow thermodynamics
1194      !!                in the case conduction flux is emulated
1195      !!
1196      !! ** Method  :   compute surface energy balance assuming neglecting heat storage
1197      !!                following the 0-layer Semtner (1976) approach
1198      !!
1199      !! ** Outputs : - ptsu    : sea-ice / snow surface temperature (K)
1200      !!              - qcn_ice : surface inner conduction flux (W/m2)
1201      !!
1202      !!---------------------------------------------------------------------
1203      LOGICAL                   , INTENT(in   ) ::   ld_virtual_itd  ! single-category option
1204      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   ptsu            ! sea ice / snow surface temperature
1205      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   ptb             ! sea ice base temperature
1206      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   phs             ! snow thickness
1207      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   phi             ! sea ice thickness
1208      !
1209      INTEGER , PARAMETER ::   nit = 10                  ! number of iterations
1210      REAL(wp), PARAMETER ::   zepsilon = 0.1_wp         ! characteristic thickness for enhanced conduction
1211      !
1212      INTEGER  ::   ji, jj, jl           ! dummy loop indices
1213      INTEGER  ::   iter                 ! local integer
1214      REAL(wp) ::   zfac, zfac2, zfac3   ! local scalars
1215      REAL(wp) ::   zkeff_h, ztsu, ztsu0 !
1216      REAL(wp) ::   zqc, zqnet           !
1217      REAL(wp) ::   zhe, zqa0            !
1218      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   zgfac   ! enhanced conduction factor
1219      !!---------------------------------------------------------------------
1220
1221      ! -------------------------------------!
1222      !      I   Enhanced conduction factor  !
1223      ! -------------------------------------!
1224      ! Emulates the enhancement of conduction by unresolved thin ice (ld_virtual_itd = T)
1225      ! Fichefet and Morales Maqueda, JGR 1997
1226      !
1227      zgfac(:,:,:) = 1._wp
1228
1229      IF( ld_virtual_itd ) THEN
1230         !
1231         zfac  = 1._wp /  ( rn_cnd_s + rcnd_i )
1232         zfac2 = EXP(1._wp) * 0.5_wp * zepsilon
1233         zfac3 = 2._wp / zepsilon
1234         !
1235         DO jl = 1, jpl
1236            DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
1237               zhe = ( rn_cnd_s * phi(ji,jj,jl) + rcnd_i * phs(ji,jj,jl) ) * zfac                            ! Effective thickness
1238               IF( zhe >=  zfac2 )   zgfac(ji,jj,jl) = MIN( 2._wp, 0.5_wp * ( 1._wp + LOG( zhe * zfac3 ) ) ) ! Enhanced conduction factor
1239            END_2D
1240         END DO
1241         !
1242      ENDIF
1243
1244      ! -------------------------------------------------------------!
1245      !      II   Surface temperature and conduction flux            !
1246      ! -------------------------------------------------------------!
1247      !
1248      zfac = rcnd_i * rn_cnd_s
1249      !
1250      DO jl = 1, jpl
1251         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
1252            !
1253            zkeff_h = zfac * zgfac(ji,jj,jl) / &                                    ! Effective conductivity of the snow-ice system divided by thickness
1254               &      ( rcnd_i * phs(ji,jj,jl) + rn_cnd_s * MAX( 0.01, phi(ji,jj,jl) ) )
1255            ztsu    = ptsu(ji,jj,jl)                                                ! Store current iteration temperature
1256            ztsu0   = ptsu(ji,jj,jl)                                                ! Store initial surface temperature
1257            zqa0    = qsr_ice(ji,jj,jl) - qtr_ice_top(ji,jj,jl) + qns_ice(ji,jj,jl) ! Net initial atmospheric heat flux
1258            !
1259            DO iter = 1, nit     ! --- Iterative loop
1260               zqc   = zkeff_h * ( ztsu - ptb(ji,jj) )                              ! Conduction heat flux through snow-ice system (>0 downwards)
1261               zqnet = zqa0 + dqns_ice(ji,jj,jl) * ( ztsu - ptsu(ji,jj,jl) ) - zqc  ! Surface energy budget
1262               ztsu  = ztsu - zqnet / ( dqns_ice(ji,jj,jl) - zkeff_h )              ! Temperature update
1263            END DO
1264            !
1265            ptsu   (ji,jj,jl) = MIN( rt0, ztsu )
1266            qcn_ice(ji,jj,jl) = zkeff_h * ( ptsu(ji,jj,jl) - ptb(ji,jj) )
1267            qns_ice(ji,jj,jl) = qns_ice(ji,jj,jl) + dqns_ice(ji,jj,jl) * ( ptsu(ji,jj,jl) - ztsu0 )
1268            qml_ice(ji,jj,jl) = ( qsr_ice(ji,jj,jl) - qtr_ice_top(ji,jj,jl) + qns_ice(ji,jj,jl) - qcn_ice(ji,jj,jl) )  &
1269               &   * MAX( 0._wp , SIGN( 1._wp, ptsu(ji,jj,jl) - rt0 ) )
1270
1271            ! --- Diagnose the heat loss due to changing non-solar flux (as in icethd_zdf_bl99) --- !
1272            hfx_err_dif(ji,jj) = hfx_err_dif(ji,jj) - ( dqns_ice(ji,jj,jl) * ( ptsu(ji,jj,jl) - ztsu0 ) ) * a_i_b(ji,jj,jl)
1273
1274         END_2D
1275         !
1276      END DO
1277      !
1278   END SUBROUTINE blk_ice_qcn
1279
1280
1281   SUBROUTINE Cdn10_Lupkes2012( pcd )
1282      !!----------------------------------------------------------------------
1283      !!                      ***  ROUTINE  Cdn10_Lupkes2012  ***
1284      !!
1285      !! ** Purpose :    Recompute the neutral air-ice drag referenced at 10m
1286      !!                 to make it dependent on edges at leads, melt ponds and flows.
1287      !!                 After some approximations, this can be resumed to a dependency
1288      !!                 on ice concentration.
1289      !!
1290      !! ** Method :     The parameterization is taken from Lupkes et al. (2012) eq.(50)
1291      !!                 with the highest level of approximation: level4, eq.(59)
1292      !!                 The generic drag over a cell partly covered by ice can be re-written as follows:
1293      !!
1294      !!                 Cd = Cdw * (1-A) + Cdi * A + Ce * (1-A)**(nu+1/(10*beta)) * A**mu
1295      !!
1296      !!                    Ce = 2.23e-3       , as suggested by Lupkes (eq. 59)
1297      !!                    nu = mu = beta = 1 , as suggested by Lupkes (eq. 59)
1298      !!                    A is the concentration of ice minus melt ponds (if any)
1299      !!
1300      !!                 This new drag has a parabolic shape (as a function of A) starting at
1301      !!                 Cdw(say 1.5e-3) for A=0, reaching 1.97e-3 for A~0.5
1302      !!                 and going down to Cdi(say 1.4e-3) for A=1
1303      !!
1304      !!                 It is theoretically applicable to all ice conditions (not only MIZ)
1305      !!                 => see Lupkes et al (2013)
1306      !!
1307      !! ** References : Lupkes et al. JGR 2012 (theory)
1308      !!                 Lupkes et al. GRL 2013 (application to GCM)
1309      !!
1310      !!----------------------------------------------------------------------
1311      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pcd
1312      REAL(wp), PARAMETER ::   zCe   = 2.23e-03_wp
1313      REAL(wp), PARAMETER ::   znu   = 1._wp
1314      REAL(wp), PARAMETER ::   zmu   = 1._wp
1315      REAL(wp), PARAMETER ::   zbeta = 1._wp
1316      REAL(wp)            ::   zcoef
1317      !!----------------------------------------------------------------------
1318      zcoef = znu + 1._wp / ( 10._wp * zbeta )
1319
1320      ! generic drag over a cell partly covered by ice
1321      !!Cd(:,:) = Cd_oce(:,:) * ( 1._wp - at_i_b(:,:) ) +  &                        ! pure ocean drag
1322      !!   &      Cd_ice      *           at_i_b(:,:)   +  &                        ! pure ice drag
1323      !!   &      zCe         * ( 1._wp - at_i_b(:,:) )**zcoef * at_i_b(:,:)**zmu   ! change due to sea-ice morphology
1324
1325      ! ice-atm drag
1326      pcd(:,:) = rCd_ice +  &                                                         ! pure ice drag
1327         &      zCe     * ( 1._wp - at_i_b(:,:) )**zcoef * at_i_b(:,:)**(zmu-1._wp)  ! change due to sea-ice morphology
1328
1329   END SUBROUTINE Cdn10_Lupkes2012
1330
1331
1332   SUBROUTINE Cdn10_Lupkes2015( ptm_su, pslp, pcd, pch )
1333      !!----------------------------------------------------------------------
1334      !!                      ***  ROUTINE  Cdn10_Lupkes2015  ***
1335      !!
1336      !! ** pUrpose :    Alternative turbulent transfert coefficients formulation
1337      !!                 between sea-ice and atmosphere with distinct momentum
1338      !!                 and heat coefficients depending on sea-ice concentration
1339      !!                 and atmospheric stability (no meltponds effect for now).
1340      !!
1341      !! ** Method :     The parameterization is adapted from Lupkes et al. (2015)
1342      !!                 and ECHAM6 atmospheric model. Compared to Lupkes2012 scheme,
1343      !!                 it considers specific skin and form drags (Andreas et al. 2010)
1344      !!                 to compute neutral transfert coefficients for both heat and
1345      !!                 momemtum fluxes. Atmospheric stability effect on transfert
1346      !!                 coefficient is also taken into account following Louis (1979).
1347      !!
1348      !! ** References : Lupkes et al. JGR 2015 (theory)
1349      !!                 Lupkes et al. ECHAM6 documentation 2015 (implementation)
1350      !!
1351      !!----------------------------------------------------------------------
1352      !
1353      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   ptm_su ! sea-ice surface temperature [K]
1354      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   pslp   ! sea-level pressure [Pa]
1355      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pcd    ! momentum transfert coefficient
1356      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   pch    ! heat transfert coefficient
1357      REAL(wp), DIMENSION(jpi,jpj)            ::   zst, zqo_sat, zqi_sat
1358      !
1359      ! ECHAM6 constants
1360      REAL(wp), PARAMETER ::   z0_skin_ice  = 0.69e-3_wp  ! Eq. 43 [m]
1361      REAL(wp), PARAMETER ::   z0_form_ice  = 0.57e-3_wp  ! Eq. 42 [m]
1362      REAL(wp), PARAMETER ::   z0_ice       = 1.00e-3_wp  ! Eq. 15 [m]
1363      REAL(wp), PARAMETER ::   zce10        = 2.80e-3_wp  ! Eq. 41
1364      REAL(wp), PARAMETER ::   zbeta        = 1.1_wp      ! Eq. 41
1365      REAL(wp), PARAMETER ::   zc           = 5._wp       ! Eq. 13
1366      REAL(wp), PARAMETER ::   zc2          = zc * zc
1367      REAL(wp), PARAMETER ::   zam          = 2. * zc     ! Eq. 14
1368      REAL(wp), PARAMETER ::   zah          = 3. * zc     ! Eq. 30
1369      REAL(wp), PARAMETER ::   z1_alpha     = 1._wp / 0.2_wp  ! Eq. 51
1370      REAL(wp), PARAMETER ::   z1_alphaf    = z1_alpha    ! Eq. 56
1371      REAL(wp), PARAMETER ::   zbetah       = 1.e-3_wp    ! Eq. 26
1372      REAL(wp), PARAMETER ::   zgamma       = 1.25_wp     ! Eq. 26
1373      REAL(wp), PARAMETER ::   z1_gamma     = 1._wp / zgamma
1374      REAL(wp), PARAMETER ::   r1_3         = 1._wp / 3._wp
1375      !
1376      INTEGER  ::   ji, jj         ! dummy loop indices
1377      REAL(wp) ::   zthetav_os, zthetav_is, zthetav_zu
1378      REAL(wp) ::   zrib_o, zrib_i
1379      REAL(wp) ::   zCdn_skin_ice, zCdn_form_ice, zCdn_ice
1380      REAL(wp) ::   zChn_skin_ice, zChn_form_ice
1381      REAL(wp) ::   z0w, z0i, zfmi, zfmw, zfhi, zfhw
1382      REAL(wp) ::   zCdn_form_tmp
1383      !!----------------------------------------------------------------------
1384
1385      ! Momentum Neutral Transfert Coefficients (should be a constant)
1386      zCdn_form_tmp = zce10 * ( LOG( 10._wp / z0_form_ice + 1._wp ) / LOG( rn_zu / z0_form_ice + 1._wp ) )**2   ! Eq. 40
1387      zCdn_skin_ice = ( vkarmn                                      / LOG( rn_zu / z0_skin_ice + 1._wp ) )**2   ! Eq. 7
1388      zCdn_ice      = zCdn_skin_ice   ! Eq. 7
1389      !zCdn_ice     = 1.89e-3         ! old ECHAM5 value (cf Eq. 32)
1390
1391      ! Heat Neutral Transfert Coefficients
1392      zChn_skin_ice = vkarmn**2 / ( LOG( rn_zu / z0_ice + 1._wp ) * LOG( rn_zu * z1_alpha / z0_skin_ice + 1._wp ) )   ! Eq. 50 + Eq. 52
1393
1394      ! Atmospheric and Surface Variables
1395      zst(:,:)     = sst_m(:,:) + rt0                                        ! convert SST from Celcius to Kelvin
1396      zqo_sat(:,:) = rdct_qsat_salt * q_sat( zst(:,:)   , pslp(:,:) )   ! saturation humidity over ocean [kg/kg]
1397      zqi_sat(:,:) =                  q_sat( ptm_su(:,:), pslp(:,:) )   ! saturation humidity over ice   [kg/kg]
1398      !
1399      DO_2D( 0, 0, 0, 0 )
1400         ! Virtual potential temperature [K]
1401         zthetav_os = zst(ji,jj)    * ( 1._wp + rctv0 * zqo_sat(ji,jj) )   ! over ocean
1402         zthetav_is = ptm_su(ji,jj) * ( 1._wp + rctv0 * zqi_sat(ji,jj) )   ! ocean ice
1403         zthetav_zu = t_zu (ji,jj)  * ( 1._wp + rctv0 * q_zu(ji,jj)    )   ! at zu
1404
1405         ! Bulk Richardson Number (could use Ri_bulk function from aerobulk instead)
1406         zrib_o = grav / zthetav_os * ( zthetav_zu - zthetav_os ) * rn_zu / MAX( 0.5, wndm(ji,jj)     )**2   ! over ocean
1407         zrib_i = grav / zthetav_is * ( zthetav_zu - zthetav_is ) * rn_zu / MAX( 0.5, wndm_ice(ji,jj) )**2   ! over ice
1408
1409         ! Momentum and Heat Neutral Transfert Coefficients
1410         zCdn_form_ice = zCdn_form_tmp * at_i_b(ji,jj) * ( 1._wp - at_i_b(ji,jj) )**zbeta  ! Eq. 40
1411         zChn_form_ice = zCdn_form_ice / ( 1._wp + ( LOG( z1_alphaf ) / vkarmn ) * SQRT( zCdn_form_ice ) )               ! Eq. 53
1412
1413         ! Momentum and Heat Stability functions (possibility to use psi_m_ecmwf instead ?)
1414         z0w = rn_zu * EXP( -1._wp * vkarmn / SQRT( Cdn_oce(ji,jj) ) ) ! over water
1415         z0i = z0_skin_ice                                             ! over ice
1416         IF( zrib_o <= 0._wp ) THEN
1417            zfmw = 1._wp - zam * zrib_o / ( 1._wp + 3._wp * zc2 * Cdn_oce(ji,jj) * SQRT( -zrib_o * ( rn_zu / z0w + 1._wp ) ) )  ! Eq. 10
1418            zfhw = ( 1._wp + ( zbetah * ( zthetav_os - zthetav_zu )**r1_3 / ( Chn_oce(ji,jj) * MAX(0.01, wndm(ji,jj)) )   &     ! Eq. 26
1419               &             )**zgamma )**z1_gamma
1420         ELSE
1421            zfmw = 1._wp / ( 1._wp + zam * zrib_o / SQRT( 1._wp + zrib_o ) )   ! Eq. 12
1422            zfhw = 1._wp / ( 1._wp + zah * zrib_o / SQRT( 1._wp + zrib_o ) )   ! Eq. 28
1423         ENDIF
1424
1425         IF( zrib_i <= 0._wp ) THEN
1426            zfmi = 1._wp - zam * zrib_i / (1._wp + 3._wp * zc2 * zCdn_ice * SQRT( -zrib_i * ( rn_zu / z0i + 1._wp)))   ! Eq.  9
1427            zfhi = 1._wp - zah * zrib_i / (1._wp + 3._wp * zc2 * zCdn_ice * SQRT( -zrib_i * ( rn_zu / z0i + 1._wp)))   ! Eq. 25
1428         ELSE
1429            zfmi = 1._wp / ( 1._wp + zam * zrib_i / SQRT( 1._wp + zrib_i ) )   ! Eq. 11
1430            zfhi = 1._wp / ( 1._wp + zah * zrib_i / SQRT( 1._wp + zrib_i ) )   ! Eq. 27
1431         ENDIF
1432
1433         ! Momentum Transfert Coefficients (Eq. 38)
1434         pcd(ji,jj) = zCdn_skin_ice *   zfmi +  &
1435            &        zCdn_form_ice * ( zfmi * at_i_b(ji,jj) + zfmw * ( 1._wp - at_i_b(ji,jj) ) ) / MAX( 1.e-06, at_i_b(ji,jj) )
1436
1437         ! Heat Transfert Coefficients (Eq. 49)
1438         pch(ji,jj) = zChn_skin_ice *   zfhi +  &
1439            &        zChn_form_ice * ( zfhi * at_i_b(ji,jj) + zfhw * ( 1._wp - at_i_b(ji,jj) ) ) / MAX( 1.e-06, at_i_b(ji,jj) )
1440         !
1441      END_2D
1442      CALL lbc_lnk_multi( 'sbcblk', pcd, 'T',  1.0_wp, pch, 'T', 1.0_wp )
1443      !
1444   END SUBROUTINE Cdn10_Lupkes2015
1445
1446#endif
1447
1448   !!======================================================================
1449END MODULE sbcblk
Note: See TracBrowser for help on using the repository browser.