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

Last change on this file since 1559 was 1482, checked in by smasson, 15 years ago

distribution of iom_put + cleaning of LIM2 outputs, see ticket:437

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