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

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 3116

Last change on this file since 3116 was 3116, checked in by cetlod, 12 years ago

dev_NEMO_MERGE_2011: add in changes dev_NOC_UKMO_MERGE developments

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