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 @ 1275

Last change on this file since 1275 was 1275, checked in by rblod, 15 years ago

First introduction off interpolation off the fly, see ticket #279

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