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

Last change on this file since 1649 was 1601, checked in by ctlod, 15 years ago

Doctor naming of OPA namelist variables , see ticket: #526

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