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

source: branches/2014/dev_CNRS_2014/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 4896

Last change on this file since 4896 was 4896, checked in by cetlod, 10 years ago

2014/dev_CNRS_2014 : merge the 1st branch onto dev_CNRS_2014, see ticket #1415

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