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

source: branches/UKMO/r6232_HZG_WAVE/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 7807

Last change on this file since 7807 was 7807, checked in by jcastill, 7 years ago

Changes as in HZG wave forcing branch, but adapted to r6232

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