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_core.F90 in branches/UKMO/dev_r5518_GO6_package_FOAMv14/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/dev_r5518_GO6_package_FOAMv14/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 9288

Last change on this file since 9288 was 9288, checked in by anaguiar, 6 years ago

Changes to read weight of sprecip from namsbc_core, branch to replace dev_r5518_GO6_package_05sprecip misleading name

File size: 50.8 KB
Line 
1MODULE sbcblk_core
2   !!======================================================================
3   !!                       ***  MODULE  sbcblk_core  ***
4   !! Ocean forcing:  momentum, heat and freshwater flux formulation
5   !!=====================================================================
6   !! History :  1.0  !  2004-08  (U. Schweckendiek)  Original code
7   !!            2.0  !  2005-04  (L. Brodeau, A.M. Treguier) additions:
8   !!                           -  new bulk routine for efficiency
9   !!                           -  WINDS ARE NOW ASSUMED TO BE AT T POINTS in input files !!!!
10   !!                           -  file names and file characteristics in namelist
11   !!                           -  Implement reading of 6-hourly fields
12   !!            3.0  !  2006-06  (G. Madec) sbc rewritting
13   !!             -   !  2006-12  (L. Brodeau) Original code for turb_core_2z
14   !!            3.2  !  2009-04  (B. Lemaire)  Introduce iom_put
15   !!            3.3  !  2010-10  (S. Masson)  add diurnal cycle
16   !!            3.4  !  2011-11  (C. Harris) Fill arrays required by CICE
17   !!            3.7  !  2014-06  (L. Brodeau) simplification and optimization of CORE bulk
18   !!----------------------------------------------------------------------
19
20   !!----------------------------------------------------------------------
21   !!   sbc_blk_core    : bulk formulation as ocean surface boundary condition (forced mode, CORE bulk formulea)
22   !!   blk_oce_core    : computes momentum, heat and freshwater fluxes over ocean
23   !!   blk_ice_core    : computes momentum, heat and freshwater fluxes over ice
24   !!   turb_core_2z    : Computes turbulent transfert coefficients
25   !!   cd_neutral_10m  : Estimate of the neutral drag coefficient at 10m
26   !!   psi_m           : universal profile stability function for momentum
27   !!   psi_h           : universal profile stability function for temperature and humidity
28   !!----------------------------------------------------------------------
29   USE oce             ! ocean dynamics and tracers
30   USE dom_oce         ! ocean space and time domain
31   USE phycst          ! physical constants
32   USE fldread         ! read input fields
33   USE sbc_oce         ! Surface boundary condition: ocean fields
34   USE cyclone         ! Cyclone 10m wind form trac of cyclone centres
35   USE sbcdcy          ! surface boundary condition: diurnal cycle
36   USE iom             ! I/O manager library
37   USE in_out_manager  ! I/O manager
38   USE lib_mpp         ! distribued memory computing library
39   USE wrk_nemo        ! work arrays
40   USE timing          ! Timing
41   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
42   USE prtctl          ! Print control
43   USE sbcwave, ONLY   :  cdn_wave ! wave module
44   USE sbc_ice         ! Surface boundary condition: ice fields
45   USE lib_fortran     ! to use key_nosignedzero
46#if defined key_lim3
47   USE ice, ONLY       : u_ice, v_ice, jpl, pfrld, a_i_b
48   USE limthd_dh       ! for CALL lim_thd_snwblow
49#elif defined key_lim2
50   USE ice_2, ONLY     : u_ice, v_ice
51   USE par_ice_2
52#endif
53
54   IMPLICIT NONE
55   PRIVATE
56
57   PUBLIC   sbc_blk_core         ! routine called in sbcmod module
58#if defined key_lim2 || defined key_lim3
59   PUBLIC   blk_ice_core_tau     ! routine called in sbc_ice_lim module
60   PUBLIC   blk_ice_core_flx     ! routine called in sbc_ice_lim module
61#endif
62   PUBLIC   turb_core_2z         ! routine calles in sbcblk_mfs module
63
64   INTEGER , PARAMETER ::   jpfld   = 9           ! maximum number of files to read
65   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point
66   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point
67   INTEGER , PARAMETER ::   jp_humi = 3           ! index of specific humidity               ( % )
68   INTEGER , PARAMETER ::   jp_qsr  = 4           ! index of solar heat                      (W/m2)
69   INTEGER , PARAMETER ::   jp_qlw  = 5           ! index of Long wave                       (W/m2)
70   INTEGER , PARAMETER ::   jp_tair = 6           ! index of 10m air temperature             (Kelvin)
71   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s)
72   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s)
73   INTEGER , PARAMETER ::   jp_tdif = 9           ! index of tau diff associated to HF tau   (N/m2)   at T-point
74
75   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read)
76
77   !                                             !!! CORE bulk parameters
78   REAL(wp), PARAMETER ::   rhoa =    1.22        ! air density
79   REAL(wp), PARAMETER ::   cpa  = 1000.5         ! specific heat of air
80   REAL(wp), PARAMETER ::   Lv   =    2.5e6       ! latent heat of vaporization
81   REAL(wp), PARAMETER ::   Ls   =    2.839e6     ! latent heat of sublimation
82   REAL(wp), PARAMETER ::   Stef =    5.67e-8     ! Stefan Boltzmann constant
83   REAL(wp), PARAMETER ::   Cice =    1.4e-3      ! iovi 1.63e-3     ! transfer coefficient over ice
84   REAL(wp), PARAMETER ::   albo =    0.066       ! ocean albedo assumed to be constant
85
86   !                                  !!* Namelist namsbc_core : CORE bulk parameters
87   LOGICAL  ::   ln_taudif   ! logical flag to use the "mean of stress module - module of mean stress" data
88   REAL(wp) ::   rn_pfac     ! multiplication factor for precipitation
89   REAL(wp) ::   rn_efac     ! multiplication factor for evaporation (clem)
90   REAL(wp) ::   rn_vfac     ! multiplication factor for ice/ocean velocity in the calculation of wind stress (clem)
91   REAL(wp) ::   rn_zqt      ! z(q,t) : height of humidity and temperature measurements
92   REAL(wp) ::   rn_zu       ! z(u)   : height of wind measurements
93   REAL(wp), PUBLIC :: rn_sfac ! multiplication factor for snow precipitation over sea-ice
94
95   !! * Substitutions
96#  include "domzgr_substitute.h90"
97#  include "vectopt_loop_substitute.h90"
98   !!----------------------------------------------------------------------
99   !! NEMO/OPA 3.7 , NEMO-consortium (2014)
100   !! $Id$
101   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
102   !!----------------------------------------------------------------------
103CONTAINS
104
105   SUBROUTINE sbc_blk_core( kt )
106      !!---------------------------------------------------------------------
107      !!                    ***  ROUTINE sbc_blk_core  ***
108      !!
109      !! ** Purpose :   provide at each time step the surface ocean fluxes
110      !!      (momentum, heat, freshwater and runoff)
111      !!
112      !! ** Method  : (1) READ each fluxes in NetCDF files:
113      !!      the 10m wind velocity (i-component) (m/s)    at T-point
114      !!      the 10m wind velocity (j-component) (m/s)    at T-point
115      !!      the 10m or 2m specific humidity     ( % )
116      !!      the solar heat                      (W/m2)
117      !!      the Long wave                       (W/m2)
118      !!      the 10m or 2m air temperature       (Kelvin)
119      !!      the total precipitation (rain+snow) (Kg/m2/s)
120      !!      the snow (solid prcipitation)       (kg/m2/s)
121      !!      the tau diff associated to HF tau   (N/m2)   at T-point   (ln_taudif=T)
122      !!              (2) CALL blk_oce_core
123      !!
124      !!      C A U T I O N : never mask the surface stress fields
125      !!                      the stress is assumed to be in the (i,j) mesh referential
126      !!
127      !! ** Action  :   defined at each time-step at the air-sea interface
128      !!              - utau, vtau  i- and j-component of the wind stress
129      !!              - taum        wind stress module at T-point
130      !!              - wndm        wind speed  module at T-point over free ocean or leads in presence of sea-ice
131      !!              - qns, qsr    non-solar and solar heat fluxes
132      !!              - emp         upward mass flux (evapo. - precip.)
133      !!              - sfx         salt flux due to freezing/melting (non-zero only if ice is present)
134      !!                            (set in limsbc(_2).F90)
135      !!
136      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
137      !!                   Brodeau et al. Ocean Modelling 2010
138      !!----------------------------------------------------------------------
139      INTEGER, INTENT(in) ::   kt   ! ocean time step
140      !
141      INTEGER  ::   ierror   ! return error code
142      INTEGER  ::   ifpr     ! dummy loop indice
143      INTEGER  ::   jfld     ! dummy loop arguments
144      INTEGER  ::   ios      ! Local integer output status for namelist read
145      !
146      CHARACTER(len=100) ::  cn_dir   !   Root directory for location of core files
147      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i     ! array of namelist informations on the fields to read
148      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read
149      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !   "                                 "
150      TYPE(FLD_N) ::   sn_tdif                                 !   "                                 "
151      NAMELIST/namsbc_core/ cn_dir , ln_taudif, rn_pfac, rn_efac, rn_vfac,  &
152         &                  sn_wndi, sn_wndj, sn_humi  , sn_qsr ,           &
153         &                  sn_qlw , sn_tair, sn_prec  , sn_snow,           &
154         &                  sn_tdif, rn_zqt,  rn_zu, rn_sfac
155      !!---------------------------------------------------------------------
156      !
157      !                                         ! ====================== !
158      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
159         !                                      ! ====================== !
160         !
161         rn_sfac = 1._wp       ! Default to one if missing from namelist
162         REWIND( numnam_ref )              ! Namelist namsbc_core in reference namelist : CORE bulk parameters
163         READ  ( numnam_ref, namsbc_core, IOSTAT = ios, ERR = 901)
164901      IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in reference namelist', lwp )
165         !
166         REWIND( numnam_cfg )              ! Namelist namsbc_core in configuration namelist : CORE bulk parameters
167         READ  ( numnam_cfg, namsbc_core, IOSTAT = ios, ERR = 902 )
168902      IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in configuration namelist', lwp )
169
170         IF(lwm) WRITE( numond, namsbc_core )
171         !                                         ! check: do we plan to use ln_dm2dc with non-daily forcing?
172         IF( ln_dm2dc .AND. sn_qsr%nfreqh /= 24 )   &
173            &   CALL ctl_stop( 'sbc_blk_core: ln_dm2dc can be activated only with daily short-wave forcing' )
174         IF( ln_dm2dc .AND. sn_qsr%ln_tint ) THEN
175            CALL ctl_warn( 'sbc_blk_core: ln_dm2dc is taking care of the temporal interpolation of daily qsr',   &
176               &         '              ==> We force time interpolation = .false. for qsr' )
177            sn_qsr%ln_tint = .false.
178         ENDIF
179         !                                         ! store namelist information in an array
180         slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
181         slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
182         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
183         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
184         slf_i(jp_tdif) = sn_tdif
185         !
186         lhftau = ln_taudif                        ! do we use HF tau information?
187         jfld = jpfld - COUNT( (/.NOT. lhftau/) )
188         !
189         ALLOCATE( sf(jfld), STAT=ierror )         ! set sf structure
190         IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_core: unable to allocate sf structure' )
191         DO ifpr= 1, jfld
192            ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
193            IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
194         END DO
195         !                                         ! fill sf with slf_i and control print
196         CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_core', 'flux formulation for ocean surface boundary condition', 'namsbc_core' )
197         !
198         sfx(:,:) = 0._wp                          ! salt flux; zero unless ice is present (computed in limsbc(_2).F90)
199         !
200      ENDIF
201
202      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
203
204      !                                            ! compute the surface ocean fluxes using CORE bulk formulea
205      IF( MOD( kt - 1, nn_fsbc ) == 0 )   CALL blk_oce_core( kt, sf, sst_m, ssu_m, ssv_m )
206
207#if defined key_cice
208      IF( MOD( kt - 1, nn_fsbc ) == 0 )   THEN
209         qlw_ice(:,:,1)   = sf(jp_qlw)%fnow(:,:,1) 
210         IF( ln_dm2dc ) THEN ; qsr_ice(:,:,1) = sbc_dcy( sf(jp_qsr)%fnow(:,:,1) )
211         ELSE                ; qsr_ice(:,:,1) =          sf(jp_qsr)%fnow(:,:,1)
212         ENDIF
213         tatm_ice(:,:)    = sf(jp_tair)%fnow(:,:,1)         
214         qatm_ice(:,:)    = sf(jp_humi)%fnow(:,:,1)
215         tprecip(:,:)     = sf(jp_prec)%fnow(:,:,1) * rn_pfac
216         sprecip(:,:)     = sf(jp_snow)%fnow(:,:,1) * rn_pfac
217         wndi_ice(:,:)    = sf(jp_wndi)%fnow(:,:,1)
218         wndj_ice(:,:)    = sf(jp_wndj)%fnow(:,:,1)
219      ENDIF
220#endif
221      !
222   END SUBROUTINE sbc_blk_core
223   
224   
225   SUBROUTINE blk_oce_core( kt, sf, pst, pu, pv )
226      !!---------------------------------------------------------------------
227      !!                     ***  ROUTINE blk_core  ***
228      !!
229      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
230      !!      the ocean surface at each time step
231      !!
232      !! ** Method  :   CORE bulk formulea for the ocean using atmospheric
233      !!      fields read in sbc_read
234      !!
235      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
236      !!              - vtau    : j-component of the stress at V-point  (N/m2)
237      !!              - taum    : Wind stress module at T-point         (N/m2)
238      !!              - wndm    : Wind speed module at T-point          (m/s)
239      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
240      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
241      !!              - emp     : evaporation minus precipitation       (kg/m2/s)
242      !!
243      !!  ** Nota  :   sf has to be a dummy argument for AGRIF on NEC
244      !!---------------------------------------------------------------------
245      INTEGER  , INTENT(in   )                 ::   kt    ! time step index
246      TYPE(fld), INTENT(inout), DIMENSION(:)   ::   sf    ! input data
247      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pst   ! surface temperature                      [Celcius]
248      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pu    ! surface current at U-point (i-component) [m/s]
249      REAL(wp) , INTENT(in)   , DIMENSION(:,:) ::   pv    ! surface current at V-point (j-component) [m/s]
250      !
251      INTEGER  ::   ji, jj               ! dummy loop indices
252      REAL(wp) ::   zcoef_qsatw, zztmp   ! local variable
253      REAL(wp), DIMENSION(:,:), POINTER ::   zwnd_i, zwnd_j    ! wind speed components at T-point
254      REAL(wp), DIMENSION(:,:), POINTER ::   zqsatw            ! specific humidity at pst
255      REAL(wp), DIMENSION(:,:), POINTER ::   zqlw, zqsb        ! long wave and sensible heat fluxes
256      REAL(wp), DIMENSION(:,:), POINTER ::   zqla, zevap       ! latent heat fluxes and evaporation
257      REAL(wp), DIMENSION(:,:), POINTER ::   Cd                ! transfer coefficient for momentum      (tau)
258      REAL(wp), DIMENSION(:,:), POINTER ::   Ch                ! transfer coefficient for sensible heat (Q_sens)
259      REAL(wp), DIMENSION(:,:), POINTER ::   Ce                ! tansfert coefficient for evaporation   (Q_lat)
260      REAL(wp), DIMENSION(:,:), POINTER ::   zst               ! surface temperature in Kelvin
261      REAL(wp), DIMENSION(:,:), POINTER ::   zt_zu             ! air temperature at wind speed height
262      REAL(wp), DIMENSION(:,:), POINTER ::   zq_zu             ! air spec. hum.  at wind speed height
263      !!---------------------------------------------------------------------
264      !
265      IF( nn_timing == 1 )  CALL timing_start('blk_oce_core')
266      !
267      CALL wrk_alloc( jpi,jpj, zwnd_i, zwnd_j, zqsatw, zqlw, zqsb, zqla, zevap )
268      CALL wrk_alloc( jpi,jpj, Cd, Ch, Ce, zst, zt_zu, zq_zu )
269      !
270      ! local scalars ( place there for vector optimisation purposes)
271      zcoef_qsatw = 0.98 * 640380. / rhoa
272     
273      zst(:,:) = pst(:,:) + rt0      ! convert SST from Celcius to Kelvin (and set minimum value far above 0 K)
274
275      ! ----------------------------------------------------------------------------- !
276      !      0   Wind components and module at T-point relative to the moving ocean   !
277      ! ----------------------------------------------------------------------------- !
278
279      ! ... components ( U10m - U_oce ) at T-point (unmasked)
280      zwnd_i(:,:) = 0.e0 
281      zwnd_j(:,:) = 0.e0
282#if defined key_cyclone
283      CALL wnd_cyc( kt, zwnd_i, zwnd_j )    ! add analytical tropical cyclone (Vincent et al. JGR 2012)
284      DO jj = 2, jpjm1
285         DO ji = fs_2, fs_jpim1   ! vect. opt.
286            sf(jp_wndi)%fnow(ji,jj,1) = sf(jp_wndi)%fnow(ji,jj,1) + zwnd_i(ji,jj)
287            sf(jp_wndj)%fnow(ji,jj,1) = sf(jp_wndj)%fnow(ji,jj,1) + zwnd_j(ji,jj)
288         END DO
289      END DO
290#endif
291      DO jj = 2, jpjm1
292         DO ji = fs_2, fs_jpim1   ! vect. opt.
293            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
294            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
295         END DO
296      END DO
297      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
298      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
299      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
300      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
301         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
302
303      ! ----------------------------------------------------------------------------- !
304      !      I   Radiative FLUXES                                                     !
305      ! ----------------------------------------------------------------------------- !
306
307      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
308      zztmp = 1. - albo
309      IF( ln_dm2dc ) THEN   ;   qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
310      ELSE                  ;   qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
311      ENDIF
312
313      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
314      ! ----------------------------------------------------------------------------- !
315      !     II    Turbulent FLUXES                                                    !
316      ! ----------------------------------------------------------------------------- !
317
318      ! ... specific humidity at SST and IST
319      zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) )
320
321      ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point :
322      CALL turb_core_2z( rn_zqt, rn_zu, zst, sf(jp_tair)%fnow, zqsatw, sf(jp_humi)%fnow, wndm,   &
323         &               Cd, Ch, Ce, zt_zu, zq_zu )
324   
325      ! ... tau module, i and j component
326      DO jj = 1, jpj
327         DO ji = 1, jpi
328            zztmp = rhoa * wndm(ji,jj) * Cd(ji,jj)
329            taum  (ji,jj) = zztmp * wndm  (ji,jj)
330            zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj)
331            zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj)
332         END DO
333      END DO
334
335      ! ... add the HF tau contribution to the wind stress module?
336      IF( lhftau ) THEN
337         taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1)
338      ENDIF
339      CALL iom_put( "taum_oce", taum )   ! output wind stress module
340
341      ! ... utau, vtau at U- and V_points, resp.
342      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
343      !     Note the use of MAX(tmask(i,j),tmask(i+1,j) is to mask tau over ice shelves
344      DO jj = 1, jpjm1
345         DO ji = 1, fs_jpim1
346            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) ) &
347               &          * MAX(tmask(ji,jj,1),tmask(ji+1,jj,1))
348            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) ) &
349               &          * MAX(tmask(ji,jj,1),tmask(ji,jj+1,1))
350         END DO
351      END DO
352      CALL lbc_lnk( utau(:,:), 'U', -1. )
353      CALL lbc_lnk( vtau(:,:), 'V', -1. )
354
355   
356      !  Turbulent fluxes over ocean
357      ! -----------------------------
358      IF( ABS( rn_zu - rn_zqt) < 0.01_wp ) THEN
359         !! q_air and t_air are (or "are almost") given at 10m (wind reference height)
360         zevap(:,:) = rn_efac*MAX( 0._wp,     rhoa*Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:,1) )*wndm(:,:) ) ! Evaporation
361         zqsb (:,:) =                     cpa*rhoa*Ch(:,:)*( zst   (:,:) - sf(jp_tair)%fnow(:,:,1) )*wndm(:,:)   ! Sensible Heat
362      ELSE
363         !! q_air and t_air are not given at 10m (wind reference height)
364         ! Values of temp. and hum. adjusted to height of wind during bulk algorithm iteration must be used!!!
365         zevap(:,:) = rn_efac*MAX( 0._wp,     rhoa*Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) )*wndm(:,:) )   ! Evaporation
366         zqsb (:,:) =                     cpa*rhoa*Ch(:,:)*( zst   (:,:) - zt_zu(:,:) )*wndm(:,:)     ! Sensible Heat
367      ENDIF
368      zqla (:,:) = Lv * zevap(:,:)                                                              ! Latent Heat
369
370      IF(ln_ctl) THEN
371         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
372         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
373         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce_core: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
374         CALL prt_ctl( tab2d_1=zqsatw, clinfo1=' blk_oce_core: zqsatw : ', tab2d_2=zst, clinfo2=' zst : ' )
375         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
376            &          tab2d_2=vtau  , clinfo2=              ' vtau : '  , mask2=vmask )
377         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce_core: wndm   : ')
378         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce_core: zst    : ')
379      ENDIF
380       
381      ! ----------------------------------------------------------------------------- !
382      !     III    Total FLUXES                                                       !
383      ! ----------------------------------------------------------------------------- !
384      !
385      emp (:,:) = (  zevap(:,:)                                          &   ! mass flux (evap. - precip.)
386         &         - sf(jp_prec)%fnow(:,:,1) * rn_pfac  ) * tmask(:,:,1)
387      !
388      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                &   ! Downward Non Solar
389         &     - sf(jp_snow)%fnow(:,:,1) * rn_pfac * lfus                         &   ! remove latent melting heat for solid precip
390         &     - zevap(:,:) * pst(:,:) * rcp                                      &   ! remove evap heat content at SST
391         &     + ( sf(jp_prec)%fnow(:,:,1) - sf(jp_snow)%fnow(:,:,1) ) * rn_pfac  &   ! add liquid precip heat content at Tair
392         &     * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp                          &
393         &     + sf(jp_snow)%fnow(:,:,1) * rn_pfac                                &   ! add solid  precip heat content at min(Tair,Tsnow)
394         &     * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1)
395      !
396#if defined key_lim3
397      qns_oce(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                ! non solar without emp (only needed by LIM3)
398      qsr_oce(:,:) = qsr(:,:)
399#endif
400      !
401      IF ( nn_ice == 0 ) THEN
402         CALL iom_put( "qlw_oce" ,   zqlw )                 ! output downward longwave heat over the ocean
403         CALL iom_put( "qsb_oce" , - zqsb )                 ! output downward sensible heat over the ocean
404         CALL iom_put( "qla_oce" , - zqla )                 ! output downward latent   heat over the ocean
405         CALL iom_put( "qemp_oce",   qns-zqlw+zqsb+zqla )   ! output downward heat content of E-P over the ocean
406         CALL iom_put( "qns_oce" ,   qns  )                 ! output downward non solar heat over the ocean
407         CALL iom_put( "qsr_oce" ,   qsr  )                 ! output downward solar heat over the ocean
408         CALL iom_put( "qt_oce"  ,   qns+qsr )              ! output total downward heat over the ocean
409         tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac   ! output total precipitation [kg/m2/s]
410         sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac   ! output solid precipitation [kg/m2/s]
411         CALL iom_put( 'snowpre', sprecip * 86400. )        ! Snow
412         CALL iom_put( 'precip' , tprecip * 86400. )        ! Total precipitation
413      ENDIF
414      !
415      IF(ln_ctl) THEN
416         CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=zqlw , clinfo2=' zqlw  : ')
417         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
418         CALL prt_ctl(tab2d_1=pst  , clinfo1=' blk_oce_core: pst    : ', tab2d_2=emp  , clinfo2=' emp   : ')
419         CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
420            &         tab2d_2=vtau , clinfo2=              ' vtau  : ' , mask2=vmask )
421      ENDIF
422      !
423      CALL wrk_dealloc( jpi,jpj, zwnd_i, zwnd_j, zqsatw, zqlw, zqsb, zqla, zevap )
424      CALL wrk_dealloc( jpi,jpj, Cd, Ch, Ce, zst, zt_zu, zq_zu )
425      !
426      IF( nn_timing == 1 )  CALL timing_stop('blk_oce_core')
427      !
428   END SUBROUTINE blk_oce_core
429 
430   
431#if defined key_lim2 || defined key_lim3
432   SUBROUTINE blk_ice_core_tau
433      !!---------------------------------------------------------------------
434      !!                     ***  ROUTINE blk_ice_core_tau  ***
435      !!
436      !! ** Purpose :   provide the surface boundary condition over sea-ice
437      !!
438      !! ** Method  :   compute momentum using CORE bulk
439      !!                formulea, ice variables and read atmospheric fields.
440      !!                NB: ice drag coefficient is assumed to be a constant
441      !!---------------------------------------------------------------------
442      INTEGER  ::   ji, jj    ! dummy loop indices
443      REAL(wp) ::   zcoef_wnorm, zcoef_wnorm2
444      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f               ! relative wind module and components at F-point
445      REAL(wp) ::             zwndi_t , zwndj_t               ! relative wind components at T-point
446      !!---------------------------------------------------------------------
447      !
448      IF( nn_timing == 1 )  CALL timing_start('blk_ice_core_tau')
449      !
450      ! local scalars ( place there for vector optimisation purposes)
451      zcoef_wnorm  = rhoa * Cice
452      zcoef_wnorm2 = rhoa * Cice * 0.5
453
454!!gm brutal....
455      utau_ice  (:,:) = 0._wp
456      vtau_ice  (:,:) = 0._wp
457      wndm_ice  (:,:) = 0._wp
458!!gm end
459
460      ! ----------------------------------------------------------------------------- !
461      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
462      ! ----------------------------------------------------------------------------- !
463      SELECT CASE( cp_ice_msh )
464      CASE( 'I' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
465         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
466         DO jj = 2, jpjm1
467            DO ji = 2, jpim1   ! B grid : NO vector opt
468               ! ... scalar wind at I-point (fld being at T-point)
469               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ,1) + sf(jp_wndi)%fnow(ji  ,jj  ,1)   &
470                  &              + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji  ,jj-1,1)  ) - rn_vfac * u_ice(ji,jj)
471               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ,1) + sf(jp_wndj)%fnow(ji  ,jj  ,1)   &
472                  &              + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji  ,jj-1,1)  ) - rn_vfac * v_ice(ji,jj)
473               zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
474               ! ... ice stress at I-point
475               utau_ice(ji,jj) = zwnorm_f * zwndi_f
476               vtau_ice(ji,jj) = zwnorm_f * zwndj_f
477               ! ... scalar wind at T-point (fld being at T-point)
478               zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  u_ice(ji,jj+1) + u_ice(ji+1,jj+1)   &
479                  &                                                    + u_ice(ji,jj  ) + u_ice(ji+1,jj  )  )
480               zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  v_ice(ji,jj+1) + v_ice(ji+1,jj+1)   &
481                  &                                                    + v_ice(ji,jj  ) + v_ice(ji+1,jj  )  )
482               wndm_ice(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
483            END DO
484         END DO
485         CALL lbc_lnk( utau_ice, 'I', -1. )
486         CALL lbc_lnk( vtau_ice, 'I', -1. )
487         CALL lbc_lnk( wndm_ice, 'T',  1. )
488         !
489      CASE( 'C' )                  ! C-grid ice dynamics :   U & V-points (same as ocean)
490         DO jj = 2, jpj
491            DO ji = fs_2, jpi   ! vect. opt.
492               zwndi_t = (  sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( u_ice(ji-1,jj  ) + u_ice(ji,jj) )  )
493               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 0.5 * ( v_ice(ji  ,jj-1) + v_ice(ji,jj) )  )
494               wndm_ice(ji,jj) = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
495            END DO
496         END DO
497         DO jj = 2, jpjm1
498            DO ji = fs_2, fs_jpim1   ! vect. opt.
499               utau_ice(ji,jj) = zcoef_wnorm2 * ( wndm_ice(ji+1,jj  ) + wndm_ice(ji,jj) )                          &
500                  &          * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj,1) + sf(jp_wndi)%fnow(ji,jj,1) ) - rn_vfac * u_ice(ji,jj) )
501               vtau_ice(ji,jj) = zcoef_wnorm2 * ( wndm_ice(ji,jj+1  ) + wndm_ice(ji,jj) )                          &
502                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - rn_vfac * v_ice(ji,jj) )
503            END DO
504         END DO
505         CALL lbc_lnk( utau_ice, 'U', -1. )
506         CALL lbc_lnk( vtau_ice, 'V', -1. )
507         CALL lbc_lnk( wndm_ice, 'T',  1. )
508         !
509      END SELECT
510
511      IF(ln_ctl) THEN
512         CALL prt_ctl(tab2d_1=utau_ice  , clinfo1=' blk_ice_core: utau_ice : ', tab2d_2=vtau_ice  , clinfo2=' vtau_ice : ')
513         CALL prt_ctl(tab2d_1=wndm_ice  , clinfo1=' blk_ice_core: wndm_ice : ')
514      ENDIF
515
516      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_core_tau')
517     
518   END SUBROUTINE blk_ice_core_tau
519
520
521   SUBROUTINE blk_ice_core_flx( ptsu, palb )
522      !!---------------------------------------------------------------------
523      !!                     ***  ROUTINE blk_ice_core_flx  ***
524      !!
525      !! ** Purpose :   provide the surface boundary condition over sea-ice
526      !!
527      !! ** Method  :   compute heat and freshwater exchanged
528      !!                between atmosphere and sea-ice using CORE bulk
529      !!                formulea, ice variables and read atmmospheric fields.
530      !!
531      !! caution : the net upward water flux has with mm/day unit
532      !!---------------------------------------------------------------------
533      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   ptsu          ! sea ice surface temperature
534      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   palb          ! ice albedo (all skies)
535      !!
536      INTEGER  ::   ji, jj, jl    ! dummy loop indices
537      REAL(wp) ::   zst2, zst3
538      REAL(wp) ::   zcoef_dqlw, zcoef_dqla, zcoef_dqsb
539      REAL(wp) ::   zztmp, z1_lsub                               ! temporary variable
540      !!
541      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qlw             ! long wave heat flux over ice
542      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qsb             ! sensible  heat flux over ice
543      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqlw            ! long wave heat sensitivity over ice
544      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqsb            ! sensible  heat sensitivity over ice
545      REAL(wp), DIMENSION(:,:)  , POINTER ::   zevap, zsnw       ! evaporation and snw distribution after wind blowing (LIM3)
546      !!---------------------------------------------------------------------
547      !
548      IF( nn_timing == 1 )  CALL timing_start('blk_ice_core_flx')
549      !
550      CALL wrk_alloc( jpi,jpj,jpl, z_qlw, z_qsb, z_dqlw, z_dqsb ) 
551
552      ! local scalars ( place there for vector optimisation purposes)
553      zcoef_dqlw   = 4.0 * 0.95 * Stef
554      zcoef_dqla   = -Ls * Cice * 11637800. * (-5897.8)
555      zcoef_dqsb   = rhoa * cpa * Cice
556
557      zztmp = 1. / ( 1. - albo )
558      !                                     ! ========================== !
559      DO jl = 1, jpl                        !  Loop over ice categories  !
560         !                                  ! ========================== !
561         DO jj = 1 , jpj
562            DO ji = 1, jpi
563               ! ----------------------------!
564               !      I   Radiative FLUXES   !
565               ! ----------------------------!
566               zst2 = ptsu(ji,jj,jl) * ptsu(ji,jj,jl)
567               zst3 = ptsu(ji,jj,jl) * zst2
568               ! Short Wave (sw)
569               qsr_ice(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
570               ! Long  Wave (lw)
571               z_qlw(ji,jj,jl) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj,1) - Stef * ptsu(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
572               ! lw sensitivity
573               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3                                               
574
575               ! ----------------------------!
576               !     II    Turbulent FLUXES  !
577               ! ----------------------------!
578
579               ! ... turbulent heat fluxes
580               ! Sensible Heat
581               z_qsb(ji,jj,jl) = rhoa * cpa * Cice * wndm_ice(ji,jj) * ( ptsu(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) )
582               ! Latent Heat
583               qla_ice(ji,jj,jl) = rn_efac * MAX( 0.e0, rhoa * Ls  * Cice * wndm_ice(ji,jj)   &                           
584                  &                         * (  11637800. * EXP( -5897.8 / ptsu(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj,1)  ) )
585              ! Latent heat sensitivity for ice (Dqla/Dt)
586               IF( qla_ice(ji,jj,jl) > 0._wp ) THEN
587                  dqla_ice(ji,jj,jl) = rn_efac * zcoef_dqla * wndm_ice(ji,jj) / ( zst2 ) * EXP( -5897.8 / ptsu(ji,jj,jl) )
588               ELSE
589                  dqla_ice(ji,jj,jl) = 0._wp
590               ENDIF
591
592               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
593               z_dqsb(ji,jj,jl) = zcoef_dqsb * wndm_ice(ji,jj)
594
595               ! ----------------------------!
596               !     III    Total FLUXES     !
597               ! ----------------------------!
598               ! Downward Non Solar flux
599               qns_ice (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - qla_ice (ji,jj,jl)
600               ! Total non solar heat flux sensitivity for ice
601               dqns_ice(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + dqla_ice(ji,jj,jl) )
602            END DO
603            !
604         END DO
605         !
606      END DO
607      !
608      tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
609      sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
610      CALL iom_put( 'snowpre', sprecip * 86400. )                  ! Snow precipitation
611      CALL iom_put( 'precip' , tprecip * 86400. )                  ! Total precipitation
612
613#if defined  key_lim3
614      CALL wrk_alloc( jpi,jpj, zevap, zsnw ) 
615
616      ! --- evaporation --- !
617      z1_lsub = 1._wp / Lsub
618      evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * z1_lsub    ! sublimation
619      devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * z1_lsub    ! d(sublimation)/dT
620      zevap    (:,:)   = rn_efac * ( emp(:,:) + tprecip(:,:) )  ! evaporation over ocean
621
622      ! --- evaporation minus precipitation --- !
623      zsnw(:,:) = 0._wp
624      CALL lim_thd_snwblow( pfrld, zsnw )  ! snow distribution over ice after wind blowing
625      emp_oce(:,:) = pfrld(:,:) * zevap(:,:) - ( tprecip(:,:) - sprecip(:,:) ) - sprecip(:,:) * (1._wp - zsnw )
626      emp_ice(:,:) = SUM( a_i_b(:,:,:) * evap_ice(:,:,:), dim=3 ) - sprecip(:,:) * zsnw
627      emp_tot(:,:) = emp_oce(:,:) + emp_ice(:,:)
628
629      ! --- heat flux associated with emp --- !
630      qemp_oce(:,:) = - pfrld(:,:) * zevap(:,:) * sst_m(:,:) * rcp                               & ! evap at sst
631         &          + ( tprecip(:,:) - sprecip(:,:) ) * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp  & ! liquid precip at Tair
632         &          +   sprecip(:,:) * ( 1._wp - zsnw ) *                                        & ! solid precip at min(Tair,Tsnow)
633         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
634      qemp_ice(:,:) =   sprecip(:,:) * zsnw *                                                    & ! solid precip (only)
635         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
636
637      ! --- total solar and non solar fluxes --- !
638      qns_tot(:,:) = pfrld(:,:) * qns_oce(:,:) + SUM( a_i_b(:,:,:) * qns_ice(:,:,:), dim=3 ) + qemp_ice(:,:) + qemp_oce(:,:)
639      qsr_tot(:,:) = pfrld(:,:) * qsr_oce(:,:) + SUM( a_i_b(:,:,:) * qsr_ice(:,:,:), dim=3 )
640
641      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
642      qprec_ice(:,:) = rhosn * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
643
644      ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) --- !
645      DO jl = 1, jpl
646         qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * cpic * tmask(:,:,1) )
647                                   ! But we do not have Tice => consider it at 0°C => evap=0
648      END DO
649
650      CALL wrk_dealloc( jpi,jpj, zevap, zsnw ) 
651#endif
652
653      !--------------------------------------------------------------------
654      ! FRACTIONs of net shortwave radiation which is not absorbed in the
655      ! thin surface layer and penetrates inside the ice cover
656      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
657      !
658      fr1_i0(:,:) = ( 0.18 * ( 1.0 - cldf_ice ) + 0.35 * cldf_ice )
659      fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice )
660      !
661      !
662      IF(ln_ctl) THEN
663         CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' blk_ice_core: qla_ice  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb    : ', kdim=jpl)
664         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice_core: z_qlw    : ', tab3d_2=dqla_ice, clinfo2=' dqla_ice : ', kdim=jpl)
665         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice_core: z_dqsb   : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw   : ', kdim=jpl)
666         CALL prt_ctl(tab3d_1=dqns_ice, clinfo1=' blk_ice_core: dqns_ice : ', tab3d_2=qsr_ice , clinfo2=' qsr_ice  : ', kdim=jpl)
667         CALL prt_ctl(tab3d_1=ptsu    , clinfo1=' blk_ice_core: ptsu     : ', tab3d_2=qns_ice , clinfo2=' qns_ice  : ', kdim=jpl)
668         CALL prt_ctl(tab2d_1=tprecip , clinfo1=' blk_ice_core: tprecip  : ', tab2d_2=sprecip , clinfo2=' sprecip  : ')
669      ENDIF
670
671      CALL wrk_dealloc( jpi,jpj,jpl, z_qlw, z_qsb, z_dqlw, z_dqsb )
672      !
673      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_core_flx')
674     
675   END SUBROUTINE blk_ice_core_flx
676#endif
677
678   SUBROUTINE turb_core_2z( zt, zu, sst, T_zt, q_sat, q_zt, dU,    &
679      &                      Cd, Ch, Ce , T_zu, q_zu )
680      !!----------------------------------------------------------------------
681      !!                      ***  ROUTINE  turb_core  ***
682      !!
683      !! ** Purpose :   Computes turbulent transfert coefficients of surface
684      !!                fluxes according to Large & Yeager (2004) and Large & Yeager (2008)
685      !!                If relevant (zt /= zu), adjust temperature and humidity from height zt to zu
686      !!
687      !! ** Method : Monin Obukhov Similarity Theory
688      !!             + Large & Yeager (2004,2008) closure: CD_n10 = f(U_n10)
689      !!
690      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
691      !!
692      !! ** Last update: Laurent Brodeau, June 2014:
693      !!    - handles both cases zt=zu and zt/=zu
694      !!    - optimized: less 2D arrays allocated and less operations
695      !!    - better first guess of stability by checking air-sea difference of virtual temperature
696      !!       rather than temperature difference only...
697      !!    - added function "cd_neutral_10m" that uses the improved parametrization of
698      !!      Large & Yeager 2008. Drag-coefficient reduction for Cyclone conditions!
699      !!    - using code-wide physical constants defined into "phycst.mod" rather than redifining them
700      !!      => 'vkarmn' and 'grav'
701      !!----------------------------------------------------------------------
702      REAL(wp), INTENT(in   )                     ::   zt       ! height for T_zt and q_zt                   [m]
703      REAL(wp), INTENT(in   )                     ::   zu       ! height for dU                              [m]
704      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   sst      ! sea surface temperature              [Kelvin]
705      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   T_zt     ! potential air temperature            [Kelvin]
706      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_sat    ! sea surface specific humidity         [kg/kg]
707      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity                 [kg/kg]
708      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   dU       ! relative wind module at zu            [m/s]
709      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
710      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
711      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
712      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   T_zu     ! air temp. shifted at zu                     [K]
713      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. hum.  shifted at zu               [kg/kg]
714      !
715      INTEGER ::   j_itt
716      INTEGER , PARAMETER ::   nb_itt = 5       ! number of itterations
717      LOGICAL ::   l_zt_equal_zu = .FALSE.      ! if q and t are given at different height than U
718      !
719      REAL(wp), DIMENSION(:,:), POINTER ::   U_zu          ! relative wind at zu                            [m/s]
720      REAL(wp), DIMENSION(:,:), POINTER ::   Ce_n10        ! 10m neutral latent coefficient
721      REAL(wp), DIMENSION(:,:), POINTER ::   Ch_n10        ! 10m neutral sensible coefficient
722      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd_n10   ! root square of Cd_n10
723      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd       ! root square of Cd
724      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_u        ! stability parameter at height zu
725      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_t        ! stability parameter at height zt
726      REAL(wp), DIMENSION(:,:), POINTER ::   zpsi_h_u, zpsi_m_u
727      REAL(wp), DIMENSION(:,:), POINTER ::   ztmp0, ztmp1, ztmp2
728      REAL(wp), DIMENSION(:,:), POINTER ::   stab          ! 1st stability test integer
729      !!----------------------------------------------------------------------
730
731      IF( nn_timing == 1 )  CALL timing_start('turb_core_2z')
732   
733      CALL wrk_alloc( jpi,jpj, U_zu, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd )
734      CALL wrk_alloc( jpi,jpj, zeta_u, stab )
735      CALL wrk_alloc( jpi,jpj, zpsi_h_u, zpsi_m_u, ztmp0, ztmp1, ztmp2 )
736
737      l_zt_equal_zu = .FALSE.
738      IF( ABS(zu - zt) < 0.01 ) l_zt_equal_zu = .TRUE.    ! testing "zu == zt" is risky with double precision
739
740      IF( .NOT. l_zt_equal_zu )   CALL wrk_alloc( jpi,jpj, zeta_t )
741
742      U_zu = MAX( 0.5 , dU )   !  relative wind speed at zu (normally 10m), we don't want to fall under 0.5 m/s
743
744      !! First guess of stability:
745      ztmp0 = T_zt*(1. + 0.608*q_zt) - sst*(1. + 0.608*q_sat) ! air-sea difference of virtual pot. temp. at zt
746      stab  = 0.5 + sign(0.5,ztmp0)                           ! stab = 1 if dTv > 0  => STABLE, 0 if unstable
747
748      !! Neutral coefficients at 10m:
749      IF( ln_cdgw ) THEN      ! wave drag case
750         cdn_wave(:,:) = cdn_wave(:,:) + rsmall * ( 1._wp - tmask(:,:,1) )
751         ztmp0   (:,:) = cdn_wave(:,:)
752      ELSE
753         ztmp0 = cd_neutral_10m( U_zu )
754      ENDIF
755      sqrt_Cd_n10 = SQRT( ztmp0 )
756      Ce_n10  = 1.e-3*( 34.6 * sqrt_Cd_n10 )
757      Ch_n10  = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab))
758   
759      !! Initializing transf. coeff. with their first guess neutral equivalents :
760      Cd = ztmp0   ;   Ce = Ce_n10   ;   Ch = Ch_n10   ;   sqrt_Cd = sqrt_Cd_n10
761
762      !! Initializing values at z_u with z_t values:
763      T_zu = T_zt   ;   q_zu = q_zt
764
765      !!  * Now starting iteration loop
766      DO j_itt=1, nb_itt
767         !
768         ztmp1 = T_zu - sst   ! Updating air/sea differences
769         ztmp2 = q_zu - q_sat 
770
771         ! Updating turbulent scales :   (L&Y 2004 eq. (7))
772         ztmp1  = Ch/sqrt_Cd*ztmp1    ! theta*
773         ztmp2  = Ce/sqrt_Cd*ztmp2    ! q*
774       
775         ztmp0 = T_zu*(1. + 0.608*q_zu) ! virtual potential temperature at zu
776
777         ! Estimate the inverse of Monin-Obukov length (1/L) at height zu:
778         ztmp0 =  (vkarmn*grav/ztmp0*(ztmp1*(1.+0.608*q_zu) + 0.608*T_zu*ztmp2)) / (Cd*U_zu*U_zu) 
779         !                                                                     ( Cd*U_zu*U_zu is U*^2 at zu)
780
781         !! Stability parameters :
782         zeta_u   = zu*ztmp0   ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
783         zpsi_h_u = psi_h( zeta_u )
784         zpsi_m_u = psi_m( zeta_u )
785       
786         !! Shifting temperature and humidity at zu (L&Y 2004 eq. (9b-9c))
787         IF ( .NOT. l_zt_equal_zu ) THEN
788            zeta_t = zt*ztmp0 ;  zeta_t = sign( min(abs(zeta_t),10.0), zeta_t )
789            stab = LOG(zu/zt) - zpsi_h_u + psi_h(zeta_t)  ! stab just used as temp array!!!
790            T_zu = T_zt + ztmp1/vkarmn*stab    ! ztmp1 is still theta*
791            q_zu = q_zt + ztmp2/vkarmn*stab    ! ztmp2 is still q*
792            q_zu = max(0., q_zu)
793         END IF
794       
795         IF( ln_cdgw ) THEN      ! surface wave case
796            sqrt_Cd = vkarmn / ( vkarmn / sqrt_Cd_n10 - zpsi_m_u ) 
797            Cd      = sqrt_Cd * sqrt_Cd
798         ELSE
799           ! Update neutral wind speed at 10m and neutral Cd at 10m (L&Y 2004 eq. 9a)...
800           !   In very rare low-wind conditions, the old way of estimating the
801           !   neutral wind speed at 10m leads to a negative value that causes the code
802           !   to crash. To prevent this a threshold of 0.25m/s is imposed.
803           ztmp0 = MAX( 0.25 , U_zu/(1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - zpsi_m_u)) ) !  U_n10
804           ztmp0 = cd_neutral_10m(ztmp0)                                                 ! Cd_n10
805           sqrt_Cd_n10 = sqrt(ztmp0)
806       
807           Ce_n10  = 1.e-3 * (34.6 * sqrt_Cd_n10)                     ! L&Y 2004 eq. (6b)
808           stab    = 0.5 + sign(0.5,zeta_u)                           ! update stability
809           Ch_n10  = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab))  ! L&Y 2004 eq. (6c-6d)
810
811           !! Update of transfer coefficients:
812           ztmp1 = 1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - zpsi_m_u)   ! L&Y 2004 eq. (10a)
813           Cd      = ztmp0 / ( ztmp1*ztmp1 )   
814           sqrt_Cd = SQRT( Cd )
815         ENDIF
816         !
817         ztmp0 = (LOG(zu/10.) - zpsi_h_u) / vkarmn / sqrt_Cd_n10
818         ztmp2 = sqrt_Cd / sqrt_Cd_n10
819         ztmp1 = 1. + Ch_n10*ztmp0               
820         Ch  = Ch_n10*ztmp2 / ztmp1  ! L&Y 2004 eq. (10b)
821         !
822         ztmp1 = 1. + Ce_n10*ztmp0               
823         Ce  = Ce_n10*ztmp2 / ztmp1  ! L&Y 2004 eq. (10c)
824         !
825      END DO
826
827      CALL wrk_dealloc( jpi,jpj, U_zu, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd )
828      CALL wrk_dealloc( jpi,jpj, zeta_u, stab )
829      CALL wrk_dealloc( jpi,jpj, zpsi_h_u, zpsi_m_u, ztmp0, ztmp1, ztmp2 )
830
831      IF( .NOT. l_zt_equal_zu ) CALL wrk_dealloc( jpi,jpj, zeta_t )
832
833      IF( nn_timing == 1 )  CALL timing_stop('turb_core_2z')
834      !
835   END SUBROUTINE turb_core_2z
836
837
838   FUNCTION cd_neutral_10m( zw10 )
839      !!----------------------------------------------------------------------
840      !! Estimate of the neutral drag coefficient at 10m as a function
841      !! of neutral wind  speed at 10m
842      !!
843      !! Origin: Large & Yeager 2008 eq.(11a) and eq.(11b)
844      !!
845      !! Author: L. Brodeau, june 2014
846      !!----------------------------------------------------------------------   
847      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   zw10           ! scalar wind speed at 10m (m/s)
848      REAL(wp), DIMENSION(jpi,jpj)             ::   cd_neutral_10m
849      !
850      REAL(wp), DIMENSION(:,:), POINTER ::   rgt33
851      !!----------------------------------------------------------------------   
852      !
853      CALL wrk_alloc( jpi,jpj, rgt33 )
854      !
855      !! When wind speed > 33 m/s => Cyclone conditions => special treatment
856      rgt33 = 0.5_wp + SIGN( 0.5_wp, (zw10 - 33._wp) )   ! If zw10 < 33. => 0, else => 1 
857      cd_neutral_10m = 1.e-3 * ( &
858         &       (1._wp - rgt33)*( 2.7_wp/zw10 + 0.142_wp + zw10/13.09_wp - 3.14807E-10*zw10**6) & ! zw10< 33.
859         &      + rgt33         *      2.34   )                                                    ! zw10 >= 33.
860      !
861      CALL wrk_dealloc( jpi,jpj, rgt33)
862      !
863   END FUNCTION cd_neutral_10m
864
865
866   FUNCTION psi_m(pta)   !! Psis, L&Y 2004 eq. (8c), (8d), (8e)
867      !-------------------------------------------------------------------------------
868      ! universal profile stability function for momentum
869      !-------------------------------------------------------------------------------
870      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pta
871      !
872      REAL(wp), DIMENSION(jpi,jpj)             :: psi_m
873      REAL(wp), DIMENSION(:,:), POINTER        :: X2, X, stabit
874      !-------------------------------------------------------------------------------
875      !
876      CALL wrk_alloc( jpi,jpj, X2, X, stabit )
877      !
878      X2 = SQRT( ABS( 1. - 16.*pta ) )  ;  X2 = MAX( X2 , 1. )   ;   X = SQRT( X2 )
879      stabit = 0.5 + SIGN( 0.5 , pta )
880      psi_m = -5.*pta*stabit  &                                                          ! Stable
881         &    + (1. - stabit)*(2.*LOG((1. + X)*0.5) + LOG((1. + X2)*0.5) - 2.*ATAN(X) + rpi*0.5)  ! Unstable
882      !
883      CALL wrk_dealloc( jpi,jpj, X2, X, stabit )
884      !
885   END FUNCTION psi_m
886
887
888   FUNCTION psi_h( pta )    !! Psis, L&Y 2004 eq. (8c), (8d), (8e)
889      !-------------------------------------------------------------------------------
890      ! universal profile stability function for temperature and humidity
891      !-------------------------------------------------------------------------------
892      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pta
893      !
894      REAL(wp), DIMENSION(jpi,jpj)             ::   psi_h
895      REAL(wp), DIMENSION(:,:), POINTER        ::   X2, X, stabit
896      !-------------------------------------------------------------------------------
897      !
898      CALL wrk_alloc( jpi,jpj, X2, X, stabit )
899      !
900      X2 = SQRT( ABS( 1. - 16.*pta ) )   ;   X2 = MAX( X2 , 1. )   ;   X = SQRT( X2 )
901      stabit = 0.5 + SIGN( 0.5 , pta )
902      psi_h = -5.*pta*stabit   &                                       ! Stable
903         &    + (1. - stabit)*(2.*LOG( (1. + X2)*0.5 ))                ! Unstable
904      !
905      CALL wrk_dealloc( jpi,jpj, X2, X, stabit )
906      !
907   END FUNCTION psi_h
908
909   !!======================================================================
910END MODULE sbcblk_core
Note: See TracBrowser for help on using the repository browser.