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

source: branches/UKMO/dev_r5518_STOPACK/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 11384

Last change on this file since 11384 was 11384, checked in by mattmartin, 5 years ago

Included Andrea Storto's STOPACK code into NEMO3.6 branch.

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