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

Last change on this file since 3963 was 3963, checked in by clem, 11 years ago

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