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

source: trunk/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 2748

Last change on this file since 2748 was 2748, checked in by cbricaud, 13 years ago

cbricaud: bugfix for sbcblk_core.F90: see ticket 814

  • Property svn:keywords set to Id
File size: 52.4 KB
Line 
1MODULE sbcblk_core
2   !!======================================================================
3   !!                       ***  MODULE  sbcblk_core  ***
4   !! Ocean forcing:  momentum, heat and freshwater flux formulation
5   !!=====================================================================
6   !! History :  1.0  !  2004-08  (U. Schweckendiek)  Original code
7   !!            2.0  !  2005-04  (L. Brodeau, A.M. Treguier) additions:
8   !!                           -  new bulk routine for efficiency
9   !!                           -  WINDS ARE NOW ASSUMED TO BE AT T POINTS in input files !!!!
10   !!                           -  file names and file characteristics in namelist
11   !!                           -  Implement reading of 6-hourly fields   
12   !!            3.0  !  2006-06  (G. Madec) sbc rewritting   
13   !!             -   !  2006-12  (L. Brodeau) Original code for TURB_CORE_2Z
14   !!            3.2  !  2009-04  (B. Lemaire)  Introduce iom_put
15   !!            3.3  !  2010-10  (S. Masson)  add diurnal cycle
16   !!----------------------------------------------------------------------
17
18   !!----------------------------------------------------------------------
19   !!   sbc_blk_core  : bulk formulation as ocean surface boundary condition
20   !!                   (forced mode, CORE bulk formulea)
21   !!   blk_oce_core  : ocean: computes momentum, heat and freshwater fluxes
22   !!   blk_ice_core  : ice  : computes momentum, heat and freshwater fluxes
23   !!   turb_core     : computes the CORE turbulent transfer coefficients
24   !!----------------------------------------------------------------------
25   USE oce             ! ocean dynamics and tracers
26   USE dom_oce         ! ocean space and time domain
27   USE phycst          ! physical constants
28   USE fldread         ! read input fields
29   USE sbc_oce         ! Surface boundary condition: ocean fields
30   USE sbcdcy          ! surface boundary condition: diurnal cycle
31   USE iom             ! I/O manager library
32   USE in_out_manager  ! I/O manager
33   USE lib_mpp         ! distribued memory computing library
34   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
35   USE prtctl          ! Print control
36#if defined key_lim3
37   USE sbc_ice         ! Surface boundary condition: ice fields
38#endif
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   sbc_blk_core         ! routine called in sbcmod module
44   PUBLIC   blk_ice_core         ! routine called in sbc_ice_lim module
45
46   INTEGER , PARAMETER ::   jpfld   = 9           ! maximum number of files to read
47   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point
48   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point
49   INTEGER , PARAMETER ::   jp_humi = 3           ! index of specific humidity               ( - )
50   INTEGER , PARAMETER ::   jp_qsr  = 4           ! index of solar heat                      (W/m2)
51   INTEGER , PARAMETER ::   jp_qlw  = 5           ! index of Long wave                       (W/m2)
52   INTEGER , PARAMETER ::   jp_tair = 6           ! index of 10m air temperature             (Kelvin)
53   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s)
54   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s)
55   INTEGER , PARAMETER ::   jp_tdif = 9           ! index of tau diff associated to HF tau   (N/m2)   at T-point
56   
57   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read)
58         
59   !                                             !!! CORE bulk parameters
60   REAL(wp), PARAMETER ::   rhoa =    1.22        ! air density
61   REAL(wp), PARAMETER ::   cpa  = 1000.5         ! specific heat of air
62   REAL(wp), PARAMETER ::   Lv   =    2.5e6       ! latent heat of vaporization
63   REAL(wp), PARAMETER ::   Ls   =    2.839e6     ! latent heat of sublimation
64   REAL(wp), PARAMETER ::   Stef =    5.67e-8     ! Stefan Boltzmann constant
65   REAL(wp), PARAMETER ::   Cice =    1.63e-3     ! transfer coefficient over ice
66   REAL(wp), PARAMETER ::   albo =    0.066       ! ocean albedo assumed to be contant
67
68   !                                  !!* Namelist namsbc_core : CORE bulk parameters
69   LOGICAL  ::   ln_2m     = .FALSE.   ! logical flag for height of air temp. and hum
70   LOGICAL  ::   ln_taudif = .FALSE.   ! logical flag to use the "mean of stress module - module of mean stress" data
71   REAL(wp) ::   rn_pfac   = 1.        ! multiplication factor for precipitation
72
73   !! * Substitutions
74#  include "domzgr_substitute.h90"
75#  include "vectopt_loop_substitute.h90"
76   !!----------------------------------------------------------------------
77   !! NEMO/OPA 3.3 , NEMO-consortium (2010)
78   !! $Id$
79   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
80   !!----------------------------------------------------------------------
81CONTAINS
82
83   SUBROUTINE sbc_blk_core( kt )
84      !!---------------------------------------------------------------------
85      !!                    ***  ROUTINE sbc_blk_core  ***
86      !!                   
87      !! ** Purpose :   provide at each time step the surface ocean fluxes
88      !!      (momentum, heat, freshwater and runoff)
89      !!
90      !! ** Method  : (1) READ each fluxes in NetCDF files:
91      !!      the 10m wind velocity (i-component) (m/s)    at T-point
92      !!      the 10m wind velocity (j-component) (m/s)    at T-point
93      !!      the specific humidity               ( - )
94      !!      the solar heat                      (W/m2)
95      !!      the Long wave                       (W/m2)
96      !!      the 10m air temperature             (Kelvin)
97      !!      the total precipitation (rain+snow) (Kg/m2/s)
98      !!      the snow (solid prcipitation)       (kg/m2/s)
99      !!   OPTIONAL parameter (see ln_taudif namelist flag):
100      !!      the tau diff associated to HF tau   (N/m2)   at T-point
101      !!              (2) CALL blk_oce_core
102      !!
103      !!      C A U T I O N : never mask the surface stress fields
104      !!                      the stress is assumed to be in the mesh referential
105      !!                      i.e. the (i,j) referential
106      !!
107      !! ** Action  :   defined at each time-step at the air-sea interface
108      !!              - utau, vtau  i- and j-component of the wind stress
109      !!              - taum        wind stress module at T-point
110      !!              - wndm        10m wind module at T-point
111      !!              - qns, qsr    non-slor and solar heat flux
112      !!              - emp, emps   evaporation minus precipitation
113      !!----------------------------------------------------------------------
114      INTEGER, INTENT(in) ::   kt   ! ocean time step
115      !!
116      INTEGER  ::   ierror   ! return error code
117      INTEGER  ::   ifpr     ! dummy loop indice
118      INTEGER  ::   jfld     ! dummy loop arguments
119      !!
120      CHARACTER(len=100) ::  cn_dir   !   Root directory for location of core files
121      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i     ! array of namelist informations on the fields to read
122      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read
123      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !   "                                 "
124      TYPE(FLD_N) ::   sn_tdif                                 !   "                                 "
125      NAMELIST/namsbc_core/ cn_dir , ln_2m  , ln_taudif, rn_pfac,           &
126         &                  sn_wndi, sn_wndj, sn_humi  , sn_qsr ,           &
127         &                  sn_qlw , sn_tair, sn_prec  , sn_snow, sn_tdif
128      !!---------------------------------------------------------------------
129
130      !                                         ! ====================== !
131      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
132         !                                      ! ====================== !
133         ! set file information (default values)
134         cn_dir = './'       ! directory in which the model is executed
135         !
136         ! (NB: frequency positive => hours, negative => months)
137         !            !    file    ! frequency ! variable ! time intep !  clim   ! 'yearly' or ! weights  ! rotation !
138         !            !    name    !  (hours)  !  name    !   (T/F)    !  (T/F)  !  'monthly'  ! filename ! pairs    !
139         sn_wndi = FLD_N( 'uwnd10m',    24     , 'u_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
140         sn_wndj = FLD_N( 'vwnd10m',    24     , 'v_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
141         sn_qsr  = FLD_N( 'qsw'    ,    24     , 'qsw'    ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
142         sn_qlw  = FLD_N( 'qlw'    ,    24     , 'qlw'    ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
143         sn_tair = FLD_N( 'tair10m',    24     , 't_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
144         sn_humi = FLD_N( 'humi10m',    24     , 'q_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
145         sn_prec = FLD_N( 'precip' ,    -1     , 'precip' ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
146         sn_snow = FLD_N( 'snow'   ,    -1     , 'snow'   ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
147         sn_tdif = FLD_N( 'taudif' ,    24     , 'taudif' ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
148         !
149         REWIND( numnam )                          ! read in namlist namsbc_core
150         READ  ( numnam, namsbc_core )
151         !                                         ! check: do we plan to use ln_dm2dc with non-daily forcing?
152         IF( ln_dm2dc .AND. sn_qsr%nfreqh /= 24 )   & 
153            &   CALL ctl_stop( 'sbc_blk_core: ln_dm2dc can be activated only with daily short-wave forcing' ) 
154         IF( ln_dm2dc .AND. sn_qsr%ln_tint ) THEN
155            CALL ctl_warn( 'sbc_blk_core: ln_dm2dc is taking care of the temporal interpolation of daily qsr',   &
156                 &         '              ==> We force time interpolation = .false. for qsr' )
157            sn_qsr%ln_tint = .false.
158         ENDIF
159         !                                         ! store namelist information in an array
160         slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
161         slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
162         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
163         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
164         slf_i(jp_tdif) = sn_tdif
165         !                 
166         lhftau = ln_taudif                        ! do we use HF tau information?
167         jfld = jpfld - COUNT( (/.NOT. lhftau/) )
168         !
169         ALLOCATE( sf(jfld), STAT=ierror )         ! set sf structure
170         IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_core: unable to allocate sf structure' )
171         DO ifpr= 1, jfld
172            ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
173            IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
174         END DO
175         !                                         ! fill sf with slf_i and control print
176         CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_core', 'flux formulation for ocean surface boundary condition', 'namsbc_core' )
177         !
178      ENDIF
179
180      CALL fld_read( kt, nn_fsbc, sf )        ! input fields provided at the current time-step
181
182#if defined key_lim3
183      tatm_ice(:,:) = sf(jp_tair)%fnow(:,:,1)                  ! LIM3: make Tair available in sea-ice
184#endif
185      !                                                        ! surface ocean fluxes computed with CLIO bulk formulea
186      IF( MOD( kt - 1, nn_fsbc ) == 0 )   CALL blk_oce_core( sf, sst_m, ssu_m, ssv_m )
187      !
188   END SUBROUTINE sbc_blk_core
189   
190   
191   SUBROUTINE blk_oce_core( sf, pst, pu, pv )
192      !!---------------------------------------------------------------------
193      !!                     ***  ROUTINE blk_core  ***
194      !!
195      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
196      !!      the ocean surface at each time step
197      !!
198      !! ** Method  :   CORE bulk formulea for the ocean using atmospheric
199      !!      fields read in sbc_read
200      !!
201      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
202      !!              - vtau    : j-component of the stress at V-point  (N/m2)
203      !!              - taum    : Wind stress module at T-point         (N/m2)
204      !!              - wndm    : Wind speed module at T-point          (m/s)
205      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
206      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
207      !!              - evap    : Evaporation over the ocean            (kg/m2/s)
208      !!              - emp(s)  : evaporation minus precipitation       (kg/m2/s)
209      !!
210      !!  ** Nota  :   sf has to be a dummy argument for AGRIF on NEC
211      !!---------------------------------------------------------------------
212      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
213      USE wrk_nemo, ONLY:   zwnd_i => wrk_2d_1  , zwnd_j => wrk_2d_2      ! wind speed components at T-point
214      USE wrk_nemo, ONLY:   zqsatw => wrk_2d_3           ! specific humidity at pst
215      USE wrk_nemo, ONLY:   zqlw   => wrk_2d_4  , zqsb   => wrk_2d_5      ! long wave and sensible heat fluxes
216      USE wrk_nemo, ONLY:   zqla   => wrk_2d_6  , zevap  => wrk_2d_7      ! latent heat fluxes and evaporation
217      USE wrk_nemo, ONLY:   Cd     => wrk_2d_8           ! transfer coefficient for momentum      (tau)
218      USE wrk_nemo, ONLY:   Ch     => wrk_2d_9           ! transfer coefficient for sensible heat (Q_sens)
219      USE wrk_nemo, ONLY:   Ce     => wrk_2d_10          ! transfer coefficient for evaporation   (Q_lat)
220      USE wrk_nemo, ONLY:   zst    => wrk_2d_11          ! surface temperature in Kelvin
221      USE wrk_nemo, ONLY:   zt_zu  => wrk_2d_12          ! air temperature at wind speed height
222      USE wrk_nemo, ONLY:   zq_zu  => wrk_2d_13          ! air spec. hum.  at wind speed height
223      !
224      TYPE(fld), INTENT(in), DIMENSION(:)   ::   sf    ! input data
225      REAL(wp) , INTENT(in), DIMENSION(:,:) ::   pst   ! surface temperature                      [Celcius]
226      REAL(wp) , INTENT(in), DIMENSION(:,:) ::   pu    ! surface current at U-point (i-component) [m/s]
227      REAL(wp) , INTENT(in), DIMENSION(:,:) ::   pv    ! surface current at V-point (j-component) [m/s]
228      !
229      INTEGER  ::   ji, jj               ! dummy loop indices
230      REAL(wp) ::   zcoef_qsatw, zztmp   ! local variable
231      !!---------------------------------------------------------------------
232
233      IF( wrk_in_use(2, 1,2,3,4,5,6,7,8,9,10,11,12,13) ) THEN
234         CALL ctl_stop('blk_oce_core: requested workspace arrays unavailable')   ;   RETURN
235      ENDIF
236      !
237      ! local scalars ( place there for vector optimisation purposes)
238      zcoef_qsatw = 0.98 * 640380. / rhoa
239     
240      zst(:,:) = pst(:,:) + rt0      ! converte Celcius to Kelvin (and set minimum value far above 0 K)
241
242      ! ----------------------------------------------------------------------------- !
243      !      0   Wind components and module at T-point relative to the moving ocean   !
244      ! ----------------------------------------------------------------------------- !
245
246      ! ... components ( U10m - U_oce ) at T-point (unmasked)
247      zwnd_i(:,:) = 0.e0 
248      zwnd_j(:,:) = 0.e0
249#if defined key_vectopt_loop
250!CDIR COLLAPSE
251#endif
252      DO jj = 2, jpjm1
253         DO ji = fs_2, fs_jpim1   ! vect. opt.
254            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj,1) - 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
255            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
256         END DO
257      END DO
258      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
259      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
260      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
261!CDIR NOVERRCHK
262!CDIR COLLAPSE
263      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
264         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
265
266      ! ----------------------------------------------------------------------------- !
267      !      I   Radiative FLUXES                                                     !
268      ! ----------------------------------------------------------------------------- !
269   
270      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
271      zztmp = 1. - albo
272      IF( ln_dm2dc ) THEN   ;   qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
273      ELSE                  ;   qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
274      ENDIF
275!CDIR COLLAPSE
276      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
277      ! ----------------------------------------------------------------------------- !
278      !     II    Turbulent FLUXES                                                    !
279      ! ----------------------------------------------------------------------------- !
280
281      ! ... specific humidity at SST and IST
282!CDIR NOVERRCHK
283!CDIR COLLAPSE
284      zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) ) 
285
286      ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point :
287      IF( ln_2m ) THEN
288         !! If air temp. and spec. hum. are given at different height (2m) than wind (10m) :
289         CALL TURB_CORE_2Z(2.,10., zst   , sf(jp_tair)%fnow,         &
290            &                      zqsatw, sf(jp_humi)%fnow, wndm,   &
291            &                      Cd    , Ch              , Ce  ,   &
292            &                      zt_zu , zq_zu                   )
293      ELSE
294         !! If air temp. and spec. hum. are given at same height than wind (10m) :
295!gm bug?  at the compiling phase, add a copy in temporary arrays...  ==> check perf
296!         CALL TURB_CORE_1Z( 10., zst   (:,:), sf(jp_tair)%fnow(:,:),              &
297!            &                    zqsatw(:,:), sf(jp_humi)%fnow(:,:), wndm(:,:),   &
298!            &                    Cd    (:,:),             Ch  (:,:), Ce  (:,:)  )
299!gm bug
300! ARPDBG - this won't compile with gfortran. Fix but check performance
301! as per comment above.
302         CALL TURB_CORE_1Z( 10., zst   , sf(jp_tair)%fnow(:,:,1),       &
303            &                    zqsatw, sf(jp_humi)%fnow(:,:,1), wndm, &
304            &                    Cd    , Ch              , Ce    )
305      ENDIF
306
307      ! ... tau module, i and j component
308      DO jj = 1, jpj
309         DO ji = 1, jpi
310            zztmp = rhoa * wndm(ji,jj) * Cd(ji,jj)
311            taum  (ji,jj) = zztmp * wndm  (ji,jj)
312            zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj)
313            zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj)
314         END DO
315      END DO
316
317      ! ... add the HF tau contribution to the wind stress module?
318      IF( lhftau ) THEN 
319!CDIR COLLAPSE
320         taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1)
321      ENDIF
322      CALL iom_put( "taum_oce", taum )   ! output wind stress module
323
324      ! ... utau, vtau at U- and V_points, resp.
325      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
326      DO jj = 1, jpjm1
327         DO ji = 1, fs_jpim1
328            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) )
329            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) )
330         END DO
331      END DO
332      CALL lbc_lnk( utau(:,:), 'U', -1. )
333      CALL lbc_lnk( vtau(:,:), 'V', -1. )
334
335      !  Turbulent fluxes over ocean
336      ! -----------------------------
337      IF( ln_2m ) THEN
338         ! Values of temp. and hum. adjusted to 10m must be used instead of 2m values
339         zevap(:,:) = MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) ) * wndm(:,:) )   ! Evaporation
340         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - zt_zu(:,:) ) * wndm(:,:)     ! Sensible Heat
341      ELSE
342!CDIR COLLAPSE
343         zevap(:,:) = MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:,1) ) * wndm(:,:) )   ! Evaporation
344!CDIR COLLAPSE
345         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - sf(jp_tair)%fnow(:,:,1) ) * wndm(:,:)     ! Sensible Heat
346      ENDIF
347!CDIR COLLAPSE
348      zqla (:,:) = Lv * zevap(:,:)                                                              ! Latent Heat
349
350      IF(ln_ctl) THEN
351         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
352         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
353         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce_core: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
354         CALL prt_ctl( tab2d_1=zqsatw, clinfo1=' blk_oce_core: zqsatw : ', tab2d_2=zst, clinfo2=' zst : ' )
355         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
356            &          tab2d_2=vtau  , clinfo2=              ' vtau : '  , mask2=vmask )
357         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce_core: wndm   : ')
358         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce_core: zst    : ')
359      ENDIF
360       
361      ! ----------------------------------------------------------------------------- !
362      !     III    Total FLUXES                                                       !
363      ! ----------------------------------------------------------------------------- !
364     
365!CDIR COLLAPSE
366      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)      ! Downward Non Solar flux
367!CDIR COLLAPSE
368      emp(:,:) = zevap(:,:) - sf(jp_prec)%fnow(:,:,1) * rn_pfac * tmask(:,:,1)
369!CDIR COLLAPSE
370      emps(:,:) = emp(:,:)
371      !
372      CALL iom_put( "qlw_oce",   zqlw )                 ! output downward longwave heat over the ocean
373      CALL iom_put( "qsb_oce", - zqsb )                 ! output downward sensible heat over the ocean
374      CALL iom_put( "qla_oce", - zqla )                 ! output downward latent   heat over the ocean
375      CALL iom_put( "qns_oce",   qns  )                 ! output downward non solar heat over the ocean
376      !
377      IF(ln_ctl) THEN
378         CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=zqlw , clinfo2=' zqlw  : ')
379         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
380         CALL prt_ctl(tab2d_1=pst  , clinfo1=' blk_oce_core: pst    : ', tab2d_2=emp  , clinfo2=' emp   : ')
381         CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
382            &         tab2d_2=vtau , clinfo2=              ' vtau  : ' , mask2=vmask )
383      ENDIF
384      !
385      IF( wrk_not_released(2, 1,2,3,4,5,6,7,8,9,10,11,12,13) )   &
386          CALL ctl_stop('blk_oce_core: failed to release workspace arrays')
387      !
388   END SUBROUTINE blk_oce_core
389   
390   
391   SUBROUTINE blk_ice_core(  pst   , pui   , pvi   , palb ,   &
392      &                      p_taui, p_tauj, p_qns , p_qsr,   &
393      &                      p_qla , p_dqns, p_dqla,          &
394      &                      p_tpr , p_spr ,                  &
395      &                      p_fr1 , p_fr2 , cd_grid, pdim  ) 
396      !!---------------------------------------------------------------------
397      !!                     ***  ROUTINE blk_ice_core  ***
398      !!
399      !! ** Purpose :   provide the surface boundary condition over sea-ice
400      !!
401      !! ** Method  :   compute momentum, heat and freshwater exchanged
402      !!                between atmosphere and sea-ice using CORE bulk
403      !!                formulea, ice variables and read atmmospheric fields.
404      !!                NB: ice drag coefficient is assumed to be a constant
405      !!
406      !! caution : the net upward water flux has with mm/day unit
407      !!---------------------------------------------------------------------
408      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
409      USE wrk_nemo, ONLY:   z_wnds_t => wrk_2d_1                ! wind speed ( = | U10m - U_ice | ) at T-point
410      USE wrk_nemo, ONLY:   wrk_3d_4 , wrk_3d_5 , wrk_3d_6 , wrk_3d_7
411      !!
412      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pst      ! ice surface temperature (>0, =rt0 over land) [Kelvin]
413      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pui      ! ice surface velocity (i- and i- components      [m/s]
414      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pvi      !    at I-point (B-grid) or U & V-point (C-grid)
415      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   palb     ! ice albedo (clear sky) (alb_ice_cs)               [%]
416      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_taui   ! i- & j-components of surface ice stress        [N/m2]
417      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_tauj   !   at I-point (B-grid) or U & V-point (C-grid)
418      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_qns    ! non solar heat flux over ice (T-point)         [W/m2]
419      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_qsr    !     solar heat flux over ice (T-point)         [W/m2]
420      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_qla    ! latent    heat flux over ice (T-point)         [W/m2]
421      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_dqns   ! non solar heat sensistivity  (T-point)         [W/m2]
422      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_dqla   ! latent    heat sensistivity  (T-point)         [W/m2]
423      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_tpr    ! total precipitation          (T-point)      [Kg/m2/s]
424      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_spr    ! solid precipitation          (T-point)      [Kg/m2/s]
425      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_fr1    ! 1sr fraction of qsr penetration in ice (T-point)  [%]
426      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_fr2    ! 2nd fraction of qsr penetration in ice (T-point)  [%]
427      CHARACTER(len=1)          , INTENT(in   ) ::   cd_grid  ! ice grid ( C or B-grid)
428      INTEGER                   , INTENT(in   ) ::   pdim     ! number of ice categories
429      !!
430      INTEGER  ::   ji, jj, jl    ! dummy loop indices
431      INTEGER  ::   ijpl          ! number of ice categories (size of 3rd dim of input arrays)
432      REAL(wp) ::   zst2, zst3
433      REAL(wp) ::   zcoef_wnorm, zcoef_wnorm2, zcoef_dqlw, zcoef_dqla, zcoef_dqsb
434      REAL(wp) ::   zztmp                                        ! temporary variable
435      REAL(wp) ::   zcoef_frca                                   ! fractional cloud amount
436      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f                  ! relative wind module and components at F-point
437      REAL(wp) ::             zwndi_t , zwndj_t                  ! relative wind components at T-point
438      !!
439      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qlw             ! long wave heat flux over ice
440      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qsb             ! sensible  heat flux over ice
441      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqlw            ! long wave heat sensitivity over ice
442      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqsb            ! sensible  heat sensitivity over ice
443      !!---------------------------------------------------------------------
444
445      ijpl  = pdim                            ! number of ice categories
446
447      ! Set-up access to workspace arrays
448      IF( wrk_in_use(2, 1) .OR. wrk_in_use(3, 4,5,6,7) ) THEN
449         CALL ctl_stop('blk_ice_core: requested workspace arrays unavailable')   ;   RETURN
450      ELSE IF(ijpl > jpk) THEN
451         CALL ctl_stop('blk_ice_core: no. of ice categories > jpk so wrk_nemo 3D workspaces cannot be used.')
452         RETURN
453      END IF
454      ! Set-up pointers to sub-arrays of workspaces
455      z_qlw  => wrk_3d_4(:,:,1:ijpl)
456      z_qsb  => wrk_3d_5(:,:,1:ijpl)
457      z_dqlw => wrk_3d_6(:,:,1:ijpl)
458      z_dqsb => wrk_3d_7(:,:,1:ijpl)
459
460      ! local scalars ( place there for vector optimisation purposes)
461      zcoef_wnorm  = rhoa * Cice
462      zcoef_wnorm2 = rhoa * Cice * 0.5
463      zcoef_dqlw   = 4.0 * 0.95 * Stef
464      zcoef_dqla   = -Ls * Cice * 11637800. * (-5897.8)
465      zcoef_dqsb   = rhoa * cpa * Cice
466      zcoef_frca   = 1.0  - 0.3
467
468!!gm brutal....
469      z_wnds_t(:,:) = 0.e0
470      p_taui  (:,:) = 0.e0
471      p_tauj  (:,:) = 0.e0
472!!gm end
473
474      ! ----------------------------------------------------------------------------- !
475      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
476      ! ----------------------------------------------------------------------------- !
477      SELECT CASE( cd_grid )
478      CASE( 'I' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
479         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
480!CDIR NOVERRCHK
481         DO jj = 2, jpjm1
482            DO ji = 2, jpim1   ! B grid : NO vector opt
483               ! ... scalar wind at I-point (fld being at T-point)
484               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ,1) + sf(jp_wndi)%fnow(ji  ,jj  ,1)   &
485                  &              + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji  ,jj-1,1)  ) - pui(ji,jj)
486               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ,1) + sf(jp_wndj)%fnow(ji  ,jj  ,1)   &
487                  &              + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji  ,jj-1,1)  ) - pvi(ji,jj)
488               zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
489               ! ... ice stress at I-point
490               p_taui(ji,jj) = zwnorm_f * zwndi_f
491               p_tauj(ji,jj) = zwnorm_f * zwndj_f
492               ! ... scalar wind at T-point (fld being at T-point)
493               zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - 0.25 * (  pui(ji,jj+1) + pui(ji+1,jj+1)   &
494                  &                                          + pui(ji,jj  ) + pui(ji+1,jj  )  )
495               zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - 0.25 * (  pvi(ji,jj+1) + pvi(ji+1,jj+1)   &
496                  &                                          + pvi(ji,jj  ) + pvi(ji+1,jj  )  )
497               z_wnds_t(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
498            END DO
499         END DO
500         CALL lbc_lnk( p_taui  , 'I', -1. )
501         CALL lbc_lnk( p_tauj  , 'I', -1. )
502         CALL lbc_lnk( z_wnds_t, 'T',  1. )
503         !
504      CASE( 'C' )                  ! C-grid ice dynamics :   U & V-points (same as ocean)
505#if defined key_vectopt_loop
506!CDIR COLLAPSE
507#endif
508         DO jj = 2, jpj
509            DO ji = fs_2, jpi   ! vect. opt.
510               zwndi_t = (  sf(jp_wndi)%fnow(ji,jj,1) - 0.5 * ( pui(ji-1,jj  ) + pui(ji,jj) )  )
511               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - 0.5 * ( pvi(ji  ,jj-1) + pvi(ji,jj) )  )
512               z_wnds_t(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
513            END DO
514         END DO
515#if defined key_vectopt_loop
516!CDIR COLLAPSE
517#endif
518         DO jj = 2, jpjm1
519            DO ji = fs_2, fs_jpim1   ! vect. opt.
520               p_taui(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji+1,jj  ) + z_wnds_t(ji,jj) )                          &
521                  &          * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj,1) + sf(jp_wndi)%fnow(ji,jj,1) ) - pui(ji,jj) )
522               p_tauj(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji,jj+1  ) + z_wnds_t(ji,jj) )                          &
523                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - pvi(ji,jj) )
524            END DO
525         END DO
526         CALL lbc_lnk( p_taui  , 'U', -1. )
527         CALL lbc_lnk( p_tauj  , 'V', -1. )
528         CALL lbc_lnk( z_wnds_t, 'T',  1. )
529         !
530      END SELECT
531
532      zztmp = 1. / ( 1. - albo )
533      !                                     ! ========================== !
534      DO jl = 1, ijpl                       !  Loop over ice categories  !
535         !                                  ! ========================== !
536!CDIR NOVERRCHK
537!CDIR COLLAPSE
538         DO jj = 1 , jpj
539!CDIR NOVERRCHK
540            DO ji = 1, jpi
541               ! ----------------------------!
542               !      I   Radiative FLUXES   !
543               ! ----------------------------!
544               zst2 = pst(ji,jj,jl) * pst(ji,jj,jl)
545               zst3 = pst(ji,jj,jl) * zst2
546               ! Short Wave (sw)
547               p_qsr(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
548               ! Long  Wave (lw)
549               z_qlw(ji,jj,jl) = 0.95 * (  sf(jp_qlw)%fnow(ji,jj,1) - Stef * pst(ji,jj,jl) * zst3  ) * tmask(ji,jj,1)
550               ! lw sensitivity
551               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3                                               
552
553               ! ----------------------------!
554               !     II    Turbulent FLUXES  !
555               ! ----------------------------!
556
557               ! ... turbulent heat fluxes
558               ! Sensible Heat
559               z_qsb(ji,jj,jl) = rhoa * cpa * Cice * z_wnds_t(ji,jj) * ( pst(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) )
560               ! Latent Heat
561               p_qla(ji,jj,jl) = MAX( 0.e0, rhoa * Ls  * Cice * z_wnds_t(ji,jj)   &                           
562                  &                    * (  11637800. * EXP( -5897.8 / pst(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj,1)  ) )
563               ! Latent heat sensitivity for ice (Dqla/Dt)
564               p_dqla(ji,jj,jl) = zcoef_dqla * z_wnds_t(ji,jj) / ( zst2 ) * EXP( -5897.8 / pst(ji,jj,jl) )
565               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
566               z_dqsb(ji,jj,jl) = zcoef_dqsb * z_wnds_t(ji,jj)
567
568               ! ----------------------------!
569               !     III    Total FLUXES     !
570               ! ----------------------------!
571               ! Downward Non Solar flux
572               p_qns (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - p_qla (ji,jj,jl)     
573               ! Total non solar heat flux sensitivity for ice
574               p_dqns(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + p_dqla(ji,jj,jl) )   
575            END DO
576            !
577         END DO
578         !
579      END DO
580      !
581      !--------------------------------------------------------------------
582      ! FRACTIONs of net shortwave radiation which is not absorbed in the
583      ! thin surface layer and penetrates inside the ice cover
584      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
585   
586!CDIR COLLAPSE
587      p_fr1(:,:) = ( 0.18 * ( 1.0 - zcoef_frca ) + 0.35 * zcoef_frca )
588!CDIR COLLAPSE
589      p_fr2(:,:) = ( 0.82 * ( 1.0 - zcoef_frca ) + 0.65 * zcoef_frca )
590       
591!CDIR COLLAPSE
592      p_tpr(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
593!CDIR COLLAPSE
594      p_spr(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
595      CALL iom_put( 'snowpre', p_spr )                  ! Snow precipitation
596      !
597      IF(ln_ctl) THEN
598         CALL prt_ctl(tab3d_1=p_qla   , clinfo1=' blk_ice_core: p_qla  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb  : ', kdim=ijpl)
599         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice_core: z_qlw  : ', tab3d_2=p_dqla  , clinfo2=' p_dqla : ', kdim=ijpl)
600         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice_core: z_dqsb : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw : ', kdim=ijpl)
601         CALL prt_ctl(tab3d_1=p_dqns  , clinfo1=' blk_ice_core: p_dqns : ', tab3d_2=p_qsr   , clinfo2=' p_qsr  : ', kdim=ijpl)
602         CALL prt_ctl(tab3d_1=pst     , clinfo1=' blk_ice_core: pst    : ', tab3d_2=p_qns   , clinfo2=' p_qns  : ', kdim=ijpl)
603         CALL prt_ctl(tab2d_1=p_tpr   , clinfo1=' blk_ice_core: p_tpr  : ', tab2d_2=p_spr   , clinfo2=' p_spr  : ')
604         CALL prt_ctl(tab2d_1=p_taui  , clinfo1=' blk_ice_core: p_taui : ', tab2d_2=p_tauj  , clinfo2=' p_tauj : ')
605         CALL prt_ctl(tab2d_1=z_wnds_t, clinfo1=' blk_ice_core: z_wnds_t : ')
606      ENDIF
607
608      IF( wrk_not_released(2, 1)       .OR.   &
609          wrk_not_released(3, 4,5,6,7) )   CALL ctl_stop('blk_ice_core: failed to release workspace arrays')
610      !
611   END SUBROUTINE blk_ice_core
612 
613
614   SUBROUTINE TURB_CORE_1Z(zu, sst, T_a, q_sat, q_a,   &
615      &                        dU , Cd , Ch   , Ce   )
616      !!----------------------------------------------------------------------
617      !!                      ***  ROUTINE  turb_core  ***
618      !!
619      !! ** Purpose :   Computes turbulent transfert coefficients of surface
620      !!                fluxes according to Large & Yeager (2004)
621      !!
622      !! ** Method  :   I N E R T I A L   D I S S I P A T I O N   M E T H O D
623      !!      Momentum, Latent and sensible heat exchange coefficients
624      !!      Caution: this procedure should only be used in cases when air
625      !!      temperature (T_air), air specific humidity (q_air) and wind (dU)
626      !!      are provided at the same height 'zzu'!
627      !!
628      !! References :   Large & Yeager, 2004 : ???
629      !!----------------------------------------------------------------------
630      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released, iwrk_in_use, iwrk_not_released
631      USE wrk_nemo, ONLY: dU10 => wrk_2d_14        ! dU                             [m/s]
632      USE wrk_nemo, ONLY: dT => wrk_2d_15          ! air/sea temperature difference   [K]
633      USE wrk_nemo, ONLY: dq => wrk_2d_16          ! air/sea humidity difference      [K]
634      USE wrk_nemo, ONLY: Cd_n10 => wrk_2d_17      ! 10m neutral drag coefficient
635      USE wrk_nemo, ONLY: Ce_n10 => wrk_2d_18      ! 10m neutral latent coefficient
636      USE wrk_nemo, ONLY: Ch_n10 => wrk_2d_19      ! 10m neutral sensible coefficient
637      USE wrk_nemo, ONLY: sqrt_Cd_n10 => wrk_2d_20 ! root square of Cd_n10
638      USE wrk_nemo, ONLY: sqrt_Cd => wrk_2d_21     ! root square of Cd
639      USE wrk_nemo, ONLY: T_vpot => wrk_2d_22      ! virtual potential temperature    [K]
640      USE wrk_nemo, ONLY: T_star => wrk_2d_23      ! turbulent scale of tem. fluct.
641      USE wrk_nemo, ONLY: q_star => wrk_2d_24      ! turbulent humidity of temp. fluct.
642      USE wrk_nemo, ONLY: U_star => wrk_2d_25      ! turb. scale of velocity fluct.
643      USE wrk_nemo, ONLY: L => wrk_2d_26           ! Monin-Obukov length              [m]
644      USE wrk_nemo, ONLY: zeta => wrk_2d_27        ! stability parameter at height zu
645      USE wrk_nemo, ONLY: U_n10 => wrk_2d_28       ! neutral wind velocity at 10m     [m]
646      USE wrk_nemo, ONLY: xlogt  => wrk_2d_29,    xct => wrk_2d_30,   &
647                          zpsi_h => wrk_2d_31, zpsi_m => wrk_2d_32
648      USE wrk_nemo, ONLY: stab => iwrk_2d_1      ! 1st guess stability test integer
649      !
650      REAL(wp)                , INTENT(in   ) ::   zu      ! altitude of wind measurement       [m]
651      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   sst     ! sea surface temperature         [Kelvin]
652      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   T_a     ! potential air temperature       [Kelvin]
653      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   q_sat   ! sea surface specific humidity   [kg/kg]
654      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   q_a     ! specific air humidity           [kg/kg]
655      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   dU      ! wind module |U(zu)-U(0)|        [m/s]
656      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Cd      ! transfert coefficient for momentum       (tau)
657      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Ch      ! transfert coefficient for temperature (Q_sens)
658      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Ce      ! transfert coefficient for evaporation  (Q_lat)
659      !!
660      INTEGER :: j_itt
661      INTEGER , PARAMETER ::   nb_itt = 3
662      REAL(wp), PARAMETER ::   grav   = 9.8   ! gravity                       
663      REAL(wp), PARAMETER ::   kappa  = 0.4   ! von Karman s constant
664      !!----------------------------------------------------------------------
665
666      IF(  wrk_in_use(2,             14,15,16,17,18,19,        &
667                         20,21,22,23,24,25,26,27,28,29,        &
668                         30,31,32)                      .OR.   &
669          iwrk_in_use(2, 1)                               ) THEN
670         CALL ctl_stop('TURB_CORE_1Z: requested workspace arrays unavailable')   ;   RETURN
671      ENDIF
672
673      !! * Start
674      !! Air/sea differences
675      dU10 = max(0.5, dU)     ! we don't want to fall under 0.5 m/s
676      dT = T_a - sst          ! assuming that T_a is allready the potential temp. at zzu
677      dq = q_a - q_sat
678      !!   
679      !! Virtual potential temperature
680      T_vpot = T_a*(1. + 0.608*q_a)
681      !!
682      !! Neutral Drag Coefficient
683      stab    = 0.5 + sign(0.5,dT)    ! stable : stab = 1 ; unstable : stab = 0
684      Cd_n10  = 1E-3 * ( 2.7/dU10 + 0.142 + dU10/13.09 )    !   L & Y eq. (6a)
685      sqrt_Cd_n10 = sqrt(Cd_n10)
686      Ce_n10  = 1E-3 * ( 34.6 * sqrt_Cd_n10 )               !   L & Y eq. (6b)
687      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1-stab)) !   L & Y eq. (6c), (6d)
688      !!
689      !! Initializing transfert coefficients with their first guess neutral equivalents :
690      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
691
692      !! * Now starting iteration loop
693      DO j_itt=1, nb_itt
694         !! Turbulent scales :
695         U_star  = sqrt_Cd*dU10                !   L & Y eq. (7a)
696         T_star  = Ch/sqrt_Cd*dT               !   L & Y eq. (7b)
697         q_star  = Ce/sqrt_Cd*dq               !   L & Y eq. (7c)
698
699         !! Estimate the Monin-Obukov length :
700         L  = (U_star**2)/( kappa*grav*(T_star/T_vpot + q_star/(q_a + 1./0.608)) )
701
702         !! Stability parameters :
703         zeta  = zu/L ;  zeta   = sign( min(abs(zeta),10.0), zeta )
704         zpsi_h  = psi_h(zeta)
705         zpsi_m  = psi_m(zeta)
706
707         !! Shifting the wind speed to 10m and neutral stability :
708         U_n10 = dU10*1./(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)) !  L & Y eq. (9a)
709
710         !! Updating the neutral 10m transfer coefficients :
711         Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)              !  L & Y eq. (6a)
712         sqrt_Cd_n10 = sqrt(Cd_n10)
713         Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                           !  L & Y eq. (6b)
714         stab    = 0.5 + sign(0.5,zeta)
715         Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab))           !  L & Y eq. (6c), (6d)
716
717         !! Shifting the neutral  10m transfer coefficients to ( zu , zeta ) :
718         !!
719         xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10) - zpsi_m)
720         Cd  = Cd_n10/(xct*xct) ;  sqrt_Cd = sqrt(Cd)
721         !!
722         xlogt = log(zu/10.) - zpsi_h
723         !!
724         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
725         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
726         !!
727         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
728         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
729         !!
730      END DO
731      !!
732      IF( wrk_not_released(2,             14,15,16,17,18,19,          &
733         &                    20,21,22,23,24,25,26,27,28,29,          &
734         &                    30,31,32                      )   .OR.  &     
735         iwrk_not_released(2, 1)                                  )   &
736         CALL ctl_stop('TURB_CORE_1Z: failed to release workspace arrays')
737      !
738    END SUBROUTINE TURB_CORE_1Z
739
740
741    SUBROUTINE TURB_CORE_2Z(zt, zu, sst, T_zt, q_sat, q_zt, dU, Cd, Ch, Ce, T_zu, q_zu)
742      !!----------------------------------------------------------------------
743      !!                      ***  ROUTINE  turb_core  ***
744      !!
745      !! ** Purpose :   Computes turbulent transfert coefficients of surface
746      !!                fluxes according to Large & Yeager (2004).
747      !!
748      !! ** Method  :   I N E R T I A L   D I S S I P A T I O N   M E T H O D
749      !!      Momentum, Latent and sensible heat exchange coefficients
750      !!      Caution: this procedure should only be used in cases when air
751      !!      temperature (T_air) and air specific humidity (q_air) are at 2m
752      !!      whereas wind (dU) is at 10m.
753      !!
754      !! References :   Large & Yeager, 2004 : ???
755      !!----------------------------------------------------------------------
756      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released, iwrk_in_use, iwrk_not_released
757      USE wrk_nemo, ONLY: dU10 => wrk_2d_14        ! dU                             [m/s]
758      USE wrk_nemo, ONLY: dT => wrk_2d_15          ! air/sea temperature difference   [K]
759      USE wrk_nemo, ONLY: dq => wrk_2d_16          ! air/sea humidity difference      [K]
760      USE wrk_nemo, ONLY: Cd_n10 => wrk_2d_17      ! 10m neutral drag coefficient
761      USE wrk_nemo, ONLY: Ce_n10 => wrk_2d_18      ! 10m neutral latent coefficient
762      USE wrk_nemo, ONLY: Ch_n10 => wrk_2d_19      ! 10m neutral sensible coefficient
763      USE wrk_nemo, ONLY: sqrt_Cd_n10 => wrk_2d_20 ! root square of Cd_n10
764      USE wrk_nemo, ONLY: sqrt_Cd => wrk_2d_21     ! root square of Cd
765      USE wrk_nemo, ONLY: T_vpot => wrk_2d_22      ! virtual potential temperature    [K]
766      USE wrk_nemo, ONLY: T_star => wrk_2d_23     ! turbulent scale of tem. fluct.
767      USE wrk_nemo, ONLY: q_star => wrk_2d_24     ! turbulent humidity of temp. fluct.
768      USE wrk_nemo, ONLY: U_star => wrk_2d_25     ! turb. scale of velocity fluct.
769      USE wrk_nemo, ONLY: L => wrk_2d_26          ! Monin-Obukov length              [m]
770      USE wrk_nemo, ONLY: zeta_u => wrk_2d_27     ! stability parameter at height zu
771      USE wrk_nemo, ONLY: zeta_t => wrk_2d_28     ! stability parameter at height zt
772      USE wrk_nemo, ONLY: U_n10 => wrk_2d_29      ! neutral wind velocity at 10m     [m]
773      USE wrk_nemo, ONLY: xlogt => wrk_2d_30, xct => wrk_2d_31, zpsi_hu => wrk_2d_32, zpsi_ht => wrk_2d_33, zpsi_m => wrk_2d_34
774      USE wrk_nemo, ONLY: stab => iwrk_2d_1      ! 1st guess stability test integer
775      !!
776      REAL(wp), INTENT(in)   :: &
777         zt,      &     ! height for T_zt and q_zt                   [m]
778         zu             ! height for dU                              [m]
779      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::  &
780         sst,      &     ! sea surface temperature              [Kelvin]
781         T_zt,     &     ! potential air temperature            [Kelvin]
782         q_sat,    &     ! sea surface specific humidity         [kg/kg]
783         q_zt,     &     ! specific air humidity                 [kg/kg]
784         dU              ! relative wind module |U(zu)-U(0)|       [m/s]
785      REAL(wp), INTENT(out), DIMENSION(jpi,jpj)  ::  &
786         Cd,       &     ! transfer coefficient for momentum         (tau)
787         Ch,       &     ! transfer coefficient for sensible heat (Q_sens)
788         Ce,       &     ! transfert coefficient for evaporation   (Q_lat)
789         T_zu,     &     ! air temp. shifted at zu                     [K]
790         q_zu            ! spec. hum.  shifted at zu               [kg/kg]
791
792      INTEGER :: j_itt
793      INTEGER, PARAMETER :: nb_itt = 3   ! number of itterations
794      REAL(wp), PARAMETER ::                        &
795         grav   = 9.8,      &  ! gravity                       
796         kappa  = 0.4          ! von Karman's constant
797      !!----------------------------------------------------------------------
798      !!  * Start
799
800      IF(  wrk_in_use(2,             14,15,16,17,18,19,        &
801                         20,21,22,23,24,25,26,27,28,29,        &         
802                         30,31,32,33,34)                .OR.   &
803          iwrk_in_use(2, 1)                               ) THEN
804         CALL ctl_stop('TURB_CORE_2Z: requested workspace arrays unavailable')   ;   RETURN
805      ENDIF
806
807      !! Initial air/sea differences
808      dU10 = max(0.5, dU)      !  we don't want to fall under 0.5 m/s
809      dT = T_zt - sst 
810      dq = q_zt - q_sat
811
812      !! Neutral Drag Coefficient :
813      stab = 0.5 + sign(0.5,dT)                 ! stab = 1  if dT > 0  -> STABLE
814      Cd_n10  = 1E-3*( 2.7/dU10 + 0.142 + dU10/13.09 ) 
815      sqrt_Cd_n10 = sqrt(Cd_n10)
816      Ce_n10  = 1E-3*( 34.6 * sqrt_Cd_n10 )
817      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1 - stab))
818
819      !! Initializing transf. coeff. with their first guess neutral equivalents :
820      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
821
822      !! Initializing z_u values with z_t values :
823      T_zu = T_zt ;  q_zu = q_zt
824
825      !!  * Now starting iteration loop
826      DO j_itt=1, nb_itt
827         dT = T_zu - sst ;  dq = q_zu - q_sat ! Updating air/sea differences
828         T_vpot = T_zu*(1. + 0.608*q_zu)    ! Updating virtual potential temperature at zu
829         U_star = sqrt_Cd*dU10                ! Updating turbulent scales :   (L & Y eq. (7))
830         T_star  = Ch/sqrt_Cd*dT              !
831         q_star  = Ce/sqrt_Cd*dq              !
832         !!
833         L = (U_star*U_star) &                ! Estimate the Monin-Obukov length at height zu
834              & / (kappa*grav/T_vpot*(T_star*(1.+0.608*q_zu) + 0.608*T_zu*q_star))
835         !! Stability parameters :
836         zeta_u  = zu/L  ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
837         zeta_t  = zt/L  ;  zeta_t = sign( min(abs(zeta_t),10.0), zeta_t )
838         zpsi_hu = psi_h(zeta_u)
839         zpsi_ht = psi_h(zeta_t)
840         zpsi_m  = psi_m(zeta_u)
841         !!
842         !! Shifting the wind speed to 10m and neutral stability : (L & Y eq.(9a))
843!        U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u)))
844         U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m))
845         !!
846         !! Shifting temperature and humidity at zu :          (L & Y eq. (9b-9c))
847!        T_zu = T_zt - T_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
848         T_zu = T_zt - T_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
849!        q_zu = q_zt - q_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
850         q_zu = q_zt - q_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
851         !!
852         !! q_zu cannot have a negative value : forcing 0
853         stab = 0.5 + sign(0.5,q_zu) ;  q_zu = stab*q_zu
854         !!
855         !! Updating the neutral 10m transfer coefficients :
856         Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)    ! L & Y eq. (6a)
857         sqrt_Cd_n10 = sqrt(Cd_n10)
858         Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                 ! L & Y eq. (6b)
859         stab    = 0.5 + sign(0.5,zeta_u)
860         Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab)) ! L & Y eq. (6c-6d)
861         !!
862         !!
863         !! Shifting the neutral 10m transfer coefficients to (zu,zeta_u) :
864!        xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u))
865         xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)
866         Cd = Cd_n10/(xct*xct) ; sqrt_Cd = sqrt(Cd)
867         !!
868!        xlogt = log(zu/10.) - psi_h(zeta_u)
869         xlogt = log(zu/10.) - zpsi_hu
870         !!
871         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
872         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
873         !!
874         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
875         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
876         !!
877         !!
878      END DO
879      !!
880     IF( wrk_not_released(2,              14,15,16,17,18,19,          &
881         &                    20,21,22,23,24,25,26,27,28,29,          &
882         &                    30,31,32,33,34                )   .OR.  & 
883         iwrk_not_released(2, 1)                                  )   &
884         CALL ctl_stop('TURB_CORE_1Z: failed to release workspace arrays')
885      !
886    END SUBROUTINE TURB_CORE_2Z
887
888
889    FUNCTION psi_m(zta)   !! Psis, L & Y eq. (8c), (8d), (8e)
890      !-------------------------------------------------------------------------------
891      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
892      USE wrk_nemo, ONLY:     X2 => wrk_2d_35
893      USE wrk_nemo, ONLY:     X  => wrk_2d_36
894      USE wrk_nemo, ONLY: stabit => wrk_2d_37
895      !!
896      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta
897
898      REAL(wp), PARAMETER :: pi = 3.141592653589793_wp
899      REAL(wp), DIMENSION(jpi,jpj)             :: psi_m
900      !-------------------------------------------------------------------------------
901
902      IF( wrk_in_use(2, 35,36,37) ) THEN
903         CALL ctl_stop('psi_m: requested workspace arrays unavailable')   ;   RETURN
904      ENDIF
905
906      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.0) ;  X  = sqrt(X2)
907      stabit    = 0.5 + sign(0.5,zta)
908      psi_m = -5.*zta*stabit  &                                                          ! Stable
909         &    + (1. - stabit)*(2*log((1. + X)/2) + log((1. + X2)/2) - 2*atan(X) + pi/2)  ! Unstable
910
911      IF( wrk_not_released(2, 35,36,37) )   CALL ctl_stop('psi_m: failed to release workspace arrays')
912      !
913    END FUNCTION psi_m
914
915
916    FUNCTION psi_h( zta )    !! Psis, L & Y eq. (8c), (8d), (8e)
917      !-------------------------------------------------------------------------------
918      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
919      USE wrk_nemo, ONLY:     X2 => wrk_2d_35
920      USE wrk_nemo, ONLY:     X  => wrk_2d_36
921      USE wrk_nemo, ONLY: stabit => wrk_2d_37
922      !
923      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   zta
924      !
925      REAL(wp), DIMENSION(jpi,jpj)             ::   psi_h
926      !-------------------------------------------------------------------------------
927
928      IF( wrk_in_use(2, 35,36,37) ) THEN
929         CALL ctl_stop('psi_h: requested workspace arrays unavailable')   ;   RETURN
930      ENDIF
931
932      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.) ;  X  = sqrt(X2)
933      stabit    = 0.5 + sign(0.5,zta)
934      psi_h = -5.*zta*stabit  &                                       ! Stable
935         &    + (1. - stabit)*(2.*log( (1. + X2)/2. ))                 ! Unstable
936
937      IF( wrk_not_released(2, 35,36,37) )   CALL ctl_stop('psi_h: failed to release workspace arrays')
938      !
939    END FUNCTION psi_h
940 
941   !!======================================================================
942END MODULE sbcblk_core
Note: See TracBrowser for help on using the repository browser.