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

source: branches/2014/dev_4728_CNRS04_coupled_interface/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 4859

Last change on this file since 4859 was 4859, checked in by smasson, 9 years ago

dev_4728_CNRS04_coupled_interface: cleaning of the coupling interface for OASIS3-MCT. 2

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