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

source: branches/UKMO/2014_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90

Last change on this file was 5275, checked in by rfurner, 9 years ago

added Charnock coefficient for wind forcing

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