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

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

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

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

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