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

source: branches/DEV_r1879_FCM/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 2292

Last change on this file since 2292 was 2292, checked in by smasson, 14 years ago

update DEV_r1879_FCM for additional tests...

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