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.
ldfslp.F90 in branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 5602

Last change on this file since 5602 was 5602, checked in by cbricaud, 9 years ago

merge change from trunk rev 5003 to 5519 ( rev where branche 3.6_stable were created )

  • Property svn:keywords set to Id
File size: 54.6 KB
Line 
1MODULE ldfslp
2   !!======================================================================
3   !!                       ***  MODULE  ldfslp  ***
4   !! Ocean physics: slopes of neutral surfaces
5   !!======================================================================
6   !! History :  OPA  ! 1994-12  (G. Madec, M. Imbard)  Original code
7   !!            8.0  ! 1997-06  (G. Madec)  optimization, lbc
8   !!            8.1  ! 1999-10  (A. Jouzeau)  NEW profile in the mixed layer
9   !!   NEMO     1.0  ! 2002-10  (G. Madec)  Free form, F90
10   !!             -   ! 2005-10  (A. Beckmann)  correction for s-coordinates
11   !!            3.3  ! 2010-10  (G. Nurser, C. Harris, G. Madec)  add Griffies operator
12   !!             -   ! 2010-11  (F. Dupond, G. Madec)  bug correction in slopes just below the ML
13   !!----------------------------------------------------------------------
14#if   defined key_ldfslp   ||   defined key_esopa
15   !!----------------------------------------------------------------------
16   !!   'key_ldfslp'                      Rotation of lateral mixing tensor
17   !!----------------------------------------------------------------------
18   !!   ldf_slp_grif  : calculates the triads of isoneutral slopes (Griffies operator)
19   !!   ldf_slp       : calculates the slopes of neutral surface   (Madec operator)
20   !!   ldf_slp_mxl   : calculates the slopes at the base of the mixed layer (Madec operator)
21   !!   ldf_slp_init  : initialization of the slopes computation
22   !!----------------------------------------------------------------------
23   USE oce            ! ocean dynamics and tracers
24   USE dom_oce        ! ocean space and time domain
25   USE ldftra_oce     ! lateral diffusion: traceur
26   USE ldfdyn_oce     ! lateral diffusion: dynamics
27   USE phycst         ! physical constants
28   USE zdfmxl         ! mixed layer depth
29   USE eosbn2         ! equation of states
30   !
31   USE in_out_manager ! I/O manager
32   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
33   USE prtctl         ! Print control
34   USE wrk_nemo       ! work arrays
35   USE timing         ! Timing
36   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
37   USE iom
38
39   IMPLICIT NONE
40   PRIVATE
41
42   PUBLIC   ldf_slp        ! routine called by step.F90
43   PUBLIC   ldf_slp_grif   ! routine called by step.F90
44   PUBLIC   ldf_slp_init   ! routine called by opa.F90
45
46   LOGICAL , PUBLIC, PARAMETER ::   lk_ldfslp = .TRUE.     !: slopes flag
47   !                                                                             !! Madec operator
48   !  Arrays allocated in ldf_slp_init() routine once we know whether we're using the Griffies or Madec operator
49   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   uslp, wslpi          !: i_slope at U- and W-points
50   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   vslp, wslpj          !: j-slope at V- and W-points
51   !                                                                !! Griffies operator
52   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   wslp2                !: wslp**2 from Griffies quarter cells
53   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi_g, triadj_g   !: skew flux  slopes relative to geopotentials
54   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi  , triadj     !: isoneutral slopes relative to model-coordinate
55
56   !                                                              !! Madec operator
57   !  Arrays allocated in ldf_slp_init() routine once we know whether we're using the Griffies or Madec operator
58   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   omlmask           ! mask of the surface mixed layer at T-pt
59   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   uslpml, wslpiml   ! i_slope at U- and W-points just below the mixed layer
60   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   vslpml, wslpjml   ! j_slope at V- and W-points just below the mixed layer
61
62   REAL(wp) ::   repsln = 1.e-25_wp       ! tiny value used as minium of di(rho), dj(rho) and dk(rho)
63
64   !! * Substitutions
65#  include "domzgr_substitute.h90"
66#  include "ldftra_substitute.h90"
67#  include "ldfeiv_substitute.h90"
68#  include "vectopt_loop_substitute.h90"
69   !!----------------------------------------------------------------------
70   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
71   !! $Id$
72   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
73   !!----------------------------------------------------------------------
74CONTAINS
75
76   SUBROUTINE ldf_slp( kt, prd, pn2 )
77      !!----------------------------------------------------------------------
78      !!                 ***  ROUTINE ldf_slp  ***
79      !!
80      !! ** Purpose :   Compute the slopes of neutral surface (slope of isopycnal
81      !!              surfaces referenced locally) (ln_traldf_iso=T).
82      !!
83      !! ** Method  :   The slope in the i-direction is computed at U- and
84      !!      W-points (uslp, wslpi) and the slope in the j-direction is
85      !!      computed at V- and W-points (vslp, wslpj).
86      !!      They are bounded by 1/100 over the whole ocean, and within the
87      !!      surface layer they are bounded by the distance to the surface
88      !!      ( slope<= depth/l  where l is the length scale of horizontal
89      !!      diffusion (here, aht=2000m2/s ==> l=20km with a typical velocity
90      !!      of 10cm/s)
91      !!        A horizontal shapiro filter is applied to the slopes
92      !!        ln_sco=T, s-coordinate, add to the previously computed slopes
93      !!      the slope of the model level surface.
94      !!        macro-tasked on horizontal slab (jk-loop)  (2, jpk-1)
95      !!      [slopes already set to zero at level 1, and to zero or the ocean
96      !!      bottom slope (ln_sco=T) at level jpk in inildf]
97      !!
98      !! ** Action : - uslp, wslpi, and vslp, wslpj, the i- and  j-slopes
99      !!               of now neutral surfaces at u-, w- and v- w-points, resp.
100      !!----------------------------------------------------------------------
101      INTEGER , INTENT(in)                   ::   kt    ! ocean time-step index
102      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   prd   ! in situ density
103      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   pn2   ! Brunt-Vaisala frequency (locally ref.)
104      !!
105      INTEGER  ::   ji , jj , jk    ! dummy loop indices
106      INTEGER  ::   ii0, ii1, iku   ! temporary integer
107      INTEGER  ::   ij0, ij1, ikv   ! temporary integer
108      REAL(wp) ::   zeps, zm1_g, zm1_2g, z1_16, zcofw ! local scalars
109      REAL(wp) ::   zci, zfi, zau, zbu, zai, zbi   !   -      -
110      REAL(wp) ::   zcj, zfj, zav, zbv, zaj, zbj   !   -      -
111      REAL(wp) ::   zck, zfk,      zbw             !   -      -
112      REAL(wp) ::   zdepv, zdepu         !   -      -
113      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zwz, zww
114      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zdzr
115      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zgru, zgrv
116      REAL(wp), POINTER, DIMENSION(:,:  ) :: zhmlpu, zhmlpv
117      !!----------------------------------------------------------------------
118      !
119      IF( nn_timing == 1 )  CALL timing_start('ldf_slp')
120      !
121      CALL wrk_alloc( jpi,jpj,jpk, zwz, zww, zdzr, zgru, zgrv )
122      CALL wrk_alloc( jpi,jpj, zhmlpu, zhmlpv )
123
124      IF ( ln_traldf_iso .OR. ln_dynldf_iso ) THEN
125     
126         zeps   =  1.e-20_wp        !==   Local constant initialization   ==!
127         z1_16  =  1.0_wp / 16._wp
128         zm1_g  = -1.0_wp / grav
129         zm1_2g = -0.5_wp / grav
130         !
131         zww(:,:,:) = 0._wp
132         zwz(:,:,:) = 0._wp
133         !
134         DO jk = 1, jpk             !==   i- & j-gradient of density   ==!
135            DO jj = 1, jpjm1
136               DO ji = 1, fs_jpim1   ! vector opt.
137                  zgru(ji,jj,jk) = umask(ji,jj,jk) * ( prd(ji+1,jj  ,jk) - prd(ji,jj,jk) )
138                  zgrv(ji,jj,jk) = vmask(ji,jj,jk) * ( prd(ji  ,jj+1,jk) - prd(ji,jj,jk) )
139               END DO
140            END DO
141         END DO
142         IF( ln_zps ) THEN                           ! partial steps correction at the bottom ocean level
143            DO jj = 1, jpjm1
144               DO ji = 1, jpim1
145                  zgru(ji,jj,mbku(ji,jj)) = gru(ji,jj)
146                  zgrv(ji,jj,mbkv(ji,jj)) = grv(ji,jj)
147               END DO
148            END DO
149         ENDIF
150         IF( ln_zps .AND. ln_isfcav ) THEN           ! partial steps correction at the bottom ocean level
151            DO jj = 1, jpjm1
152               DO ji = 1, jpim1
153               IF ( miku(ji,jj) > 1 ) zgru(ji,jj,miku(ji,jj)) = grui(ji,jj) 
154               IF ( mikv(ji,jj) > 1 ) zgrv(ji,jj,mikv(ji,jj)) = grvi(ji,jj)
155               END DO
156            END DO
157         ENDIF
158         !
159         !==   Local vertical density gradient at T-point   == !   (evaluated from N^2)
160         ! interior value
161         DO jk = 2, jpkm1
162            !                                ! zdzr = d/dz(prd)= - ( prd ) / grav * mk(pn2) -- at t point
163            !                                !   trick: tmask(ik  )  = 0   =>   all pn2   = 0   =>   zdzr = 0
164            !                                !    else  tmask(ik+1)  = 0   =>   pn2(ik+1) = 0   =>   zdzr divides by 1
165            !                                !          umask(ik+1) /= 0   =>   all pn2  /= 0   =>   zdzr divides by 2
166            !                                ! NB: 1/(tmask+1) = (1-.5*tmask)  substitute a / by a *  ==> faster
167            zdzr(:,:,jk) = zm1_g * ( prd(:,:,jk) + 1._wp )              &
168               &                 * ( pn2(:,:,jk) + pn2(:,:,jk+1) ) * ( 1._wp - 0.5_wp * tmask(:,:,jk+1) )
169         END DO
170         ! surface initialisation
171         zdzr(:,:,1) = 0._wp 
172         IF ( ln_isfcav ) THEN
173            ! if isf need to overwrite the interior value at at the first ocean point
174            DO jj = 1, jpjm1
175               DO ji = 1, jpim1
176                  zdzr(ji,jj,1:mikt(ji,jj)) = 0._wp
177               END DO
178            END DO
179         END IF
180         !
181         !                          !==   Slopes just below the mixed layer   ==!
182         CALL ldf_slp_mxl( prd, pn2, zgru, zgrv, zdzr )        ! output: uslpml, vslpml, wslpiml, wslpjml
183
184
185         ! I.  slopes at u and v point      | uslp = d/di( prd ) / d/dz( prd )
186         ! ===========================      | vslp = d/dj( prd ) / d/dz( prd )
187         !
188         IF ( ln_isfcav ) THEN
189            DO jj = 2, jpjm1
190               DO ji = fs_2, fs_jpim1   ! vector opt.
191                  IF (miku(ji,jj) .GT. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji  ,jj  ),                   5._wp)
192                  IF (miku(ji,jj) .LT. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji+1,jj  ),                   5._wp)
193                  IF (miku(ji,jj) .EQ. miku(ji+1,jj)) zhmlpu(ji,jj) = MAX(hmlpt(ji  ,jj  ), hmlpt(ji+1,jj  ), 5._wp)
194                  IF (mikv(ji,jj) .GT. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji  ,jj  ),                   5._wp)
195                  IF (mikv(ji,jj) .LT. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji  ,jj+1),                   5._wp)
196                  IF (mikv(ji,jj) .EQ. miku(ji,jj+1)) zhmlpv(ji,jj) = MAX(hmlpt(ji  ,jj  ), hmlpt(ji  ,jj+1), 5._wp)
197               ENDDO
198            ENDDO
199         ELSE
200            DO jj = 2, jpjm1
201               DO ji = fs_2, fs_jpim1   ! vector opt.
202                  zhmlpu(ji,jj) = MAX(hmlpt(ji,jj), hmlpt(ji+1,jj  ), 5._wp)
203                  zhmlpv(ji,jj) = MAX(hmlpt(ji,jj), hmlpt(ji  ,jj+1), 5._wp)
204               ENDDO
205            ENDDO
206         END IF
207         DO jk = 2, jpkm1                            !* Slopes at u and v points
208            DO jj = 2, jpjm1
209               DO ji = fs_2, fs_jpim1   ! vector opt.
210                  !                                      ! horizontal and vertical density gradient at u- and v-points
211                  zau = zgru(ji,jj,jk) / e1u(ji,jj)
212                  zav = zgrv(ji,jj,jk) / e2v(ji,jj)
213                  zbu = 0.5_wp * ( zdzr(ji,jj,jk) + zdzr(ji+1,jj  ,jk) )
214                  zbv = 0.5_wp * ( zdzr(ji,jj,jk) + zdzr(ji  ,jj+1,jk) )
215                  !                                      ! bound the slopes: abs(zw.)<= 1/100 and zb..<0
216                  !                                      ! + kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
217                  zbu = MIN(  zbu, -100._wp* ABS( zau ) , -7.e+3_wp/fse3u(ji,jj,jk)* ABS( zau )  )
218                  zbv = MIN(  zbv, -100._wp* ABS( zav ) , -7.e+3_wp/fse3v(ji,jj,jk)* ABS( zav )  )
219                  !                                      ! uslp and vslp output in zwz and zww, resp.
220                  zfi = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj  ,jk) )
221                  zfj = MAX( omlmask(ji,jj,jk), omlmask(ji  ,jj+1,jk) )
222                  ! thickness of water column between surface and level k at u/v point
223                  zdepu = 0.5_wp * ( ( fsdept(ji,jj,jk) + fsdept(ji+1,jj  ,jk) )                              &
224                                   - 2 * MAX( risfdep(ji,jj), risfdep(ji+1,jj  ) ) - fse3u(ji,jj,miku(ji,jj)) )
225                  zdepv = 0.5_wp * ( ( fsdept(ji,jj,jk) + fsdept(ji  ,jj+1,jk) )                              &
226                                   - 2 * MAX( risfdep(ji,jj), risfdep(ji  ,jj+1) ) - fse3v(ji,jj,mikv(ji,jj)) )
227                  !
228                  zwz(ji,jj,jk) = ( 1. - zfi) * zau / ( zbu - zeps )                                          &
229                     &                 + zfi  * uslpml(ji,jj) * zdepu / zhmlpu(ji,jj)
230                  zwz(ji,jj,jk) = zwz(ji,jj,jk) * wumask(ji,jj,jk)
231                  zww(ji,jj,jk) = ( 1. - zfj) * zav / ( zbv - zeps )                                          &
232                     &                 + zfj  * vslpml(ji,jj) * zdepv / zhmlpv(ji,jj) 
233                  zww(ji,jj,jk) = zww(ji,jj,jk) * wvmask(ji,jj,jk)
234                 
235                 
236!!gm  modif to suppress omlmask.... (as in Griffies case)
237!                  !                                         ! jk must be >= ML level for zf=1. otherwise  zf=0.
238!                  zfi = REAL( 1 - 1/(1 + jk / MAX( nmln(ji+1,jj), nmln(ji,jj) ) ), wp )
239!                  zfj = REAL( 1 - 1/(1 + jk / MAX( nmln(ji,jj+1), nmln(ji,jj) ) ), wp )
240!                  zci = 0.5 * ( fsdept(ji+1,jj,jk)+fsdept(ji,jj,jk) ) / MAX( hmlpt(ji,jj), hmlpt(ji+1,jj), 10. ) )
241!                  zcj = 0.5 * ( fsdept(ji,jj+1,jk)+fsdept(ji,jj,jk) ) / MAX( hmlpt(ji,jj), hmlpt(ji,jj+1), 10. ) )
242!                  zwz(ji,jj,jk) = ( zfi * zai / ( zbi - zeps ) + ( 1._wp - zfi ) * wslpiml(ji,jj) * zci ) * tmask(ji,jj,jk)
243!                  zww(ji,jj,jk) = ( zfj * zaj / ( zbj - zeps ) + ( 1._wp - zfj ) * wslpjml(ji,jj) * zcj ) * tmask(ji,jj,jk)
244!!gm end modif
245               END DO
246            END DO
247         END DO
248         CALL lbc_lnk( zwz, 'U', -1. )   ;   CALL lbc_lnk( zww, 'V', -1. )      ! lateral boundary conditions
249
250         !
251         !                                            !* horizontal Shapiro filter
252         DO jk = 2, jpkm1
253            DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only
254               DO ji = 2, jpim1
255                  uslp(ji,jj,jk) = z1_16 * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      &
256                     &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      &
257                     &                       + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)      &
258                     &                       +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )    &
259                     &                       + 4.*  zwz(ji  ,jj  ,jk)                       )
260                  vslp(ji,jj,jk) = z1_16 * (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)      &
261                     &                       +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)      &
262                     &                       + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)      &
263                     &                       +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )    &
264                     &                       + 4.*  zww(ji,jj    ,jk)                       )
265               END DO
266            END DO
267            DO jj = 3, jpj-2                               ! other rows
268               DO ji = fs_2, fs_jpim1   ! vector opt.
269                  uslp(ji,jj,jk) = z1_16 * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      &
270                     &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      &
271                     &                       + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)      &
272                     &                       +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )    &
273                     &                       + 4.*  zwz(ji  ,jj  ,jk)                       )
274                  vslp(ji,jj,jk) = z1_16 * (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)      &
275                     &                       +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)      &
276                     &                       + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)         &
277                     &                       +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )    &
278                     &                       + 4.*  zww(ji,jj    ,jk)                       )
279               END DO
280            END DO
281            !                                        !* decrease along coastal boundaries
282            DO jj = 2, jpjm1
283               DO ji = fs_2, fs_jpim1   ! vector opt.
284                  uslp(ji,jj,jk) = uslp(ji,jj,jk) * ( umask(ji,jj+1,jk) + umask(ji,jj-1,jk  ) ) * 0.5_wp   &
285                     &                            * ( umask(ji,jj  ,jk) + umask(ji,jj  ,jk+1) ) * 0.5_wp   &
286                     &                            *   umask(ji,jj,jk-1)
287                  vslp(ji,jj,jk) = vslp(ji,jj,jk) * ( vmask(ji+1,jj,jk) + vmask(ji-1,jj,jk  ) ) * 0.5_wp   &
288                     &                            * ( vmask(ji  ,jj,jk) + vmask(ji  ,jj,jk+1) ) * 0.5_wp   &
289                     &                            *   vmask(ji,jj,jk-1)
290               END DO
291            END DO
292         END DO
293
294
295         ! II.  slopes at w point           | wslpi = mij( d/di( prd ) / d/dz( prd )
296         ! ===========================      | wslpj = mij( d/dj( prd ) / d/dz( prd )
297         !
298         DO jk = 2, jpkm1
299            DO jj = 2, jpjm1
300               DO ji = fs_2, fs_jpim1   ! vector opt.
301                  !                                  !* Local vertical density gradient evaluated from N^2
302                  zbw = zm1_2g * pn2 (ji,jj,jk) * ( prd (ji,jj,jk) + prd (ji,jj,jk-1) + 2. ) * wmask(ji,jj,jk)
303                  !                                  !* Slopes at w point
304                  !                                        ! i- & j-gradient of density at w-points
305                  zci = MAX(  umask(ji-1,jj,jk  ) + umask(ji,jj,jk  )           &
306                     &      + umask(ji-1,jj,jk-1) + umask(ji,jj,jk-1) , zeps  ) * e1t(ji,jj)
307                  zcj = MAX(  vmask(ji,jj-1,jk  ) + vmask(ji,jj,jk-1)           &
308                     &      + vmask(ji,jj-1,jk-1) + vmask(ji,jj,jk  ) , zeps  ) *  e2t(ji,jj)
309                  zai =    (  zgru (ji-1,jj,jk  ) + zgru (ji,jj,jk-1)           &
310                     &      + zgru (ji-1,jj,jk-1) + zgru (ji,jj,jk  )   ) / zci
311                  zaj =    (  zgrv (ji,jj-1,jk  ) + zgrv (ji,jj,jk-1)           &
312                     &      + zgrv (ji,jj-1,jk-1) + zgrv (ji,jj,jk  )   ) / zcj
313                  !                                        ! bound the slopes: abs(zw.)<= 1/100 and zb..<0.
314                  !                                        ! + kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
315                  zbi = MIN( zbw ,- 100._wp* ABS( zai ) , -7.e+3_wp/fse3w(ji,jj,jk)* ABS( zai )  )
316                  zbj = MIN( zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/fse3w(ji,jj,jk)* ABS( zaj )  )
317                  !                                        ! wslpi and wslpj with ML flattening (output in zwz and zww, resp.)
318                  zfk = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) )   ! zfk=1 in the ML otherwise zfk=0
319                  zck = ( fsdepw(ji,jj,jk) - fsdepw(ji,jj,mikt(ji,jj) ) ) / MAX( hmlp(ji,jj), 10._wp )
320                  zwz(ji,jj,jk) = (  zai / ( zbi - zeps ) * ( 1._wp - zfk ) &
321                     &            + zck * wslpiml(ji,jj) * zfk  ) * wmask(ji,jj,jk)
322                  zww(ji,jj,jk) = (  zaj / ( zbj - zeps ) * ( 1._wp - zfk ) &
323                     &            + zck * wslpjml(ji,jj) * zfk  ) * wmask(ji,jj,jk)
324
325!!gm  modif to suppress omlmask....  (as in Griffies operator)
326!                  !                                         ! jk must be >= ML level for zfk=1. otherwise  zfk=0.
327!                  zfk = REAL( 1 - 1/(1 + jk / nmln(ji+1,jj)), wp )
328!                  zck = fsdepw(ji,jj,jk)    / MAX( hmlp(ji,jj), 10. )
329!                  zwz(ji,jj,jk) = ( zfk * zai / ( zbi - zeps ) + ( 1._wp - zfk ) * wslpiml(ji,jj) * zck ) * tmask(ji,jj,jk)
330!                  zww(ji,jj,jk) = ( zfk * zaj / ( zbj - zeps ) + ( 1._wp - zfk ) * wslpjml(ji,jj) * zck ) * tmask(ji,jj,jk)
331!!gm end modif
332               END DO
333            END DO
334         END DO
335         CALL lbc_lnk( zwz, 'T', -1. )   ;    CALL lbc_lnk( zww, 'T', -1. )      ! lateral boundary conditions
336         !
337         !                                           !* horizontal Shapiro filter
338         DO jk = 2, jpkm1
339            DO jj = 2, jpjm1, MAX(1, jpj-3)                        ! rows jj=2 and =jpjm1 only
340               DO ji = 2, jpim1
341                  zcofw = tmask(ji,jj,jk) * z1_16
342                  wslpi(ji,jj,jk) = (          zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     &
343                       &                +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     &
344                       &                + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)     &
345                       &                +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )   &
346                       &                + 4.*  zwz(ji  ,jj  ,jk)                         ) * zcofw
347
348                  wslpj(ji,jj,jk) = (          zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)     &
349                       &                +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)     &
350                       &                + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)     &
351                       &                +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )   &
352                       &                + 4.*  zww(ji  ,jj  ,jk)                         ) * zcofw
353               END DO
354            END DO
355            DO jj = 3, jpj-2                               ! other rows
356               DO ji = fs_2, fs_jpim1   ! vector opt.
357                  zcofw = tmask(ji,jj,jk) * z1_16
358                  wslpi(ji,jj,jk) = (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     &
359                       &                +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     &
360                       &                + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)     &
361                       &                +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )   &
362                       &                + 4.*  zwz(ji  ,jj  ,jk)                         ) * zcofw
363
364                  wslpj(ji,jj,jk) = (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)     &
365                       &                +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)     &
366                       &                + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)     &
367                       &                +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )   &
368                       &                + 4.*  zww(ji  ,jj  ,jk)                         ) * zcofw
369               END DO
370            END DO
371            !                                        !* decrease along coastal boundaries
372            DO jj = 2, jpjm1
373               DO ji = fs_2, fs_jpim1   ! vector opt.
374                  zck =   ( umask(ji,jj,jk) + umask(ji-1,jj,jk) )   &
375                     &  * ( vmask(ji,jj,jk) + vmask(ji,jj-1,jk) ) * 0.25
376                  wslpi(ji,jj,jk) = wslpi(ji,jj,jk) * zck * wmask(ji,jj,jk)
377                  wslpj(ji,jj,jk) = wslpj(ji,jj,jk) * zck * wmask(ji,jj,jk)
378               END DO
379            END DO
380         END DO
381
382         ! III.  Specific grid points
383         ! ===========================
384         !
385         IF( cp_cfg == "orca" .AND. jp_cfg == 4 ) THEN     !  ORCA_R4 configuration: horizontal diffusion in specific area
386            !                                                    ! Gibraltar Strait
387            ij0 =  50   ;   ij1 =  53
388            ii0 =  69   ;   ii1 =  71   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
389            ij0 =  51   ;   ij1 =  53
390            ii0 =  68   ;   ii1 =  71   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
391            ii0 =  69   ;   ii1 =  71   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
392            ii0 =  69   ;   ii1 =  71   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
393            !
394            !                                                    ! Mediterrannean Sea
395            ij0 =  49   ;   ij1 =  56
396            ii0 =  71   ;   ii1 =  90   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
397            ij0 =  50   ;   ij1 =  56
398            ii0 =  70   ;   ii1 =  90   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
399            ii0 =  71   ;   ii1 =  90   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
400            ii0 =  71   ;   ii1 =  90   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0._wp
401         ENDIF
402
403
404         ! IV. Lateral boundary conditions
405         ! ===============================
406         CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. )
407         CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. )
408
409
410         IF(ln_ctl) THEN
411            CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp  - u : ', tab3d_2=vslp,  clinfo2=' v : ', kdim=jpk)
412            CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk)
413         ENDIF
414         !
415
416      ELSEIF ( lk_vvl ) THEN
417 
418         IF(lwp) THEN
419            WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces' 
420         ENDIF 
421
422         ! geopotential diffusion in s-coordinates on tracers and/or momentum
423         ! The slopes of s-surfaces are computed at each time step due to vvl
424         ! The slopes for momentum diffusion are i- or j- averaged of those on tracers
425
426         ! set the slope of diffusion to the slope of s-surfaces
427         !      ( c a u t i o n : minus sign as fsdep has positive value )
428         DO jj = 2, jpjm1 
429            DO ji = fs_2, fs_jpim1   ! vector opt.
430               uslp(ji,jj,1) = -1./e1u(ji,jj) * ( fsdept_b(ji+1,jj,1) - fsdept_b(ji ,jj ,1) )  * umask(ji,jj,1) 
431               vslp(ji,jj,1) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,1) - fsdept_b(ji ,jj ,1) )  * vmask(ji,jj,1) 
432               wslpi(ji,jj,1) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,1) - fsdepw_b(ji-1,jj,1) ) * tmask(ji,jj,1) * 0.5 
433               wslpj(ji,jj,1) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,1) - fsdepw_b(ji,jj-1,1) ) * tmask(ji,jj,1) * 0.5 
434            END DO
435         END DO
436
437         DO jk = 2, jpk 
438            DO jj = 2, jpjm1 
439               DO ji = fs_2, fs_jpim1   ! vector opt.
440                  uslp(ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept_b(ji+1,jj,jk) - fsdept_b(ji ,jj ,jk) ) * umask(ji,jj,jk) 
441                  vslp(ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,jk) - fsdept_b(ji ,jj ,jk) ) * vmask(ji,jj,jk) 
442                  wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,jk) - fsdepw_b(ji-1,jj,jk) ) &
443                    &                              * wmask(ji,jj,jk) * 0.5 
444                  wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,jk) - fsdepw_b(ji,jj-1,jk) ) &
445                    &                              * wmask(ji,jj,jk) * 0.5 
446               END DO
447            END DO
448         END DO 
449
450         ! Lateral boundary conditions on the slopes
451         CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. ) 
452         CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. ) 
453 
454         if( kt == nit000 ) then
455            IF(lwp) WRITE(numout,*) ' max slop: u',SQRT( MAXVAL(uslp*uslp)), ' v ', SQRT(MAXVAL(vslp)),  & 
456               &                             ' wi', sqrt(MAXVAL(wslpi)), ' wj', sqrt(MAXVAL(wslpj)) 
457         endif
458 
459         IF(ln_ctl) THEN
460            CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp  - u : ', tab3d_2=vslp,  clinfo2=' v : ', kdim=jpk) 
461            CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk) 
462         ENDIF
463
464      ENDIF
465
466      CALL iom_put("zgru",zgru)
467      CALL iom_put("zgrv",zgrv)
468      CALL iom_put("zdzr",zdzr)
469      CALL iom_put("zwz",zwz)
470      CALL iom_put("zww",zww)
471      CALL iom_put("uslp",uslp)
472      CALL iom_put("vslp",vslp)
473     
474      CALL wrk_dealloc( jpi,jpj,jpk, zwz, zww, zdzr, zgru, zgrv )
475      CALL wrk_dealloc( jpi,jpj,     zhmlpu, zhmlpv)
476      !
477      IF( nn_timing == 1 )  CALL timing_stop('ldf_slp')
478      !
479   END SUBROUTINE ldf_slp
480
481
482   SUBROUTINE ldf_slp_grif ( kt )
483      !!----------------------------------------------------------------------
484      !!                 ***  ROUTINE ldf_slp_grif  ***
485      !!
486      !! ** Purpose :   Compute the squared slopes of neutral surfaces (slope
487      !!      of iso-pycnal surfaces referenced locally) (ln_traldf_grif=T)
488      !!      at W-points using the Griffies quarter-cells.
489      !!
490      !! ** Method  :   calculates alpha and beta at T-points
491      !!
492      !! ** Action : - triadi_g, triadj_g   T-pts i- and j-slope triads relative to geopot. (used for eiv)
493      !!             - triadi , triadj    T-pts i- and j-slope triads relative to model-coordinate
494      !!             - wslp2              squared slope of neutral surfaces at w-points.
495      !!----------------------------------------------------------------------
496      INTEGER, INTENT( in ) ::   kt             ! ocean time-step index
497      !!
498      INTEGER  ::   ji, jj, jk, jl, ip, jp, kp  ! dummy loop indices
499      INTEGER  ::   iku, ikv                    ! local integer
500      REAL(wp) ::   zfacti, zfactj              ! local scalars
501      REAL(wp) ::   znot_thru_surface           ! local scalars
502      REAL(wp) ::   zdit, zdis, zdjt, zdjs, zdkt, zdks, zbu, zbv, zbti, zbtj
503      REAL(wp) ::   zdxrho_raw, zti_coord, zti_raw, zti_lim, zti_g_raw, zti_g_lim
504      REAL(wp) ::   zdyrho_raw, ztj_coord, ztj_raw, ztj_lim, ztj_g_raw, ztj_g_lim
505      REAL(wp) ::   zdzrho_raw
506      REAL(wp), POINTER, DIMENSION(:,:)     ::   z1_mlbw
507      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zdxrho , zdyrho, zdzrho     ! Horizontal and vertical density gradients
508      REAL(wp), POINTER, DIMENSION(:,:,:,:) ::   zti_mlb, ztj_mlb            ! for Griffies operator only
509      !!----------------------------------------------------------------------
510      !
511      IF( nn_timing == 1 )  CALL timing_start('ldf_slp_grif')
512      !
513      CALL wrk_alloc( jpi,jpj, z1_mlbw )
514      CALL wrk_alloc( jpi,jpj,jpk,2, zdxrho , zdyrho, zdzrho,              klstart = 0  )
515      CALL wrk_alloc( jpi,jpj,  2,2, zti_mlb, ztj_mlb,        kkstart = 0, klstart = 0  )
516      !
517      !--------------------------------!
518      !  Some preliminary calculation  !
519      !--------------------------------!
520      !
521      DO jl = 0, 1                            !==  unmasked before density i- j-, k-gradients  ==!
522         !
523         ip = jl   ;   jp = jl                ! guaranteed nonzero gradients ( absolute value larger than repsln)
524         DO jk = 1, jpkm1                     ! done each pair of triad
525            DO jj = 1, jpjm1                  ! NB: not masked ==>  a minimum value is set
526               DO ji = 1, fs_jpim1            ! vector opt.
527                  zdit = ( tsb(ji+1,jj,jk,jp_tem) - tsb(ji,jj,jk,jp_tem) )    ! i-gradient of T & S at u-point
528                  zdis = ( tsb(ji+1,jj,jk,jp_sal) - tsb(ji,jj,jk,jp_sal) )
529                  zdjt = ( tsb(ji,jj+1,jk,jp_tem) - tsb(ji,jj,jk,jp_tem) )    ! j-gradient of T & S at v-point
530                  zdjs = ( tsb(ji,jj+1,jk,jp_sal) - tsb(ji,jj,jk,jp_sal) )
531                  zdxrho_raw = ( - rab_b(ji+ip,jj   ,jk,jp_tem) * zdit + rab_b(ji+ip,jj   ,jk,jp_sal) * zdis ) / e1u(ji,jj)
532                  zdyrho_raw = ( - rab_b(ji   ,jj+jp,jk,jp_tem) * zdjt + rab_b(ji   ,jj+jp,jk,jp_sal) * zdjs ) / e2v(ji,jj)
533                  zdxrho(ji+ip,jj   ,jk,1-ip) = SIGN( MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw )   ! keep the sign
534                  zdyrho(ji   ,jj+jp,jk,1-jp) = SIGN( MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw )
535               END DO
536            END DO
537         END DO
538         !
539         IF( ln_zps .AND. l_grad_zps ) THEN     ! partial steps: correction of i- & j-grad on bottom
540            DO jj = 1, jpjm1
541               DO ji = 1, jpim1
542                  iku  = mbku(ji,jj)          ;   ikv  = mbkv(ji,jj)             ! last ocean level (u- & v-points)
543                  zdit = gtsu(ji,jj,jp_tem)   ;   zdjt = gtsv(ji,jj,jp_tem)      ! i- & j-gradient of Temperature
544                  zdis = gtsu(ji,jj,jp_sal)   ;   zdjs = gtsv(ji,jj,jp_sal)      ! i- & j-gradient of Salinity
545                  zdxrho_raw = ( - rab_b(ji+ip,jj   ,iku,jp_tem) * zdit + rab_b(ji+ip,jj   ,iku,jp_sal) * zdis ) / e1u(ji,jj)
546                  zdyrho_raw = ( - rab_b(ji   ,jj+jp,ikv,jp_tem) * zdjt + rab_b(ji   ,jj+jp,ikv,jp_sal) * zdjs ) / e2v(ji,jj)
547                  zdxrho(ji+ip,jj   ,iku,1-ip) = SIGN( MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw )   ! keep the sign
548                  zdyrho(ji   ,jj+jp,ikv,1-jp) = SIGN( MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw )
549               END DO
550            END DO
551         ENDIF
552         !
553      END DO
554
555      DO kp = 0, 1                            !==  unmasked before density i- j-, k-gradients  ==!
556         DO jk = 1, jpkm1                     ! done each pair of triad
557            DO jj = 1, jpj                    ! NB: not masked ==>  a minimum value is set
558               DO ji = 1, jpi                 ! vector opt.
559                  IF( jk+kp > 1 ) THEN        ! k-gradient of T & S a jk+kp
560                     zdkt = ( tsb(ji,jj,jk+kp-1,jp_tem) - tsb(ji,jj,jk+kp,jp_tem) )
561                     zdks = ( tsb(ji,jj,jk+kp-1,jp_sal) - tsb(ji,jj,jk+kp,jp_sal) )
562                  ELSE
563                     zdkt = 0._wp                                             ! 1st level gradient set to zero
564                     zdks = 0._wp
565                  ENDIF
566                  zdzrho_raw = ( - rab_b(ji,jj,jk,jp_tem) * zdkt + rab_b(ji,jj,jk,jp_sal) * zdks ) / fse3w(ji,jj,jk+kp)
567                  zdzrho(ji,jj,jk,kp) = - MIN( - repsln, zdzrho_raw )    ! force zdzrho >= repsln
568                 END DO
569            END DO
570         END DO
571      END DO
572      !
573      DO jj = 1, jpj                          !==  Reciprocal depth of the w-point below ML base  ==!
574         DO ji = 1, jpi
575            jk = MIN( nmln(ji,jj), mbkt(ji,jj) ) + 1     ! MIN in case ML depth is the ocean depth
576            z1_mlbw(ji,jj) = 1._wp / fsdepw(ji,jj,jk)
577         END DO
578      END DO
579      !
580      !                                       !==  intialisations to zero  ==!
581      !
582      wslp2  (:,:,:)     = 0._wp              ! wslp2 will be cumulated 3D field set to zero
583      triadi_g(:,:,1,:,:) = 0._wp   ;   triadi_g(:,:,jpk,:,:) = 0._wp   ! set surface and bottom slope to zero
584      triadj_g(:,:,1,:,:) = 0._wp   ;   triadj_g(:,:,jpk,:,:) = 0._wp
585      !!gm _iso set to zero missing
586      triadi  (:,:,1,:,:) = 0._wp   ;   triadj  (:,:,jpk,:,:) = 0._wp   ! set surface and bottom slope to zero
587      triadj  (:,:,1,:,:) = 0._wp   ;   triadj  (:,:,jpk,:,:) = 0._wp
588
589      !-------------------------------------!
590      !  Triads just below the Mixed Layer  !
591      !-------------------------------------!
592      !
593      DO jl = 0, 1                            ! calculate slope of the 4 triads immediately ONE level below mixed-layer base
594         DO kp = 0, 1                         ! with only the slope-max limit   and   MASKED
595            DO jj = 1, jpjm1
596               DO ji = 1, fs_jpim1
597                  ip = jl   ;   jp = jl
598                  !
599                  jk = nmln(ji+ip,jj) + 1
600                  IF( jk .GT. mbkt(ji+ip,jj) ) THEN  !ML reaches bottom
601                    zti_mlb(ji+ip,jj   ,1-ip,kp) = 0.0_wp
602                  ELSE
603                    ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth)
604                    zti_g_raw = (  zdxrho(ji+ip,jj,jk-kp,1-ip) / zdzrho(ji+ip,jj,jk-kp,kp)      &
605                       &      - ( fsdept(ji+1,jj,jk-kp) - fsdept(ji,jj,jk-kp) ) / e1u(ji,jj)  ) * umask(ji,jj,jk)
606                    zti_mlb(ji+ip,jj   ,1-ip,kp) = SIGN( MIN( rn_slpmax, ABS( zti_g_raw ) ), zti_g_raw )
607                  ENDIF
608                  !
609                  jk = nmln(ji,jj+jp) + 1
610                  IF( jk .GT. mbkt(ji,jj+jp) ) THEN  !ML reaches bottom
611                    ztj_mlb(ji   ,jj+jp,1-jp,kp) = 0.0_wp
612                  ELSE
613                    ztj_g_raw = (  zdyrho(ji,jj+jp,jk-kp,1-jp) / zdzrho(ji,jj+jp,jk-kp,kp)      &
614                       &      - ( fsdept(ji,jj+1,jk-kp) - fsdept(ji,jj,jk-kp) ) / e2v(ji,jj)  ) * vmask(ji,jj,jk)
615                    ztj_mlb(ji   ,jj+jp,1-jp,kp) = SIGN( MIN( rn_slpmax, ABS( ztj_g_raw ) ), ztj_g_raw )
616                  ENDIF
617               END DO
618            END DO
619         END DO
620      END DO
621
622      !-------------------------------------!
623      !  Triads with surface limits         !
624      !-------------------------------------!
625      !
626      DO kp = 0, 1                            ! k-index of triads
627         DO jl = 0, 1
628            ip = jl   ;   jp = jl             ! i- and j-indices of triads (i-k and j-k planes)
629            DO jk = 1, jpkm1
630               ! Must mask contribution to slope from dz/dx at constant s for triads jk=1,kp=0 that poke up though ocean surface
631               znot_thru_surface = REAL( 1-1/(jk+kp), wp )  !jk+kp=1,=0.; otherwise=1.0
632               DO jj = 1, jpjm1
633                  DO ji = 1, fs_jpim1         ! vector opt.
634                     !
635                     ! Calculate slope relative to geopotentials used for GM skew fluxes
636                     ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth)
637                     ! Limit by slope *relative to geopotentials* by rn_slpmax, and mask by psi-point
638                     ! masked by umask taken at the level of dz(rho)
639                     !
640                     ! raw slopes: unmasked unbounded slopes (relative to geopotential (zti_g) and model surface (zti)
641                     !
642                     zti_raw   = zdxrho(ji+ip,jj   ,jk,1-ip) / zdzrho(ji+ip,jj   ,jk,kp)                   ! unmasked
643                     ztj_raw   = zdyrho(ji   ,jj+jp,jk,1-jp) / zdzrho(ji   ,jj+jp,jk,kp)
644
645                     ! Must mask contribution to slope for triad jk=1,kp=0 that poke up though ocean surface
646                     zti_coord = znot_thru_surface * ( fsdept(ji+1,jj  ,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
647                     ztj_coord = znot_thru_surface * ( fsdept(ji  ,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)                  ! unmasked
648                     zti_g_raw = zti_raw - zti_coord      ! ref to geopot surfaces
649                     ztj_g_raw = ztj_raw - ztj_coord
650                     zti_g_lim = SIGN( MIN( rn_slpmax, ABS( zti_g_raw ) ), zti_g_raw )
651                     ztj_g_lim = SIGN( MIN( rn_slpmax, ABS( ztj_g_raw ) ), ztj_g_raw )
652                     !
653                     ! Below  ML use limited zti_g as is & mask
654                     ! Inside ML replace by linearly reducing sx_mlb towards surface & mask
655                     !
656                     zfacti = REAL( 1 - 1/(1 + (jk+kp-1)/nmln(ji+ip,jj)), wp )  ! k index of uppermost point(s) of triad is jk+kp-1
657                     zfactj = REAL( 1 - 1/(1 + (jk+kp-1)/nmln(ji,jj+jp)), wp )  ! must be .ge. nmln(ji,jj) for zfact=1
658                     !                                                          !                   otherwise  zfact=0
659                     zti_g_lim =          ( zfacti   * zti_g_lim                       &
660                        &      + ( 1._wp - zfacti ) * zti_mlb(ji+ip,jj,1-ip,kp)   &
661                        &                           * fsdepw(ji+ip,jj,jk+kp) * z1_mlbw(ji+ip,jj) ) * umask(ji,jj,jk+kp)
662                     ztj_g_lim =          ( zfactj   * ztj_g_lim                       &
663                        &      + ( 1._wp - zfactj ) * ztj_mlb(ji,jj+jp,1-jp,kp)   &
664                        &                           * fsdepw(ji,jj+jp,jk+kp) * z1_mlbw(ji,jj+jp) ) * vmask(ji,jj,jk+kp)
665                     !
666                     triadi_g(ji+ip,jj   ,jk,1-ip,kp) = zti_g_lim
667                     triadj_g(ji   ,jj+jp,jk,1-jp,kp) = ztj_g_lim
668                     !
669                     ! Get coefficients of isoneutral diffusion tensor
670                     ! 1. Utilise gradients *relative* to s-coordinate, so add t-point slopes (*subtract* depth gradients)
671                     ! 2. We require that isoneutral diffusion  gives no vertical buoyancy flux
672                     !     i.e. 33 term = (real slope* 31, 13 terms)
673                     ! To do this, retain limited sx**2  in vertical flux, but divide by real slope for 13/31 terms
674                     ! Equivalent to tapering A_iso = sx_limited**2/(real slope)**2
675                     !
676                     zti_lim  = ( zti_g_lim + zti_coord ) * umask(ji,jj,jk+kp)    ! remove coordinate slope => relative to coordinate surfaces
677                     ztj_lim  = ( ztj_g_lim + ztj_coord ) * vmask(ji,jj,jk+kp)
678                     !
679                     IF( ln_triad_iso ) THEN
680                        zti_raw = zti_lim**2 / zti_raw
681                        ztj_raw = ztj_lim**2 / ztj_raw
682                        zti_raw = SIGN( MIN( ABS(zti_lim), ABS( zti_raw ) ), zti_raw )
683                        ztj_raw = SIGN( MIN( ABS(ztj_lim), ABS( ztj_raw ) ), ztj_raw )
684                        zti_lim =           zfacti   * zti_lim                       &
685                        &      + ( 1._wp - zfacti ) * zti_raw
686                        ztj_lim =           zfactj   * ztj_lim                       &
687                        &      + ( 1._wp - zfactj ) * ztj_raw
688                     ENDIF
689                     triadi(ji+ip,jj   ,jk,1-ip,kp) = zti_lim
690                     triadj(ji   ,jj+jp,jk,1-jp,kp) = ztj_lim
691                    !
692                     zbu = e1u(ji    ,jj) * e2u(ji   ,jj) * fse3u(ji   ,jj,jk   )
693                     zbv = e1v(ji    ,jj) * e2v(ji   ,jj) * fse3v(ji   ,jj,jk   )
694                     zbti = e1t(ji+ip,jj) * e2t(ji+ip,jj) * fse3w(ji+ip,jj,jk+kp)
695                     zbtj = e1t(ji,jj+jp) * e2t(ji,jj+jp) * fse3w(ji,jj+jp,jk+kp)
696                     !
697                     !!gm this may inhibit vectorization on Vect Computers, and even on scalar computers....  ==> to be checked
698                     wslp2 (ji+ip,jj,jk+kp) = wslp2(ji+ip,jj,jk+kp) + 0.25_wp * zbu / zbti * zti_g_lim**2      ! masked
699                     wslp2 (ji,jj+jp,jk+kp) = wslp2(ji,jj+jp,jk+kp) + 0.25_wp * zbv / zbtj * ztj_g_lim**2
700                  END DO
701               END DO
702            END DO
703         END DO
704      END DO
705      !
706      wslp2(:,:,1) = 0._wp                ! force the surface wslp to zero
707
708      CALL lbc_lnk( wslp2, 'W', 1. )      ! lateral boundary confition on wslp2 only   ==>>> gm : necessary ? to be checked
709      !
710      CALL wrk_dealloc( jpi,jpj, z1_mlbw )
711      CALL wrk_dealloc( jpi,jpj,jpk,2, zdxrho , zdyrho, zdzrho,              klstart = 0  )
712      CALL wrk_dealloc( jpi,jpj,  2,2, zti_mlb, ztj_mlb,        kkstart = 0, klstart = 0  )
713      !
714      IF( nn_timing == 1 )  CALL timing_stop('ldf_slp_grif')
715      !
716   END SUBROUTINE ldf_slp_grif
717
718
719   SUBROUTINE ldf_slp_mxl( prd, pn2, p_gru, p_grv, p_dzr )
720      !!----------------------------------------------------------------------
721      !!                  ***  ROUTINE ldf_slp_mxl  ***
722      !!
723      !! ** Purpose :   Compute the slopes of iso-neutral surface just below
724      !!              the mixed layer.
725      !!
726      !! ** Method  :   The slope in the i-direction is computed at u- & w-points
727      !!              (uslpml, wslpiml) and the slope in the j-direction is computed
728      !!              at v- and w-points (vslpml, wslpjml) with the same bounds as
729      !!              in ldf_slp.
730      !!
731      !! ** Action  :   uslpml, wslpiml :  i- &  j-slopes of neutral surfaces
732      !!                vslpml, wslpjml    just below the mixed layer
733      !!                omlmask         :  mixed layer mask
734      !!----------------------------------------------------------------------
735      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   prd            ! in situ density
736      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   pn2            ! Brunt-Vaisala frequency (locally ref.)
737      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_gru, p_grv   ! i- & j-gradient of density (u- & v-pts)
738      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_dzr          ! z-gradient of density      (T-point)
739      !!
740      INTEGER  ::   ji , jj , jk                   ! dummy loop indices
741      INTEGER  ::   iku, ikv, ik, ikm1             ! local integers
742      REAL(wp) ::   zeps, zm1_g, zm1_2g            ! local scalars
743      REAL(wp) ::   zci, zfi, zau, zbu, zai, zbi   !   -      -
744      REAL(wp) ::   zcj, zfj, zav, zbv, zaj, zbj   !   -      -
745      REAL(wp) ::   zck, zfk,      zbw             !   -      -
746      !!----------------------------------------------------------------------
747      !
748      IF( nn_timing == 1 )  CALL timing_start('ldf_slp_mxl')
749      !
750      zeps   =  1.e-20_wp        !==   Local constant initialization   ==!
751      zm1_g  = -1.0_wp / grav
752      zm1_2g = -0.5_wp / grav
753      !
754      uslpml (1,:) = 0._wp      ;      uslpml (jpi,:) = 0._wp
755      vslpml (1,:) = 0._wp      ;      vslpml (jpi,:) = 0._wp
756      wslpiml(1,:) = 0._wp      ;      wslpiml(jpi,:) = 0._wp
757      wslpjml(1,:) = 0._wp      ;      wslpjml(jpi,:) = 0._wp
758      !
759      !                                            !==   surface mixed layer mask   !
760      DO jk = 1, jpk                               ! =1 inside the mixed layer, =0 otherwise
761         DO jj = 1, jpj
762            DO ji = 1, jpi
763               ik = nmln(ji,jj) - 1
764               IF( jk <= ik .AND. jk >= mikt(ji,jj) ) THEN
765                  omlmask(ji,jj,jk) = 1._wp
766               ELSE
767                  omlmask(ji,jj,jk) = 0._wp
768               ENDIF
769            END DO
770         END DO
771      END DO
772
773
774      ! Slopes of isopycnal surfaces just before bottom of mixed layer
775      ! --------------------------------------------------------------
776      ! The slope are computed as in the 3D case.
777      ! A key point here is the definition of the mixed layer at u- and v-points.
778      ! It is assumed to be the maximum of the two neighbouring T-point mixed layer depth.
779      ! Otherwise, a n2 value inside the mixed layer can be involved in the computation
780      ! of the slope, resulting in a too steep diagnosed slope and thus a spurious eddy
781      ! induce velocity field near the base of the mixed layer.
782      !-----------------------------------------------------------------------
783      !
784      DO jj = 2, jpjm1
785         DO ji = 2, jpim1
786            !                        !==   Slope at u- & v-points just below the Mixed Layer   ==!
787            !
788            !                        !- vertical density gradient for u- and v-slopes (from dzr at T-point)
789            iku = MIN(  MAX( miku(ji,jj)+1, nmln(ji,jj) , nmln(ji+1,jj) ) , jpkm1  )   ! ML (MAX of T-pts, bound by jpkm1)
790            ikv = MIN(  MAX( mikv(ji,jj)+1, nmln(ji,jj) , nmln(ji,jj+1) ) , jpkm1  )   !
791            zbu = 0.5_wp * ( p_dzr(ji,jj,iku) + p_dzr(ji+1,jj  ,iku) )
792            zbv = 0.5_wp * ( p_dzr(ji,jj,ikv) + p_dzr(ji  ,jj+1,ikv) )
793            !                        !- horizontal density gradient at u- & v-points
794            zau = p_gru(ji,jj,iku) / e1u(ji,jj)
795            zav = p_grv(ji,jj,ikv) / e2v(ji,jj)
796            !                        !- bound the slopes: abs(zw.)<= 1/100 and zb..<0
797            !                           kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
798            zbu = MIN(  zbu , -100._wp* ABS( zau ) , -7.e+3_wp/fse3u(ji,jj,iku)* ABS( zau )  )
799            zbv = MIN(  zbv , -100._wp* ABS( zav ) , -7.e+3_wp/fse3v(ji,jj,ikv)* ABS( zav )  )
800            !                        !- Slope at u- & v-points (uslpml, vslpml)
801            uslpml(ji,jj) = zau / ( zbu - zeps ) * umask(ji,jj,iku)
802            vslpml(ji,jj) = zav / ( zbv - zeps ) * vmask(ji,jj,ikv)
803            !
804            !                        !==   i- & j-slopes at w-points just below the Mixed Layer   ==!
805            !
806            ik   = MIN( nmln(ji,jj) + 1, jpk )
807            ikm1 = MAX( 1, ik-1 )
808            !                        !- vertical density gradient for w-slope (from N^2)
809            zbw = zm1_2g * pn2 (ji,jj,ik) * ( prd (ji,jj,ik) + prd (ji,jj,ikm1) + 2. )
810            !                        !- horizontal density i- & j-gradient at w-points
811            zci = MAX(   umask(ji-1,jj,ik  ) + umask(ji,jj,ik  )           &
812               &       + umask(ji-1,jj,ikm1) + umask(ji,jj,ikm1) , zeps  ) * e1t(ji,jj)
813            zcj = MAX(   vmask(ji,jj-1,ik  ) + vmask(ji,jj,ik  )           &
814               &       + vmask(ji,jj-1,ikm1) + vmask(ji,jj,ikm1) , zeps  ) * e2t(ji,jj)
815            zai =    (   p_gru(ji-1,jj,ik  ) + p_gru(ji,jj,ik)           &
816               &       + p_gru(ji-1,jj,ikm1) + p_gru(ji,jj,ikm1  )  ) / zci  * tmask(ji,jj,ik)
817            zaj =    (   p_grv(ji,jj-1,ik  ) + p_grv(ji,jj,ik  )           &
818               &       + p_grv(ji,jj-1,ikm1) + p_grv(ji,jj,ikm1)  ) / zcj  * tmask(ji,jj,ik)
819            !                        !- bound the slopes: abs(zw.)<= 1/100 and zb..<0.
820            !                           kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
821            zbi = MIN(  zbw , -100._wp* ABS( zai ) , -7.e+3_wp/fse3w(ji,jj,ik)* ABS( zai )  )
822            zbj = MIN(  zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/fse3w(ji,jj,ik)* ABS( zaj )  )
823            !                        !- i- & j-slope at w-points (wslpiml, wslpjml)
824            wslpiml(ji,jj) = zai / ( zbi - zeps ) * wmask (ji,jj,ik)
825            wslpjml(ji,jj) = zaj / ( zbj - zeps ) * wmask (ji,jj,ik)
826         END DO
827      END DO
828      !!gm this lbc_lnk should be useless....
829      CALL lbc_lnk( uslpml , 'U', -1. )   ;   CALL lbc_lnk( vslpml , 'V', -1. )   ! lateral boundary cond. (sign change)
830      CALL lbc_lnk( wslpiml, 'W', -1. )   ;   CALL lbc_lnk( wslpjml, 'W', -1. )   ! lateral boundary conditions
831      !
832      IF( nn_timing == 1 )  CALL timing_stop('ldf_slp_mxl')
833      !
834   END SUBROUTINE ldf_slp_mxl
835
836
837   SUBROUTINE ldf_slp_init
838      !!----------------------------------------------------------------------
839      !!                  ***  ROUTINE ldf_slp_init  ***
840      !!
841      !! ** Purpose :   Initialization for the isopycnal slopes computation
842      !!
843      !! ** Method  :   read the nammbf namelist and check the parameter
844      !!      values called by tra_dmp at the first timestep (nit000)
845      !!----------------------------------------------------------------------
846      INTEGER ::   ji, jj, jk   ! dummy loop indices
847      INTEGER ::   ierr         ! local integer
848      !!----------------------------------------------------------------------
849      !
850      IF( nn_timing == 1 )  CALL timing_start('ldf_slp_init')
851      !
852      IF(lwp) THEN
853         WRITE(numout,*)
854         WRITE(numout,*) 'ldf_slp_init : direction of lateral mixing'
855         WRITE(numout,*) '~~~~~~~~~~~~'
856      ENDIF
857
858      IF( ln_traldf_grif ) THEN        ! Griffies operator : triad of slopes
859         ALLOCATE( triadi_g(jpi,jpj,jpk,0:1,0:1) , triadj_g(jpi,jpj,jpk,0:1,0:1) , wslp2(jpi,jpj,jpk) , STAT=ierr )
860         ALLOCATE( triadi  (jpi,jpj,jpk,0:1,0:1) , triadj  (jpi,jpj,jpk,0:1,0:1)                      , STAT=ierr )
861         IF( ierr > 0             )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Griffies operator slope' )
862         !
863         IF( ln_dynldf_iso )   CALL ctl_stop( 'ldf_slp_init: Griffies operator on momentum not supported' )
864         !
865      ELSE                             ! Madec operator : slopes at u-, v-, and w-points
866         ALLOCATE( uslp(jpi,jpj,jpk) , vslp(jpi,jpj,jpk) , wslpi(jpi,jpj,jpk) , wslpj(jpi,jpj,jpk) ,                &
867            &   omlmask(jpi,jpj,jpk) , uslpml(jpi,jpj)   , vslpml(jpi,jpj)    , wslpiml(jpi,jpj)   , wslpjml(jpi,jpj) , STAT=ierr )
868         IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Madec operator slope ' )
869
870         ! Direction of lateral diffusion (tracers and/or momentum)
871         ! ------------------------------
872         uslp (:,:,:) = 0._wp   ;   uslpml (:,:) = 0._wp      ! set the slope to zero (even in s-coordinates)
873         vslp (:,:,:) = 0._wp   ;   vslpml (:,:) = 0._wp
874         wslpi(:,:,:) = 0._wp   ;   wslpiml(:,:) = 0._wp
875         wslpj(:,:,:) = 0._wp   ;   wslpjml(:,:) = 0._wp
876
877         IF(ln_sco .AND.  (ln_traldf_hor .OR. ln_dynldf_hor )) THEN
878            IF(lwp)   WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces'
879
880            ! geopotential diffusion in s-coordinates on tracers and/or momentum
881            ! The slopes of s-surfaces are computed once (no call to ldfslp in step)
882            ! The slopes for momentum diffusion are i- or j- averaged of those on tracers
883
884            ! set the slope of diffusion to the slope of s-surfaces
885            !      ( c a u t i o n : minus sign as fsdep has positive value )
886            DO jk = 1, jpk
887               DO jj = 2, jpjm1
888                  DO ji = fs_2, fs_jpim1   ! vector opt.
889                     uslp (ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept_b(ji+1,jj,jk) - fsdept_b(ji ,jj ,jk) ) * umask(ji,jj,jk)
890                     vslp (ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept_b(ji,jj+1,jk) - fsdept_b(ji ,jj ,jk) ) * vmask(ji,jj,jk)
891                     wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw_b(ji+1,jj,jk) - fsdepw_b(ji-1,jj,jk) ) * tmask(ji,jj,jk) * 0.5
892                     wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw_b(ji,jj+1,jk) - fsdepw_b(ji,jj-1,jk) ) * tmask(ji,jj,jk) * 0.5
893                  END DO
894               END DO
895            END DO
896            CALL lbc_lnk( uslp , 'U', -1. )   ;   CALL lbc_lnk( vslp , 'V', -1. )      ! Lateral boundary conditions
897            CALL lbc_lnk( wslpi, 'W', -1. )   ;   CALL lbc_lnk( wslpj, 'W', -1. )
898         ENDIF
899      ENDIF
900      !
901      IF( nn_timing == 1 )  CALL timing_stop('ldf_slp_init')
902      !
903   END SUBROUTINE ldf_slp_init
904
905#else
906   !!------------------------------------------------------------------------
907   !!   Dummy module :                 NO Rotation of lateral mixing tensor
908   !!------------------------------------------------------------------------
909   LOGICAL, PUBLIC, PARAMETER ::   lk_ldfslp = .FALSE.    !: slopes flag
910CONTAINS
911   SUBROUTINE ldf_slp( kt, prd, pn2 )   ! Dummy routine
912      INTEGER, INTENT(in) :: kt
913      REAL, DIMENSION(:,:,:), INTENT(in) :: prd, pn2
914      WRITE(*,*) 'ldf_slp: You should not have seen this print! error?', kt, prd(1,1,1), pn2(1,1,1)
915   END SUBROUTINE ldf_slp
916   SUBROUTINE ldf_slp_grif( kt )        ! Dummy routine
917      INTEGER, INTENT(in) :: kt
918      WRITE(*,*) 'ldf_slp_grif: You should not have seen this print! error?', kt
919   END SUBROUTINE ldf_slp_grif
920   SUBROUTINE ldf_slp_init              ! Dummy routine
921   END SUBROUTINE ldf_slp_init
922#endif
923
924   !!======================================================================
925END MODULE ldfslp
Note: See TracBrowser for help on using the repository browser.