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

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

Branch 2013/dev_ECMWF_waves ticket #1216. Corrections to first set of code imports to solve compilation issues

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