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

source: trunk/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 1156

Last change on this file since 1156 was 1156, checked in by rblod, 16 years ago

Update Id and licence information, see ticket #210

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