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/2017/dev_r8126_ROBUST08_no_ghost/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2017/dev_r8126_ROBUST08_no_ghost/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk.F90 @ 8809

Last change on this file since 8809 was 7753, checked in by mocavero, 7 years ago

Reverting trunk to remove OpenMP

File size: 58.0 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      IF ( ln_wave ) THEN
226      !Activated wave module but neither drag nor stokes drift activated
227         IF ( .NOT.(ln_cdgw .OR. ln_sdw .OR. ln_tauoc .OR. ln_stcor ) )   THEN
228            CALL ctl_warn( 'Ask for wave coupling but ln_cdgw=F, ln_sdw=F, ln_tauoc=F, ln_stcor=F')
229      !drag coefficient read from wave model definable only with mfs bulk formulae and core
230         ELSEIF (ln_cdgw .AND. .NOT. ln_NCAR )       THEN       
231             CALL ctl_stop( 'drag coefficient read from wave model definable only with mfs bulk formulae and core')
232         ELSEIF (ln_stcor .AND. .NOT. ln_sdw)                             THEN
233             CALL ctl_stop( 'Stokes-Coriolis term calculated only if activated Stokes Drift ln_sdw=T')
234         ENDIF
235      ELSE
236      IF ( ln_cdgw .OR. ln_sdw .OR. ln_tauoc .OR. ln_stcor )                & 
237         &   CALL ctl_stop( 'Not Activated Wave Module (ln_wave=F) but asked coupling ',    &
238         &                  'with drag coefficient (ln_cdgw =T) '  ,                        &
239         &                  'or Stokes Drift (ln_sdw=T) ' ,                                 &
240         &                  'or ocean stress modification due to waves (ln_tauoc=T) ',      & 
241         &                  'or Stokes-Coriolis term (ln_stcori=T)'  )
242      ENDIF 
243      !
244      !           
245      IF(lwp) THEN                     !** Control print
246         !
247         WRITE(numout,*)                  !* namelist
248         WRITE(numout,*) '   Namelist namsbc_blk (other than data information):'
249         WRITE(numout,*) '      "NCAR"      algorithm   (Large and Yeager 2008)     ln_NCAR      = ', ln_NCAR
250         WRITE(numout,*) '      "COARE 3.0" algorithm   (Fairall et al. 2003)       ln_COARE_3p0 = ', ln_COARE_3p0
251         WRITE(numout,*) '      "COARE 3.5" algorithm   (Edson et al. 2013)         ln_COARE_3p5 = ', ln_COARE_3p0
252         WRITE(numout,*) '      "ECMWF"     algorithm   (IFS cycle 31)              ln_ECMWF     = ', ln_ECMWF
253         WRITE(numout,*) '      add High freq.contribution to the stress module     ln_taudif    = ', ln_taudif
254         WRITE(numout,*) '      Air temperature and humidity reference height (m)   rn_zqt       = ', rn_zqt
255         WRITE(numout,*) '      Wind vector reference height (m)                    rn_zu        = ', rn_zu
256         WRITE(numout,*) '      factor applied on precipitation (total & snow)      rn_pfac      = ', rn_pfac
257         WRITE(numout,*) '      factor applied on evaporation                       rn_efac      = ', rn_efac
258         WRITE(numout,*) '      factor applied on ocean/ice velocity                rn_vfac      = ', rn_vfac
259         WRITE(numout,*) '         (form absolute (=0) to relative winds(=1))'
260         !
261         WRITE(numout,*)
262         SELECT CASE( nblk )              !* Print the choice of bulk algorithm
263         CASE( np_NCAR      )   ;   WRITE(numout,*) '      ===>>   "NCAR" algorithm        (Large and Yeager 2008)'
264         CASE( np_COARE_3p0 )   ;   WRITE(numout,*) '      ===>>   "COARE 3.0" algorithm   (Fairall et al. 2003)'
265         CASE( np_COARE_3p5 )   ;   WRITE(numout,*) '      ===>>   "COARE 3.5" algorithm   (Edson et al. 2013)'
266         CASE( np_ECMWF     )   ;   WRITE(numout,*) '      ===>>   "ECMWF" algorithm       (IFS cycle 31)'
267         END SELECT
268         !
269      ENDIF
270      !
271   END SUBROUTINE sbc_blk_init
272
273
274   SUBROUTINE sbc_blk( kt )
275      !!---------------------------------------------------------------------
276      !!                    ***  ROUTINE sbc_blk  ***
277      !!
278      !! ** Purpose :   provide at each time step the surface ocean fluxes
279      !!      (momentum, heat, freshwater and runoff)
280      !!
281      !! ** Method  : (1) READ each fluxes in NetCDF files:
282      !!      the 10m wind velocity (i-component) (m/s)    at T-point
283      !!      the 10m wind velocity (j-component) (m/s)    at T-point
284      !!      the 10m or 2m specific humidity     ( % )
285      !!      the solar heat                      (W/m2)
286      !!      the Long wave                       (W/m2)
287      !!      the 10m or 2m air temperature       (Kelvin)
288      !!      the total precipitation (rain+snow) (Kg/m2/s)
289      !!      the snow (solid prcipitation)       (kg/m2/s)
290      !!      the tau diff associated to HF tau   (N/m2)   at T-point   (ln_taudif=T)
291      !!              (2) CALL blk_oce
292      !!
293      !!      C A U T I O N : never mask the surface stress fields
294      !!                      the stress is assumed to be in the (i,j) mesh referential
295      !!
296      !! ** Action  :   defined at each time-step at the air-sea interface
297      !!              - utau, vtau  i- and j-component of the wind stress
298      !!              - taum        wind stress module at T-point
299      !!              - wndm        wind speed  module at T-point over free ocean or leads in presence of sea-ice
300      !!              - qns, qsr    non-solar and solar heat fluxes
301      !!              - emp         upward mass flux (evapo. - precip.)
302      !!              - sfx         salt flux due to freezing/melting (non-zero only if ice is present)
303      !!                            (set in limsbc(_2).F90)
304      !!
305      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
306      !!                   Brodeau et al. Ocean Modelling 2010
307      !!----------------------------------------------------------------------
308      INTEGER, INTENT(in) ::   kt   ! ocean time step
309      !!---------------------------------------------------------------------
310      !
311      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
312      !
313      !                                            ! compute the surface ocean fluxes using bulk formulea
314      IF( MOD( kt - 1, nn_fsbc ) == 0 )   CALL blk_oce( kt, sf, sst_m, ssu_m, ssv_m )
315
316#if defined key_cice
317      IF( MOD( kt - 1, nn_fsbc ) == 0 )   THEN
318         qlw_ice(:,:,1)   = sf(jp_qlw )%fnow(:,:,1)
319         IF( ln_dm2dc ) THEN ; qsr_ice(:,:,1) = sbc_dcy( sf(jp_qsr)%fnow(:,:,1) )
320         ELSE                ; qsr_ice(:,:,1) =          sf(jp_qsr)%fnow(:,:,1) 
321         ENDIF
322         tatm_ice(:,:)    = sf(jp_tair)%fnow(:,:,1)
323         qatm_ice(:,:)    = sf(jp_humi)%fnow(:,:,1)
324         tprecip(:,:)     = sf(jp_prec)%fnow(:,:,1) * rn_pfac
325         sprecip(:,:)     = sf(jp_snow)%fnow(:,:,1) * rn_pfac
326         wndi_ice(:,:)    = sf(jp_wndi)%fnow(:,:,1)
327         wndj_ice(:,:)    = sf(jp_wndj)%fnow(:,:,1)
328      ENDIF
329#endif
330      !
331   END SUBROUTINE sbc_blk
332
333
334   SUBROUTINE blk_oce( kt, sf, pst, pu, pv )
335      !!---------------------------------------------------------------------
336      !!                     ***  ROUTINE blk_oce  ***
337      !!
338      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
339      !!      the ocean surface at each time step
340      !!
341      !! ** Method  :   bulk formulea for the ocean using atmospheric
342      !!      fields read in sbc_read
343      !!
344      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
345      !!              - vtau    : j-component of the stress at V-point  (N/m2)
346      !!              - taum    : Wind stress module at T-point         (N/m2)
347      !!              - wndm    : Wind speed module at T-point          (m/s)
348      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
349      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
350      !!              - emp     : evaporation minus precipitation       (kg/m2/s)
351      !!
352      !!  ** Nota  :   sf has to be a dummy argument for AGRIF on NEC
353      !!---------------------------------------------------------------------
354      INTEGER  , INTENT(in   )                 ::   kt    ! time step index
355      TYPE(fld), INTENT(inout), DIMENSION(:)   ::   sf    ! input data
356      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pst   ! surface temperature                      [Celcius]
357      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pu    ! surface current at U-point (i-component) [m/s]
358      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pv    ! surface current at V-point (j-component) [m/s]
359      !
360      INTEGER  ::   ji, jj               ! dummy loop indices
361      REAL(wp) ::   zztmp                ! local variable
362      REAL(wp), DIMENSION(:,:), POINTER ::   zwnd_i, zwnd_j    ! wind speed components at T-point
363      REAL(wp), DIMENSION(:,:), POINTER ::   zsq               ! specific humidity at pst
364      REAL(wp), DIMENSION(:,:), POINTER ::   zqlw, zqsb        ! long wave and sensible heat fluxes
365      REAL(wp), DIMENSION(:,:), POINTER ::   zqla, zevap       ! latent heat fluxes and evaporation
366      REAL(wp), DIMENSION(:,:), POINTER ::   Cd                ! transfer coefficient for momentum      (tau)
367      REAL(wp), DIMENSION(:,:), POINTER ::   Ch                ! transfer coefficient for sensible heat (Q_sens)
368      REAL(wp), DIMENSION(:,:), POINTER ::   Ce                ! tansfert coefficient for evaporation   (Q_lat)
369      REAL(wp), DIMENSION(:,:), POINTER ::   zst               ! surface temperature in Kelvin
370      REAL(wp), DIMENSION(:,:), POINTER ::   zt_zu             ! air temperature at wind speed height
371      REAL(wp), DIMENSION(:,:), POINTER ::   zq_zu             ! air spec. hum.  at wind speed height
372      REAL(wp), DIMENSION(:,:), POINTER ::   zU_zu             ! bulk wind speed at height zu  [m/s]
373      REAL(wp), DIMENSION(:,:), POINTER ::   ztpot             ! potential temperature of air at z=rn_zqt [K]
374      REAL(wp), DIMENSION(:,:), POINTER ::   zrhoa             ! density of air   [kg/m^3]
375      !!---------------------------------------------------------------------
376      !
377      IF( nn_timing == 1 )  CALL timing_start('blk_oce')
378      !
379      CALL wrk_alloc( jpi,jpj,   zwnd_i, zwnd_j, zsq, zqlw, zqsb, zqla, zevap )
380      CALL wrk_alloc( jpi,jpj,   Cd, Ch, Ce, zst, zt_zu, zq_zu )
381      CALL wrk_alloc( jpi,jpj,   zU_zu, ztpot, zrhoa )
382      !
383
384      ! local scalars ( place there for vector optimisation purposes)
385      zst(:,:) = pst(:,:) + rt0      ! convert SST from Celcius to Kelvin (and set minimum value far above 0 K)
386
387      ! ----------------------------------------------------------------------------- !
388      !      0   Wind components and module at T-point relative to the moving ocean   !
389      ! ----------------------------------------------------------------------------- !
390
391      ! ... components ( U10m - U_oce ) at T-point (unmasked)
392!!gm    move zwnd_i (_j) set to zero  inside the key_cyclone ???
393      zwnd_i(:,:) = 0._wp
394      zwnd_j(:,:) = 0._wp
395#if defined key_cyclone
396      CALL wnd_cyc( kt, zwnd_i, zwnd_j )    ! add analytical tropical cyclone (Vincent et al. JGR 2012)
397      DO jj = 2, jpjm1
398         DO ji = fs_2, fs_jpim1   ! vect. opt.
399            sf(jp_wndi)%fnow(ji,jj,1) = sf(jp_wndi)%fnow(ji,jj,1) + zwnd_i(ji,jj)
400            sf(jp_wndj)%fnow(ji,jj,1) = sf(jp_wndj)%fnow(ji,jj,1) + zwnd_j(ji,jj)
401         END DO
402      END DO
403#endif
404      DO jj = 2, jpjm1
405         DO ji = fs_2, fs_jpim1   ! vect. opt.
406            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
407            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
408         END DO
409      END DO
410      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
411      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
412      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
413      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
414         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
415
416      ! ----------------------------------------------------------------------------- !
417      !      I   Radiative FLUXES                                                     !
418      ! ----------------------------------------------------------------------------- !
419
420      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
421      zztmp = 1. - albo
422      IF( ln_dm2dc ) THEN   ;   qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
423      ELSE                  ;   qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
424      ENDIF
425
426      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
427
428
429
430      ! ----------------------------------------------------------------------------- !
431      !     II    Turbulent FLUXES                                                    !
432      ! ----------------------------------------------------------------------------- !
433
434      ! ... specific humidity at SST and IST tmask(
435      zsq(:,:) = 0.98 * q_sat( zst(:,:), sf(jp_slp)%fnow(:,:,1) )
436      !!
437      !! Estimate of potential temperature at z=rn_zqt, based on adiabatic lapse-rate
438      !!    (see Josey, Gulev & Yu, 2013) / doi=10.1016/B978-0-12-391851-2.00005-2
439      !!    (since reanalysis products provide T at z, not theta !)
440      ztpot = sf(jp_tair)%fnow(:,:,1) + gamma_moist( sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1) ) * rn_zqt
441
442      SELECT CASE( nblk )        !==  transfer coefficients  ==!   Cd, Ch, Ce at T-point
443      !
444      CASE( np_NCAR      )   ;   CALL turb_ncar    ( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! NCAR-COREv2
445         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
446      CASE( np_COARE_3p0 )   ;   CALL turb_coare   ( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! COARE v3.0
447         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
448      CASE( np_COARE_3p5 )   ;   CALL turb_coare3p5( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! COARE v3.5
449         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
450      CASE( np_ECMWF     )   ;   CALL turb_ecmwf   ( rn_zqt, rn_zu, zst, ztpot, zsq, sf(jp_humi)%fnow, wndm,   &  ! ECMWF
451         &                                               Cd, Ch, Ce, zt_zu, zq_zu, zU_zu )
452      CASE DEFAULT
453         CALL ctl_stop( 'STOP', 'sbc_oce: non-existing bulk formula selected' )
454      END SELECT
455
456      !                          ! Compute true air density :
457      IF( ABS(rn_zu - rn_zqt) > 0.01 ) THEN     ! At zu: (probably useless to remove zrho*grav*rn_zu from SLP...)
458         zrhoa(:,:) = rho_air( zt_zu(:,:)             , zq_zu(:,:)             , sf(jp_slp)%fnow(:,:,1) )
459      ELSE                                      ! At zt:
460         zrhoa(:,:) = rho_air( sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1) )
461      END IF
462
463      Cd_oce(:,:) = Cd(:,:)  ! record value of pure ocean-atm. drag (clem)
464
465      DO jj = 1, jpj             ! tau module, i and j component
466         DO ji = 1, jpi
467            zztmp = zrhoa(ji,jj)  * zU_zu(ji,jj) * Cd(ji,jj)   ! using bulk wind speed
468            taum  (ji,jj) = zztmp * wndm  (ji,jj)
469            zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj)
470            zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj)
471         END DO
472      END DO
473
474      !                          ! add the HF tau contribution to the wind stress module
475      IF( lhftau )   taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1)
476
477      CALL iom_put( "taum_oce", taum )   ! output wind stress module
478
479      ! ... utau, vtau at U- and V_points, resp.
480      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
481      !     Note the use of MAX(tmask(i,j),tmask(i+1,j) is to mask tau over ice shelves
482      DO jj = 1, jpjm1
483         DO ji = 1, fs_jpim1
484            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) ) &
485               &          * MAX(tmask(ji,jj,1),tmask(ji+1,jj,1))
486            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) ) &
487               &          * MAX(tmask(ji,jj,1),tmask(ji,jj+1,1))
488         END DO
489      END DO
490      CALL lbc_lnk( utau(:,:), 'U', -1. )
491      CALL lbc_lnk( vtau(:,:), 'V', -1. )
492
493
494      !  Turbulent fluxes over ocean
495      ! -----------------------------
496
497      ! zqla used as temporary array, for rho*U (common term of bulk formulae):
498      zqla(:,:) = zrhoa(:,:) * zU_zu(:,:)
499
500      IF( ABS( rn_zu - rn_zqt) < 0.01_wp ) THEN
501         !! q_air and t_air are given at 10m (wind reference height)
502         zevap(:,:) = rn_efac*MAX( 0._wp,             zqla(:,:)*Ce(:,:)*(zsq(:,:) - sf(jp_humi)%fnow(:,:,1)) ) ! Evaporation, using bulk wind speed
503         zqsb (:,:) = cp_air(sf(jp_humi)%fnow(:,:,1))*zqla(:,:)*Ch(:,:)*(zst(:,:) - ztpot(:,:)             )   ! Sensible Heat, using bulk wind speed
504      ELSE
505         !! q_air and t_air are not given at 10m (wind reference height)
506         ! Values of temp. and hum. adjusted to height of wind during bulk algorithm iteration must be used!!!
507         zevap(:,:) = rn_efac*MAX( 0._wp,             zqla(:,:)*Ce(:,:)*(zsq(:,:) - zq_zu(:,:) ) ) ! Evaporation ! using bulk wind speed
508         zqsb (:,:) = cp_air(sf(jp_humi)%fnow(:,:,1))*zqla(:,:)*Ch(:,:)*(zst(:,:) - zt_zu(:,:) )   ! Sensible Heat ! using bulk wind speed
509      ENDIF
510
511      zqla(:,:) = L_vap(zst(:,:)) * zevap(:,:)     ! Latent Heat flux
512
513
514      IF(ln_ctl) THEN
515         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
516         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
517         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
518         CALL prt_ctl( tab2d_1=zsq   , clinfo1=' blk_oce: zsq    : ', tab2d_2=zst, clinfo2=' zst : ' )
519         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce: utau   : ', mask1=umask,   &
520            &          tab2d_2=vtau  , clinfo2=           ' vtau : ', mask2=vmask )
521         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce: wndm   : ')
522         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce: zst    : ')
523      ENDIF
524
525      ! ----------------------------------------------------------------------------- !
526      !     III    Total FLUXES                                                       !
527      ! ----------------------------------------------------------------------------- !
528      !
529      emp (:,:) = (  zevap(:,:)                                          &   ! mass flux (evap. - precip.)
530         &         - sf(jp_prec)%fnow(:,:,1) * rn_pfac  ) * tmask(:,:,1)
531      !
532      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                &   ! Downward Non Solar
533         &     - sf(jp_snow)%fnow(:,:,1) * rn_pfac * lfus                         &   ! remove latent melting heat for solid precip
534         &     - zevap(:,:) * pst(:,:) * rcp                                      &   ! remove evap heat content at SST
535         &     + ( sf(jp_prec)%fnow(:,:,1) - sf(jp_snow)%fnow(:,:,1) ) * rn_pfac  &   ! add liquid precip heat content at Tair
536         &     * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp                          &
537         &     + sf(jp_snow)%fnow(:,:,1) * rn_pfac                                &   ! add solid  precip heat content at min(Tair,Tsnow)
538         &     * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1)
539      !
540#if defined key_lim3
541      qns_oce(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                ! non solar without emp (only needed by LIM3)
542      qsr_oce(:,:) = qsr(:,:)
543#endif
544      !
545      IF ( nn_ice == 0 ) THEN
546         CALL iom_put( "qlw_oce" ,   zqlw )                 ! output downward longwave heat over the ocean
547         CALL iom_put( "qsb_oce" , - zqsb )                 ! output downward sensible heat over the ocean
548         CALL iom_put( "qla_oce" , - zqla )                 ! output downward latent   heat over the ocean
549         CALL iom_put( "qemp_oce",   qns-zqlw+zqsb+zqla )   ! output downward heat content of E-P over the ocean
550         CALL iom_put( "qns_oce" ,   qns  )                 ! output downward non solar heat over the ocean
551         CALL iom_put( "qsr_oce" ,   qsr  )                 ! output downward solar heat over the ocean
552         CALL iom_put( "qt_oce"  ,   qns+qsr )              ! output total downward heat over the ocean
553         tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac   ! output total precipitation [kg/m2/s]
554         sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac   ! output solid precipitation [kg/m2/s]
555         CALL iom_put( 'snowpre', sprecip * 86400. )        ! Snow
556         CALL iom_put( 'precip' , tprecip * 86400. )        ! Total precipitation
557      ENDIF
558      !
559      IF(ln_ctl) THEN
560         CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce: zqsb   : ', tab2d_2=zqlw , clinfo2=' zqlw  : ')
561         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce: zqla   : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
562         CALL prt_ctl(tab2d_1=pst  , clinfo1=' blk_oce: pst    : ', tab2d_2=emp  , clinfo2=' emp   : ')
563         CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce: utau   : ', mask1=umask,   &
564            &         tab2d_2=vtau , clinfo2=              ' vtau  : ' , mask2=vmask )
565      ENDIF
566      !
567      CALL wrk_dealloc( jpi,jpj,   zwnd_i, zwnd_j, zsq, zqlw, zqsb, zqla, zevap )
568      CALL wrk_dealloc( jpi,jpj,   Cd, Ch, Ce, zst, zt_zu, zq_zu )
569      CALL wrk_dealloc( jpi,jpj,   zU_zu, ztpot, zrhoa )
570      !
571      IF( nn_timing == 1 )  CALL timing_stop('blk_oce')
572      !
573   END SUBROUTINE blk_oce
574
575#if defined key_lim2 || defined key_lim3
576
577   SUBROUTINE blk_ice_tau
578      !!---------------------------------------------------------------------
579      !!                     ***  ROUTINE blk_ice_tau  ***
580      !!
581      !! ** Purpose :   provide the surface boundary condition over sea-ice
582      !!
583      !! ** Method  :   compute momentum using bulk formulation
584      !!                formulea, ice variables and read atmospheric fields.
585      !!                NB: ice drag coefficient is assumed to be a constant
586      !!---------------------------------------------------------------------
587      INTEGER  ::   ji, jj    ! dummy loop indices
588      !
589      REAL(wp), DIMENSION(:,:)  , POINTER :: zrhoa
590      !
591      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f               ! relative wind module and components at F-point
592      REAL(wp) ::             zwndi_t , zwndj_t               ! relative wind components at T-point
593      REAL(wp), DIMENSION(:,:), POINTER ::   Cd               ! transfer coefficient for momentum      (tau)
594      !!---------------------------------------------------------------------
595      !
596      IF( nn_timing == 1 )  CALL timing_start('blk_ice_tau')
597      !
598      CALL wrk_alloc( jpi,jpj, zrhoa )
599      CALL wrk_alloc( jpi,jpj, Cd )
600
601      Cd(:,:) = Cd_ice
602
603      ! Make ice-atm. drag dependent on ice concentration (see Lupkes et al. 2012) (clem)
604#if defined key_lim3
605      IF( ln_Cd_L12 ) THEN
606         CALL Cdn10_Lupkes2012( Cd ) ! calculate new drag from Lupkes(2012) equations
607      ENDIF
608#endif
609
610      ! local scalars ( place there for vector optimisation purposes)
611      ! Computing density of air! Way denser that 1.2 over sea-ice !!!
612      !!
613      zrhoa (:,:) =  rho_air(sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1))
614
615      !!gm brutal....
616      utau_ice  (:,:) = 0._wp
617      vtau_ice  (:,:) = 0._wp
618      wndm_ice  (:,:) = 0._wp
619      !!gm end
620
621      ! ----------------------------------------------------------------------------- !
622      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
623      ! ----------------------------------------------------------------------------- !
624      SELECT CASE( cp_ice_msh )
625      CASE( 'I' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
626         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
627         DO jj = 2, jpjm1
628            DO ji = 2, jpim1   ! B grid : NO vector opt
629               ! ... scalar wind at I-point (fld being at T-point)
630               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ,1) + sf(jp_wndi)%fnow(ji  ,jj  ,1)   &
631                  &              + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji  ,jj-1,1)  ) - rn_vfac * u_ice(ji,jj)
632               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ,1) + sf(jp_wndj)%fnow(ji  ,jj  ,1)   &
633                  &              + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji  ,jj-1,1)  ) - rn_vfac * v_ice(ji,jj)
634               zwnorm_f = zrhoa(ji,jj) * Cd(ji,jj) * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
635               ! ... ice stress at I-point
636               utau_ice(ji,jj) = zwnorm_f * zwndi_f
637               vtau_ice(ji,jj) = zwnorm_f * zwndj_f
638               ! ... scalar wind at T-point (fld being at T-point)
639               zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  u_ice(ji,jj+1) + u_ice(ji+1,jj+1)   &
640                  &                                                    + u_ice(ji,jj  ) + u_ice(ji+1,jj  )  )
641               zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  v_ice(ji,jj+1) + v_ice(ji+1,jj+1)   &
642                  &                                                    + v_ice(ji,jj  ) + v_ice(ji+1,jj  )  )
643               wndm_ice(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
644            END DO
645         END DO
646         CALL lbc_lnk( utau_ice, 'I', -1. )
647         CALL lbc_lnk( vtau_ice, 'I', -1. )
648         CALL lbc_lnk( wndm_ice, 'T',  1. )
649         !
650      CASE( 'C' )                  ! C-grid ice dynamics :   U & V-points (same as ocean)
651         DO jj = 2, jpj
652            DO ji = fs_2, jpi   ! vect. opt.
653               zwndi_t = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( u_ice(ji-1,jj  ) + u_ice(ji,jj) )  )
654               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( v_ice(ji  ,jj-1) + v_ice(ji,jj) )  )
655               wndm_ice(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
656            END DO
657         END DO
658         DO jj = 2, jpjm1
659            DO ji = fs_2, fs_jpim1   ! vect. opt.
660               utau_ice(ji,jj) = 0.5 * zrhoa(ji,jj) * Cd(ji,jj) * ( wndm_ice(ji+1,jj  ) + wndm_ice(ji,jj) )                          &
661                  &          * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj,1) + sf(jp_wndi)%fnow(ji,jj,1) ) - rn_vfac * u_ice(ji,jj) )
662               vtau_ice(ji,jj) = 0.5 * zrhoa(ji,jj) * Cd(ji,jj) * ( wndm_ice(ji,jj+1  ) + wndm_ice(ji,jj) )                          &
663                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - rn_vfac * v_ice(ji,jj) )
664            END DO
665         END DO
666         CALL lbc_lnk( utau_ice, 'U', -1. )
667         CALL lbc_lnk( vtau_ice, 'V', -1. )
668         CALL lbc_lnk( wndm_ice, 'T',  1. )
669         !
670      END SELECT
671
672      IF(ln_ctl) THEN
673         CALL prt_ctl(tab2d_1=utau_ice  , clinfo1=' blk_ice: utau_ice : ', tab2d_2=vtau_ice  , clinfo2=' vtau_ice : ')
674         CALL prt_ctl(tab2d_1=wndm_ice  , clinfo1=' blk_ice: wndm_ice : ')
675      ENDIF
676
677      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_tau')
678
679   END SUBROUTINE blk_ice_tau
680
681
682   SUBROUTINE blk_ice_flx( ptsu, palb )
683      !!---------------------------------------------------------------------
684      !!                     ***  ROUTINE blk_ice_flx  ***
685      !!
686      !! ** Purpose :   provide the surface boundary condition over sea-ice
687      !!
688      !! ** Method  :   compute heat and freshwater exchanged
689      !!                between atmosphere and sea-ice using bulk formulation
690      !!                formulea, ice variables and read atmmospheric fields.
691      !!
692      !! caution : the net upward water flux has with mm/day unit
693      !!---------------------------------------------------------------------
694      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   ptsu   ! sea ice surface temperature
695      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   palb   ! ice albedo (all skies)
696      !!
697      INTEGER  ::   ji, jj, jl               ! dummy loop indices
698      REAL(wp) ::   zst2, zst3               ! local variable
699      REAL(wp) ::   zcoef_dqlw, zcoef_dqla   !   -      -
700      REAL(wp) ::   zztmp, z1_lsub           !   -      -
701      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qlw         ! long wave heat flux over ice
702      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qsb         ! sensible  heat flux over ice
703      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqlw        ! long wave heat sensitivity over ice
704      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqsb        ! sensible  heat sensitivity over ice
705      REAL(wp), DIMENSION(:,:)  , POINTER ::   zevap, zsnw   ! evaporation and snw distribution after wind blowing (LIM3)
706      REAL(wp), DIMENSION(:,:)  , POINTER ::   zrhoa
707      REAL(wp), DIMENSION(:,:)  , POINTER ::   Cd            ! transfer coefficient for momentum      (tau)
708      !!---------------------------------------------------------------------
709      !
710      IF( nn_timing == 1 )  CALL timing_start('blk_ice_flx')
711      !
712      CALL wrk_alloc( jpi,jpj,jpl,   z_qlw, z_qsb, z_dqlw, z_dqsb )
713      CALL wrk_alloc( jpi,jpj,       zrhoa)
714      CALL wrk_alloc( jpi,jpj, Cd )
715
716      Cd(:,:) = Cd_ice
717
718      ! Make ice-atm. drag dependent on ice concentration (see Lupkes et al.  2012) (clem)
719#if defined key_lim3
720      IF( ln_Cd_L12 ) THEN
721         CALL Cdn10_Lupkes2012( Cd ) ! calculate new drag from Lupkes(2012) equations
722      ENDIF
723#endif
724
725      !
726      ! local scalars ( place there for vector optimisation purposes)
727      zcoef_dqlw   = 4.0 * 0.95 * Stef
728      zcoef_dqla   = -Ls * 11637800. * (-5897.8)
729      !
730      zrhoa(:,:) = rho_air( sf(jp_tair)%fnow(:,:,1), sf(jp_humi)%fnow(:,:,1), sf(jp_slp)%fnow(:,:,1) )
731      !
732      zztmp = 1. / ( 1. - albo )
733      !                                     ! ========================== !
734      DO jl = 1, jpl                        !  Loop over ice categories  !
735         !                                  ! ========================== !
736         DO jj = 1 , jpj
737            DO ji = 1, jpi
738               ! ----------------------------!
739               !      I   Radiative FLUXES   !
740               ! ----------------------------!
741               zst2 = ptsu(ji,jj,jl) * ptsu(ji,jj,jl)
742               zst3 = ptsu(ji,jj,jl) * zst2
743               ! Short Wave (sw)
744               qsr_ice(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
745               ! Long  Wave (lw)
746               z_qlw(ji,jj,jl) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj,1) - Stef * ptsu(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
747               ! lw sensitivity
748               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3
749
750               ! ----------------------------!
751               !     II    Turbulent FLUXES  !
752               ! ----------------------------!
753
754               ! ... turbulent heat fluxes
755               ! Sensible Heat
756               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) )
757               ! Latent Heat
758               qla_ice(ji,jj,jl) = rn_efac * MAX( 0.e0, zrhoa(ji,jj) * Ls  * Cd(ji,jj) * wndm_ice(ji,jj)   &
759                  &                         * (  11637800. * EXP( -5897.8 / ptsu(ji,jj,jl) ) / zrhoa(ji,jj) - sf(jp_humi)%fnow(ji,jj,1)  ) )
760               ! Latent heat sensitivity for ice (Dqla/Dt)
761               IF( qla_ice(ji,jj,jl) > 0._wp ) THEN
762                  dqla_ice(ji,jj,jl) = rn_efac * zcoef_dqla * Cd(ji,jj) * wndm_ice(ji,jj) / ( zst2 ) * EXP( -5897.8 / ptsu(ji,jj,jl) )
763               ELSE
764                  dqla_ice(ji,jj,jl) = 0._wp
765               ENDIF
766
767               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
768               z_dqsb(ji,jj,jl) = zrhoa(ji,jj) * cpa * Cd(ji,jj) * wndm_ice(ji,jj)
769
770               ! ----------------------------!
771               !     III    Total FLUXES     !
772               ! ----------------------------!
773               ! Downward Non Solar flux
774               qns_ice (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - qla_ice (ji,jj,jl)
775               ! Total non solar heat flux sensitivity for ice
776               dqns_ice(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + dqla_ice(ji,jj,jl) )
777            END DO
778            !
779         END DO
780         !
781      END DO
782      !
783      tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
784      sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
785      CALL iom_put( 'snowpre', sprecip * 86400. )                  ! Snow precipitation
786      CALL iom_put( 'precip' , tprecip * 86400. )                  ! Total precipitation
787
788#if defined  key_lim3
789      CALL wrk_alloc( jpi,jpj,   zevap, zsnw )
790
791      ! --- evaporation --- !
792      z1_lsub = 1._wp / Lsub
793      evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * z1_lsub    ! sublimation
794      devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * z1_lsub    ! d(sublimation)/dT
795      zevap    (:,:)   = rn_efac * ( emp(:,:) + tprecip(:,:) )  ! evaporation over ocean
796
797      ! --- evaporation minus precipitation --- !
798      zsnw(:,:) = 0._wp
799      CALL lim_thd_snwblow( pfrld, zsnw )  ! snow distribution over ice after wind blowing
800      emp_oce(:,:) = pfrld(:,:) * zevap(:,:) - ( tprecip(:,:) - sprecip(:,:) ) - sprecip(:,:) * (1._wp - zsnw )
801      emp_ice(:,:) = SUM( a_i_b(:,:,:) * evap_ice(:,:,:), dim=3 ) - sprecip(:,:) * zsnw
802      emp_tot(:,:) = emp_oce(:,:) + emp_ice(:,:)
803
804      ! --- heat flux associated with emp --- !
805      qemp_oce(:,:) = - pfrld(:,:) * zevap(:,:) * sst_m(:,:) * rcp                               & ! evap at sst
806         &          + ( tprecip(:,:) - sprecip(:,:) ) * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp  & ! liquid precip at Tair
807         &          +   sprecip(:,:) * ( 1._wp - zsnw ) *                                        & ! solid precip at min(Tair,Tsnow)
808         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
809      qemp_ice(:,:) =   sprecip(:,:) * zsnw *                                                    & ! solid precip (only)
810         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
811
812      ! --- total solar and non solar fluxes --- !
813      qns_tot(:,:) = pfrld(:,:) * qns_oce(:,:) + SUM( a_i_b(:,:,:) * qns_ice(:,:,:), dim=3 ) + qemp_ice(:,:) + qemp_oce(:,:)
814      qsr_tot(:,:) = pfrld(:,:) * qsr_oce(:,:) + SUM( a_i_b(:,:,:) * qsr_ice(:,:,:), dim=3 )
815
816      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
817      qprec_ice(:,:) = rhosn * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
818
819      ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) ---
820      DO jl = 1, jpl
821         qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * cpic * tmask(:,:,1) )
822                                   ! But we do not have Tice => consider it at 0degC => evap=0
823      END DO
824
825      CALL wrk_dealloc( jpi,jpj,   zevap, zsnw )
826#endif
827
828      !--------------------------------------------------------------------
829      ! FRACTIONs of net shortwave radiation which is not absorbed in the
830      ! thin surface layer and penetrates inside the ice cover
831      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
832      !
833      fr1_i0(:,:) = ( 0.18 * ( 1.0 - cldf_ice ) + 0.35 * cldf_ice )
834      fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice )
835      !
836      !
837      IF(ln_ctl) THEN
838         CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' blk_ice: qla_ice  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb    : ', kdim=jpl)
839         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice: z_qlw    : ', tab3d_2=dqla_ice, clinfo2=' dqla_ice : ', kdim=jpl)
840         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice: z_dqsb   : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw   : ', kdim=jpl)
841         CALL prt_ctl(tab3d_1=dqns_ice, clinfo1=' blk_ice: dqns_ice : ', tab3d_2=qsr_ice , clinfo2=' qsr_ice  : ', kdim=jpl)
842         CALL prt_ctl(tab3d_1=ptsu    , clinfo1=' blk_ice: ptsu     : ', tab3d_2=qns_ice , clinfo2=' qns_ice  : ', kdim=jpl)
843         CALL prt_ctl(tab2d_1=tprecip , clinfo1=' blk_ice: tprecip  : ', tab2d_2=sprecip , clinfo2=' sprecip  : ')
844      ENDIF
845
846      CALL wrk_dealloc( jpi,jpj,jpl,   z_qlw, z_qsb, z_dqlw, z_dqsb )
847      CALL wrk_dealloc( jpi,jpj,       zrhoa )
848      CALL wrk_dealloc( jpi,jpj, Cd )
849      !
850      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_flx')
851
852   END SUBROUTINE blk_ice_flx
853   
854#endif
855
856   FUNCTION rho_air( ptak, pqa, pslp )
857      !!-------------------------------------------------------------------------------
858      !!                           ***  FUNCTION rho_air  ***
859      !!
860      !! ** Purpose : compute density of (moist) air using the eq. of state of the atmosphere
861      !!
862      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
863      !!-------------------------------------------------------------------------------
864      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   ptak      ! air temperature             [K]
865      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pqa       ! air specific humidity   [kg/kg]
866      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pslp      ! pressure in                [Pa]
867      REAL(wp), DIMENSION(jpi,jpj)             ::   rho_air   ! density of moist air   [kg/m^3]
868      !!-------------------------------------------------------------------------------
869      !
870      rho_air = pslp / (  R_dry*ptak * ( 1._wp + rctv0*pqa )  )
871      !
872   END FUNCTION rho_air
873
874
875   FUNCTION cp_air( pqa )
876      !!-------------------------------------------------------------------------------
877      !!                           ***  FUNCTION cp_air  ***
878      !!
879      !! ** Purpose : Compute specific heat (Cp) of moist air
880      !!
881      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
882      !!-------------------------------------------------------------------------------
883      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pqa      ! air specific humidity         [kg/kg]
884      REAL(wp), DIMENSION(jpi,jpj)             ::   cp_air   ! specific heat of moist air   [J/K/kg]
885      !!-------------------------------------------------------------------------------
886      !
887      Cp_air = Cp_dry + Cp_vap * pqa
888      !
889   END FUNCTION cp_air
890
891
892   FUNCTION q_sat( ptak, pslp )
893      !!----------------------------------------------------------------------------------
894      !!                           ***  FUNCTION q_sat  ***
895      !!
896      !! ** Purpose : Specific humidity at saturation in [kg/kg]
897      !!              Based on accurate estimate of "e_sat"
898      !!              aka saturation water vapor (Goff, 1957)
899      !!
900      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
901      !!----------------------------------------------------------------------------------
902      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   ptak    ! air temperature                       [K]
903      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pslp    ! sea level atmospheric pressure       [Pa]
904      REAL(wp), DIMENSION(jpi,jpj)             ::   q_sat   ! Specific humidity at saturation   [kg/kg]
905      !
906      INTEGER  ::   ji, jj         ! dummy loop indices
907      REAL(wp) ::   ze_sat, ztmp   ! local scalar
908      !!----------------------------------------------------------------------------------
909      !
910      DO jj = 1, jpj
911         DO ji = 1, jpi
912            !
913            ztmp = rt0 / ptak(ji,jj)
914            !
915            ! Vapour pressure at saturation [hPa] : WMO, (Goff, 1957)
916            ze_sat = 10.**( 10.79574*(1. - ztmp) - 5.028*LOG10(ptak(ji,jj)/rt0)        &
917               &    + 1.50475*10.**(-4)*(1. - 10.**(-8.2969*(ptak(ji,jj)/rt0 - 1.)) )  &
918               &    + 0.42873*10.**(-3)*(10.**(4.76955*(1. - ztmp)) - 1.) + 0.78614  )
919               !
920            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]
921            !
922         END DO
923      END DO
924      !
925   END FUNCTION q_sat
926
927
928   FUNCTION gamma_moist( ptak, pqa )
929      !!----------------------------------------------------------------------------------
930      !!                           ***  FUNCTION gamma_moist  ***
931      !!
932      !! ** Purpose : Compute the moist adiabatic lapse-rate.
933      !!     => http://glossary.ametsoc.org/wiki/Moist-adiabatic_lapse_rate
934      !!     => http://www.geog.ucsb.edu/~joel/g266_s10/lecture_notes/chapt03/oh10_3_01/oh10_3_01.html
935      !!
936      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
937      !!----------------------------------------------------------------------------------
938      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   ptak          ! air temperature       [K]
939      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pqa           ! specific humidity [kg/kg]
940      REAL(wp), DIMENSION(jpi,jpj)             ::   gamma_moist   ! moist adiabatic lapse-rate
941      !
942      INTEGER  ::   ji, jj         ! dummy loop indices
943      REAL(wp) :: zrv, ziRT        ! local scalar
944      !!----------------------------------------------------------------------------------
945      !
946      DO jj = 1, jpj
947         DO ji = 1, jpi
948            zrv = pqa(ji,jj) / (1. - pqa(ji,jj))
949            ziRT = 1. / (R_dry*ptak(ji,jj))    ! 1/RT
950            gamma_moist(ji,jj) = grav * ( 1. + cevap*zrv*ziRT ) / ( Cp_dry + cevap*cevap*zrv*reps0*ziRT/ptak(ji,jj) )
951         END DO
952      END DO
953      !
954   END FUNCTION gamma_moist
955
956
957   FUNCTION L_vap( psst )
958      !!---------------------------------------------------------------------------------
959      !!                           ***  FUNCTION L_vap  ***
960      !!
961      !! ** Purpose : Compute the latent heat of vaporization of water from temperature
962      !!
963      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk)
964      !!----------------------------------------------------------------------------------
965      REAL(wp), DIMENSION(jpi,jpj)             ::   L_vap   ! latent heat of vaporization   [J/kg]
966      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   psst   ! water temperature                [K]
967      !!----------------------------------------------------------------------------------
968      !
969      L_vap = (  2.501 - 0.00237 * ( psst(:,:) - rt0)  ) * 1.e6
970      !
971   END FUNCTION L_vap
972
973
974#if defined key_lim3
975   SUBROUTINE Cdn10_Lupkes2012( Cd )
976      !!----------------------------------------------------------------------
977      !!                      ***  ROUTINE  Cdn10_Lupkes2012  ***
978      !!
979      !! ** Purpose :    Recompute the ice-atm drag at 10m height to make
980      !!                 it dependent on edges at leads, melt ponds and flows.
981      !!                 After some approximations, this can be resumed to a dependency
982      !!                 on ice concentration.
983      !!               
984      !! ** Method :     The parameterization is taken from Lupkes et al. (2012) eq.(50)
985      !!                 with the highest level of approximation: level4, eq.(59)
986      !!                 The generic drag over a cell partly covered by ice can be re-written as follows:
987      !!
988      !!                 Cd = Cdw * (1-A) + Cdi * A + Ce * (1-A)**(nu+1/(10*beta)) * A**mu
989      !!
990      !!                    Ce = 2.23e-3       , as suggested by Lupkes (eq. 59)
991      !!                    nu = mu = beta = 1 , as suggested by Lupkes (eq. 59)
992      !!                    A is the concentration of ice minus melt ponds (if any)
993      !!
994      !!                 This new drag has a parabolic shape (as a function of A) starting at
995      !!                 Cdw(say 1.5e-3) for A=0, reaching 1.97e-3 for A~0.5
996      !!                 and going down to Cdi(say 1.4e-3) for A=1
997      !!
998      !!                 It is theoretically applicable to all ice conditions (not only MIZ)
999      !!                 => see Lupkes et al (2013)
1000      !!
1001      !! ** References : Lupkes et al. JGR 2012 (theory)
1002      !!                 Lupkes et al. GRL 2013 (application to GCM)
1003      !!
1004      !!----------------------------------------------------------------------
1005      REAL(wp), DIMENSION(:,:), INTENT(inout) ::   Cd
1006      REAL(wp), PARAMETER ::   zCe   = 2.23e-03_wp
1007      REAL(wp), PARAMETER ::   znu   = 1._wp
1008      REAL(wp), PARAMETER ::   zmu   = 1._wp
1009      REAL(wp), PARAMETER ::   zbeta = 1._wp
1010      REAL(wp)            ::   zcoef
1011      !!----------------------------------------------------------------------
1012      zcoef = znu + 1._wp / ( 10._wp * zbeta )
1013
1014      ! generic drag over a cell partly covered by ice
1015      !!Cd(:,:) = Cd_oce(:,:) * ( 1._wp - at_i_b(:,:) ) +  &                        ! pure ocean drag
1016      !!   &      Cd_ice      *           at_i_b(:,:)   +  &                        ! pure ice drag
1017      !!   &      zCe         * ( 1._wp - at_i_b(:,:) )**zcoef * at_i_b(:,:)**zmu   ! change due to sea-ice morphology
1018
1019      ! ice-atm drag
1020      Cd(:,:) = Cd_ice +  &                                                         ! pure ice drag
1021         &      zCe    * ( 1._wp - at_i_b(:,:) )**zcoef * at_i_b(:,:)**(zmu-1._wp)  ! change due to sea-ice morphology
1022     
1023   END SUBROUTINE Cdn10_Lupkes2012
1024#endif
1025   
1026
1027   !!======================================================================
1028END MODULE sbcblk
Note: See TracBrowser for help on using the repository browser.