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

source: branches/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 2236

Last change on this file since 2236 was 2236, checked in by cetlod, 14 years ago

First guess of NEMO_v3.3

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