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

source: branches/UKMO/smagorinsky_diagnostics/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 5699

Last change on this file since 5699 was 5699, checked in by davestorkey, 9 years ago

Remove SVN keyword updating from UKMO/smagorinsky_diagnostics branch.

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