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

source: branches/2013/dev_ECMWF_waves/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 4351

Last change on this file since 4351 was 4351, checked in by acc, 10 years ago

Branch 2013/dev_ECMWF_waves ticket #1216. First set of code imports from ECMWF. Untested

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