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

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 2412

Last change on this file since 2412 was 2412, checked in by gm, 13 years ago

v3.3beta: #762 Bug correction in CORE bulk on NEC

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