New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
sbcblk_core.F90 in branches/UKMO/test_moci_test_suite_namelist_read/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

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

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

#2050 fixes and changes

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