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 @ 6727

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

#1751 - branch SIMPLIF_6_aerobulk: minor correction

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