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 branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 3186

Last change on this file since 3186 was 3186, checked in by smasson, 12 years ago

dev_NEMO_MERGE_2011: replace the old wrk_nemo with the new wrk_nemo

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