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/2013/dev_r3406_CNRS_LIM3/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2013/dev_r3406_CNRS_LIM3/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 3938

Last change on this file since 3938 was 3938, checked in by flavoni, 11 years ago

dev_r3406_CNRS_LIM3: update LIM3, see ticket #1116

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