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

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 2633

Last change on this file since 2633 was 2633, checked in by trackstand2, 13 years ago

Renamed wrk_use => wrk_in_use and wrk_release => wrk_not_released

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