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 branches/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2016/dev_merge_2016/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk.F90 @ 7421

Last change on this file since 7421 was 7421, checked in by flavoni, 8 years ago

#1811 merge dev_CNRS_MERATOR_2016 with dev_merge_2016 branch

File size: 56.3 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 (http://aerobulk.sourceforge.net/)
18   !!            4.0  !  2016-10  (G. Madec)  introduce a sbc_blk_init routine
19   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   sbc_blk_init  : initialisation of the chosen bulk formulation as ocean surface boundary condition
23   !!   sbc_blk       : bulk formulation as ocean surface boundary condition
24   !!   blk_oce       : computes momentum, heat and freshwater fluxes over ocean
25   !!   blk_ice       : computes momentum, heat and freshwater fluxes over sea ice
26   !!   rho_air       : density of (moist) air (depends on T_air, q_air and SLP
27   !!   cp_air        : specific heat of (moist) air (depends spec. hum. q_air)
28   !!   q_sat         : saturation humidity as a function of SLP and temperature
29   !!   L_vap         : latent heat of vaporization of water as a function of temperature
30   !!----------------------------------------------------------------------
31   USE oce            ! ocean dynamics and tracers
32   USE dom_oce        ! ocean space and time domain
33   USE phycst         ! physical constants
34   USE fldread        ! read input fields
35   USE sbc_oce        ! Surface boundary condition: ocean fields
36   USE cyclone        ! Cyclone 10m wind form trac of cyclone centres
37   USE sbcdcy         ! surface boundary condition: diurnal cycle
38   USE sbcwave , ONLY :   cdn_wave ! wave module
39   USE sbc_ice        ! Surface boundary condition: ice fields
40   USE lib_fortran    ! to use key_nosignedzero
41#if defined key_lim3
42   USE ice     , ONLY :   u_ice, v_ice, jpl, pfrld, a_i_b, at_i_b
43   USE limthd_dh      ! for CALL lim_thd_snwblow
44#elif defined key_lim2
45   USE ice_2   , ONLY :   u_ice, v_ice
46   USE par_ice_2      ! LIM-2 parameters
47#endif
48   USE sbcblk_algo_ncar     ! => turb_ncar     : NCAR - CORE (Large & Yeager, 2009)
49   USE sbcblk_algo_coare    ! => turb_coare    : COAREv3.0 (Fairall et al. 2003)
50   USE sbcblk_algo_coare3p5 ! => turb_coare3p5 : COAREv3.5 (Edson et al. 2013)
51   USE sbcblk_algo_ecmwf    ! => turb_ecmwf    : ECMWF (IFS cycle 31)
52   !
53   USE iom            ! I/O manager library
54   USE in_out_manager ! I/O manager
55   USE lib_mpp        ! distribued memory computing library
56   USE wrk_nemo       ! work arrays
57   USE timing         ! Timing
58   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
59   USE prtctl         ! Print control
60
61   IMPLICIT NONE
62   PRIVATE
63
64   PUBLIC   sbc_blk_init  ! called in sbcmod
65   PUBLIC   sbc_blk       ! called in sbcmod
66#if defined key_lim2 || defined key_lim3
67   PUBLIC   blk_ice_tau   ! routine called in sbc_ice_lim module
68   PUBLIC   blk_ice_flx   ! routine called in sbc_ice_lim module
69#endif
70
71!!Lolo: should ultimately be moved in the module with all physical constants ?
72!!gm  : In principle, yes.
73   REAL(wp), PARAMETER ::   Cp_dry = 1005.0       !: Specic heat of dry air, constant pressure      [J/K/kg]
74   REAL(wp), PARAMETER ::   Cp_vap = 1860.0       !: Specic heat of water vapor, constant pressure  [J/K/kg]
75   REAL(wp), PARAMETER ::   R_dry = 287.05_wp     !: Specific gas constant for dry air              [J/K/kg]
76   REAL(wp), PARAMETER ::   R_vap = 461.495_wp    !: Specific gas constant for water vapor          [J/K/kg]
77   REAL(wp), PARAMETER ::   reps0 = R_dry/R_vap   !: ratio of gas constant for dry air and water vapor => ~ 0.622
78   REAL(wp), PARAMETER ::   rctv0 = R_vap/R_dry   !: for virtual temperature (== (1-eps)/eps) => ~ 0.608
79
80   INTEGER , PARAMETER ::   jpfld   =10           ! maximum number of files to read
81   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point
82   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point
83   INTEGER , PARAMETER ::   jp_tair = 3           ! index of 10m air temperature             (Kelvin)
84   INTEGER , PARAMETER ::   jp_humi = 4           ! index of specific humidity               ( % )
85   INTEGER , PARAMETER ::   jp_qsr  = 5           ! index of solar heat                      (W/m2)
86   INTEGER , PARAMETER ::   jp_qlw  = 6           ! index of Long wave                       (W/m2)
87   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s)
88   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s)
89   INTEGER , PARAMETER ::   jp_slp  = 9           ! index of sea level pressure              (Pa)
90   INTEGER , PARAMETER ::   jp_tdif =10           ! index of tau diff associated to HF tau   (N/m2)   at T-point
91
92   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read)
93
94   !                                             !!! Bulk parameters
95   REAL(wp), PARAMETER ::   cpa    = 1000.5         ! specific heat of air (only used for ice fluxes now...)
96   REAL(wp), PARAMETER ::   Ls     =    2.839e6     ! latent heat of sublimation
97   REAL(wp), PARAMETER ::   Stef   =    5.67e-8     ! Stefan Boltzmann constant
98   REAL(wp), PARAMETER ::   Cd_ice =    1.4e-3      ! iovi 1.63e-3     ! transfer coefficient over ice
99   REAL(wp), PARAMETER ::   albo   =    0.066       ! ocean albedo assumed to be constant
100   !
101   !                           !!* Namelist namsbc_blk : bulk parameters
102   LOGICAL  ::   ln_NCAR        ! "NCAR"      algorithm   (Large and Yeager 2008)
103   LOGICAL  ::   ln_COARE_3p0   ! "COARE 3.0" algorithm   (Fairall et al. 2003)
104   LOGICAL  ::   ln_COARE_3p5   ! "COARE 3.5" algorithm   (Edson et al. 2013)
105   LOGICAL  ::   ln_ECMWF       ! "ECMWF"     algorithm   (IFS cycle 31)
106   !
107   LOGICAL  ::   ln_taudif      ! logical flag to use the "mean of stress module - module of mean stress" data
108   REAL(wp) ::   rn_pfac        ! multiplication factor for precipitation
109   REAL(wp) ::   rn_efac        ! multiplication factor for evaporation (clem)
110   REAL(wp) ::   rn_vfac        ! multiplication factor for ice/ocean velocity in the calculation of wind stress (clem)
111   REAL(wp) ::   rn_zqt         ! z(q,t) : height of humidity and temperature measurements
112   REAL(wp) ::   rn_zu          ! z(u)   : height of wind measurements
113   LOGICAL  ::   ln_Cd_L12 = .FALSE. !  Modify the drag ice-atm and oce-atm depending on ice concentration (from Lupkes et al. JGR2012)
114   !
115   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   Cd_oce   ! air-ocean drag (clem)
116
117   INTEGER  ::   nblk           ! choice of the bulk algorithm
118   !                            ! associated indices:
119   INTEGER, PARAMETER ::   np_NCAR      = 1   ! "NCAR" algorithm        (Large and Yeager 2008)
120   INTEGER, PARAMETER ::   np_COARE_3p0 = 2   ! "COARE 3.0" algorithm   (Fairall et al. 2003)
121   INTEGER, PARAMETER ::   np_COARE_3p5 = 3   ! "COARE 3.5" algorithm   (Edson et al. 2013)
122   INTEGER, PARAMETER ::   np_ECMWF     = 4   ! "ECMWF" algorithm       (IFS cycle 31)
123
124   !! * Substitutions
125#  include "vectopt_loop_substitute.h90"
126   !!----------------------------------------------------------------------
127   !! NEMO/OPA 3.7 , NEMO-consortium (2014)
128   !! $Id: sbcblk.F90 6416 2016-04-01 12:22:17Z clem $
129   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
130   !!----------------------------------------------------------------------
131CONTAINS
132
133   INTEGER FUNCTION sbc_blk_alloc()
134      !!-------------------------------------------------------------------
135      !!             ***  ROUTINE sbc_blk_alloc ***
136      !!-------------------------------------------------------------------
137      ALLOCATE( Cd_oce(jpi,jpj) , STAT=sbc_blk_alloc )
138      !
139      IF( lk_mpp             )   CALL mpp_sum ( sbc_blk_alloc )
140      IF( sbc_blk_alloc /= 0 )   CALL ctl_warn('sbc_blk_alloc: failed to allocate arrays')
141   END FUNCTION sbc_blk_alloc
142
143   SUBROUTINE sbc_blk_init
144      !!---------------------------------------------------------------------
145      !!                    ***  ROUTINE sbc_blk_init  ***
146      !!
147      !! ** Purpose :   choose and initialize a bulk formulae formulation
148      !!
149      !! ** Method  :
150      !!
151      !!      C A U T I O N : never mask the surface stress fields
152      !!                      the stress is assumed to be in the (i,j) mesh referential
153      !!
154      !! ** Action  :   
155      !!
156      !!----------------------------------------------------------------------
157      INTEGER  ::   ifpr, jfld            ! dummy loop indice and argument
158      INTEGER  ::   ios, ierror, ioptio   ! Local integer
159      !!
160      CHARACTER(len=100)            ::   cn_dir                ! Root directory for location of atmospheric forcing files
161      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i                 ! array of namelist informations on the fields to read
162      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read
163      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !       "                        "
164      TYPE(FLD_N) ::   sn_slp , sn_tdif                        !       "                        "
165      NAMELIST/namsbc_blk/ sn_wndi, sn_wndj, sn_humi, sn_qsr, sn_qlw ,                &   ! input fields
166         &                 sn_tair, sn_prec, sn_snow, sn_slp, sn_tdif,                &
167         &                 ln_NCAR, ln_COARE_3p0, ln_COARE_3p5, ln_ECMWF,             &   ! bulk algorithm
168         &                 cn_dir , ln_taudif, rn_zqt, rn_zu,                         & 
169         &                 rn_pfac, rn_efac, rn_vfac, ln_Cd_L12
170      !!---------------------------------------------------------------------
171      !
172      !                                      ! allocate sbc_blk_core array
173      IF( sbc_blk_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'sbc_blk : unable to allocate standard arrays' )
174      !
175      !                             !** read bulk namelist 
176      REWIND( numnam_ref )                !* Namelist namsbc_blk in reference namelist : bulk parameters
177      READ  ( numnam_ref, namsbc_blk, IOSTAT = ios, ERR = 901)
178901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_blk in reference namelist', lwp )
179      !
180      REWIND( numnam_cfg )                !* Namelist namsbc_blk in configuration namelist : bulk parameters
181      READ  ( numnam_cfg, namsbc_blk, IOSTAT = ios, ERR = 902 )
182902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_blk in configuration namelist', lwp )
183      !
184      IF(lwm) WRITE( numond, namsbc_blk )
185      !
186      !                             !** initialization of the chosen bulk formulae (+ check)
187      !                                   !* select the bulk chosen in the namelist and check the choice
188      ;                                                        ioptio = 0
189      IF( ln_NCAR      ) THEN   ;   nblk =  np_NCAR        ;   ioptio = ioptio + 1   ;   ENDIF
190      IF( ln_COARE_3p0 ) THEN   ;   nblk =  np_COARE_3p0   ;   ioptio = ioptio + 1   ;   ENDIF
191      IF( ln_COARE_3p5 ) THEN   ;   nblk =  np_COARE_3p5   ;   ioptio = ioptio + 1   ;   ENDIF
192      IF( ln_ECMWF     ) THEN   ;   nblk =  np_ECMWF       ;   ioptio = ioptio + 1   ;   ENDIF
193      !
194      IF( ioptio /= 1 )   CALL ctl_stop( 'sbc_blk_init: Choose one and only one bulk algorithm' )
195      !
196      IF( ln_dm2dc ) THEN                 !* check: diurnal cycle on Qsr
197         IF( sn_qsr%nfreqh /= 24 )   CALL ctl_stop( 'sbc_blk_init: ln_dm2dc=T only with daily short-wave input' )
198         IF( sn_qsr%ln_tint ) THEN
199            CALL ctl_warn( 'sbc_blk_init: ln_dm2dc=T daily qsr time interpolation done by sbcdcy module',   &
200               &           '              ==> We force time interpolation = .false. for qsr' )
201            sn_qsr%ln_tint = .false.
202         ENDIF
203      ENDIF
204      !                                   !* set the bulk structure
205      !                                      !- store namelist information in an array
206      slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
207      slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
208      slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
209      slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
210      slf_i(jp_slp)  = sn_slp    ;   slf_i(jp_tdif) = sn_tdif
211      !
212      lhftau = ln_taudif                     !- add an extra field if HF stress is used
213      jfld = jpfld - COUNT( (/.NOT.lhftau/) )
214      !
215      !                                      !- allocate the bulk structure
216      ALLOCATE( sf(jfld), STAT=ierror )
217      IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_init: unable to allocate sf structure' )
218      DO ifpr= 1, jfld
219         ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
220         IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
221      END DO
222      !                                      !- fill the bulk structure with namelist informations
223      CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_init', 'surface boundary condition -- bulk formulae', 'namsbc_blk' )
224      !
225      !           
226      IF(lwp) THEN                     !** Control print
227         !
228         WRITE(numout,*)                  !* namelist
229         WRITE(numout,*) '   Namelist namsbc_blk (other than data information):'
230         WRITE(numout,*) '      "NCAR"      algorithm   (Large and Yeager 2008)     ln_NCAR      = ', ln_NCAR
231         WRITE(numout,*) '      "COARE 3.0" algorithm   (Fairall et al. 2003)       ln_COARE_3p0 = ', ln_COARE_3p0
232         WRITE(numout,*) '      "COARE 3.5" algorithm   (Edson et al. 2013)         ln_COARE_3p5 = ', ln_COARE_3p0
233         WRITE(numout,*) '      "ECMWF"     algorithm   (IFS cycle 31)              ln_ECMWF     = ', ln_ECMWF
234         WRITE(numout,*) '      add High freq.contribution to the stress module     ln_taudif    = ', ln_taudif
235         WRITE(numout,*) '      Air temperature and humidity reference height (m)   rn_zqt       = ', rn_zqt
236         WRITE(numout,*) '      Wind vector reference height (m)                    rn_zu        = ', rn_zu
237         WRITE(numout,*) '      factor applied on precipitation (total & snow)      rn_pfac      = ', rn_pfac
238         WRITE(numout,*) '      factor applied on evaporation                       rn_efac      = ', rn_efac
239         WRITE(numout,*) '      factor applied on ocean/ice velocity                rn_vfac      = ', rn_vfac
240         WRITE(numout,*) '         (form absolute (=0) to relative winds(=1))'
241         !
242         WRITE(numout,*)
243         SELECT CASE( nblk )              !* Print the choice of bulk algorithm
244         CASE( np_NCAR      )   ;   WRITE(numout,*) '      ===>>   "NCAR" algorithm        (Large and Yeager 2008)'
245         CASE( np_COARE_3p0 )   ;   WRITE(numout,*) '      ===>>   "COARE 3.0" algorithm   (Fairall et al. 2003)'
246         CASE( np_COARE_3p5 )   ;   WRITE(numout,*) '      ===>>   "COARE 3.5" algorithm   (Edson et al. 2013)'
247         CASE( np_ECMWF     )   ;   WRITE(numout,*) '      ===>>   "ECMWF" algorithm       (IFS cycle 31)'
248         END SELECT
249         !
250      ENDIF
251      !
252   END SUBROUTINE sbc_blk_init
253
254
255   SUBROUTINE sbc_blk( kt )
256      !!---------------------------------------------------------------------
257      !!                    ***  ROUTINE sbc_blk  ***
258      !!
259      !! ** Purpose :   provide at each time step the surface ocean fluxes
260      !!      (momentum, heat, freshwater and runoff)
261      !!
262      !! ** Method  : (1) READ each fluxes in NetCDF files:
263      !!      the 10m wind velocity (i-component) (m/s)    at T-point
264      !!      the 10m wind velocity (j-component) (m/s)    at T-point
265      !!      the 10m or 2m specific humidity     ( % )
266      !!      the solar heat                      (W/m2)
267      !!      the Long wave                       (W/m2)
268      !!      the 10m or 2m air temperature       (Kelvin)
269      !!      the total precipitation (rain+snow) (Kg/m2/s)
270      !!      the snow (solid prcipitation)       (kg/m2/s)
271      !!      the tau diff associated to HF tau   (N/m2)   at T-point   (ln_taudif=T)
272      !!              (2) CALL blk_oce
273      !!
274      !!      C A U T I O N : never mask the surface stress fields
275      !!                      the stress is assumed to be in the (i,j) mesh referential
276      !!
277      !! ** Action  :   defined at each time-step at the air-sea interface
278      !!              - utau, vtau  i- and j-component of the wind stress
279      !!              - taum        wind stress module at T-point
280      !!              - wndm        wind speed  module at T-point over free ocean or leads in presence of sea-ice
281      !!              - qns, qsr    non-solar and solar heat fluxes
282      !!              - emp         upward mass flux (evapo. - precip.)
283      !!              - sfx         salt flux due to freezing/melting (non-zero only if ice is present)
284      !!                            (set in limsbc(_2).F90)
285      !!
286      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
287      !!                   Brodeau et al. Ocean Modelling 2010
288      !!----------------------------------------------------------------------
289      INTEGER, INTENT(in) ::   kt   ! ocean time step
290      !!---------------------------------------------------------------------
291      !
292      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
293      !
294      !                                            ! compute the surface ocean fluxes using bulk formulea
295      IF( MOD( kt - 1, nn_fsbc ) == 0 )   CALL blk_oce( kt, sf, sst_m, ssu_m, ssv_m )
296
297#if defined key_cice
298      IF( MOD( kt - 1, nn_fsbc ) == 0 )   THEN
299         qlw_ice(:,:,1)   = sf(jp_qlw )%fnow(:,:,1)
300         IF( ln_dm2dc ) THEN ; qsr_ice(:,:,1) = sbc_dcy( sf(jp_qsr)%fnow(:,:,1) )
301         ELSE                ; qsr_ice(:,:,1) =          sf(jp_qsr)%fnow(:,:,1) 
302         ENDIF
303         tatm_ice(:,:)    = sf(jp_tair)%fnow(:,:,1)
304         qatm_ice(:,:)    = sf(jp_humi)%fnow(:,:,1)
305         tprecip(:,:)     = sf(jp_prec)%fnow(:,:,1) * rn_pfac
306         sprecip(:,:)     = sf(jp_snow)%fnow(:,:,1) * rn_pfac
307         wndi_ice(:,:)    = sf(jp_wndi)%fnow(:,:,1)
308         wndj_ice(:,:)    = sf(jp_wndj)%fnow(:,:,1)
309      ENDIF
310#endif
311      !
312   END SUBROUTINE sbc_blk
313
314
315   SUBROUTINE blk_oce( kt, sf, pst, pu, pv )
316      !!---------------------------------------------------------------------
317      !!                     ***  ROUTINE blk_oce  ***
318      !!
319      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
320      !!      the ocean surface at each time step
321      !!
322      !! ** Method  :   bulk formulea for the ocean using atmospheric
323      !!      fields read in sbc_read
324      !!
325      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
326      !!              - vtau    : j-component of the stress at V-point  (N/m2)
327      !!              - taum    : Wind stress module at T-point         (N/m2)
328      !!              - wndm    : Wind speed module at T-point          (m/s)
329      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
330      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
331      !!              - emp     : evaporation minus precipitation       (kg/m2/s)
332      !!
333      !!  ** Nota  :   sf has to be a dummy argument for AGRIF on NEC
334      !!---------------------------------------------------------------------
335      INTEGER  , INTENT(in   )                 ::   kt    ! time step index
336      TYPE(fld), INTENT(inout), DIMENSION(:)   ::   sf    ! input data
337      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pst   ! surface temperature                      [Celcius]
338      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pu    ! surface current at U-point (i-component) [m/s]
339      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pv    ! surface current at V-point (j-component) [m/s]
340      !
341      INTEGER  ::   ji, jj               ! dummy loop indices
342      REAL(wp) ::   zztmp                ! local variable
343      REAL(wp), DIMENSION(:,:), POINTER ::   zwnd_i, zwnd_j    ! wind speed components at T-point
344      REAL(wp), DIMENSION(:,:), POINTER ::   zsq               ! specific humidity at pst
345      REAL(wp), DIMENSION(:,:), POINTER ::   zqlw, zqsb        ! long wave and sensible heat fluxes
346      REAL(wp), DIMENSION(:,:), POINTER ::   zqla, zevap       ! latent heat fluxes and evaporation
347      REAL(wp), DIMENSION(:,:), POINTER ::   Cd                ! transfer coefficient for momentum      (tau)
348      REAL(wp), DIMENSION(:,:), POINTER ::   Ch                ! transfer coefficient for sensible heat (Q_sens)
349      REAL(wp), DIMENSION(:,:), POINTER ::   Ce                ! tansfert coefficient for evaporation   (Q_lat)
350      REAL(wp), DIMENSION(:,:), POINTER ::   zst               ! surface temperature in Kelvin
351      REAL(wp), DIMENSION(:,:), POINTER ::   zt_zu             ! air temperature at wind speed height
352      REAL(wp), DIMENSION(:,:), POINTER ::   zq_zu             ! air spec. hum.  at wind speed height
353      REAL(wp), DIMENSION(:,:), POINTER ::   zU_zu             ! bulk wind speed at height zu  [m/s]
354      REAL(wp), DIMENSION(:,:), POINTER ::   ztpot             ! potential temperature of air at z=rn_zqt [K]
355      REAL(wp), DIMENSION(:,:), POINTER ::   zrhoa             ! density of air   [kg/m^3]
356      !!---------------------------------------------------------------------
357      !
358      IF( nn_timing == 1 )  CALL timing_start('blk_oce')
359      !
360      CALL wrk_alloc( jpi,jpj,   zwnd_i, zwnd_j, zsq, zqlw, zqsb, zqla, zevap )
361      CALL wrk_alloc( jpi,jpj,   Cd, Ch, Ce, zst, zt_zu, zq_zu )
362      CALL wrk_alloc( jpi,jpj,   zU_zu, ztpot, zrhoa )
363      !
364
365      ! local scalars ( place there for vector optimisation purposes)
366      zst(:,:) = pst(:,:) + rt0      ! convert SST from Celcius to Kelvin (and set minimum value far above 0 K)
367
368      ! ----------------------------------------------------------------------------- !
369      !      0   Wind components and module at T-point relative to the moving ocean   !
370      ! ----------------------------------------------------------------------------- !
371
372      ! ... components ( U10m - U_oce ) at T-point (unmasked)
373!!gm    move zwnd_i (_j) set to zero  inside the key_cyclone ???
374      zwnd_i(:,:) = 0._wp
375      zwnd_j(:,:) = 0._wp
376#if defined key_cyclone
377      CALL wnd_cyc( kt, zwnd_i, zwnd_j )    ! add analytical tropical cyclone (Vincent et al. JGR 2012)
378      DO jj = 2, jpjm1
379         DO ji = fs_2, fs_jpim1   ! vect. opt.
380            sf(jp_wndi)%fnow(ji,jj,1) = sf(jp_wndi)%fnow(ji,jj,1) + zwnd_i(ji,jj)
381            sf(jp_wndj)%fnow(ji,jj,1) = sf(jp_wndj)%fnow(ji,jj,1) + zwnd_j(ji,jj)
382         END DO
383      END DO
384#endif
385      DO jj = 2, jpjm1
386         DO ji = fs_2, fs_jpim1   ! vect. opt.
387            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
388            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
389         END DO
390      END DO
391      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
392      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
393      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
394      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
395         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
396
397      ! ----------------------------------------------------------------------------- !
398      !      I   Radiative FLUXES                                                     !
399      ! ----------------------------------------------------------------------------- !
400
401      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
402      zztmp = 1. - albo
403      IF( ln_dm2dc ) THEN   ;   qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
404      ELSE                  ;   qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
405      ENDIF
406
407      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
408
409
410
411      ! ----------------------------------------------------------------------------- !
412      !     II    Turbulent FLUXES                                                    !
413      ! ----------------------------------------------------------------------------- !
414
415      ! ... specific humidity at SST and IST tmask(
416      zsq(:,:) = 0.98 * q_sat( zst(:,:), sf(jp_slp)%fnow(:,:,1) )
417      !!
418      !! Estimate of potential temperature at z=rn_zqt, based on adiabatic lapse-rate
419      !!    (see Josey, Gulev & Yu, 2013) / doi=10.1016/B978-0-12-391851-2.00005-2
420      !!    (since reanalysis products provide T at z, not theta !)
421      ztpot = sf(jp_tair)%fnow(:,:,1) + gamma_moist( sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1) ) * rn_zqt
422
423      SELECT CASE( nblk )        !==  transfer coefficients  ==!   Cd, Ch, Ce at T-point
424      !
425      CASE( np_NCAR      )   ;   CALL turb_ncar    ( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! NCAR-COREv2
426         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
427      CASE( np_COARE_3p0 )   ;   CALL turb_coare   ( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! COARE v3.0
428         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
429      CASE( np_COARE_3p5 )   ;   CALL turb_coare3p5( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! COARE v3.5
430         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
431      CASE( np_ECMWF     )   ;   CALL turb_ecmwf   ( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! ECMWF
432         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
433      CASE DEFAULT
434         CALL ctl_stop( 'STOP', 'sbc_oce: non-existing bulk formula selected' )
435      END SELECT
436
437      !                          ! Compute true air density :
438      IF( ABS(rn_zu - rn_zqt) > 0.01 ) THEN     ! At zu: (probably useless to remove zrho*grav*rn_zu from SLP...)
439         zrhoa(:,:) = rho_air( zt_zu(:,:)             , zq_zu(:,:)             , sf(jp_slp)%fnow(:,:,1) )
440      ELSE                                      ! At zt:
441         zrhoa(:,:) = rho_air( sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1) )
442      END IF
443
444      Cd_oce(:,:) = Cd(:,:)  ! record value of pure ocean-atm. drag (clem)
445
446      DO jj = 1, jpj             ! tau module, i and j component
447         DO ji = 1, jpi
448            zztmp = zrhoa(ji,jj)  * zU_zu(ji,jj) * Cd(ji,jj)   ! using bulk wind speed
449            taum  (ji,jj) = zztmp * wndm  (ji,jj)
450            zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj)
451            zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj)
452         END DO
453      END DO
454
455      !                          ! add the HF tau contribution to the wind stress module
456      IF( lhftau )   taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1)
457
458      CALL iom_put( "taum_oce", taum )   ! output wind stress module
459
460      ! ... utau, vtau at U- and V_points, resp.
461      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
462      !     Note the use of MAX(tmask(i,j),tmask(i+1,j) is to mask tau over ice shelves
463      DO jj = 1, jpjm1
464         DO ji = 1, fs_jpim1
465            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) ) &
466               &          * MAX(tmask(ji,jj,1),tmask(ji+1,jj,1))
467            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) ) &
468               &          * MAX(tmask(ji,jj,1),tmask(ji,jj+1,1))
469         END DO
470      END DO
471      CALL lbc_lnk( utau(:,:), 'U', -1. )
472      CALL lbc_lnk( vtau(:,:), 'V', -1. )
473
474
475      !  Turbulent fluxes over ocean
476      ! -----------------------------
477
478      ! zqla used as temporary array, for rho*U (common term of bulk formulae):
479      zqla(:,:) = zrhoa(:,:) * zU_zu(:,:)
480
481      IF( ABS( rn_zu - rn_zqt) < 0.01_wp ) THEN
482         !! q_air and t_air are given at 10m (wind reference height)
483         zevap(:,:) = rn_efac*MAX( 0._wp,             zqla(:,:)*Ce(:,:)*(zsq(:,:) - sf(jp_humi)%fnow(:,:,1)) ) ! Evaporation, using bulk wind speed
484         zqsb (:,:) = cp_air(sf(jp_humi)%fnow(:,:,1))*zqla(:,:)*Ch(:,:)*(zst(:,:) - ztpot(:,:)             )   ! Sensible Heat, using bulk wind speed
485      ELSE
486         !! q_air and t_air are not given at 10m (wind reference height)
487         ! Values of temp. and hum. adjusted to height of wind during bulk algorithm iteration must be used!!!
488         zevap(:,:) = rn_efac*MAX( 0._wp,             zqla(:,:)*Ce(:,:)*(zsq(:,:) - zq_zu(:,:) ) ) ! Evaporation ! using bulk wind speed
489         zqsb (:,:) = cp_air(sf(jp_humi)%fnow(:,:,1))*zqla(:,:)*Ch(:,:)*(zst(:,:) - zt_zu(:,:) )   ! Sensible Heat ! using bulk wind speed
490      ENDIF
491
492      zqla(:,:) = L_vap(zst(:,:)) * zevap(:,:)     ! Latent Heat flux
493
494
495      IF(ln_ctl) THEN
496         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
497         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
498         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
499         CALL prt_ctl( tab2d_1=zsq   , clinfo1=' blk_oce: zsq    : ', tab2d_2=zst, clinfo2=' zst : ' )
500         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce: utau   : ', mask1=umask,   &
501            &          tab2d_2=vtau  , clinfo2=           ' vtau : ', mask2=vmask )
502         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce: wndm   : ')
503         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce: zst    : ')
504      ENDIF
505
506      ! ----------------------------------------------------------------------------- !
507      !     III    Total FLUXES                                                       !
508      ! ----------------------------------------------------------------------------- !
509      !
510      emp (:,:) = (  zevap(:,:)                                          &   ! mass flux (evap. - precip.)
511         &         - sf(jp_prec)%fnow(:,:,1) * rn_pfac  ) * tmask(:,:,1)
512      !
513      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                &   ! Downward Non Solar
514         &     - sf(jp_snow)%fnow(:,:,1) * rn_pfac * lfus                         &   ! remove latent melting heat for solid precip
515         &     - zevap(:,:) * pst(:,:) * rcp                                      &   ! remove evap heat content at SST
516         &     + ( sf(jp_prec)%fnow(:,:,1) - sf(jp_snow)%fnow(:,:,1) ) * rn_pfac  &   ! add liquid precip heat content at Tair
517         &     * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp                          &
518         &     + sf(jp_snow)%fnow(:,:,1) * rn_pfac                                &   ! add solid  precip heat content at min(Tair,Tsnow)
519         &     * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1)
520      !
521#if defined key_lim3
522      qns_oce(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                ! non solar without emp (only needed by LIM3)
523      qsr_oce(:,:) = qsr(:,:)
524#endif
525      !
526      IF ( nn_ice == 0 ) THEN
527         CALL iom_put( "qlw_oce" ,   zqlw )                 ! output downward longwave heat over the ocean
528         CALL iom_put( "qsb_oce" , - zqsb )                 ! output downward sensible heat over the ocean
529         CALL iom_put( "qla_oce" , - zqla )                 ! output downward latent   heat over the ocean
530         CALL iom_put( "qemp_oce",   qns-zqlw+zqsb+zqla )   ! output downward heat content of E-P over the ocean
531         CALL iom_put( "qns_oce" ,   qns  )                 ! output downward non solar heat over the ocean
532         CALL iom_put( "qsr_oce" ,   qsr  )                 ! output downward solar heat over the ocean
533         CALL iom_put( "qt_oce"  ,   qns+qsr )              ! output total downward heat over the ocean
534         tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac   ! output total precipitation [kg/m2/s]
535         sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac   ! output solid precipitation [kg/m2/s]
536         CALL iom_put( 'snowpre', sprecip * 86400. )        ! Snow
537         CALL iom_put( 'precip' , tprecip * 86400. )        ! Total precipitation
538      ENDIF
539      !
540      IF(ln_ctl) THEN
541         CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce: zqsb   : ', tab2d_2=zqlw , clinfo2=' zqlw  : ')
542         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce: zqla   : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
543         CALL prt_ctl(tab2d_1=pst  , clinfo1=' blk_oce: pst    : ', tab2d_2=emp  , clinfo2=' emp   : ')
544         CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce: utau   : ', mask1=umask,   &
545            &         tab2d_2=vtau , clinfo2=              ' vtau  : ' , mask2=vmask )
546      ENDIF
547      !
548      CALL wrk_dealloc( jpi,jpj,   zwnd_i, zwnd_j, zsq, zqlw, zqsb, zqla, zevap )
549      CALL wrk_dealloc( jpi,jpj,   Cd, Ch, Ce, zst, zt_zu, zq_zu )
550      CALL wrk_dealloc( jpi,jpj,   zU_zu, ztpot, zrhoa )
551      !
552      IF( nn_timing == 1 )  CALL timing_stop('blk_oce')
553      !
554   END SUBROUTINE blk_oce
555
556#if defined key_lim2 || defined key_lim3
557
558   SUBROUTINE blk_ice_tau
559      !!---------------------------------------------------------------------
560      !!                     ***  ROUTINE blk_ice_tau  ***
561      !!
562      !! ** Purpose :   provide the surface boundary condition over sea-ice
563      !!
564      !! ** Method  :   compute momentum using bulk formulation
565      !!                formulea, ice variables and read atmospheric fields.
566      !!                NB: ice drag coefficient is assumed to be a constant
567      !!---------------------------------------------------------------------
568      INTEGER  ::   ji, jj    ! dummy loop indices
569      !
570      REAL(wp), DIMENSION(:,:)  , POINTER :: zrhoa
571      !
572      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f               ! relative wind module and components at F-point
573      REAL(wp) ::             zwndi_t , zwndj_t               ! relative wind components at T-point
574      REAL(wp), DIMENSION(:,:), POINTER ::   Cd               ! transfer coefficient for momentum      (tau)
575      !!---------------------------------------------------------------------
576      !
577      IF( nn_timing == 1 )  CALL timing_start('blk_ice_tau')
578      !
579      CALL wrk_alloc( jpi,jpj, zrhoa )
580      CALL wrk_alloc( jpi,jpj, Cd )
581
582      Cd(:,:) = Cd_ice
583
584      ! Make ice-atm. drag dependent on ice concentration (see Lupkes et al. 2012) (clem)
585#if defined key_lim3
586      IF( ln_Cd_L12 ) THEN
587         CALL Cdn10_Lupkes2012( Cd ) ! calculate new drag from Lupkes(2012) equations
588      ENDIF
589#endif
590
591      ! local scalars ( place there for vector optimisation purposes)
592      ! Computing density of air! Way denser that 1.2 over sea-ice !!!
593      !!
594      zrhoa (:,:) =  rho_air(sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1))
595
596      !!gm brutal....
597      utau_ice  (:,:) = 0._wp
598      vtau_ice  (:,:) = 0._wp
599      wndm_ice  (:,:) = 0._wp
600      !!gm end
601
602      ! ----------------------------------------------------------------------------- !
603      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
604      ! ----------------------------------------------------------------------------- !
605      SELECT CASE( cp_ice_msh )
606      CASE( 'I' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
607         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
608         DO jj = 2, jpjm1
609            DO ji = 2, jpim1   ! B grid : NO vector opt
610               ! ... scalar wind at I-point (fld being at T-point)
611               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ,1) + sf(jp_wndi)%fnow(ji  ,jj  ,1)   &
612                  &              + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji  ,jj-1,1)  ) - rn_vfac * u_ice(ji,jj)
613               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ,1) + sf(jp_wndj)%fnow(ji  ,jj  ,1)   &
614                  &              + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji  ,jj-1,1)  ) - rn_vfac * v_ice(ji,jj)
615               zwnorm_f = zrhoa(ji,jj) * Cd(ji,jj) * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
616               ! ... ice stress at I-point
617               utau_ice(ji,jj) = zwnorm_f * zwndi_f
618               vtau_ice(ji,jj) = zwnorm_f * zwndj_f
619               ! ... scalar wind at T-point (fld being at T-point)
620               zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  u_ice(ji,jj+1) + u_ice(ji+1,jj+1)   &
621                  &                                                    + u_ice(ji,jj  ) + u_ice(ji+1,jj  )  )
622               zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  v_ice(ji,jj+1) + v_ice(ji+1,jj+1)   &
623                  &                                                    + v_ice(ji,jj  ) + v_ice(ji+1,jj  )  )
624               wndm_ice(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
625            END DO
626         END DO
627         CALL lbc_lnk( utau_ice, 'I', -1. )
628         CALL lbc_lnk( vtau_ice, 'I', -1. )
629         CALL lbc_lnk( wndm_ice, 'T',  1. )
630         !
631      CASE( 'C' )                  ! C-grid ice dynamics :   U & V-points (same as ocean)
632         DO jj = 2, jpj
633            DO ji = fs_2, jpi   ! vect. opt.
634               zwndi_t = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( u_ice(ji-1,jj  ) + u_ice(ji,jj) )  )
635               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( v_ice(ji  ,jj-1) + v_ice(ji,jj) )  )
636               wndm_ice(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
637            END DO
638         END DO
639         DO jj = 2, jpjm1
640            DO ji = fs_2, fs_jpim1   ! vect. opt.
641               utau_ice(ji,jj) = 0.5 * zrhoa(ji,jj) * Cd(ji,jj) * ( wndm_ice(ji+1,jj  ) + wndm_ice(ji,jj) )                          &
642                  &          * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj,1) + sf(jp_wndi)%fnow(ji,jj,1) ) - rn_vfac * u_ice(ji,jj) )
643               vtau_ice(ji,jj) = 0.5 * zrhoa(ji,jj) * Cd(ji,jj) * ( wndm_ice(ji,jj+1  ) + wndm_ice(ji,jj) )                          &
644                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - rn_vfac * v_ice(ji,jj) )
645            END DO
646         END DO
647         CALL lbc_lnk( utau_ice, 'U', -1. )
648         CALL lbc_lnk( vtau_ice, 'V', -1. )
649         CALL lbc_lnk( wndm_ice, 'T',  1. )
650         !
651      END SELECT
652
653      IF(ln_ctl) THEN
654         CALL prt_ctl(tab2d_1=utau_ice  , clinfo1=' blk_ice: utau_ice : ', tab2d_2=vtau_ice  , clinfo2=' vtau_ice : ')
655         CALL prt_ctl(tab2d_1=wndm_ice  , clinfo1=' blk_ice: wndm_ice : ')
656      ENDIF
657
658      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_tau')
659
660   END SUBROUTINE blk_ice_tau
661
662
663   SUBROUTINE blk_ice_flx( ptsu, palb )
664      !!---------------------------------------------------------------------
665      !!                     ***  ROUTINE blk_ice_flx  ***
666      !!
667      !! ** Purpose :   provide the surface boundary condition over sea-ice
668      !!
669      !! ** Method  :   compute heat and freshwater exchanged
670      !!                between atmosphere and sea-ice using bulk formulation
671      !!                formulea, ice variables and read atmmospheric fields.
672      !!
673      !! caution : the net upward water flux has with mm/day unit
674      !!---------------------------------------------------------------------
675      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   ptsu   ! sea ice surface temperature
676      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   palb   ! ice albedo (all skies)
677      !!
678      INTEGER  ::   ji, jj, jl               ! dummy loop indices
679      REAL(wp) ::   zst2, zst3               ! local variable
680      REAL(wp) ::   zcoef_dqlw, zcoef_dqla   !   -      -
681      REAL(wp) ::   zztmp, z1_lsub           !   -      -
682      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qlw         ! long wave heat flux over ice
683      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qsb         ! sensible  heat flux over ice
684      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqlw        ! long wave heat sensitivity over ice
685      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqsb        ! sensible  heat sensitivity over ice
686      REAL(wp), DIMENSION(:,:)  , POINTER ::   zevap, zsnw   ! evaporation and snw distribution after wind blowing (LIM3)
687      REAL(wp), DIMENSION(:,:)  , POINTER ::   zrhoa
688      REAL(wp), DIMENSION(:,:)  , POINTER ::   Cd            ! transfer coefficient for momentum      (tau)
689      !!---------------------------------------------------------------------
690      !
691      IF( nn_timing == 1 )  CALL timing_start('blk_ice_flx')
692      !
693      CALL wrk_alloc( jpi,jpj,jpl,   z_qlw, z_qsb, z_dqlw, z_dqsb )
694      CALL wrk_alloc( jpi,jpj,       zrhoa)
695      CALL wrk_alloc( jpi,jpj, Cd )
696
697      Cd(:,:) = Cd_ice
698
699      ! Make ice-atm. drag dependent on ice concentration (see Lupkes et al.  2012) (clem)
700#if defined key_lim3
701      IF( ln_Cd_L12 ) THEN
702         CALL Cdn10_Lupkes2012( Cd ) ! calculate new drag from Lupkes(2012) equations
703      ENDIF
704#endif
705
706      !
707      ! local scalars ( place there for vector optimisation purposes)
708      zcoef_dqlw   = 4.0 * 0.95 * Stef
709      zcoef_dqla   = -Ls * 11637800. * (-5897.8)
710      !
711      zrhoa(:,:) = rho_air( sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1) )
712      !
713      zztmp = 1. / ( 1. - albo )
714      !                                     ! ========================== !
715      DO jl = 1, jpl                        !  Loop over ice categories  !
716         !                                  ! ========================== !
717         DO jj = 1 , jpj
718            DO ji = 1, jpi
719               ! ----------------------------!
720               !      I   Radiative FLUXES   !
721               ! ----------------------------!
722               zst2 = ptsu(ji,jj,jl) * ptsu(ji,jj,jl)
723               zst3 = ptsu(ji,jj,jl) * zst2
724               ! Short Wave (sw)
725               qsr_ice(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
726               ! Long  Wave (lw)
727               z_qlw(ji,jj,jl) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj,1) - Stef * ptsu(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
728               ! lw sensitivity
729               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3
730
731               ! ----------------------------!
732               !     II    Turbulent FLUXES  !
733               ! ----------------------------!
734
735               ! ... turbulent heat fluxes
736               ! Sensible Heat
737               z_qsb(ji,jj,jl) = zrhoa(ji,jj) * cpa * Cd(ji,jj) * wndm_ice(ji,jj) * ( ptsu(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) )
738               ! Latent Heat
739               qla_ice(ji,jj,jl) = rn_efac * MAX( 0.e0, zrhoa(ji,jj) * Ls  * Cd(ji,jj) * wndm_ice(ji,jj)   &
740                  &                         * (  11637800. * EXP( -5897.8 / ptsu(ji,jj,jl) ) / zrhoa(ji,jj) - sf(jp_humi)%fnow(ji,jj,1)  ) )
741               ! Latent heat sensitivity for ice (Dqla/Dt)
742               IF( qla_ice(ji,jj,jl) > 0._wp ) THEN
743                  dqla_ice(ji,jj,jl) = rn_efac * zcoef_dqla * wndm_ice(ji,jj) / ( zst2 ) * EXP( -5897.8 / ptsu(ji,jj,jl) )
744               ELSE
745                  dqla_ice(ji,jj,jl) = 0._wp
746               ENDIF
747
748               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
749               z_dqsb(ji,jj,jl) = zrhoa(ji,jj) * cpa * Cd(ji,jj) * wndm_ice(ji,jj)
750
751               ! ----------------------------!
752               !     III    Total FLUXES     !
753               ! ----------------------------!
754               ! Downward Non Solar flux
755               qns_ice (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - qla_ice (ji,jj,jl)
756               ! Total non solar heat flux sensitivity for ice
757               dqns_ice(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + dqla_ice(ji,jj,jl) )
758            END DO
759            !
760         END DO
761         !
762      END DO
763      !
764      tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
765      sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
766      CALL iom_put( 'snowpre', sprecip * 86400. )                  ! Snow precipitation
767      CALL iom_put( 'precip' , tprecip * 86400. )                  ! Total precipitation
768
769#if defined  key_lim3
770      CALL wrk_alloc( jpi,jpj,   zevap, zsnw )
771
772      ! --- evaporation --- !
773      z1_lsub = 1._wp / Lsub
774      evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * z1_lsub    ! sublimation
775      devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * z1_lsub    ! d(sublimation)/dT
776      zevap    (:,:)   = rn_efac * ( emp(:,:) + tprecip(:,:) )  ! evaporation over ocean
777
778      ! --- evaporation minus precipitation --- !
779      zsnw(:,:) = 0._wp
780      CALL lim_thd_snwblow( pfrld, zsnw )  ! snow distribution over ice after wind blowing
781      emp_oce(:,:) = pfrld(:,:) * zevap(:,:) - ( tprecip(:,:) - sprecip(:,:) ) - sprecip(:,:) * (1._wp - zsnw )
782      emp_ice(:,:) = SUM( a_i_b(:,:,:) * evap_ice(:,:,:), dim=3 ) - sprecip(:,:) * zsnw
783      emp_tot(:,:) = emp_oce(:,:) + emp_ice(:,:)
784
785      ! --- heat flux associated with emp --- !
786      qemp_oce(:,:) = - pfrld(:,:) * zevap(:,:) * sst_m(:,:) * rcp                               & ! evap at sst
787         &          + ( tprecip(:,:) - sprecip(:,:) ) * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp  & ! liquid precip at Tair
788         &          +   sprecip(:,:) * ( 1._wp - zsnw ) *                                        & ! solid precip at min(Tair,Tsnow)
789         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
790      qemp_ice(:,:) =   sprecip(:,:) * zsnw *                                                    & ! solid precip (only)
791         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
792
793      ! --- total solar and non solar fluxes --- !
794      qns_tot(:,:) = pfrld(:,:) * qns_oce(:,:) + SUM( a_i_b(:,:,:) * qns_ice(:,:,:), dim=3 ) + qemp_ice(:,:) + qemp_oce(:,:)
795      qsr_tot(:,:) = pfrld(:,:) * qsr_oce(:,:) + SUM( a_i_b(:,:,:) * qsr_ice(:,:,:), dim=3 )
796
797      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
798      qprec_ice(:,:) = rhosn * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
799
800      CALL wrk_dealloc( jpi,jpj,   zevap, zsnw )
801#endif
802
803      !--------------------------------------------------------------------
804      ! FRACTIONs of net shortwave radiation which is not absorbed in the
805      ! thin surface layer and penetrates inside the ice cover
806      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
807      !
808      fr1_i0(:,:) = ( 0.18 * ( 1.0 - cldf_ice ) + 0.35 * cldf_ice )
809      fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice )
810      !
811      !
812      IF(ln_ctl) THEN
813         CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' blk_ice: qla_ice  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb    : ', kdim=jpl)
814         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice: z_qlw    : ', tab3d_2=dqla_ice, clinfo2=' dqla_ice : ', kdim=jpl)
815         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice: z_dqsb   : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw   : ', kdim=jpl)
816         CALL prt_ctl(tab3d_1=dqns_ice, clinfo1=' blk_ice: dqns_ice : ', tab3d_2=qsr_ice , clinfo2=' qsr_ice  : ', kdim=jpl)
817         CALL prt_ctl(tab3d_1=ptsu    , clinfo1=' blk_ice: ptsu     : ', tab3d_2=qns_ice , clinfo2=' qns_ice  : ', kdim=jpl)
818         CALL prt_ctl(tab2d_1=tprecip , clinfo1=' blk_ice: tprecip  : ', tab2d_2=sprecip , clinfo2=' sprecip  : ')
819      ENDIF
820
821      CALL wrk_dealloc( jpi,jpj,jpl,   z_qlw, z_qsb, z_dqlw, z_dqsb )
822      CALL wrk_dealloc( jpi,jpj,       zrhoa )
823      CALL wrk_dealloc( jpi,jpj, Cd )
824      !
825      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_flx')
826
827   END SUBROUTINE blk_ice_flx
828   
829#endif
830
831   FUNCTION rho_air( ptak, pqa, pslp )
832      !!-------------------------------------------------------------------------------
833      !!                           ***  FUNCTION rho_air  ***
834      !!
835      !! ** Purpose : compute density of (moist) air using the eq. of state of the atmosphere
836      !!
837      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
838      !!-------------------------------------------------------------------------------
839      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   ptak      ! air temperature             [K]
840      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pqa       ! air specific humidity   [kg/kg]
841      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pslp      ! pressure in                [Pa]
842      REAL(wp), DIMENSION(jpi,jpj)             ::   rho_air   ! density of moist air   [kg/m^3]
843      !!-------------------------------------------------------------------------------
844      !
845      rho_air = pslp / (  R_dry*ptak * ( 1._wp + rctv0*pqa )  )
846      !
847   END FUNCTION rho_air
848
849
850   FUNCTION cp_air( pqa )
851      !!-------------------------------------------------------------------------------
852      !!                           ***  FUNCTION cp_air  ***
853      !!
854      !! ** Purpose : Compute specific heat (Cp) of moist air
855      !!
856      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
857      !!-------------------------------------------------------------------------------
858      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pqa      ! air specific humidity         [kg/kg]
859      REAL(wp), DIMENSION(jpi,jpj)             ::   cp_air   ! specific heat of moist air   [J/K/kg]
860      !!-------------------------------------------------------------------------------
861      !
862      Cp_air = Cp_dry + Cp_vap * pqa
863      !
864   END FUNCTION cp_air
865
866
867   FUNCTION q_sat( ptak, pslp )
868      !!----------------------------------------------------------------------------------
869      !!                           ***  FUNCTION q_sat  ***
870      !!
871      !! ** Purpose : Specific humidity at saturation in [kg/kg]
872      !!              Based on accurate estimate of "e_sat"
873      !!              aka saturation water vapor (Goff, 1957)
874      !!
875      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
876      !!----------------------------------------------------------------------------------
877      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   ptak    ! air temperature                       [K]
878      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pslp    ! sea level atmospheric pressure       [Pa]
879      REAL(wp), DIMENSION(jpi,jpj)             ::   q_sat   ! Specific humidity at saturation   [kg/kg]
880      !
881      INTEGER  ::   ji, jj         ! dummy loop indices
882      REAL(wp) ::   ze_sat, ztmp   ! local scalar
883      !!----------------------------------------------------------------------------------
884      !
885      DO jj = 1, jpj
886         DO ji = 1, jpi
887            !
888            ztmp = rt0 / ptak(ji,jj)
889            !
890            ! Vapour pressure at saturation [hPa] : WMO, (Goff, 1957)
891            ze_sat = 10.**( 10.79574*(1. - ztmp) - 5.028*LOG10(ptak(ji,jj)/rt0)        &
892               &    + 1.50475*10.**(-4)*(1. - 10.**(-8.2969*(ptak(ji,jj)/rt0 - 1.)) )  &
893               &    + 0.42873*10.**(-3)*(10.**(4.76955*(1. - ztmp)) - 1.) + 0.78614  )
894               !
895            q_sat(ji,jj) = reps0 * ze_sat/( 0.01_wp*pslp(ji,jj) - (1._wp - reps0)*ze_sat )   ! 0.01 because SLP is in [Pa]
896            !
897         END DO
898      END DO
899      !
900   END FUNCTION q_sat
901
902
903   FUNCTION gamma_moist( ptak, pqa )
904      !!----------------------------------------------------------------------------------
905      !!                           ***  FUNCTION gamma_moist  ***
906      !!
907      !! ** Purpose : Compute the moist adiabatic lapse-rate.
908      !!     => http://glossary.ametsoc.org/wiki/Moist-adiabatic_lapse_rate
909      !!     => http://www.geog.ucsb.edu/~joel/g266_s10/lecture_notes/chapt03/oh10_3_01/oh10_3_01.html
910      !!
911      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
912      !!----------------------------------------------------------------------------------
913      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   ptak          ! air temperature       [K]
914      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pqa           ! specific humidity [kg/kg]
915      REAL(wp), DIMENSION(jpi,jpj)             ::   gamma_moist   ! moist adiabatic lapse-rate
916      !
917      INTEGER  ::   ji, jj         ! dummy loop indices
918      REAL(wp) :: zrv, ziRT        ! local scalar
919      !!----------------------------------------------------------------------------------
920      !
921      DO jj = 1, jpj
922         DO ji = 1, jpi
923            zrv = pqa(ji,jj) / (1. - pqa(ji,jj))
924            ziRT = 1. / (R_dry*ptak(ji,jj))    ! 1/RT
925            gamma_moist(ji,jj) = grav * ( 1. + cevap*zrv*ziRT ) / ( Cp_dry + cevap*cevap*zrv*reps0*ziRT/ptak(ji,jj) )
926         END DO
927      END DO
928      !
929   END FUNCTION gamma_moist
930
931
932   FUNCTION L_vap( psst )
933      !!---------------------------------------------------------------------------------
934      !!                           ***  FUNCTION L_vap  ***
935      !!
936      !! ** Purpose : Compute the latent heat of vaporization of water from temperature
937      !!
938      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
939      !!----------------------------------------------------------------------------------
940      REAL(wp), DIMENSION(jpi,jpj)             ::   L_vap   ! latent heat of vaporization   [J/kg]
941      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   psst   ! water temperature                [K]
942      !!----------------------------------------------------------------------------------
943      !
944      L_vap = (  2.501 - 0.00237 * ( psst(:,:) - rt0)  ) * 1.e6
945      !
946   END FUNCTION L_vap
947
948
949#if defined key_lim3
950   SUBROUTINE Cdn10_Lupkes2012( Cd )
951      !!----------------------------------------------------------------------
952      !!                      ***  ROUTINE  Cdn10_Lupkes2012  ***
953      !!
954      !! ** Purpose :    Recompute the ice-atm drag at 10m height to make
955      !!                 it dependent on edges at leads, melt ponds and flows.
956      !!                 After some approximations, this can be resumed to a dependency
957      !!                 on ice concentration.
958      !!               
959      !! ** Method :     The parameterization is taken from Lupkes et al. (2012) eq.(50)
960      !!                 with the highest level of approximation: level4, eq.(59)
961      !!                 The generic drag over a cell partly covered by ice can be re-written as follows:
962      !!
963      !!                 Cd = Cdw * (1-A) + Cdi * A + Ce * (1-A)**(nu+1/(10*beta)) * A**mu
964      !!
965      !!                    Ce = 2.23e-3       , as suggested by Lupkes (eq. 59)
966      !!                    nu = mu = beta = 1 , as suggested by Lupkes (eq. 59)
967      !!                    A is the concentration of ice minus melt ponds (if any)
968      !!
969      !!                 This new drag has a parabolic shape (as a function of A) starting at
970      !!                 Cdw(say 1.5e-3) for A=0, reaching 1.97e-3 for A~0.5
971      !!                 and going down to Cdi(say 1.4e-3) for A=1
972      !!
973      !!                 It is theoretically applicable to all ice conditions !(not only MIZ)
974      !!                 => see Lupkes et al (2013)
975      !!
976      !! ** References : Lupkes et al. JGR 2012 (theory)
977      !!                 Lupkes et al. GRL 2013 (application to GCM)
978      !!
979      !!----------------------------------------------------------------------
980      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   Cd
981      REAL(wp), PARAMETER ::   zCe   = 2.23e-03_wp
982      REAL(wp), PARAMETER ::   znu   = 1._wp
983      REAL(wp), PARAMETER ::   zmu   = 1._wp
984      REAL(wp), PARAMETER ::   zbeta = 1._wp
985      REAL(wp)            ::   zcoef
986      !!----------------------------------------------------------------------
987      zcoef = znu + 1._wp / ( 10._wp * zbeta )
988
989      ! generic drag over a cell partly covered by ice
990      !!Cd(:,:) = Cd_oce(:,:) * ( 1._wp - at_i_b(:,:) ) +  & !! pure ocean drag
991      !!   &      Cd_ice      *           at_i_b(:,:)   +  & !! pure ice drag
992      !!   &      zCe         * ( 1._wp - at_i_b(:,:) )**zcoef * !at_i_b(:,:)**zmu   ! change due to sea-ice morphology
993
994      ! ice-atm drag
995      Cd(:,:) = Cd_ice +  & ! pure ice drag
996         &      zCe    * ( 1._wp - at_i_b(:,:) )**zcoef * at_i_b(:,:)**(zmu-1._wp)   ! change due to sea-ice morphology
997     
998   END SUBROUTINE Cdn10_Lupkes2012
999#endif
1000   
1001
1002   !!======================================================================
1003END MODULE sbcblk
Note: See TracBrowser for help on using the repository browser.