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

source: branches/UKMO/test_moci_test_suite/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 8243

Last change on this file since 8243 was 8243, checked in by andmirek, 7 years ago

#1914 working XIOS read, XIOS write and single processor read

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