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

source: branches/UKMO/test_moci_test_suite_namelist_read/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 9366

Last change on this file since 9366 was 9366, checked in by andmirek, 6 years ago

#2050 first version. Compiled OK in moci test suite

File size: 52.4 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   PRIVATE  core_namelist
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) ::   rn_zqt      ! z(q,t) : height of humidity and temperature measurements
93   REAL(wp) ::   rn_zu       ! z(u)   : height of wind measurements
94   LOGICAL  ::   ln_core_sio ! single processor read flag
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, ln_core_sio
156      !!---------------------------------------------------------------------
157      !
158      !                                         ! ====================== !
159      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
160         !                                      ! ====================== !
161         !
162         ln_core_sio = .FALSE.
163         IF(lwm) THEN
164            REWIND( numnam_ref )              ! Namelist namsbc_core in reference namelist : CORE bulk parameters
165            READ  ( numnam_ref, namsbc_core, IOSTAT = ios, ERR = 901)
166901         IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in reference namelist', lwm )
167            REWIND( numnam_cfg )              ! Namelist namsbc_core in configuration namelist : CORE bulk parameters
168            READ  ( numnam_cfg, namsbc_core, IOSTAT = ios, ERR = 902 )
169902         IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_core in configuration namelist', lwm )
170         ENDIF
171
172         IF(lwm) WRITE( numond, namsbc_core )
173         !                                         ! check: do we plan to use ln_dm2dc with non-daily forcing?
174         CALL core_namelist(cn_dir, sn_wndi, sn_wndj, sn_humi, sn_qsr, sn_qlw, &
175         &                  sn_tair, sn_prec, sn_snow, sn_tdif)
176         !
177         IF( ln_dm2dc .AND. sn_qsr%nfreqh /= 24 )   &
178            &   CALL ctl_stop( 'sbc_blk_core: ln_dm2dc can be activated only with daily short-wave forcing' )
179         IF( ln_dm2dc .AND. sn_qsr%ln_tint ) THEN
180            CALL ctl_warn( 'sbc_blk_core: ln_dm2dc is taking care of the temporal interpolation of daily qsr',   &
181               &         '              ==> We force time interpolation = .false. for qsr' )
182            sn_qsr%ln_tint = .false.
183         ENDIF
184         !                                         ! store namelist information in an array
185         slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
186         slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
187         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
188         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
189         slf_i(jp_tdif) = sn_tdif
190         !
191         lhftau = ln_taudif                        ! do we use HF tau information?
192         jfld = jpfld - COUNT( (/.NOT. lhftau/) )
193         !
194         ALLOCATE( sf(jfld), STAT=ierror )         ! set sf structure
195         IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_core: unable to allocate sf structure' )
196         DO ifpr= 1, jfld
197            ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
198            IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
199         END DO
200         !                                         ! fill sf with slf_i and control print
201         CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_core', 'flux formulation for ocean surface boundary condition', 'namsbc_core' )
202         !
203         sfx(:,:) = 0._wp                          ! salt flux; zero unless ice is present (computed in limsbc(_2).F90)
204         !
205      ENDIF
206      lspr = ln_core_sio
207      CALL fld_read( kt, nn_fsbc, sf )             ! input fields provided at the current time-step
208      lspr = .false.
209
210      !                                            ! compute the surface ocean fluxes using CORE bulk formulea
211      IF( MOD( kt - 1, nn_fsbc ) == 0 )   CALL blk_oce_core( kt, sf, sst_m, ssu_m, ssv_m )
212
213#if defined key_cice
214      IF( MOD( kt - 1, nn_fsbc ) == 0 )   THEN
215         qlw_ice(:,:,1)   = sf(jp_qlw)%fnow(:,:,1) 
216         IF( ln_dm2dc ) THEN ; qsr_ice(:,:,1) = sbc_dcy( sf(jp_qsr)%fnow(:,:,1) )
217         ELSE                ; qsr_ice(:,:,1) =          sf(jp_qsr)%fnow(:,:,1)
218         ENDIF
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_vfac * 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
300            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 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         tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac   ! output total precipitation [kg/m2/s]
416         sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac   ! output solid precipitation [kg/m2/s]
417         CALL iom_put( 'snowpre', sprecip * 86400. )        ! Snow
418         CALL iom_put( 'precip' , tprecip * 86400. )        ! Total precipitation
419      ENDIF
420      !
421      IF(ln_ctl) THEN
422         CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=zqlw , clinfo2=' zqlw  : ')
423         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
424         CALL prt_ctl(tab2d_1=pst  , clinfo1=' blk_oce_core: pst    : ', tab2d_2=emp  , clinfo2=' emp   : ')
425         CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
426            &         tab2d_2=vtau , clinfo2=              ' vtau  : ' , mask2=vmask )
427      ENDIF
428      !
429      CALL wrk_dealloc( jpi,jpj, zwnd_i, zwnd_j, zqsatw, zqlw, zqsb, zqla, zevap )
430      CALL wrk_dealloc( jpi,jpj, Cd, Ch, Ce, zst, zt_zu, zq_zu )
431      !
432      IF( nn_timing == 1 )  CALL timing_stop('blk_oce_core')
433      !
434   END SUBROUTINE blk_oce_core
435 
436   
437#if defined key_lim2 || defined key_lim3
438   SUBROUTINE blk_ice_core_tau
439      !!---------------------------------------------------------------------
440      !!                     ***  ROUTINE blk_ice_core_tau  ***
441      !!
442      !! ** Purpose :   provide the surface boundary condition over sea-ice
443      !!
444      !! ** Method  :   compute momentum using CORE bulk
445      !!                formulea, ice variables and read atmospheric fields.
446      !!                NB: ice drag coefficient is assumed to be a constant
447      !!---------------------------------------------------------------------
448      INTEGER  ::   ji, jj    ! dummy loop indices
449      REAL(wp) ::   zcoef_wnorm, zcoef_wnorm2
450      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f               ! relative wind module and components at F-point
451      REAL(wp) ::             zwndi_t , zwndj_t               ! relative wind components at T-point
452      !!---------------------------------------------------------------------
453      !
454      IF( nn_timing == 1 )  CALL timing_start('blk_ice_core_tau')
455      !
456      ! local scalars ( place there for vector optimisation purposes)
457      zcoef_wnorm  = rhoa * Cice
458      zcoef_wnorm2 = rhoa * Cice * 0.5
459
460!!gm brutal....
461      utau_ice  (:,:) = 0._wp
462      vtau_ice  (:,:) = 0._wp
463      wndm_ice  (:,:) = 0._wp
464!!gm end
465
466      ! ----------------------------------------------------------------------------- !
467      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
468      ! ----------------------------------------------------------------------------- !
469      SELECT CASE( cp_ice_msh )
470      CASE( 'I' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
471         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
472         DO jj = 2, jpjm1
473            DO ji = 2, jpim1   ! B grid : NO vector opt
474               ! ... scalar wind at I-point (fld being at T-point)
475               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ,1) + sf(jp_wndi)%fnow(ji  ,jj  ,1)   &
476                  &              + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji  ,jj-1,1)  ) - rn_vfac * u_ice(ji,jj)
477               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ,1) + sf(jp_wndj)%fnow(ji  ,jj  ,1)   &
478                  &              + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji  ,jj-1,1)  ) - rn_vfac * v_ice(ji,jj)
479               zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
480               ! ... ice stress at I-point
481               utau_ice(ji,jj) = zwnorm_f * zwndi_f
482               vtau_ice(ji,jj) = zwnorm_f * zwndj_f
483               ! ... scalar wind at T-point (fld being at T-point)
484               zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - rn_vfac * 0.25 * (  u_ice(ji,jj+1) + u_ice(ji+1,jj+1)   &
485                  &                                                    + u_ice(ji,jj  ) + u_ice(ji+1,jj  )  )
486               zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 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_vfac * 0.5 * ( u_ice(ji-1,jj  ) + u_ice(ji,jj) )  )
499               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - rn_vfac * 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) ) - rn_vfac * u_ice(ji,jj) )
507               vtau_ice(ji,jj) = zcoef_wnorm2 * ( wndm_ice(ji,jj+1  ) + wndm_ice(ji,jj) )                          &
508                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - rn_vfac * v_ice(ji,jj) )
509            END DO
510         END DO
511         CALL lbc_lnk( utau_ice, 'U', -1. )
512         CALL lbc_lnk( vtau_ice, 'V', -1. )
513         CALL lbc_lnk( wndm_ice, 'T',  1. )
514         !
515      END SELECT
516
517      IF(ln_ctl) THEN
518         CALL prt_ctl(tab2d_1=utau_ice  , clinfo1=' blk_ice_core: utau_ice : ', tab2d_2=vtau_ice  , clinfo2=' vtau_ice : ')
519         CALL prt_ctl(tab2d_1=wndm_ice  , clinfo1=' blk_ice_core: wndm_ice : ')
520      ENDIF
521
522      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_core_tau')
523     
524   END SUBROUTINE blk_ice_core_tau
525
526
527   SUBROUTINE blk_ice_core_flx( ptsu, palb )
528      !!---------------------------------------------------------------------
529      !!                     ***  ROUTINE blk_ice_core_flx  ***
530      !!
531      !! ** Purpose :   provide the surface boundary condition over sea-ice
532      !!
533      !! ** Method  :   compute heat and freshwater exchanged
534      !!                between atmosphere and sea-ice using CORE bulk
535      !!                formulea, ice variables and read atmmospheric fields.
536      !!
537      !! caution : the net upward water flux has with mm/day unit
538      !!---------------------------------------------------------------------
539      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   ptsu          ! sea ice surface temperature
540      REAL(wp), DIMENSION(:,:,:), INTENT(in)  ::   palb          ! ice albedo (all skies)
541      !!
542      INTEGER  ::   ji, jj, jl    ! dummy loop indices
543      REAL(wp) ::   zst2, zst3
544      REAL(wp) ::   zcoef_dqlw, zcoef_dqla, zcoef_dqsb
545      REAL(wp) ::   zztmp, z1_lsub                               ! temporary variable
546      !!
547      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qlw             ! long wave heat flux over ice
548      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qsb             ! sensible  heat flux over ice
549      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqlw            ! long wave heat sensitivity over ice
550      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqsb            ! sensible  heat sensitivity over ice
551      REAL(wp), DIMENSION(:,:)  , POINTER ::   zevap, zsnw       ! evaporation and snw distribution after wind blowing (LIM3)
552      !!---------------------------------------------------------------------
553      !
554      IF( nn_timing == 1 )  CALL timing_start('blk_ice_core_flx')
555      !
556      CALL wrk_alloc( jpi,jpj,jpl, z_qlw, z_qsb, z_dqlw, z_dqsb ) 
557
558      ! local scalars ( place there for vector optimisation purposes)
559      zcoef_dqlw   = 4.0 * 0.95 * Stef
560      zcoef_dqla   = -Ls * Cice * 11637800. * (-5897.8)
561      zcoef_dqsb   = rhoa * cpa * Cice
562
563      zztmp = 1. / ( 1. - albo )
564      !                                     ! ========================== !
565      DO jl = 1, jpl                        !  Loop over ice categories  !
566         !                                  ! ========================== !
567         DO jj = 1 , jpj
568            DO ji = 1, jpi
569               ! ----------------------------!
570               !      I   Radiative FLUXES   !
571               ! ----------------------------!
572               zst2 = ptsu(ji,jj,jl) * ptsu(ji,jj,jl)
573               zst3 = ptsu(ji,jj,jl) * zst2
574               ! Short Wave (sw)
575               qsr_ice(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
576               ! Long  Wave (lw)
577               z_qlw(ji,jj,jl) = 0.95 * ( sf(jp_qlw)%fnow(ji,jj,1) - Stef * ptsu(ji,jj,jl) * zst3 ) * tmask(ji,jj,1)
578               ! lw sensitivity
579               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3                                               
580
581               ! ----------------------------!
582               !     II    Turbulent FLUXES  !
583               ! ----------------------------!
584
585               ! ... turbulent heat fluxes
586               ! Sensible Heat
587               z_qsb(ji,jj,jl) = rhoa * cpa * Cice * wndm_ice(ji,jj) * ( ptsu(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) )
588               ! Latent Heat
589               qla_ice(ji,jj,jl) = rn_efac * MAX( 0.e0, rhoa * Ls  * Cice * wndm_ice(ji,jj)   &                           
590                  &                         * (  11637800. * EXP( -5897.8 / ptsu(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj,1)  ) )
591              ! Latent heat sensitivity for ice (Dqla/Dt)
592               IF( qla_ice(ji,jj,jl) > 0._wp ) THEN
593                  dqla_ice(ji,jj,jl) = rn_efac * zcoef_dqla * wndm_ice(ji,jj) / ( zst2 ) * EXP( -5897.8 / ptsu(ji,jj,jl) )
594               ELSE
595                  dqla_ice(ji,jj,jl) = 0._wp
596               ENDIF
597
598               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
599               z_dqsb(ji,jj,jl) = zcoef_dqsb * wndm_ice(ji,jj)
600
601               ! ----------------------------!
602               !     III    Total FLUXES     !
603               ! ----------------------------!
604               ! Downward Non Solar flux
605               qns_ice (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - qla_ice (ji,jj,jl)
606               ! Total non solar heat flux sensitivity for ice
607               dqns_ice(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + dqla_ice(ji,jj,jl) )
608            END DO
609            !
610         END DO
611         !
612      END DO
613      !
614      tprecip(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
615      sprecip(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
616      CALL iom_put( 'snowpre', sprecip * 86400. )                  ! Snow precipitation
617      CALL iom_put( 'precip' , tprecip * 86400. )                  ! Total precipitation
618
619#if defined  key_lim3
620      CALL wrk_alloc( jpi,jpj, zevap, zsnw ) 
621
622      ! --- evaporation --- !
623      z1_lsub = 1._wp / Lsub
624      evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * z1_lsub    ! sublimation
625      devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * z1_lsub    ! d(sublimation)/dT
626      zevap    (:,:)   = rn_efac * ( emp(:,:) + tprecip(:,:) )  ! evaporation over ocean
627
628      ! --- evaporation minus precipitation --- !
629      zsnw(:,:) = 0._wp
630      CALL lim_thd_snwblow( pfrld, zsnw )  ! snow distribution over ice after wind blowing
631      emp_oce(:,:) = pfrld(:,:) * zevap(:,:) - ( tprecip(:,:) - sprecip(:,:) ) - sprecip(:,:) * (1._wp - zsnw )
632      emp_ice(:,:) = SUM( a_i_b(:,:,:) * evap_ice(:,:,:), dim=3 ) - sprecip(:,:) * zsnw
633      emp_tot(:,:) = emp_oce(:,:) + emp_ice(:,:)
634
635      ! --- heat flux associated with emp --- !
636      qemp_oce(:,:) = - pfrld(:,:) * zevap(:,:) * sst_m(:,:) * rcp                               & ! evap at sst
637         &          + ( tprecip(:,:) - sprecip(:,:) ) * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp  & ! liquid precip at Tair
638         &          +   sprecip(:,:) * ( 1._wp - zsnw ) *                                        & ! solid precip at min(Tair,Tsnow)
639         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
640      qemp_ice(:,:) =   sprecip(:,:) * zsnw *                                                    & ! solid precip (only)
641         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
642
643      ! --- total solar and non solar fluxes --- !
644      qns_tot(:,:) = pfrld(:,:) * qns_oce(:,:) + SUM( a_i_b(:,:,:) * qns_ice(:,:,:), dim=3 ) + qemp_ice(:,:) + qemp_oce(:,:)
645      qsr_tot(:,:) = pfrld(:,:) * qsr_oce(:,:) + SUM( a_i_b(:,:,:) * qsr_ice(:,:,:), dim=3 )
646
647      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- !
648      qprec_ice(:,:) = rhosn * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus )
649
650      ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) --- !
651      DO jl = 1, jpl
652         qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * cpic * tmask(:,:,1) )
653                                   ! But we do not have Tice => consider it at 0°C => evap=0
654      END DO
655
656      CALL wrk_dealloc( jpi,jpj, zevap, zsnw ) 
657#endif
658
659      !--------------------------------------------------------------------
660      ! FRACTIONs of net shortwave radiation which is not absorbed in the
661      ! thin surface layer and penetrates inside the ice cover
662      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
663      !
664      fr1_i0(:,:) = ( 0.18 * ( 1.0 - cldf_ice ) + 0.35 * cldf_ice )
665      fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice )
666      !
667      !
668      IF(ln_ctl) THEN
669         CALL prt_ctl(tab3d_1=qla_ice , clinfo1=' blk_ice_core: qla_ice  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb    : ', kdim=jpl)
670         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice_core: z_qlw    : ', tab3d_2=dqla_ice, clinfo2=' dqla_ice : ', kdim=jpl)
671         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice_core: z_dqsb   : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw   : ', kdim=jpl)
672         CALL prt_ctl(tab3d_1=dqns_ice, clinfo1=' blk_ice_core: dqns_ice : ', tab3d_2=qsr_ice , clinfo2=' qsr_ice  : ', kdim=jpl)
673         CALL prt_ctl(tab3d_1=ptsu    , clinfo1=' blk_ice_core: ptsu     : ', tab3d_2=qns_ice , clinfo2=' qns_ice  : ', kdim=jpl)
674         CALL prt_ctl(tab2d_1=tprecip , clinfo1=' blk_ice_core: tprecip  : ', tab2d_2=sprecip , clinfo2=' sprecip  : ')
675      ENDIF
676
677      CALL wrk_dealloc( jpi,jpj,jpl, z_qlw, z_qsb, z_dqlw, z_dqsb )
678      !
679      IF( nn_timing == 1 )  CALL timing_stop('blk_ice_core_flx')
680     
681   END SUBROUTINE blk_ice_core_flx
682#endif
683
684   SUBROUTINE turb_core_2z( zt, zu, sst, T_zt, q_sat, q_zt, dU,    &
685      &                      Cd, Ch, Ce , T_zu, q_zu )
686      !!----------------------------------------------------------------------
687      !!                      ***  ROUTINE  turb_core  ***
688      !!
689      !! ** Purpose :   Computes turbulent transfert coefficients of surface
690      !!                fluxes according to Large & Yeager (2004) and Large & Yeager (2008)
691      !!                If relevant (zt /= zu), adjust temperature and humidity from height zt to zu
692      !!
693      !! ** Method : Monin Obukhov Similarity Theory
694      !!             + Large & Yeager (2004,2008) closure: CD_n10 = f(U_n10)
695      !!
696      !! ** References :   Large & Yeager, 2004 / Large & Yeager, 2008
697      !!
698      !! ** Last update: Laurent Brodeau, June 2014:
699      !!    - handles both cases zt=zu and zt/=zu
700      !!    - optimized: less 2D arrays allocated and less operations
701      !!    - better first guess of stability by checking air-sea difference of virtual temperature
702      !!       rather than temperature difference only...
703      !!    - added function "cd_neutral_10m" that uses the improved parametrization of
704      !!      Large & Yeager 2008. Drag-coefficient reduction for Cyclone conditions!
705      !!    - using code-wide physical constants defined into "phycst.mod" rather than redifining them
706      !!      => 'vkarmn' and 'grav'
707      !!----------------------------------------------------------------------
708      REAL(wp), INTENT(in   )                     ::   zt       ! height for T_zt and q_zt                   [m]
709      REAL(wp), INTENT(in   )                     ::   zu       ! height for dU                              [m]
710      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   sst      ! sea surface temperature              [Kelvin]
711      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   T_zt     ! potential air temperature            [Kelvin]
712      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_sat    ! sea surface specific humidity         [kg/kg]
713      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   q_zt     ! specific air humidity                 [kg/kg]
714      REAL(wp), INTENT(in   ), DIMENSION(jpi,jpj) ::   dU       ! relative wind module at zu            [m/s]
715      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Cd       ! transfer coefficient for momentum         (tau)
716      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ch       ! transfer coefficient for sensible heat (Q_sens)
717      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   Ce       ! transfert coefficient for evaporation   (Q_lat)
718      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   T_zu     ! air temp. shifted at zu                     [K]
719      REAL(wp), INTENT(  out), DIMENSION(jpi,jpj) ::   q_zu     ! spec. hum.  shifted at zu               [kg/kg]
720      !
721      INTEGER ::   j_itt
722      INTEGER , PARAMETER ::   nb_itt = 5       ! number of itterations
723      LOGICAL ::   l_zt_equal_zu = .FALSE.      ! if q and t are given at different height than U
724      !
725      REAL(wp), DIMENSION(:,:), POINTER ::   U_zu          ! relative wind at zu                            [m/s]
726      REAL(wp), DIMENSION(:,:), POINTER ::   Ce_n10        ! 10m neutral latent coefficient
727      REAL(wp), DIMENSION(:,:), POINTER ::   Ch_n10        ! 10m neutral sensible coefficient
728      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd_n10   ! root square of Cd_n10
729      REAL(wp), DIMENSION(:,:), POINTER ::   sqrt_Cd       ! root square of Cd
730      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_u        ! stability parameter at height zu
731      REAL(wp), DIMENSION(:,:), POINTER ::   zeta_t        ! stability parameter at height zt
732      REAL(wp), DIMENSION(:,:), POINTER ::   zpsi_h_u, zpsi_m_u
733      REAL(wp), DIMENSION(:,:), POINTER ::   ztmp0, ztmp1, ztmp2
734      REAL(wp), DIMENSION(:,:), POINTER ::   stab          ! 1st stability test integer
735      !!----------------------------------------------------------------------
736
737      IF( nn_timing == 1 )  CALL timing_start('turb_core_2z')
738   
739      CALL wrk_alloc( jpi,jpj, U_zu, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd )
740      CALL wrk_alloc( jpi,jpj, zeta_u, stab )
741      CALL wrk_alloc( jpi,jpj, zpsi_h_u, zpsi_m_u, ztmp0, ztmp1, ztmp2 )
742
743      l_zt_equal_zu = .FALSE.
744      IF( ABS(zu - zt) < 0.01 ) l_zt_equal_zu = .TRUE.    ! testing "zu == zt" is risky with double precision
745
746      IF( .NOT. l_zt_equal_zu )   CALL wrk_alloc( jpi,jpj, zeta_t )
747
748      U_zu = MAX( 0.5 , dU )   !  relative wind speed at zu (normally 10m), we don't want to fall under 0.5 m/s
749
750      !! First guess of stability:
751      ztmp0 = T_zt*(1. + 0.608*q_zt) - sst*(1. + 0.608*q_sat) ! air-sea difference of virtual pot. temp. at zt
752      stab  = 0.5 + sign(0.5,ztmp0)                           ! stab = 1 if dTv > 0  => STABLE, 0 if unstable
753
754      !! Neutral coefficients at 10m:
755      IF( ln_cdgw ) THEN      ! wave drag case
756         cdn_wave(:,:) = cdn_wave(:,:) + rsmall * ( 1._wp - tmask(:,:,1) )
757         ztmp0   (:,:) = cdn_wave(:,:)
758      ELSE
759         ztmp0 = cd_neutral_10m( U_zu )
760      ENDIF
761      sqrt_Cd_n10 = SQRT( ztmp0 )
762      Ce_n10  = 1.e-3*( 34.6 * sqrt_Cd_n10 )
763      Ch_n10  = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab))
764   
765      !! Initializing transf. coeff. with their first guess neutral equivalents :
766      Cd = ztmp0   ;   Ce = Ce_n10   ;   Ch = Ch_n10   ;   sqrt_Cd = sqrt_Cd_n10
767
768      !! Initializing values at z_u with z_t values:
769      T_zu = T_zt   ;   q_zu = q_zt
770
771      !!  * Now starting iteration loop
772      DO j_itt=1, nb_itt
773         !
774         ztmp1 = T_zu - sst   ! Updating air/sea differences
775         ztmp2 = q_zu - q_sat 
776
777         ! Updating turbulent scales :   (L&Y 2004 eq. (7))
778         ztmp1  = Ch/sqrt_Cd*ztmp1    ! theta*
779         ztmp2  = Ce/sqrt_Cd*ztmp2    ! q*
780       
781         ztmp0 = T_zu*(1. + 0.608*q_zu) ! virtual potential temperature at zu
782
783         ! Estimate the inverse of Monin-Obukov length (1/L) at height zu:
784         ztmp0 =  (vkarmn*grav/ztmp0*(ztmp1*(1.+0.608*q_zu) + 0.608*T_zu*ztmp2)) / (Cd*U_zu*U_zu) 
785         !                                                                     ( Cd*U_zu*U_zu is U*^2 at zu)
786
787         !! Stability parameters :
788         zeta_u   = zu*ztmp0   ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
789         zpsi_h_u = psi_h( zeta_u )
790         zpsi_m_u = psi_m( zeta_u )
791       
792         !! Shifting temperature and humidity at zu (L&Y 2004 eq. (9b-9c))
793         IF ( .NOT. l_zt_equal_zu ) THEN
794            zeta_t = zt*ztmp0 ;  zeta_t = sign( min(abs(zeta_t),10.0), zeta_t )
795            stab = LOG(zu/zt) - zpsi_h_u + psi_h(zeta_t)  ! stab just used as temp array!!!
796            T_zu = T_zt + ztmp1/vkarmn*stab    ! ztmp1 is still theta*
797            q_zu = q_zt + ztmp2/vkarmn*stab    ! ztmp2 is still q*
798            q_zu = max(0., q_zu)
799         END IF
800       
801         IF( ln_cdgw ) THEN      ! surface wave case
802            sqrt_Cd = vkarmn / ( vkarmn / sqrt_Cd_n10 - zpsi_m_u ) 
803            Cd      = sqrt_Cd * sqrt_Cd
804         ELSE
805           ! Update neutral wind speed at 10m and neutral Cd at 10m (L&Y 2004 eq. 9a)...
806           !   In very rare low-wind conditions, the old way of estimating the
807           !   neutral wind speed at 10m leads to a negative value that causes the code
808           !   to crash. To prevent this a threshold of 0.25m/s is imposed.
809           ztmp0 = MAX( 0.25 , U_zu/(1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - zpsi_m_u)) ) !  U_n10
810           ztmp0 = cd_neutral_10m(ztmp0)                                                 ! Cd_n10
811           sqrt_Cd_n10 = sqrt(ztmp0)
812       
813           Ce_n10  = 1.e-3 * (34.6 * sqrt_Cd_n10)                     ! L&Y 2004 eq. (6b)
814           stab    = 0.5 + sign(0.5,zeta_u)                           ! update stability
815           Ch_n10  = 1.e-3*sqrt_Cd_n10*(18.*stab + 32.7*(1. - stab))  ! L&Y 2004 eq. (6c-6d)
816
817           !! Update of transfer coefficients:
818           ztmp1 = 1. + sqrt_Cd_n10/vkarmn*(LOG(zu/10.) - zpsi_m_u)   ! L&Y 2004 eq. (10a)
819           Cd      = ztmp0 / ( ztmp1*ztmp1 )   
820           sqrt_Cd = SQRT( Cd )
821         ENDIF
822         !
823         ztmp0 = (LOG(zu/10.) - zpsi_h_u) / vkarmn / sqrt_Cd_n10
824         ztmp2 = sqrt_Cd / sqrt_Cd_n10
825         ztmp1 = 1. + Ch_n10*ztmp0               
826         Ch  = Ch_n10*ztmp2 / ztmp1  ! L&Y 2004 eq. (10b)
827         !
828         ztmp1 = 1. + Ce_n10*ztmp0               
829         Ce  = Ce_n10*ztmp2 / ztmp1  ! L&Y 2004 eq. (10c)
830         !
831      END DO
832
833      CALL wrk_dealloc( jpi,jpj, U_zu, Ce_n10, Ch_n10, sqrt_Cd_n10, sqrt_Cd )
834      CALL wrk_dealloc( jpi,jpj, zeta_u, stab )
835      CALL wrk_dealloc( jpi,jpj, zpsi_h_u, zpsi_m_u, ztmp0, ztmp1, ztmp2 )
836
837      IF( .NOT. l_zt_equal_zu ) CALL wrk_dealloc( jpi,jpj, zeta_t )
838
839      IF( nn_timing == 1 )  CALL timing_stop('turb_core_2z')
840      !
841   END SUBROUTINE turb_core_2z
842
843   SUBROUTINE core_namelist(cd_dir, sd_wndi, sd_wndj, sd_humi, sd_qsr, sd_qlw, &
844         &                  sd_tair, sd_prec, sd_snow, sd_tdif)
845     !!---------------------------------------------------------------------
846     !!                   ***  ROUTINE core_namelist  ***
847     !!                     
848     !! ** Purpose :   Broadcast namelist variables read by procesor lwm
849     !!
850     !! ** Method  :   use lib_mpp
851     !!----------------------------------------------------------------------
852     CHARACTER(len=100), INTENT(INOUT)  :: cd_dir   !   Root directory for location of core files
853     TYPE(FLD_N), INTENT(INOUT)         :: sd_wndi, sd_wndj, sd_humi, sd_qsr       ! informations about the fields to be read
854     TYPE(FLD_N), INTENT(INOUT)         :: sd_qlw , sd_tair, sd_prec, sd_snow      !
855     TYPE(FLD_N), INTENT(INOUT)         :: sd_tdif
856#if defined key_mpp_mpi
857      CALL mpp_bcast(cd_dir, 100)
858      CALL mpp_bcast(ln_taudif)
859      CALL mpp_bcast(rn_pfac)
860      CALL mpp_bcast(rn_efac)
861      CALL mpp_bcast(rn_vfac)
862      CALL fld_n_bcast(sd_wndi)
863      CALL fld_n_bcast(sd_wndj)
864      CALL fld_n_bcast(sd_humi)
865      CALL fld_n_bcast(sd_qsr)
866      CALL fld_n_bcast(sd_qlw)
867      CALL fld_n_bcast(sd_tair)
868      CALL fld_n_bcast(sd_prec)
869      CALL fld_n_bcast(sd_snow)
870      CALL fld_n_bcast(sd_tdif)
871      CALL mpp_bcast(rn_zqt)
872      CALL mpp_bcast(rn_zu)
873      CALL mpp_bcast(ln_core_sio)
874#endif
875   END SUBROUTINE core_namelist
876
877   FUNCTION cd_neutral_10m( zw10 )
878      !!----------------------------------------------------------------------
879      !! Estimate of the neutral drag coefficient at 10m as a function
880      !! of neutral wind  speed at 10m
881      !!
882      !! Origin: Large & Yeager 2008 eq.(11a) and eq.(11b)
883      !!
884      !! Author: L. Brodeau, june 2014
885      !!----------------------------------------------------------------------   
886      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   zw10           ! scalar wind speed at 10m (m/s)
887      REAL(wp), DIMENSION(jpi,jpj)             ::   cd_neutral_10m
888      !
889      REAL(wp), DIMENSION(:,:), POINTER ::   rgt33
890      !!----------------------------------------------------------------------   
891      !
892      CALL wrk_alloc( jpi,jpj, rgt33 )
893      !
894      !! When wind speed > 33 m/s => Cyclone conditions => special treatment
895      rgt33 = 0.5_wp + SIGN( 0.5_wp, (zw10 - 33._wp) )   ! If zw10 < 33. => 0, else => 1 
896      cd_neutral_10m = 1.e-3 * ( &
897         &       (1._wp - rgt33)*( 2.7_wp/zw10 + 0.142_wp + zw10/13.09_wp - 3.14807E-10*zw10**6) & ! zw10< 33.
898         &      + rgt33         *      2.34   )                                                    ! zw10 >= 33.
899      !
900      CALL wrk_dealloc( jpi,jpj, rgt33)
901      !
902   END FUNCTION cd_neutral_10m
903
904
905   FUNCTION psi_m(pta)   !! Psis, L&Y 2004 eq. (8c), (8d), (8e)
906      !-------------------------------------------------------------------------------
907      ! universal profile stability function for momentum
908      !-------------------------------------------------------------------------------
909      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pta
910      !
911      REAL(wp), DIMENSION(jpi,jpj)             :: psi_m
912      REAL(wp), DIMENSION(:,:), POINTER        :: X2, X, stabit
913      !-------------------------------------------------------------------------------
914      !
915      CALL wrk_alloc( jpi,jpj, X2, X, stabit )
916      !
917      X2 = SQRT( ABS( 1. - 16.*pta ) )  ;  X2 = MAX( X2 , 1. )   ;   X = SQRT( X2 )
918      stabit = 0.5 + SIGN( 0.5 , pta )
919      psi_m = -5.*pta*stabit  &                                                          ! Stable
920         &    + (1. - stabit)*(2.*LOG((1. + X)*0.5) + LOG((1. + X2)*0.5) - 2.*ATAN(X) + rpi*0.5)  ! Unstable
921      !
922      CALL wrk_dealloc( jpi,jpj, X2, X, stabit )
923      !
924   END FUNCTION psi_m
925
926
927   FUNCTION psi_h( pta )    !! Psis, L&Y 2004 eq. (8c), (8d), (8e)
928      !-------------------------------------------------------------------------------
929      ! universal profile stability function for temperature and humidity
930      !-------------------------------------------------------------------------------
931      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pta
932      !
933      REAL(wp), DIMENSION(jpi,jpj)             ::   psi_h
934      REAL(wp), DIMENSION(:,:), POINTER        ::   X2, X, stabit
935      !-------------------------------------------------------------------------------
936      !
937      CALL wrk_alloc( jpi,jpj, X2, X, stabit )
938      !
939      X2 = SQRT( ABS( 1. - 16.*pta ) )   ;   X2 = MAX( X2 , 1. )   ;   X = SQRT( X2 )
940      stabit = 0.5 + SIGN( 0.5 , pta )
941      psi_h = -5.*pta*stabit   &                                       ! Stable
942         &    + (1. - stabit)*(2.*LOG( (1. + X2)*0.5 ))                ! Unstable
943      !
944      CALL wrk_dealloc( jpi,jpj, X2, X, stabit )
945      !
946   END FUNCTION psi_h
947
948   !!======================================================================
949END MODULE sbcblk_core
Note: See TracBrowser for help on using the repository browser.