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

source: branches/2011/dev_r2784_CMCC1_topbfm/NEMOGCM/NEMO/OPA_SRC/SBC/sbcblk_core.F90 @ 2801

Last change on this file since 2801 was 2801, checked in by vichi, 13 years ago

Bug corrected in sbcblk_core.F90

  • Property svn:keywords set to Id
File size: 64.5 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   !!----------------------------------------------------------------------
17
18   !!----------------------------------------------------------------------
19   !!   sbc_blk_core  : bulk formulation as ocean surface boundary condition
20   !!                   (forced mode, CORE bulk formulea)
21   !!   blk_oce_core  : ocean: computes momentum, heat and freshwater fluxes
22   !!   blk_ice_core  : ice  : computes momentum, heat and freshwater fluxes
23   !!   turb_core     : computes the CORE turbulent transfer coefficients
24   !!----------------------------------------------------------------------
25   USE oce             ! ocean dynamics and tracers
26   USE dom_oce         ! ocean space and time domain
27   USE phycst          ! physical constants
28   USE fldread         ! read input fields
29   USE sbc_oce         ! Surface boundary condition: ocean fields
30   USE sbcdcy          ! surface boundary condition: diurnal cycle
31   USE iom             ! I/O manager library
32   USE in_out_manager  ! I/O manager
33   USE lib_mpp         ! distribued memory computing library
34   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
35   USE prtctl          ! Print control
36#if defined key_lim3
37   USE sbc_ice         ! Surface boundary condition: ice fields
38#endif
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   sbc_blk_core         ! routine called in sbcmod module
44   PUBLIC   blk_ice_core         ! routine called in sbc_ice_lim module
45
46   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point
47   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point
48   INTEGER , PARAMETER ::   jp_humi = 3           ! index of specific humidity               ( - )
49   INTEGER , PARAMETER ::   jp_qsr  = 4           ! index of solar heat                      (W/m2)
50   INTEGER , PARAMETER ::   jp_qlw  = 5           ! index of Long wave                       (W/m2)
51   INTEGER , PARAMETER ::   jp_tair = 6           ! index of 10m air temperature             (Kelvin)
52   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s)
53   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s)
54   INTEGER , PARAMETER ::   jp_tdif = 9           ! index of tau diff associated to HF tau   (N/m2)   at T-point
55#if defined key_orca_r025
56   INTEGER , PARAMETER ::   jp_swc  = 10          ! index of GEWEX correction for SW radiation  at T-point
57   INTEGER , PARAMETER ::   jp_lwc  = 11          ! index of GEWEX correction for LW radiation  at T-point
58   INTEGER , PARAMETER ::   jp_prc  = 12          ! index of PMWC correction forat T-point
59   INTEGER , PARAMETER ::   jpfld   = 12          ! maximum number of files to read
60#else
61   INTEGER , PARAMETER ::   jpfld   = 9           ! maximum number of files to read
62#endif
63   
64   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read)
65         
66   !                                             !!! CORE bulk parameters
67   REAL(wp), PARAMETER ::   rhoa =    1.22        ! air density
68   REAL(wp), PARAMETER ::   cpa  = 1000.5         ! specific heat of air
69   REAL(wp), PARAMETER ::   Lv   =    2.5e6       ! latent heat of vaporization
70   REAL(wp), PARAMETER ::   Ls   =    2.839e6     ! latent heat of sublimation
71   REAL(wp), PARAMETER ::   Stef =    5.67e-8     ! Stefan Boltzmann constant
72   REAL(wp), PARAMETER ::   Cice =    1.63e-3     ! transfer coefficient over ice
73   REAL(wp), PARAMETER ::   albo =    0.066       ! ocean albedo assumed to be contant
74
75   !                                  !!* Namelist namsbc_core : CORE bulk parameters
76   LOGICAL  ::   ln_2m     = .FALSE.   ! logical flag for height of air temp. and hum
77   LOGICAL  ::   ln_taudif = .FALSE.   ! logical flag to use the "mean of stress module - module of mean stress" data
78   REAL(wp) ::   rn_pfac   = 1.        ! multiplication factor for precipitation
79#if defined key_orca_r025
80   LOGICAL  ::   ln_printdia= .TRUE.     ! logical flag for height of air temp. and hum
81   LOGICAL  ::   ln_netsw   = .TRUE.     ! logical flag for height of air temp. and hum
82   LOGICAL  ::   ln_core_graceopt=.FALSE., ln_core_spinup=.FALSE.
83   LOGICAL  ::   ln_gwxc = .TRUE.
84   LOGICAL  ::   ln_corad_antar =.FALSE., ln_corad_arc =.FALSE. , ln_cotair_arc = .FALSE.
85   LOGICAL  ::   ln_coprecip =.FALSE.
86   REAL(wp) ::   rn_qns_bias = 0._wp     ! heat flux bias
87
88   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::  area 
89   REAL(wp) :: araux
90   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::  zqlw, zqsb  ! long wave and sensible heat fluxes
91   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::  zqla, zevap ! latent heat fluxes and evaporation
92
93   REAL(wp), PARAMETER :: zalph = 2.408724e-06_wp, &
94   & zbet = -0.006936579_wp, zgam = 449.9094_wp ! GRACE regression coefficients
95#endif
96
97   !! * Substitutions
98#  include "domzgr_substitute.h90"
99#  include "vectopt_loop_substitute.h90"
100   !!----------------------------------------------------------------------
101   !! NEMO/OPA 3.3 , NEMO-consortium (2010)
102   !! $Id$
103   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
104   !!----------------------------------------------------------------------
105CONTAINS
106   
107   SUBROUTINE sbc_blk_core( kt )
108      !!---------------------------------------------------------------------
109      !!                    ***  ROUTINE sbc_blk_core  ***
110      !!                   
111      !! ** Purpose :   provide at each time step the surface ocean fluxes
112      !!      (momentum, heat, freshwater and runoff)
113      !!
114      !! ** Method  : (1) READ each fluxes in NetCDF files:
115      !!      the 10m wind velocity (i-component) (m/s)    at T-point
116      !!      the 10m wind velocity (j-component) (m/s)    at T-point
117      !!      the specific humidity               ( - )
118      !!      the solar heat                      (W/m2)
119      !!      the Long wave                       (W/m2)
120      !!      the 10m air temperature             (Kelvin)
121      !!      the total precipitation (rain+snow) (Kg/m2/s)
122      !!      the snow (solid prcipitation)       (kg/m2/s)
123      !!   OPTIONAL parameter (see ln_taudif namelist flag):
124      !!      the tau diff associated to HF tau   (N/m2)   at T-point
125      !!              (2) CALL blk_oce_core
126      !!
127      !!      C A U T I O N : never mask the surface stress fields
128      !!                      the stress is assumed to be in the mesh referential
129      !!                      i.e. the (i,j) referential
130      !!
131      !! ** Action  :   defined at each time-step at the air-sea interface
132      !!              - utau, vtau  i- and j-component of the wind stress
133      !!              - taum        wind stress module at T-point
134      !!              - wndm        10m wind module at T-point
135      !!              - qns, qsr    non-slor and solar heat flux
136      !!              - emp, emps   evaporation minus precipitation
137      !!----------------------------------------------------------------------
138#if defined key_orca_r025
139      USE ice_2
140#endif
141      INTEGER, INTENT(in) ::   kt   ! ocean time step
142      !!
143      INTEGER  ::   ierror   ! return error code
144      INTEGER  ::   ifpr     ! dummy loop indice
145      INTEGER  ::   jfld     ! dummy loop arguments
146      INTEGER  ::   ji, jj
147      !!
148      CHARACTER(len=100) ::  cn_dir   !   Root directory for location of core files
149      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i     ! array of namelist informations on the fields to read
150      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read
151      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !   "                                 "
152      TYPE(FLD_N) ::   sn_tdif                                 !   "                                 "
153      TYPE(FLD_N) ::   sn_swc, sn_lwc                          !   "                                 "
154      TYPE(FLD_N) ::   sn_prc
155      INTEGER  ::   iter_shapiro = 250
156      REAL :: zzlat, zzlat1, zzlat2, zfm, zfrld, ztmp
157      REAL(wp), DIMENSION(jpi,jpj):: xyt,z_qsr,z_qlw,z_qsr1,z_qlw1, z_hum, z_tair
158      REAL(wp), DIMENSION(jpi,jpj):: zqsr_lr, zqsr_hr, zqlw_lr, zqlw_hr, zprec_hr, zprec_lr
159      CHARACTER(len=20)  ::  c_kind='ORCA_GLOB'
160#if defined key_orca_r025
161      NAMELIST/namsbc_core/ cn_dir , ln_2m  , ln_taudif, rn_pfac,           &
162         &                  sn_wndi, sn_wndj, sn_humi  , sn_qsr ,           &
163         &                  sn_qlw , sn_tair, sn_prec  , sn_snow, sn_tdif,  &
164         &                  sn_swc , sn_lwc , sn_prc   , ln_gwxc,           &
165         &                  ln_corad_antar, ln_corad_arc, ln_cotair_arc, ln_coprecip ,  &
166         &                  rn_qns_bias, ln_printdia, ln_netsw, ln_core_graceopt,ln_core_spinup
167      !!---------------------------------------------------------------------
168#else
169      NAMELIST/namsbc_core/ cn_dir , ln_2m  , ln_taudif, rn_pfac,           &
170         &                  sn_wndi, sn_wndj, sn_humi  , sn_qsr ,           &
171         &                  sn_qlw , sn_tair, sn_prec  , sn_snow, sn_tdif
172      !!---------------------------------------------------------------------
173#endif
174
175      !                                         ! ====================== !
176      IF( kt == nit000 ) THEN                   !  First call kt=nit000  !
177         !                                      ! ====================== !
178#if defined key_orca_r025
179         !                                      !==  allocate sbc arrays
180            IF( sbc_blk_core_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_core_alloc : unable to allocate arrays' )
181#endif
182
183         ! set file information (default values)
184         cn_dir = './'       ! directory in which the model is executed
185         !
186         ! (NB: frequency positive => hours, negative => months)
187         !            !    file    ! frequency ! variable ! time intep !  clim   ! 'yearly' or ! weights  ! rotation !
188         !            !    name    !  (hours)  !  name    !   (T/F)    !  (T/F)  !  'monthly'  ! filename ! pairs    !
189         sn_wndi = FLD_N( 'uwnd10m',    24     , 'u_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
190         sn_wndj = FLD_N( 'vwnd10m',    24     , 'v_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
191         sn_qsr  = FLD_N( 'qsw'    ,    24     , 'qsw'    ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
192         sn_qlw  = FLD_N( 'qlw'    ,    24     , 'qlw'    ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
193         sn_tair = FLD_N( 'tair10m',    24     , 't_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
194         sn_humi = FLD_N( 'humi10m',    24     , 'q_10'   ,  .false.   , .false. ,   'yearly'  , ''       , ''       )
195         sn_prec = FLD_N( 'precip' ,    -1     , 'precip' ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
196         sn_snow = FLD_N( 'snow'   ,    -1     , 'snow'   ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
197         sn_tdif = FLD_N( 'taudif' ,    24     , 'taudif' ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
198#if defined key_orca_r025
199         sn_swc  = FLD_N( 'swc'    ,    24     ,  'swc'   ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
200         sn_lwc  = FLD_N( 'lwc'    ,    24     ,  'lwc'   ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
201         sn_prc  = FLD_N( 'prc'    ,    24     ,  'prc'   ,  .true.    , .false. ,   'yearly'  , ''       , ''       )
202#endif
203         !
204         REWIND( numnam )                          ! read in namlist namsbc_core
205         READ  ( numnam, namsbc_core )
206         !                                         ! check: do we plan to use ln_dm2dc with non-daily forcing?
207         IF( ln_dm2dc .AND. sn_qsr%nfreqh /= 24 )   & 
208            &   CALL ctl_stop( 'sbc_blk_core: ln_dm2dc can be activated only with daily short-wave forcing' ) 
209         IF( ln_dm2dc .AND. sn_qsr%ln_tint ) THEN
210            CALL ctl_warn( 'sbc_blk_core: ln_dm2dc is taking care of the temporal interpolation of daily qsr',   &
211                 &         '              ==> We force time interpolation = .false. for qsr' )
212            sn_qsr%ln_tint = .false.
213         ENDIF
214         !                                         ! store namelist information in an array
215         slf_i(jp_wndi) = sn_wndi   ;   slf_i(jp_wndj) = sn_wndj
216         slf_i(jp_qsr ) = sn_qsr    ;   slf_i(jp_qlw ) = sn_qlw
217         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi
218         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow
219         slf_i(jp_tdif) = sn_tdif
220#if defined key_orca_r025
221            slf_i(jp_swc ) = sn_swc
222            slf_i(jp_lwc ) = sn_lwc
223            slf_i(jp_prc ) = sn_prc
224#endif
225         !                 
226         lhftau = ln_taudif                        ! do we use HF tau information?
227         jfld = jpfld - COUNT( (/.NOT. lhftau/) )
228#if defined key_orca_r025
229         IF( .NOT. ln_gwxc )     jfld = jfld - 2
230         IF( .NOT. ln_coprecip ) jfld = jfld - 1
231#endif
232         !
233         ALLOCATE( sf(jfld), STAT=ierror )         ! set sf structure
234         IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_blk_core: unable to allocate sf structure' )
235         DO ifpr= 1, jfld
236            ALLOCATE( sf(ifpr)%fnow(jpi,jpj,1) )
237            IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf(ifpr)%fdta(jpi,jpj,1,2) )
238         END DO
239         !                                         ! fill sf with slf_i and control print
240         CALL fld_fill( sf, slf_i, cn_dir, 'sbc_blk_core', 'flux formulation for ocean surface boundary condition', 'namsbc_core' )
241         !
242#if defined key_orca_r025
243         IF( ln_printdia .OR. ln_core_graceopt ) THEN
244           area = (e1t * e2t)
245           araux = sum ( area * tmask(:,:,1) )
246           IF(lk_mpp) CALL mpp_sum ( araux )
247         ENDIF
248#endif
249      ENDIF
250
251      CALL fld_read( kt, nn_fsbc, sf )        ! input fields provided at the current time-step
252
253      IF( MOD( kt - 1, nn_fsbc ) == 0 ) THEN
254
255#if defined key_orca_r025
256      ! Introduce ERA-Interim filtering and correction
257
258         IF( ln_gwxc ) THEN
259
260           call Shapiro_1D(sf(jp_qsr)%fnow(:,:,1),iter_shapiro, c_kind, zqsr_lr)
261           zqsr_hr(:,:)=sf(jp_qsr)%fnow(:,:,1)-zqsr_lr(:,:)          ! We get large scale and small scale
262
263           call Shapiro_1D(sf(jp_qlw)%fnow(:,:,1),iter_shapiro, c_kind, zqlw_lr)
264           zqlw_hr(:,:)=sf(jp_qlw)%fnow(:,:,1)-zqlw_lr(:,:)          ! We get large scale and small scale
265
266           z_qsr1(:,:)=zqsr_lr(:,:)*sf(jp_swc)%fnow(:,:,1) + zqsr_hr(:,:)
267           z_qlw1(:,:)=zqlw_lr(:,:)*sf(jp_lwc)%fnow(:,:,1) + zqlw_hr(:,:)
268
269           DO jj=1,jpj
270             DO ji=1,jpi
271               z_qsr1(ji,jj)=max(z_qsr1(ji,jj),0.0)
272               z_qlw1(ji,jj)=max(z_qlw1(ji,jj),0.0)
273             END DO
274           END DO
275
276         ENDIF
277
278         IF( ln_coprecip ) THEN
279
280           call Shapiro_1D(sf(jp_prec)%fnow(:,:,1),iter_shapiro,c_kind,zprec_lr)
281           zprec_hr(:,:)=sf(jp_prec)%fnow(:,:,1)-zprec_lr(:,:)       ! We get large scale and small scale
282
283           DO jj=1,jpj
284             DO ji=1,jpi
285               IF( zprec_lr(ji,jj) .GT. 0._wp ) THEN
286                  ztmp = LOG( ( 1000._wp + sf(jp_prc)%fnow(ji,jj,1) ) * EXP( zprec_lr(ji,jj) ) / 1000._wp )
287                  sf(jp_prec)%fnow(ji,jj,1) = max(ztmp+zprec_hr(ji,jj),0.0)
288               ENDIF
289             END DO
290           END DO
291
292         ENDIF
293
294         IF ( ln_corad_antar ) THEN           ! correction of SW and LW in the Southern Ocean
295 
296           z_qsr(:,:)=0.8*z_qsr1(:,:)
297           z_qlw(:,:)=1.1*z_qlw1(:,:)
298           xyt(:,:) = 0.e0
299           zzlat1 = -65.
300           zzlat2 = -60.
301           DO jj = 1, jpj
302             DO ji = 1, jpi
303               zzlat = gphit(ji,jj)
304               IF ( zzlat >= zzlat1 .AND. zzlat <= zzlat2 ) THEN
305                  xyt(ji,jj) = (zzlat2-zzlat)/(zzlat2-zzlat1)
306               ELSE IF ( zzlat < zzlat1 ) THEN
307                  xyt(ji,jj) = 1
308               ENDIF
309             END DO
310           END DO
311           z_qsr1(:,:)=z_qsr(:,:)*xyt(:,:)+(1.0-xyt(:,:))*z_qsr1(:,:)
312           z_qlw1(:,:)=z_qlw(:,:)*xyt(:,:)+(1.0-xyt(:,:))*z_qlw1(:,:)
313
314         ENDIF
315
316         IF ( ln_corad_arc ) THEN         ! correction of SW in the Arctic Ocean
317
318           z_qsr(:,:)=0.7*z_qsr1(:,:)
319           xyt(:,:) = 0.e0
320           zzlat1 = 78.
321           zzlat2 = 82.
322           DO jj = 1, jpj
323             DO ji = 1, jpi
324               zzlat = gphit(ji,jj)
325               IF ( zzlat >= zzlat1 .AND. zzlat <= zzlat2 ) THEN
326                  xyt(ji,jj) = (zzlat-zzlat1)/(zzlat2-zzlat1)
327               ELSE IF ( zzlat > zzlat2 ) THEN
328                  xyt(ji,jj) = 1
329               ENDIF
330             END DO
331           END DO
332           z_qsr1(:,:)=z_qsr(:,:)*xyt(:,:)+(1.0-xyt(:,:))*z_qsr1(:,:)
333
334         ENDIF
335
336         sf(jp_qsr)%fnow(:,:,1)=z_qsr1(:,:)
337         sf(jp_qlw)%fnow(:,:,1)=z_qlw1(:,:)
338
339         IF ( ln_cotair_arc ) THEN           ! correction of Air Temperature in the Arctic Ocean
340
341           z_tair(:,:)=sf(jp_tair)%fnow(:,:,1) - 2.0
342           xyt(:,:) = 0.e0 ; zzlat1 = 78. ; zzlat2 = 82.
343           DO jj = 1, jpj
344             DO ji = 1, jpi
345               zzlat = gphit(ji,jj) ; zfrld=frld(ji,jj)
346               IF ( zzlat >= zzlat1 .AND. zzlat <= zzlat2 .AND. zfrld < 0.85 ) THEN
347                  xyt(ji,jj) = (zzlat-zzlat1)/(zzlat2-zzlat1)
348               ELSE IF ( zzlat > zzlat2 .AND. zfrld < 0.85 ) THEN
349                  xyt(ji,jj) = 1
350               ENDIF
351             END DO
352           END DO
353           sf(jp_tair)%fnow(:,:,1)=z_tair(:,:)*xyt(:,:)+(1.0-xyt(:,:))*sf(jp_tair)%fnow(:,:,1)
354
355         ENDIF
356
357#endif
358
359         CALL blk_oce_core( sf, sst_m, ssu_m, ssv_m )   ! compute the surface ocean fluxes using CLIO bulk formulea
360
361      ENDIF
362      !                                                  ! using CORE bulk formulea
363   END SUBROUTINE sbc_blk_core
364   
365   
366   SUBROUTINE blk_oce_core( sf, pst, pu, pv )
367      !!---------------------------------------------------------------------
368      !!                     ***  ROUTINE blk_core  ***
369      !!
370      !! ** Purpose :   provide the momentum, heat and freshwater fluxes at
371      !!      the ocean surface at each time step
372      !!
373      !! ** Method  :   CORE bulk formulea for the ocean using atmospheric
374      !!      fields read in sbc_read
375      !!
376      !! ** Outputs : - utau    : i-component of the stress at U-point  (N/m2)
377      !!              - vtau    : j-component of the stress at V-point  (N/m2)
378      !!              - taum    : Wind stress module at T-point         (N/m2)
379      !!              - wndm    : Wind speed module at T-point          (m/s)
380      !!              - qsr     : Solar heat flux over the ocean        (W/m2)
381      !!              - qns     : Non Solar heat flux over the ocean    (W/m2)
382      !!              - evap    : Evaporation over the ocean            (kg/m2/s)
383      !!              - emp(s)  : evaporation minus precipitation       (kg/m2/s)
384      !!
385      !!  ** Nota  :   sf has to be a dummy argument for AGRIF on NEC
386      !!---------------------------------------------------------------------
387      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
388      USE wrk_nemo, ONLY:   zwnd_i => wrk_2d_1  , zwnd_j => wrk_2d_2      ! wind speed components at T-point
389      USE wrk_nemo, ONLY:   zqsatw => wrk_2d_3           ! specific humidity at pst
390      USE wrk_nemo, ONLY:   zqlw   => wrk_2d_4  , zqsb   => wrk_2d_5      ! long wave and sensible heat fluxes
391      USE wrk_nemo, ONLY:   zqla   => wrk_2d_6  , zevap  => wrk_2d_7      ! latent heat fluxes and evaporation
392      USE wrk_nemo, ONLY:   Cd     => wrk_2d_8           ! transfer coefficient for momentum      (tau)
393      USE wrk_nemo, ONLY:   Ch     => wrk_2d_9           ! transfer coefficient for sensible heat (Q_sens)
394      USE wrk_nemo, ONLY:   Ce     => wrk_2d_10          ! transfer coefficient for evaporation   (Q_lat)
395      USE wrk_nemo, ONLY:   zst    => wrk_2d_11          ! surface temperature in Kelvin
396      USE wrk_nemo, ONLY:   zt_zu  => wrk_2d_12          ! air temperature at wind speed height
397      USE wrk_nemo, ONLY:   zq_zu  => wrk_2d_13          ! air spec. hum.  at wind speed height
398      !
399      TYPE(fld), INTENT(in), DIMENSION(:)   ::   sf    ! input data
400      REAL(wp) , INTENT(in), DIMENSION(:,:) ::   pst   ! surface temperature                      [Celcius]
401      REAL(wp) , INTENT(in), DIMENSION(:,:) ::   pu    ! surface current at U-point (i-component) [m/s]
402      REAL(wp) , INTENT(in), DIMENSION(:,:) ::   pv    ! surface current at V-point (j-component) [m/s]
403      !
404      INTEGER  ::   ji, jj               ! dummy loop indices
405      REAL(wp) ::   zcoef_qsatw, zztmp   ! local variable
406      !!---------------------------------------------------------------------
407
408      IF( wrk_in_use(2, 1,2,3,4,5,6,7,8,9,10,11,12,13) ) THEN
409         CALL ctl_stop('blk_oce_core: requested workspace arrays unavailable')   ;   RETURN
410      ENDIF
411      !
412      ! local scalars ( place there for vector optimisation purposes)
413      zcoef_qsatw = 0.98 * 640380. / rhoa
414     
415      zst(:,:) = pst(:,:) + rt0      ! converte Celcius to Kelvin (and set minimum value far above 0 K)
416
417      ! ----------------------------------------------------------------------------- !
418      !      0   Wind components and module at T-point relative to the moving ocean   !
419      ! ----------------------------------------------------------------------------- !
420
421      ! ... components ( U10m - U_oce ) at T-point (unmasked)
422      zwnd_i(:,:) = 0.e0 
423      zwnd_j(:,:) = 0.e0
424#if defined key_vectopt_loop
425!CDIR COLLAPSE
426#endif
427      DO jj = 2, jpjm1
428         DO ji = fs_2, fs_jpim1   ! vect. opt.
429            zwnd_i(ji,jj) = (  sf(jp_wndi)%fnow(ji,jj,1) - 0.5 * ( pu(ji-1,jj  ) + pu(ji,jj) )  )
430            zwnd_j(ji,jj) = (  sf(jp_wndj)%fnow(ji,jj,1) - 0.5 * ( pv(ji  ,jj-1) + pv(ji,jj) )  )
431         END DO
432      END DO
433      CALL lbc_lnk( zwnd_i(:,:) , 'T', -1. )
434      CALL lbc_lnk( zwnd_j(:,:) , 'T', -1. )
435      ! ... scalar wind ( = | U10m - U_oce | ) at T-point (masked)
436!CDIR NOVERRCHK
437!CDIR COLLAPSE
438      wndm(:,:) = SQRT(  zwnd_i(:,:) * zwnd_i(:,:)   &
439         &             + zwnd_j(:,:) * zwnd_j(:,:)  ) * tmask(:,:,1)
440
441      ! ----------------------------------------------------------------------------- !
442      !      I   Radiative FLUXES                                                     !
443      ! ----------------------------------------------------------------------------- !
444   
445      ! ocean albedo assumed to be constant + modify now Qsr to include the diurnal cycle                    ! Short Wave
446      zztmp = 1. - albo
447      IF( ln_dm2dc ) THEN   ;   qsr(:,:) = zztmp * sbc_dcy( sf(jp_qsr)%fnow(:,:,1) ) * tmask(:,:,1)
448      ELSE                  ;   qsr(:,:) = zztmp *          sf(jp_qsr)%fnow(:,:,1)   * tmask(:,:,1)
449      ENDIF
450!CDIR COLLAPSE
451      zqlw(:,:) = (  sf(jp_qlw)%fnow(:,:,1) - Stef * zst(:,:)*zst(:,:)*zst(:,:)*zst(:,:)  ) * tmask(:,:,1)   ! Long  Wave
452      ! ----------------------------------------------------------------------------- !
453      !     II    Turbulent FLUXES                                                    !
454      ! ----------------------------------------------------------------------------- !
455
456      ! ... specific humidity at SST and IST
457!CDIR NOVERRCHK
458!CDIR COLLAPSE
459      zqsatw(:,:) = zcoef_qsatw * EXP( -5107.4 / zst(:,:) ) 
460
461      ! ... NCAR Bulk formulae, computation of Cd, Ch, Ce at T-point :
462      IF( ln_2m ) THEN
463         !! If air temp. and spec. hum. are given at different height (2m) than wind (10m) :
464         CALL TURB_CORE_2Z(2.,10., zst   , sf(jp_tair)%fnow,         &
465            &                      zqsatw, sf(jp_humi)%fnow, wndm,   &
466            &                      Cd    , Ch              , Ce  ,   &
467            &                      zt_zu , zq_zu                   )
468      ELSE
469         !! If air temp. and spec. hum. are given at same height than wind (10m) :
470!gm bug?  at the compiling phase, add a copy in temporary arrays...  ==> check perf
471!         CALL TURB_CORE_1Z( 10., zst   (:,:), sf(jp_tair)%fnow(:,:),              &
472!            &                    zqsatw(:,:), sf(jp_humi)%fnow(:,:), wndm(:,:),   &
473!            &                    Cd    (:,:),             Ch  (:,:), Ce  (:,:)  )
474!gm bug
475! ARPDBG - this won't compile with gfortran. Fix but check performance
476! as per comment above.
477         CALL TURB_CORE_1Z( 10., zst   , sf(jp_tair)%fnow(:,:,1),       &
478            &                    zqsatw, sf(jp_humi)%fnow(:,:,1), wndm, &
479            &                    Cd    , Ch              , Ce    )
480      ENDIF
481
482      ! ... tau module, i and j component
483      DO jj = 1, jpj
484         DO ji = 1, jpi
485            zztmp = rhoa * wndm(ji,jj) * Cd(ji,jj)
486            taum  (ji,jj) = zztmp * wndm  (ji,jj)
487            zwnd_i(ji,jj) = zztmp * zwnd_i(ji,jj)
488            zwnd_j(ji,jj) = zztmp * zwnd_j(ji,jj)
489         END DO
490      END DO
491
492      ! ... add the HF tau contribution to the wind stress module?
493      IF( lhftau ) THEN 
494!CDIR COLLAPSE
495#if defined key_orca_r025
496         ! Changed!!! Multiply by QSCAT correction
497           zwnd_i(:,:) = zwnd_i(:,:) * sf(jp_tdif)%fnow(:,:,1)
498           zwnd_j(:,:) = zwnd_j(:,:) * sf(jp_tdif)%fnow(:,:,1)
499#endif
500         taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:,1)
501      ENDIF
502      CALL iom_put( "taum_oce", taum )   ! output wind stress module
503
504      ! ... utau, vtau at U- and V_points, resp.
505      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines
506      DO jj = 1, jpjm1
507         DO ji = 1, fs_jpim1
508            utau(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zwnd_i(ji,jj) + zwnd_i(ji+1,jj  ) )
509            vtau(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zwnd_j(ji,jj) + zwnd_j(ji  ,jj+1) )
510         END DO
511      END DO
512      CALL lbc_lnk( utau(:,:), 'U', -1. )
513      CALL lbc_lnk( vtau(:,:), 'V', -1. )
514
515      !  Turbulent fluxes over ocean
516      ! -----------------------------
517      IF( ln_2m ) THEN
518         ! Values of temp. and hum. adjusted to 10m must be used instead of 2m values
519         zevap(:,:) = MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - zq_zu(:,:) ) * wndm(:,:) )   ! Evaporation
520         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - zt_zu(:,:) ) * wndm(:,:)     ! Sensible Heat
521      ELSE
522!CDIR COLLAPSE
523         zevap(:,:) = MAX( 0.e0, rhoa    *Ce(:,:)*( zqsatw(:,:) - sf(jp_humi)%fnow(:,:,1) ) * wndm(:,:) )   ! Evaporation
524!CDIR COLLAPSE
525         zqsb (:,:) =            rhoa*cpa*Ch(:,:)*( zst   (:,:) - sf(jp_tair)%fnow(:,:,1) ) * wndm(:,:)     ! Sensible Heat
526      ENDIF
527!CDIR COLLAPSE
528      zqla (:,:) = Lv * zevap(:,:)                                                              ! Latent Heat
529
530      IF(ln_ctl) THEN
531         CALL prt_ctl( tab2d_1=zqla  , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=Ce , clinfo2=' Ce  : ' )
532         CALL prt_ctl( tab2d_1=zqsb  , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=Ch , clinfo2=' Ch  : ' )
533         CALL prt_ctl( tab2d_1=zqlw  , clinfo1=' blk_oce_core: zqlw   : ', tab2d_2=qsr, clinfo2=' qsr : ' )
534         CALL prt_ctl( tab2d_1=zqsatw, clinfo1=' blk_oce_core: zqsatw : ', tab2d_2=zst, clinfo2=' zst : ' )
535         CALL prt_ctl( tab2d_1=utau  , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
536            &          tab2d_2=vtau  , clinfo2=              ' vtau : '  , mask2=vmask )
537         CALL prt_ctl( tab2d_1=wndm  , clinfo1=' blk_oce_core: wndm   : ')
538         CALL prt_ctl( tab2d_1=zst   , clinfo1=' blk_oce_core: zst    : ')
539      ENDIF
540       
541      ! ----------------------------------------------------------------------------- !
542      !     III    Total FLUXES                                                       !
543      ! ----------------------------------------------------------------------------- !
544     
545!CDIR COLLAPSE
546      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)      ! Downward Non Solar flux
547!CDIR COLLAPSE
548      emp(:,:) = zevap(:,:) - sf(jp_prec)%fnow(:,:,1) * rn_pfac * tmask(:,:,1)
549!CDIR COLLAPSE
550      emps(:,:) = emp(:,:)
551      !
552      CALL iom_put( "qlw_oce",   zqlw )                 ! output downward longwave heat over the ocean
553      CALL iom_put( "qsb_oce", - zqsb )                 ! output downward sensible heat over the ocean
554      CALL iom_put( "qla_oce", - zqla )                 ! output downward latent   heat over the ocean
555      CALL iom_put( "qns_oce",   qns  )                 ! output downward non solar heat over the ocean
556      !
557      IF(ln_ctl) THEN
558         CALL prt_ctl(tab2d_1=zqsb , clinfo1=' blk_oce_core: zqsb   : ', tab2d_2=zqlw , clinfo2=' zqlw  : ')
559         CALL prt_ctl(tab2d_1=zqla , clinfo1=' blk_oce_core: zqla   : ', tab2d_2=qsr  , clinfo2=' qsr   : ')
560         CALL prt_ctl(tab2d_1=pst  , clinfo1=' blk_oce_core: pst    : ', tab2d_2=emp  , clinfo2=' emp   : ')
561         CALL prt_ctl(tab2d_1=utau , clinfo1=' blk_oce_core: utau   : ', mask1=umask,   &
562            &         tab2d_2=vtau , clinfo2=              ' vtau  : ' , mask2=vmask )
563      ENDIF
564      !
565      IF( wrk_not_released(2, 1,2,3,4,5,6,7,8,9,10,11,12,13) )   &
566          CALL ctl_stop('blk_oce_core: failed to release workspace arrays')
567      !
568   END SUBROUTINE blk_oce_core
569   
570   
571   SUBROUTINE blk_ice_core(  pst   , pui   , pvi   , palb ,   &
572      &                      p_taui, p_tauj, p_qns , p_qsr,   &
573      &                      p_qla , p_dqns, p_dqla,          &
574      &                      p_tpr , p_spr ,                  &
575      &                      p_fr1 , p_fr2 , cd_grid, pdim  ) 
576      !!---------------------------------------------------------------------
577      !!                     ***  ROUTINE blk_ice_core  ***
578      !!
579      !! ** Purpose :   provide the surface boundary condition over sea-ice
580      !!
581      !! ** Method  :   compute momentum, heat and freshwater exchanged
582      !!                between atmosphere and sea-ice using CORE bulk
583      !!                formulea, ice variables and read atmmospheric fields.
584      !!                NB: ice drag coefficient is assumed to be a constant
585      !!
586      !! caution : the net upward water flux has with mm/day unit
587      !!---------------------------------------------------------------------
588      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
589      USE wrk_nemo, ONLY:   z_wnds_t => wrk_2d_1                ! wind speed ( = | U10m - U_ice | ) at T-point
590      USE wrk_nemo, ONLY:   wrk_3d_4 , wrk_3d_5 , wrk_3d_6 , wrk_3d_7
591      !!
592      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   pst      ! ice surface temperature (>0, =rt0 over land) [Kelvin]
593      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pui      ! ice surface velocity (i- and i- components      [m/s]
594      REAL(wp), DIMENSION(:,:)  , INTENT(in   ) ::   pvi      !    at I-point (B-grid) or U & V-point (C-grid)
595      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   palb     ! ice albedo (clear sky) (alb_ice_cs)               [%]
596      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_taui   ! i- & j-components of surface ice stress        [N/m2]
597      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_tauj   !   at I-point (B-grid) or U & V-point (C-grid)
598      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_qns    ! non solar heat flux over ice (T-point)         [W/m2]
599      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_qsr    !     solar heat flux over ice (T-point)         [W/m2]
600      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_qla    ! latent    heat flux over ice (T-point)         [W/m2]
601      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_dqns   ! non solar heat sensistivity  (T-point)         [W/m2]
602      REAL(wp), DIMENSION(:,:,:), INTENT(  out) ::   p_dqla   ! latent    heat sensistivity  (T-point)         [W/m2]
603      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_tpr    ! total precipitation          (T-point)      [Kg/m2/s]
604      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_spr    ! solid precipitation          (T-point)      [Kg/m2/s]
605      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_fr1    ! 1sr fraction of qsr penetration in ice (T-point)  [%]
606      REAL(wp), DIMENSION(:,:)  , INTENT(  out) ::   p_fr2    ! 2nd fraction of qsr penetration in ice (T-point)  [%]
607      CHARACTER(len=1)          , INTENT(in   ) ::   cd_grid  ! ice grid ( C or B-grid)
608      INTEGER                   , INTENT(in   ) ::   pdim     ! number of ice categories
609      !!
610      INTEGER  ::   ji, jj, jl    ! dummy loop indices
611      INTEGER  ::   ijpl          ! number of ice categories (size of 3rd dim of input arrays)
612      REAL(wp) ::   zst2, zst3
613      REAL(wp) ::   zcoef_wnorm, zcoef_wnorm2, zcoef_dqlw, zcoef_dqla, zcoef_dqsb
614      REAL(wp) ::   zztmp                                        ! temporary variable
615      REAL(wp) ::   zcoef_frca                                   ! fractional cloud amount
616      REAL(wp) ::   zwnorm_f, zwndi_f , zwndj_f                  ! relative wind module and components at F-point
617      REAL(wp) ::             zwndi_t , zwndj_t                  ! relative wind components at T-point
618      !!
619      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qlw             ! long wave heat flux over ice
620      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_qsb             ! sensible  heat flux over ice
621      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqlw            ! long wave heat sensitivity over ice
622      REAL(wp), DIMENSION(:,:,:), POINTER ::   z_dqsb            ! sensible  heat sensitivity over ice
623      !!---------------------------------------------------------------------
624
625      ijpl  = pdim                            ! number of ice categories
626
627      ! Set-up access to workspace arrays
628      IF( wrk_in_use(2, 1) .OR. wrk_in_use(3, 4,5,6,7) ) THEN
629         CALL ctl_stop('blk_ice_core: requested workspace arrays unavailable')   ;   RETURN
630      ELSE IF(ijpl > jpk) THEN
631         CALL ctl_stop('blk_ice_core: no. of ice categories > jpk so wrk_nemo 3D workspaces cannot be used.')
632         RETURN
633      END IF
634      ! Set-up pointers to sub-arrays of workspaces
635      z_qlw  => wrk_3d_4(:,:,1:ijpl)
636      z_qsb  => wrk_3d_5(:,:,1:ijpl)
637      z_dqlw => wrk_3d_6(:,:,1:ijpl)
638      z_dqsb => wrk_3d_7(:,:,1:ijpl)
639
640      ! local scalars ( place there for vector optimisation purposes)
641      zcoef_wnorm  = rhoa * Cice
642      zcoef_wnorm2 = rhoa * Cice * 0.5
643      zcoef_dqlw   = 4.0 * 0.95 * Stef
644      zcoef_dqla   = -Ls * Cice * 11637800. * (-5897.8)
645      zcoef_dqsb   = rhoa * cpa * Cice
646      zcoef_frca   = 1.0  - 0.3
647
648!!gm brutal....
649      z_wnds_t(:,:) = 0.e0
650      p_taui  (:,:) = 0.e0
651      p_tauj  (:,:) = 0.e0
652!!gm end
653
654#if defined key_lim3
655      tatm_ice(:,:) = sf(jp_tair)%fnow(:,:,1)   ! LIM3: make Tair available in sea-ice. WARNING allocated after call to ice_init
656#endif
657      ! ----------------------------------------------------------------------------- !
658      !    Wind components and module relative to the moving ocean ( U10m - U_ice )   !
659      ! ----------------------------------------------------------------------------- !
660      SELECT CASE( cd_grid )
661      CASE( 'I' )                  ! B-grid ice dynamics :   I-point (i.e. F-point with sea-ice indexation)
662         !                           and scalar wind at T-point ( = | U10m - U_ice | ) (masked)
663!CDIR NOVERRCHK
664         DO jj = 2, jpjm1
665            DO ji = 2, jpim1   ! B grid : NO vector opt
666               ! ... scalar wind at I-point (fld being at T-point)
667               zwndi_f = 0.25 * (  sf(jp_wndi)%fnow(ji-1,jj  ,1) + sf(jp_wndi)%fnow(ji  ,jj  ,1)   &
668                  &              + sf(jp_wndi)%fnow(ji-1,jj-1,1) + sf(jp_wndi)%fnow(ji  ,jj-1,1)  ) - pui(ji,jj)
669               zwndj_f = 0.25 * (  sf(jp_wndj)%fnow(ji-1,jj  ,1) + sf(jp_wndj)%fnow(ji  ,jj  ,1)   &
670                  &              + sf(jp_wndj)%fnow(ji-1,jj-1,1) + sf(jp_wndj)%fnow(ji  ,jj-1,1)  ) - pvi(ji,jj)
671               zwnorm_f = zcoef_wnorm * SQRT( zwndi_f * zwndi_f + zwndj_f * zwndj_f )
672               ! ... ice stress at I-point
673               p_taui(ji,jj) = zwnorm_f * zwndi_f
674               p_tauj(ji,jj) = zwnorm_f * zwndj_f
675               ! ... scalar wind at T-point (fld being at T-point)
676               zwndi_t = sf(jp_wndi)%fnow(ji,jj,1) - 0.25 * (  pui(ji,jj+1) + pui(ji+1,jj+1)   &
677                  &                                          + pui(ji,jj  ) + pui(ji+1,jj  )  )
678               zwndj_t = sf(jp_wndj)%fnow(ji,jj,1) - 0.25 * (  pvi(ji,jj+1) + pvi(ji+1,jj+1)   &
679                  &                                          + pvi(ji,jj  ) + pvi(ji+1,jj  )  )
680               z_wnds_t(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
681            END DO
682         END DO
683         CALL lbc_lnk( p_taui  , 'I', -1. )
684         CALL lbc_lnk( p_tauj  , 'I', -1. )
685         CALL lbc_lnk( z_wnds_t, 'T',  1. )
686         !
687      CASE( 'C' )                  ! C-grid ice dynamics :   U & V-points (same as ocean)
688#if defined key_vectopt_loop
689!CDIR COLLAPSE
690#endif
691         DO jj = 2, jpj
692            DO ji = fs_2, jpi   ! vect. opt.
693               zwndi_t = (  sf(jp_wndi)%fnow(ji,jj,1) - 0.5 * ( pui(ji-1,jj  ) + pui(ji,jj) )  )
694               zwndj_t = (  sf(jp_wndj)%fnow(ji,jj,1) - 0.5 * ( pvi(ji  ,jj-1) + pvi(ji,jj) )  )
695               z_wnds_t(ji,jj)  = SQRT( zwndi_t * zwndi_t + zwndj_t * zwndj_t ) * tmask(ji,jj,1)
696            END DO
697         END DO
698#if defined key_vectopt_loop
699!CDIR COLLAPSE
700#endif
701         DO jj = 2, jpjm1
702            DO ji = fs_2, fs_jpim1   ! vect. opt.
703               p_taui(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji+1,jj  ) + z_wnds_t(ji,jj) )                          &
704                  &          * ( 0.5 * (sf(jp_wndi)%fnow(ji+1,jj,1) + sf(jp_wndi)%fnow(ji,jj,1) ) - pui(ji,jj) )
705               p_tauj(ji,jj) = zcoef_wnorm2 * ( z_wnds_t(ji,jj+1  ) + z_wnds_t(ji,jj) )                          &
706                  &          * ( 0.5 * (sf(jp_wndj)%fnow(ji,jj+1,1) + sf(jp_wndj)%fnow(ji,jj,1) ) - pvi(ji,jj) )
707            END DO
708         END DO
709         CALL lbc_lnk( p_taui  , 'U', -1. )
710         CALL lbc_lnk( p_tauj  , 'V', -1. )
711         CALL lbc_lnk( z_wnds_t, 'T',  1. )
712         !
713      END SELECT
714
715      zztmp = 1. / ( 1. - albo )
716      !                                     ! ========================== !
717      DO jl = 1, ijpl                       !  Loop over ice categories  !
718         !                                  ! ========================== !
719!CDIR NOVERRCHK
720!CDIR COLLAPSE
721         DO jj = 1 , jpj
722!CDIR NOVERRCHK
723            DO ji = 1, jpi
724               ! ----------------------------!
725               !      I   Radiative FLUXES   !
726               ! ----------------------------!
727               zst2 = pst(ji,jj,jl) * pst(ji,jj,jl)
728               zst3 = pst(ji,jj,jl) * zst2
729               ! Short Wave (sw)
730               p_qsr(ji,jj,jl) = zztmp * ( 1. - palb(ji,jj,jl) ) * qsr(ji,jj)
731               ! Long  Wave (lw)
732               z_qlw(ji,jj,jl) = 0.95 * (  sf(jp_qlw)%fnow(ji,jj,1) - Stef * pst(ji,jj,jl) * zst3  ) * tmask(ji,jj,1)
733               ! lw sensitivity
734               z_dqlw(ji,jj,jl) = zcoef_dqlw * zst3                                               
735
736               ! ----------------------------!
737               !     II    Turbulent FLUXES  !
738               ! ----------------------------!
739
740               ! ... turbulent heat fluxes
741               ! Sensible Heat
742               z_qsb(ji,jj,jl) = rhoa * cpa * Cice * z_wnds_t(ji,jj) * ( pst(ji,jj,jl) - sf(jp_tair)%fnow(ji,jj,1) )
743               ! Latent Heat
744               p_qla(ji,jj,jl) = MAX( 0.e0, rhoa * Ls  * Cice * z_wnds_t(ji,jj)   &                           
745                  &                    * (  11637800. * EXP( -5897.8 / pst(ji,jj,jl) ) / rhoa - sf(jp_humi)%fnow(ji,jj,1)  ) )
746               ! Latent heat sensitivity for ice (Dqla/Dt)
747               p_dqla(ji,jj,jl) = zcoef_dqla * z_wnds_t(ji,jj) / ( zst2 ) * EXP( -5897.8 / pst(ji,jj,jl) )
748               ! Sensible heat sensitivity (Dqsb_ice/Dtn_ice)
749               z_dqsb(ji,jj,jl) = zcoef_dqsb * z_wnds_t(ji,jj)
750
751               ! ----------------------------!
752               !     III    Total FLUXES     !
753               ! ----------------------------!
754               ! Downward Non Solar flux
755               p_qns (ji,jj,jl) =     z_qlw (ji,jj,jl) - z_qsb (ji,jj,jl) - p_qla (ji,jj,jl)     
756               ! Total non solar heat flux sensitivity for ice
757               p_dqns(ji,jj,jl) = - ( z_dqlw(ji,jj,jl) + z_dqsb(ji,jj,jl) + p_dqla(ji,jj,jl) )   
758            END DO
759            !
760         END DO
761         !
762      END DO
763      !
764      !--------------------------------------------------------------------
765      ! FRACTIONs of net shortwave radiation which is not absorbed in the
766      ! thin surface layer and penetrates inside the ice cover
767      ! ( Maykut and Untersteiner, 1971 ; Ebert and Curry, 1993 )
768   
769!CDIR COLLAPSE
770      p_fr1(:,:) = ( 0.18 * ( 1.0 - zcoef_frca ) + 0.35 * zcoef_frca )
771!CDIR COLLAPSE
772      p_fr2(:,:) = ( 0.82 * ( 1.0 - zcoef_frca ) + 0.65 * zcoef_frca )
773       
774!CDIR COLLAPSE
775      p_tpr(:,:) = sf(jp_prec)%fnow(:,:,1) * rn_pfac      ! total precipitation [kg/m2/s]
776!CDIR COLLAPSE
777      p_spr(:,:) = sf(jp_snow)%fnow(:,:,1) * rn_pfac      ! solid precipitation [kg/m2/s]
778      CALL iom_put( 'snowpre', p_spr )                  ! Snow precipitation
779      !
780      IF(ln_ctl) THEN
781         CALL prt_ctl(tab3d_1=p_qla   , clinfo1=' blk_ice_core: p_qla  : ', tab3d_2=z_qsb   , clinfo2=' z_qsb  : ', kdim=ijpl)
782         CALL prt_ctl(tab3d_1=z_qlw   , clinfo1=' blk_ice_core: z_qlw  : ', tab3d_2=p_dqla  , clinfo2=' p_dqla : ', kdim=ijpl)
783         CALL prt_ctl(tab3d_1=z_dqsb  , clinfo1=' blk_ice_core: z_dqsb : ', tab3d_2=z_dqlw  , clinfo2=' z_dqlw : ', kdim=ijpl)
784         CALL prt_ctl(tab3d_1=p_dqns  , clinfo1=' blk_ice_core: p_dqns : ', tab3d_2=p_qsr   , clinfo2=' p_qsr  : ', kdim=ijpl)
785         CALL prt_ctl(tab3d_1=pst     , clinfo1=' blk_ice_core: pst    : ', tab3d_2=p_qns   , clinfo2=' p_qns  : ', kdim=ijpl)
786         CALL prt_ctl(tab2d_1=p_tpr   , clinfo1=' blk_ice_core: p_tpr  : ', tab2d_2=p_spr   , clinfo2=' p_spr  : ')
787         CALL prt_ctl(tab2d_1=p_taui  , clinfo1=' blk_ice_core: p_taui : ', tab2d_2=p_tauj  , clinfo2=' p_tauj : ')
788         CALL prt_ctl(tab2d_1=z_wnds_t, clinfo1=' blk_ice_core: z_wnds_t : ')
789      ENDIF
790
791      IF( wrk_not_released(2, 1)       .OR.   &
792          wrk_not_released(3, 4,5,6,7) )   CALL ctl_stop('blk_ice_core: failed to release workspace arrays')
793      !
794   END SUBROUTINE blk_ice_core
795 
796
797   SUBROUTINE TURB_CORE_1Z(zu, sst, T_a, q_sat, q_a,   &
798      &                        dU , Cd , Ch   , Ce   )
799      !!----------------------------------------------------------------------
800      !!                      ***  ROUTINE  turb_core  ***
801      !!
802      !! ** Purpose :   Computes turbulent transfert coefficients of surface
803      !!                fluxes according to Large & Yeager (2004)
804      !!
805      !! ** 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
806      !!      Momentum, Latent and sensible heat exchange coefficients
807      !!      Caution: this procedure should only be used in cases when air
808      !!      temperature (T_air), air specific humidity (q_air) and wind (dU)
809      !!      are provided at the same height 'zzu'!
810      !!
811      !! References :   Large & Yeager, 2004 : ???
812      !!----------------------------------------------------------------------
813      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released, iwrk_in_use, iwrk_not_released
814      USE wrk_nemo, ONLY: dU10 => wrk_2d_14        ! dU                             [m/s]
815      USE wrk_nemo, ONLY: dT => wrk_2d_15          ! air/sea temperature difference   [K]
816      USE wrk_nemo, ONLY: dq => wrk_2d_16          ! air/sea humidity difference      [K]
817      USE wrk_nemo, ONLY: Cd_n10 => wrk_2d_17      ! 10m neutral drag coefficient
818      USE wrk_nemo, ONLY: Ce_n10 => wrk_2d_18      ! 10m neutral latent coefficient
819      USE wrk_nemo, ONLY: Ch_n10 => wrk_2d_19      ! 10m neutral sensible coefficient
820      USE wrk_nemo, ONLY: sqrt_Cd_n10 => wrk_2d_20 ! root square of Cd_n10
821      USE wrk_nemo, ONLY: sqrt_Cd => wrk_2d_21     ! root square of Cd
822      USE wrk_nemo, ONLY: T_vpot => wrk_2d_22      ! virtual potential temperature    [K]
823      USE wrk_nemo, ONLY: T_star => wrk_2d_23      ! turbulent scale of tem. fluct.
824      USE wrk_nemo, ONLY: q_star => wrk_2d_24      ! turbulent humidity of temp. fluct.
825      USE wrk_nemo, ONLY: U_star => wrk_2d_25      ! turb. scale of velocity fluct.
826      USE wrk_nemo, ONLY: L => wrk_2d_26           ! Monin-Obukov length              [m]
827      USE wrk_nemo, ONLY: zeta => wrk_2d_27        ! stability parameter at height zu
828      USE wrk_nemo, ONLY: U_n10 => wrk_2d_28       ! neutral wind velocity at 10m     [m]
829      USE wrk_nemo, ONLY: xlogt  => wrk_2d_29,    xct => wrk_2d_30,   &
830                          zpsi_h => wrk_2d_31, zpsi_m => wrk_2d_32
831      USE wrk_nemo, ONLY: stab => iwrk_2d_1      ! 1st guess stability test integer
832      !
833      REAL(wp)                , INTENT(in   ) ::   zu      ! altitude of wind measurement       [m]
834      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   sst     ! sea surface temperature         [Kelvin]
835      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   T_a     ! potential air temperature       [Kelvin]
836      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   q_sat   ! sea surface specific humidity   [kg/kg]
837      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   q_a     ! specific air humidity           [kg/kg]
838      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   dU      ! wind module |U(zu)-U(0)|        [m/s]
839      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Cd      ! transfert coefficient for momentum       (tau)
840      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Ch      ! transfert coefficient for temperature (Q_sens)
841      REAL(wp), DIMENSION(:,:), INTENT(  out) ::   Ce      ! transfert coefficient for evaporation  (Q_lat)
842      !!
843      INTEGER :: j_itt
844      INTEGER , PARAMETER ::   nb_itt = 3
845      REAL(wp), PARAMETER ::   grav   = 9.8   ! gravity                       
846      REAL(wp), PARAMETER ::   kappa  = 0.4   ! von Karman s constant
847      !!----------------------------------------------------------------------
848
849      IF(  wrk_in_use(2,             14,15,16,17,18,19,        &
850                         20,21,22,23,24,25,26,27,28,29,        &
851                         30,31,32)                      .OR.   &
852          iwrk_in_use(2, 1)                               ) THEN
853         CALL ctl_stop('TURB_CORE_1Z: requested workspace arrays unavailable')   ;   RETURN
854      ENDIF
855
856      !! * Start
857      !! Air/sea differences
858      dU10 = max(0.5, dU)     ! we don't want to fall under 0.5 m/s
859      dT = T_a - sst          ! assuming that T_a is allready the potential temp. at zzu
860      dq = q_a - q_sat
861      !!   
862      !! Virtual potential temperature
863      T_vpot = T_a*(1. + 0.608*q_a)
864      !!
865      !! Neutral Drag Coefficient
866      stab    = 0.5 + sign(0.5,dT)    ! stable : stab = 1 ; unstable : stab = 0
867      Cd_n10  = 1E-3 * ( 2.7/dU10 + 0.142 + dU10/13.09 )    !   L & Y eq. (6a)
868      sqrt_Cd_n10 = sqrt(Cd_n10)
869      Ce_n10  = 1E-3 * ( 34.6 * sqrt_Cd_n10 )               !   L & Y eq. (6b)
870      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1-stab)) !   L & Y eq. (6c), (6d)
871      !!
872      !! Initializing transfert coefficients with their first guess neutral equivalents :
873      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
874
875      !! * Now starting iteration loop
876      DO j_itt=1, nb_itt
877         !! Turbulent scales :
878         U_star  = sqrt_Cd*dU10                !   L & Y eq. (7a)
879         T_star  = Ch/sqrt_Cd*dT               !   L & Y eq. (7b)
880         q_star  = Ce/sqrt_Cd*dq               !   L & Y eq. (7c)
881
882         !! Estimate the Monin-Obukov length :
883         L  = (U_star**2)/( kappa*grav*(T_star/T_vpot + q_star/(q_a + 1./0.608)) )
884
885         !! Stability parameters :
886         zeta  = zu/L ;  zeta   = sign( min(abs(zeta),10.0), zeta )
887         zpsi_h  = psi_h(zeta)
888         zpsi_m  = psi_m(zeta)
889
890         !! Shifting the wind speed to 10m and neutral stability :
891         U_n10 = dU10*1./(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)) !  L & Y eq. (9a)
892
893         !! Updating the neutral 10m transfer coefficients :
894         Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)              !  L & Y eq. (6a)
895         sqrt_Cd_n10 = sqrt(Cd_n10)
896         Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                           !  L & Y eq. (6b)
897         stab    = 0.5 + sign(0.5,zeta)
898         Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab))           !  L & Y eq. (6c), (6d)
899
900         !! Shifting the neutral  10m transfer coefficients to ( zu , zeta ) :
901         !!
902         xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10) - zpsi_m)
903         Cd  = Cd_n10/(xct*xct) ;  sqrt_Cd = sqrt(Cd)
904         !!
905         xlogt = log(zu/10.) - zpsi_h
906         !!
907         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
908         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
909         !!
910         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
911         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
912         !!
913      END DO
914      !!
915      IF( wrk_not_released(2,             14,15,16,17,18,19,          &
916         &                    20,21,22,23,24,25,26,27,28,29,          &
917         &                    30,31,32                      )   .OR.  &     
918         iwrk_not_released(2, 1)                                  )   &
919         CALL ctl_stop('TURB_CORE_1Z: failed to release workspace arrays')
920      !
921    END SUBROUTINE TURB_CORE_1Z
922
923
924    SUBROUTINE TURB_CORE_2Z(zt, zu, sst, T_zt, q_sat, q_zt, dU, Cd, Ch, Ce, T_zu, q_zu)
925      !!----------------------------------------------------------------------
926      !!                      ***  ROUTINE  turb_core  ***
927      !!
928      !! ** Purpose :   Computes turbulent transfert coefficients of surface
929      !!                fluxes according to Large & Yeager (2004).
930      !!
931      !! ** 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
932      !!      Momentum, Latent and sensible heat exchange coefficients
933      !!      Caution: this procedure should only be used in cases when air
934      !!      temperature (T_air) and air specific humidity (q_air) are at 2m
935      !!      whereas wind (dU) is at 10m.
936      !!
937      !! References :   Large & Yeager, 2004 : ???
938      !!----------------------------------------------------------------------
939      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released, iwrk_in_use, iwrk_not_released
940      USE wrk_nemo, ONLY: dU10 => wrk_2d_14        ! dU                             [m/s]
941      USE wrk_nemo, ONLY: dT => wrk_2d_15          ! air/sea temperature difference   [K]
942      USE wrk_nemo, ONLY: dq => wrk_2d_16          ! air/sea humidity difference      [K]
943      USE wrk_nemo, ONLY: Cd_n10 => wrk_2d_17      ! 10m neutral drag coefficient
944      USE wrk_nemo, ONLY: Ce_n10 => wrk_2d_18      ! 10m neutral latent coefficient
945      USE wrk_nemo, ONLY: Ch_n10 => wrk_2d_19      ! 10m neutral sensible coefficient
946      USE wrk_nemo, ONLY: sqrt_Cd_n10 => wrk_2d_20 ! root square of Cd_n10
947      USE wrk_nemo, ONLY: sqrt_Cd => wrk_2d_21     ! root square of Cd
948      USE wrk_nemo, ONLY: T_vpot => wrk_2d_22      ! virtual potential temperature    [K]
949      USE wrk_nemo, ONLY: T_star => wrk_2d_23     ! turbulent scale of tem. fluct.
950      USE wrk_nemo, ONLY: q_star => wrk_2d_24     ! turbulent humidity of temp. fluct.
951      USE wrk_nemo, ONLY: U_star => wrk_2d_25     ! turb. scale of velocity fluct.
952      USE wrk_nemo, ONLY: L => wrk_2d_26          ! Monin-Obukov length              [m]
953      USE wrk_nemo, ONLY: zeta_u => wrk_2d_27     ! stability parameter at height zu
954      USE wrk_nemo, ONLY: zeta_t => wrk_2d_28     ! stability parameter at height zt
955      USE wrk_nemo, ONLY: U_n10 => wrk_2d_29      ! neutral wind velocity at 10m     [m]
956      USE wrk_nemo, ONLY: xlogt => wrk_2d_30, xct => wrk_2d_31, zpsi_hu => wrk_2d_32, zpsi_ht => wrk_2d_33, zpsi_m => wrk_2d_34
957      USE wrk_nemo, ONLY: stab => iwrk_2d_1      ! 1st guess stability test integer
958      !!
959      REAL(wp), INTENT(in)   :: &
960         zt,      &     ! height for T_zt and q_zt                   [m]
961         zu             ! height for dU                              [m]
962      REAL(wp), INTENT(in), DIMENSION(jpi,jpj) ::  &
963         sst,      &     ! sea surface temperature              [Kelvin]
964         T_zt,     &     ! potential air temperature            [Kelvin]
965         q_sat,    &     ! sea surface specific humidity         [kg/kg]
966         q_zt,     &     ! specific air humidity                 [kg/kg]
967         dU              ! relative wind module |U(zu)-U(0)|       [m/s]
968      REAL(wp), INTENT(out), DIMENSION(jpi,jpj)  ::  &
969         Cd,       &     ! transfer coefficient for momentum         (tau)
970         Ch,       &     ! transfer coefficient for sensible heat (Q_sens)
971         Ce,       &     ! transfert coefficient for evaporation   (Q_lat)
972         T_zu,     &     ! air temp. shifted at zu                     [K]
973         q_zu            ! spec. hum.  shifted at zu               [kg/kg]
974
975      INTEGER :: j_itt
976      INTEGER, PARAMETER :: nb_itt = 3   ! number of itterations
977      REAL(wp), PARAMETER ::                        &
978         grav   = 9.8,      &  ! gravity                       
979         kappa  = 0.4          ! von Karman's constant
980      !!----------------------------------------------------------------------
981      !!  * Start
982
983      IF(  wrk_in_use(2,             14,15,16,17,18,19,        &
984                         20,21,22,23,24,25,26,27,28,29,        &         
985                         30,31,32,33,34)                .OR.   &
986          iwrk_in_use(2, 1)                               ) THEN
987         CALL ctl_stop('TURB_CORE_2Z: requested workspace arrays unavailable')   ;   RETURN
988      ENDIF
989
990      !! Initial air/sea differences
991      dU10 = max(0.5, dU)      !  we don't want to fall under 0.5 m/s
992      dT = T_zt - sst 
993      dq = q_zt - q_sat
994
995      !! Neutral Drag Coefficient :
996      stab = 0.5 + sign(0.5,dT)                 ! stab = 1  if dT > 0  -> STABLE
997      Cd_n10  = 1E-3*( 2.7/dU10 + 0.142 + dU10/13.09 ) 
998      sqrt_Cd_n10 = sqrt(Cd_n10)
999      Ce_n10  = 1E-3*( 34.6 * sqrt_Cd_n10 )
1000      Ch_n10  = 1E-3*sqrt_Cd_n10*(18*stab + 32.7*(1 - stab))
1001
1002      !! Initializing transf. coeff. with their first guess neutral equivalents :
1003      Cd = Cd_n10 ;  Ce = Ce_n10 ;  Ch = Ch_n10 ;  sqrt_Cd = sqrt(Cd)
1004
1005      !! Initializing z_u values with z_t values :
1006      T_zu = T_zt ;  q_zu = q_zt
1007
1008      !!  * Now starting iteration loop
1009      DO j_itt=1, nb_itt
1010         dT = T_zu - sst ;  dq = q_zu - q_sat ! Updating air/sea differences
1011         T_vpot = T_zu*(1. + 0.608*q_zu)    ! Updating virtual potential temperature at zu
1012         U_star = sqrt_Cd*dU10                ! Updating turbulent scales :   (L & Y eq. (7))
1013         T_star  = Ch/sqrt_Cd*dT              !
1014         q_star  = Ce/sqrt_Cd*dq              !
1015         !!
1016         L = (U_star*U_star) &                ! Estimate the Monin-Obukov length at height zu
1017              & / (kappa*grav/T_vpot*(T_star*(1.+0.608*q_zu) + 0.608*T_zu*q_star))
1018         !! Stability parameters :
1019         zeta_u  = zu/L  ;  zeta_u = sign( min(abs(zeta_u),10.0), zeta_u )
1020         zeta_t  = zt/L  ;  zeta_t = sign( min(abs(zeta_t),10.0), zeta_t )
1021         zpsi_hu = psi_h(zeta_u)
1022         zpsi_ht = psi_h(zeta_t)
1023         zpsi_m  = psi_m(zeta_u)
1024         !!
1025         !! Shifting the wind speed to 10m and neutral stability : (L & Y eq.(9a))
1026!        U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u)))
1027         U_n10 = dU10/(1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m))
1028         !!
1029         !! Shifting temperature and humidity at zu :          (L & Y eq. (9b-9c))
1030!        T_zu = T_zt - T_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
1031         T_zu = T_zt - T_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
1032!        q_zu = q_zt - q_star/kappa*(log(zt/zu) + psi_h(zeta_u) - psi_h(zeta_t))
1033         q_zu = q_zt - q_star/kappa*(log(zt/zu) + zpsi_hu - zpsi_ht)
1034         !!
1035         !! q_zu cannot have a negative value : forcing 0
1036         stab = 0.5 + sign(0.5,q_zu) ;  q_zu = stab*q_zu
1037         !!
1038         !! Updating the neutral 10m transfer coefficients :
1039         Cd_n10  = 1E-3 * (2.7/U_n10 + 0.142 + U_n10/13.09)    ! L & Y eq. (6a)
1040         sqrt_Cd_n10 = sqrt(Cd_n10)
1041         Ce_n10  = 1E-3 * (34.6 * sqrt_Cd_n10)                 ! L & Y eq. (6b)
1042         stab    = 0.5 + sign(0.5,zeta_u)
1043         Ch_n10  = 1E-3*sqrt_Cd_n10*(18.*stab + 32.7*(1-stab)) ! L & Y eq. (6c-6d)
1044         !!
1045         !!
1046         !! Shifting the neutral 10m transfer coefficients to (zu,zeta_u) :
1047!        xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - psi_m(zeta_u))
1048         xct = 1. + sqrt_Cd_n10/kappa*(log(zu/10.) - zpsi_m)
1049         Cd = Cd_n10/(xct*xct) ; sqrt_Cd = sqrt(Cd)
1050         !!
1051!        xlogt = log(zu/10.) - psi_h(zeta_u)
1052         xlogt = log(zu/10.) - zpsi_hu
1053         !!
1054         xct = 1. + Ch_n10*xlogt/kappa/sqrt_Cd_n10
1055         Ch  = Ch_n10*sqrt_Cd/sqrt_Cd_n10/xct
1056         !!
1057         xct = 1. + Ce_n10*xlogt/kappa/sqrt_Cd_n10
1058         Ce  = Ce_n10*sqrt_Cd/sqrt_Cd_n10/xct
1059         !!
1060         !!
1061      END DO
1062      !!
1063     IF( wrk_not_released(2,              14,15,16,17,18,19,          &
1064         &                    20,21,22,23,24,25,26,27,28,29,          &
1065         &                    30,31,32,33,34                )   .OR.  & 
1066         iwrk_not_released(2, 1)                                  )   &
1067         CALL ctl_stop('TURB_CORE_2Z: failed to release workspace arrays')
1068      !
1069    END SUBROUTINE TURB_CORE_2Z
1070
1071
1072    FUNCTION psi_m(zta)   !! Psis, L & Y eq. (8c), (8d), (8e)
1073      !-------------------------------------------------------------------------------
1074      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
1075      USE wrk_nemo, ONLY:     X2 => wrk_2d_35
1076      USE wrk_nemo, ONLY:     X  => wrk_2d_36
1077      USE wrk_nemo, ONLY: stabit => wrk_2d_37
1078      !!
1079      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: zta
1080
1081      REAL(wp), PARAMETER :: pi = 3.141592653589793_wp
1082      REAL(wp), DIMENSION(jpi,jpj)             :: psi_m
1083      !-------------------------------------------------------------------------------
1084
1085      IF( wrk_in_use(2, 35,36,37) ) THEN
1086         CALL ctl_stop('psi_m: requested workspace arrays unavailable')   ;   RETURN
1087      ENDIF
1088
1089      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.0) ;  X  = sqrt(X2)
1090      stabit    = 0.5 + sign(0.5,zta)
1091      psi_m = -5.*zta*stabit  &                                                          ! Stable
1092         &    + (1. - stabit)*(2*log((1. + X)/2) + log((1. + X2)/2) - 2*atan(X) + pi/2)  ! Unstable
1093
1094      IF( wrk_not_released(2, 35,36,37) )   CALL ctl_stop('psi_m: failed to release workspace arrays')
1095      !
1096    END FUNCTION psi_m
1097
1098
1099    FUNCTION psi_h( zta )    !! Psis, L & Y eq. (8c), (8d), (8e)
1100      !-------------------------------------------------------------------------------
1101      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
1102      USE wrk_nemo, ONLY:     X2 => wrk_2d_35
1103      USE wrk_nemo, ONLY:     X  => wrk_2d_36
1104      USE wrk_nemo, ONLY: stabit => wrk_2d_37
1105      !
1106      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   zta
1107      !
1108      REAL(wp), DIMENSION(jpi,jpj)             ::   psi_h
1109      !-------------------------------------------------------------------------------
1110
1111      IF( wrk_in_use(2, 35,36,37) ) THEN
1112         CALL ctl_stop('psi_h: requested workspace arrays unavailable')   ;   RETURN
1113      ENDIF
1114
1115      X2 = sqrt(abs(1. - 16.*zta))  ;  X2 = max(X2 , 1.) ;  X  = sqrt(X2)
1116      stabit    = 0.5 + sign(0.5,zta)
1117      psi_h = -5.*zta*stabit  &                                       ! Stable
1118         &    + (1. - stabit)*(2.*log( (1. + X2)/2. ))                 ! Unstable
1119
1120      IF( wrk_not_released(2, 35,36,37) )   CALL ctl_stop('psi_h: failed to release workspace arrays')
1121      !
1122    END FUNCTION psi_h
1123
1124#if defined key_orca_r025
1125    INTEGER FUNCTION sbc_blk_core_alloc()
1126      !!----------------------------------------------------------------------
1127      !!                ***  ROUTINE sbc_blk_core_alloc  ***
1128      !!----------------------------------------------------------------------
1129      ALLOCATE( area(jpi,jpj)         , zqlw(jpi,jpj)      ,     &
1130         &      zqsb(jpi,jpj)         , zqla(jpi,jpj)      ,     &
1131         &      zevap(jpi,jpj)        , STAT=sbc_blk_core_alloc )
1132      zqlw=0._wp
1133      zqsb=0._wp
1134   !
1135      IF( lk_mpp            )   CALL mpp_sum ( sbc_blk_core_alloc )
1136      IF( sbc_blk_core_alloc > 0 )   CALL ctl_warn('sbc_blk_core_alloc: allocation of arrays failed')
1137    END FUNCTION sbc_blk_core_alloc
1138#endif
1139
1140    SUBROUTINE Shapiro_1D(rla_varin,id_np, cd_overlap, rlpa_varout) !GIG
1141      !!=====================================================================
1142      !!
1143      !! Description: This function applies a 1D Shapiro filter
1144      !!              (3 points filter) horizontally to a 2D field
1145      !!              in regular grid
1146      !! Arguments :
1147      !!            rla_varin   : Input variable to filter
1148      !!            zla_mask    : Input mask variable
1149      !!            id_np       : Number of Shapiro filter iterations
1150      !!            cd_overlap  : Logical argument for periodical condition
1151      !!                          (global ocean case)
1152      !!            rlpa_varout : Output filtered variable
1153      !!
1154      !! History : 08/2009  S. CAILLEAU : from 1st version of N. FERRY
1155      !!           09/2009  C. REGNIER  : Corrections
1156      !!
1157      !!=====================================================================
1158      IMPLICIT NONE
1159      INTEGER, INTENT(IN)                       :: id_np
1160      REAL(wp), DIMENSION(jpi,jpj), INTENT(IN)  :: rla_varin !GIG
1161      CHARACTER(len=20), INTENT(IN)             :: cd_overlap !GIG
1162      REAL(wp), DIMENSION(jpi,jpj), INTENT(OUT) :: rlpa_varout !GIG
1163
1164      REAL(wp), DIMENSION(jpi,jpj)              :: rlpa_varout_tmp
1165      REAL, PARAMETER                           :: rl_alpha = 1./2.    ! fixed stability coefficient (isotrope case)
1166      REAL, parameter                           :: rap_aniso_diff_XY=2.25 ! anisotrope case
1167      REAL                                      :: alphax,alphay, znum, zden,test
1168      INTEGER                                   :: ji, jj, jn, nn
1169!
1170!! rap_aniso_diff_XY=2.25 : valeur trouvée empiriquement pour 140 itération pour le filtre de Shapiro et
1171!! pour un rapport d'anisotopie de 1.5 : on filtre de plus rapidement en x qu'eny.
1172!------------------------------------------------------------------------------
1173!
1174! Loop on several filter iterations
1175
1176!     Global ocean case
1177      IF (( cd_overlap == 'MERCA_GLOB' )   .OR.   &
1178          ( cd_overlap == 'REGULAR_GLOB' ) .OR.   &
1179          ( cd_overlap == 'ORCA_GLOB' )) THEN
1180             rlpa_varout(:,:) = rla_varin(:,:)
1181             rlpa_varout_tmp(:,:) = rlpa_varout(:,:)
1182!
1183
1184       alphax=1./2.
1185       alphay=1./2.
1186!  Dx/Dy=rap_aniso_diff_XY  , D_ = vitesse de diffusion
1187!  140 passes du fitre, Lx/Ly=1.5, le rap_aniso_diff_XY correspondant est:
1188       IF ( rap_aniso_diff_XY .GE. 1. ) alphay=alphay/rap_aniso_diff_XY
1189       IF ( rap_aniso_diff_XY .LT. 1. ) alphax=alphax*rap_aniso_diff_XY
1190
1191        DO jn = 1,id_np   ! number of passes of the filter
1192            DO ji = 2,jpim1
1193               DO jj = 2,jpjm1
1194                  ! We crop on the coast
1195                   znum = rlpa_varout_tmp(ji,jj)   &
1196                          + 0.25*alphax*(rlpa_varout_tmp(ji-1,jj  )-rlpa_varout_tmp(ji,jj))*tmask(ji-1,jj  ,1)  &
1197                          + 0.25*alphax*(rlpa_varout_tmp(ji+1,jj  )-rlpa_varout_tmp(ji,jj))*tmask(ji+1,jj  ,1)  &
1198                          + 0.25*alphay*(rlpa_varout_tmp(ji  ,jj-1)-rlpa_varout_tmp(ji,jj))*tmask(ji  ,jj-1,1)  &
1199                          + 0.25*alphay*(rlpa_varout_tmp(ji  ,jj+1)-rlpa_varout_tmp(ji,jj))*tmask(ji  ,jj+1,1)
1200                   rlpa_varout(ji,jj)=znum*tmask(ji,jj,1)+rla_varin(ji,jj)*(1.-tmask(ji,jj,1))
1201                ENDDO  ! end loop ji
1202            ENDDO  ! end loop jj
1203!
1204!
1205!           Periodical condition in case of cd_overlap (global ocean)
1206!           - on a mercator projection grid we consider that singular point at poles
1207!             are a mean of the values at points of the previous latitude
1208!           - on ORCA and regular grid we copy the values at points of the previous latitude
1209            IF ( cd_overlap == 'MERCAT_GLOB' ) THEN
1210!GIG case unchecked
1211               rlpa_varout(1,1) = SUM(rlpa_varout(:,2)) / jpi
1212               rlpa_varout(jpi,jpj) = SUM(rlpa_varout(:,jpj-1)) / jpi
1213            ELSE
1214               call lbc_lnk(rlpa_varout, 'T', 1.) ! Boundary condition
1215            ENDIF
1216            rlpa_varout_tmp(:,:) = rlpa_varout(:,:)
1217         ENDDO  ! end loop jn
1218      ENDIF
1219
1220!
1221    END SUBROUTINE Shapiro_1D
1222
1223 
1224   !!======================================================================
1225END MODULE sbcblk_core
Note: See TracBrowser for help on using the repository browser.