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 on Ticket #583 – Attachment – NEMO

Ticket #583: sbcblk_core.F90

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