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 NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/LDF – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/LDF/ldfslp.F90 @ 12340

Last change on this file since 12340 was 12340, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. This commit introduces basic do loop macro
substitution to the 2019 option 1, merge branch. These changes have been SETTE
tested. The only addition is the do_loop_substitute.h90 file in the OCE directory but
the macros defined therein are used throughout the code to replace identifiable, 2D-
and 3D- nested loop opening and closing statements with single-line alternatives. Code
indents are also adjusted accordingly.

The following explanation is taken from comments in the new header file:

This header file contains preprocessor definitions and macros used in the do-loop
substitutions introduced between version 4.0 and 4.2. The primary aim of these macros
is to assist in future applications of tiling to improve performance. This is expected
to be achieved by alternative versions of these macros in selected locations. The
initial introduction of these macros simply replaces all identifiable nested 2D- and
3D-loops with single line statements (and adjusts indenting accordingly). Do loops
are identifiable if they comform to either:

DO jk = ....

DO jj = .... DO jj = ...

DO ji = .... DO ji = ...
. OR .
. .

END DO END DO

END DO END DO

END DO

and white-space variants thereof.

Additionally, only loops with recognised jj and ji loops limits are treated; these are:
Lower limits of 1, 2 or fs_2
Upper limits of jpi, jpim1 or fs_jpim1 (for ji) or jpj, jpjm1 or fs_jpjm1 (for jj)

The macro naming convention takes the form: DO_2D_BT_LR where:

B is the Bottom offset from the PE's inner domain;
T is the Top offset from the PE's inner domain;
L is the Left offset from the PE's inner domain;
R is the Right offset from the PE's inner domain

So, given an inner domain of 2,jpim1 and 2,jpjm1, a typical example would replace:

DO jj = 2, jpj

DO ji = 1, jpim1
.
.

END DO

END DO

with:

DO_2D_01_10
.
.
END_2D

similar conventions apply to the 3D loops macros. jk loop limits are retained
through macro arguments and are not restricted. This includes the possibility of
strides for which an extra set of DO_3DS macros are defined.

In the example definition below the inner PE domain is defined by start indices of
(kIs, kJs) and end indices of (kIe, KJe)

#define DO_2D_00_00 DO jj = kJs, kJe ; DO ji = kIs, kIe
#define END_2D END DO ; END DO

TO DO:


Only conventional nested loops have been identified and replaced by this step. There are constructs such as:

DO jk = 2, jpkm1

z2d(:,:) = z2d(:,:) + e3w(:,:,jk,Kmm) * z3d(:,:,jk) * wmask(:,:,jk)

END DO

which may need to be considered.

  • Property svn:keywords set to Id
File size: 46.9 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   !!            3.7  ! 2013-12  (F. Lemarie, G. Madec)  add limiter on triad slopes
14   !!----------------------------------------------------------------------
15
16   !!----------------------------------------------------------------------
17   !!   ldf_slp       : calculates the slopes of neutral surface   (Madec operator)
18   !!   ldf_slp_triad : calculates the triads of isoneutral slopes (Griffies operator)
19   !!   ldf_slp_mxl   : calculates the slopes at the base of the mixed layer (Madec operator)
20   !!   ldf_slp_init  : initialization of the slopes computation
21   !!----------------------------------------------------------------------
22   USE oce            ! ocean dynamics and tracers
23   USE isf_oce        ! ice shelf
24   USE dom_oce        ! ocean space and time domain
25!   USE ldfdyn         ! lateral diffusion: eddy viscosity coef.
26   USE phycst         ! physical constants
27   USE zdfmxl         ! mixed layer depth
28   USE eosbn2         ! equation of states
29   !
30   USE in_out_manager ! I/O manager
31   USE prtctl         ! Print control
32   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
33   USE lib_mpp        ! distribued memory computing library
34   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
35   USE timing         ! Timing
36
37   IMPLICIT NONE
38   PRIVATE
39
40   PUBLIC   ldf_slp         ! routine called by step.F90
41   PUBLIC   ldf_slp_triad   ! routine called by step.F90
42   PUBLIC   ldf_slp_init    ! routine called by nemogcm.F90
43
44   LOGICAL , PUBLIC ::   l_ldfslp = .FALSE.     !: slopes flag
45
46   LOGICAL , PUBLIC ::   ln_traldf_iso   = .TRUE.       !: iso-neutral direction                           (nam_traldf namelist)
47   LOGICAL , PUBLIC ::   ln_traldf_triad = .FALSE.      !: griffies triad scheme                           (nam_traldf namelist)
48   LOGICAL , PUBLIC ::   ln_dynldf_iso                  !: iso-neutral direction                           (nam_dynldf namelist)
49
50   LOGICAL , PUBLIC ::   ln_triad_iso    = .FALSE.      !: pure horizontal mixing in ML                    (nam_traldf namelist)
51   LOGICAL , PUBLIC ::   ln_botmix_triad = .FALSE.      !: mixing on bottom                                (nam_traldf namelist)
52   REAL(wp), PUBLIC ::   rn_sw_triad     = 1._wp        !: =1 switching triads ; =0 all four triads used   (nam_traldf namelist)
53   REAL(wp), PUBLIC ::   rn_slpmax       = 0.01_wp      !: slope limit                                     (nam_traldf namelist)
54
55   LOGICAL , PUBLIC ::   l_grad_zps = .FALSE.           !: special treatment for Horz Tgradients w partial steps (triad operator)
56   
57   !                                                     !! Classic operator (Madec)
58   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   uslp, wslpi          !: i_slope at U- and W-points
59   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   vslp, wslpj          !: j-slope at V- and W-points
60   !                                                     !! triad operator (Griffies)
61   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   wslp2                !: wslp**2 from Griffies quarter cells
62   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi_g, triadj_g   !: skew flux  slopes relative to geopotentials
63   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:,:) ::   triadi  , triadj     !: isoneutral slopes relative to model-coordinate
64   !                                                     !! both operators
65   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   ah_wslp2             !: ah * slope^2 at w-point
66   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     ::   akz                  !: stabilizing vertical diffusivity
67   
68   !                                                     !! Madec operator
69   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   omlmask           ! mask of the surface mixed layer at T-pt
70   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   uslpml, wslpiml   ! i_slope at U- and W-points just below the mixed layer
71   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   vslpml, wslpjml   ! j_slope at V- and W-points just below the mixed layer
72
73   REAL(wp) ::   repsln = 1.e-25_wp       ! tiny value used as minium of di(rho), dj(rho) and dk(rho)
74
75   !! * Substitutions
76#  include "vectopt_loop_substitute.h90"
77#  include "do_loop_substitute.h90"
78   !!----------------------------------------------------------------------
79   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
80   !! $Id$
81   !! Software governed by the CeCILL license (see ./LICENSE)
82   !!----------------------------------------------------------------------
83CONTAINS
84
85   SUBROUTINE ldf_slp( kt, prd, pn2, Kbb, Kmm )
86      !!----------------------------------------------------------------------
87      !!                 ***  ROUTINE ldf_slp  ***
88      !!
89      !! ** Purpose :   Compute the slopes of neutral surface (slope of isopycnal
90      !!              surfaces referenced locally) (ln_traldf_iso=T).
91      !!
92      !! ** Method  :   The slope in the i-direction is computed at U- and
93      !!      W-points (uslp, wslpi) and the slope in the j-direction is
94      !!      computed at V- and W-points (vslp, wslpj).
95      !!      They are bounded by 1/100 over the whole ocean, and within the
96      !!      surface layer they are bounded by the distance to the surface
97      !!      ( slope<= depth/l  where l is the length scale of horizontal
98      !!      diffusion (here, aht=2000m2/s ==> l=20km with a typical velocity
99      !!      of 10cm/s)
100      !!        A horizontal shapiro filter is applied to the slopes
101      !!        ln_sco=T, s-coordinate, add to the previously computed slopes
102      !!      the slope of the model level surface.
103      !!        macro-tasked on horizontal slab (jk-loop)  (2, jpk-1)
104      !!      [slopes already set to zero at level 1, and to zero or the ocean
105      !!      bottom slope (ln_sco=T) at level jpk in inildf]
106      !!
107      !! ** Action : - uslp, wslpi, and vslp, wslpj, the i- and  j-slopes
108      !!               of now neutral surfaces at u-, w- and v- w-points, resp.
109      !!----------------------------------------------------------------------
110      INTEGER , INTENT(in)                   ::   kt    ! ocean time-step index
111      INTEGER , INTENT(in)                   ::   Kbb, Kmm   ! ocean time level indices
112      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   prd   ! in situ density
113      REAL(wp), INTENT(in), DIMENSION(:,:,:) ::   pn2   ! Brunt-Vaisala frequency (locally ref.)
114      !!
115      INTEGER  ::   ji , jj , jk    ! dummy loop indices
116      INTEGER  ::   ii0, ii1        ! temporary integer
117      INTEGER  ::   ij0, ij1        ! temporary integer
118      REAL(wp) ::   zeps, zm1_g, zm1_2g, z1_16, zcofw, z1_slpmax ! local scalars
119      REAL(wp) ::   zci, zfi, zau, zbu, zai, zbi   !   -      -
120      REAL(wp) ::   zcj, zfj, zav, zbv, zaj, zbj   !   -      -
121      REAL(wp) ::   zck, zfk,      zbw             !   -      -
122      REAL(wp) ::   zdepu, zdepv                   !   -      -
123      REAL(wp), DIMENSION(jpi,jpj)     ::  zslpml_hmlpu, zslpml_hmlpv
124      REAL(wp), DIMENSION(jpi,jpj,jpk) ::  zgru, zwz, zdzr
125      REAL(wp), DIMENSION(jpi,jpj,jpk) ::  zgrv, zww
126      !!----------------------------------------------------------------------
127      !
128      IF( ln_timing )   CALL timing_start('ldf_slp')
129      !
130      zeps   =  1.e-20_wp        !==   Local constant initialization   ==!
131      z1_16  =  1.0_wp / 16._wp
132      zm1_g  = -1.0_wp / grav
133      zm1_2g = -0.5_wp / grav
134      z1_slpmax = 1._wp / rn_slpmax
135      !
136      zww(:,:,:) = 0._wp
137      zwz(:,:,:) = 0._wp
138      !
139      DO_3D_10_10( 1, jpk )
140         zgru(ji,jj,jk) = umask(ji,jj,jk) * ( prd(ji+1,jj  ,jk) - prd(ji,jj,jk) )
141         zgrv(ji,jj,jk) = vmask(ji,jj,jk) * ( prd(ji  ,jj+1,jk) - prd(ji,jj,jk) )
142      END_3D
143      IF( ln_zps ) THEN                           ! partial steps correction at the bottom ocean level
144         DO_2D_10_10
145            zgru(ji,jj,mbku(ji,jj)) = gru(ji,jj)
146            zgrv(ji,jj,mbkv(ji,jj)) = grv(ji,jj)
147         END_2D
148      ENDIF
149      IF( ln_zps .AND. ln_isfcav ) THEN           ! partial steps correction at the bottom ocean level
150         DO_2D_10_10
151            IF( miku(ji,jj) > 1 )   zgru(ji,jj,miku(ji,jj)) = grui(ji,jj) 
152            IF( mikv(ji,jj) > 1 )   zgrv(ji,jj,mikv(ji,jj)) = grvi(ji,jj)
153         END_2D
154      ENDIF
155      !
156      zdzr(:,:,1) = 0._wp        !==   Local vertical density gradient at T-point   == !   (evaluated from N^2)
157      DO jk = 2, jpkm1
158         !                                ! zdzr = d/dz(prd)= - ( prd ) / grav * mk(pn2) -- at t point
159         !                                !   trick: tmask(ik  )  = 0   =>   all pn2   = 0   =>   zdzr = 0
160         !                                !    else  tmask(ik+1)  = 0   =>   pn2(ik+1) = 0   =>   zdzr divides by 1
161         !                                !          umask(ik+1) /= 0   =>   all pn2  /= 0   =>   zdzr divides by 2
162         !                                ! NB: 1/(tmask+1) = (1-.5*tmask)  substitute a / by a *  ==> faster
163         zdzr(:,:,jk) = zm1_g * ( prd(:,:,jk) + 1._wp )              &
164            &                 * ( pn2(:,:,jk) + pn2(:,:,jk+1) ) * ( 1._wp - 0.5_wp * tmask(:,:,jk+1) )
165      END DO
166      !
167      !                          !==   Slopes just below the mixed layer   ==!
168      CALL ldf_slp_mxl( prd, pn2, zgru, zgrv, zdzr, Kmm )        ! output: uslpml, vslpml, wslpiml, wslpjml
169
170
171      ! I.  slopes at u and v point      | uslp = d/di( prd ) / d/dz( prd )
172      ! ===========================      | vslp = d/dj( prd ) / d/dz( prd )
173      !
174      IF ( ln_isfcav ) THEN
175         DO_2D_00_00
176            zslpml_hmlpu(ji,jj) = uslpml(ji,jj) / ( MAX(hmlpt  (ji,jj), hmlpt  (ji+1,jj  ), 5._wp) &
177               &                                  - MAX(risfdep(ji,jj), risfdep(ji+1,jj  )       ) ) 
178            zslpml_hmlpv(ji,jj) = vslpml(ji,jj) / ( MAX(hmlpt  (ji,jj), hmlpt  (ji  ,jj+1), 5._wp) &
179               &                                  - MAX(risfdep(ji,jj), risfdep(ji  ,jj+1)       ) )
180         END_2D
181      ELSE
182         DO_2D_00_00
183            zslpml_hmlpu(ji,jj) = uslpml(ji,jj) / MAX(hmlpt(ji,jj), hmlpt(ji+1,jj  ), 5._wp)
184            zslpml_hmlpv(ji,jj) = vslpml(ji,jj) / MAX(hmlpt(ji,jj), hmlpt(ji  ,jj+1), 5._wp)
185         END_2D
186      END IF
187
188      DO_3D_00_00( 2, jpkm1 )
189         !                                      ! horizontal and vertical density gradient at u- and v-points
190         zau = zgru(ji,jj,jk) * r1_e1u(ji,jj)
191         zav = zgrv(ji,jj,jk) * r1_e2v(ji,jj)
192         zbu = 0.5_wp * ( zdzr(ji,jj,jk) + zdzr(ji+1,jj  ,jk) )
193         zbv = 0.5_wp * ( zdzr(ji,jj,jk) + zdzr(ji  ,jj+1,jk) )
194         !                                      ! bound the slopes: abs(zw.)<= 1/100 and zb..<0
195         !                                      ! + kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
196         zbu = MIN(  zbu, - z1_slpmax * ABS( zau ) , -7.e+3_wp/e3u(ji,jj,jk,Kmm)* ABS( zau )  )
197         zbv = MIN(  zbv, - z1_slpmax * ABS( zav ) , -7.e+3_wp/e3v(ji,jj,jk,Kmm)* ABS( zav )  )
198         !                                      ! Fred Dupont: add a correction for bottom partial steps:
199         !                                      !              max slope = 1/2 * e3 / e1
200         IF (ln_zps .AND. jk==mbku(ji,jj)) &
201            zbu = MIN(  zbu, - z1_slpmax * ABS( zau ) , - 2._wp * e1u(ji,jj) / e3u(ji,jj,jk,Kmm)* ABS( zau )  )
202         IF (ln_zps .AND. jk==mbkv(ji,jj)) &
203            zbv = MIN(  zbv, - z1_slpmax * ABS( zav ) , - 2._wp * e2v(ji,jj) / e3v(ji,jj,jk,Kmm)* ABS( zav )  )
204         !                                      ! uslp and vslp output in zwz and zww, resp.
205         zfi = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj,jk) )
206         zfj = MAX( omlmask(ji,jj,jk), omlmask(ji,jj+1,jk) )
207         ! thickness of water column between surface and level k at u/v point
208         zdepu = 0.5_wp * ( ( gdept (ji,jj,jk,Kmm) + gdept (ji+1,jj,jk,Kmm) )                            &
209                          - 2 * MAX( risfdep(ji,jj), risfdep(ji+1,jj) ) - e3u(ji,jj,miku(ji,jj),Kmm)   )
210         zdepv = 0.5_wp * ( ( gdept (ji,jj,jk,Kmm) + gdept (ji,jj+1,jk,Kmm) )                            &
211                          - 2 * MAX( risfdep(ji,jj), risfdep(ji,jj+1) ) - e3v(ji,jj,mikv(ji,jj),Kmm)   )
212         !
213         zwz(ji,jj,jk) = ( ( 1._wp - zfi) * zau / ( zbu - zeps )                                     &
214            &                      + zfi  * zdepu * zslpml_hmlpu(ji,jj) ) * umask(ji,jj,jk)
215         zww(ji,jj,jk) = ( ( 1._wp - zfj) * zav / ( zbv - zeps )                                     &
216            &                      + zfj  * zdepv * zslpml_hmlpv(ji,jj) ) * vmask(ji,jj,jk)
217!!gm  modif to suppress omlmask.... (as in Griffies case)
218!               !                                         ! jk must be >= ML level for zf=1. otherwise  zf=0.
219!               zfi = REAL( 1 - 1/(1 + jk / MAX( nmln(ji+1,jj), nmln(ji,jj) ) ), wp )
220!               zfj = REAL( 1 - 1/(1 + jk / MAX( nmln(ji,jj+1), nmln(ji,jj) ) ), wp )
221!               zci = 0.5 * ( gdept(ji+1,jj,jk,Kmm)+gdept(ji,jj,jk,Kmm) ) / MAX( hmlpt(ji,jj), hmlpt(ji+1,jj), 10. ) )
222!               zcj = 0.5 * ( gdept(ji,jj+1,jk,Kmm)+gdept(ji,jj,jk,Kmm) ) / MAX( hmlpt(ji,jj), hmlpt(ji,jj+1), 10. ) )
223!               zwz(ji,jj,jk) = ( zfi * zai / ( zbi - zeps ) + ( 1._wp - zfi ) * wslpiml(ji,jj) * zci ) * tmask(ji,jj,jk)
224!               zww(ji,jj,jk) = ( zfj * zaj / ( zbj - zeps ) + ( 1._wp - zfj ) * wslpjml(ji,jj) * zcj ) * tmask(ji,jj,jk)
225!!gm end modif
226      END_3D
227      CALL lbc_lnk_multi( 'ldfslp', zwz, 'U', -1.,  zww, 'V', -1. )      ! lateral boundary conditions
228      !
229      !                                            !* horizontal Shapiro filter
230      DO jk = 2, jpkm1
231         DO_2D_00_00
232            uslp(ji,jj,jk) = z1_16 * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      &
233               &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      &
234               &                       + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)      &
235               &                       +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )    &
236               &                       + 4.*  zwz(ji  ,jj  ,jk)                       )
237            vslp(ji,jj,jk) = z1_16 * (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)      &
238               &                       +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)      &
239               &                       + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)      &
240               &                       +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )    &
241               &                       + 4.*  zww(ji,jj    ,jk)                       )
242         END_2D
243         DO jj = 3, jpj-2                               ! other rows
244            DO ji = fs_2, fs_jpim1   ! vector opt.
245               uslp(ji,jj,jk) = z1_16 * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      &
246                  &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      &
247                  &                       + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)      &
248                  &                       +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )    &
249                  &                       + 4.*  zwz(ji  ,jj  ,jk)                       )
250               vslp(ji,jj,jk) = z1_16 * (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)      &
251                  &                       +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)      &
252                  &                       + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)      &
253                  &                       +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )    &
254                  &                       + 4.*  zww(ji,jj    ,jk)                       )
255            END DO
256         END DO
257         !                                        !* decrease along coastal boundaries
258         DO_2D_00_00
259            uslp(ji,jj,jk) = uslp(ji,jj,jk) * ( umask(ji,jj+1,jk) + umask(ji,jj-1,jk  ) ) * 0.5_wp   &
260               &                            * ( umask(ji,jj  ,jk) + umask(ji,jj  ,jk+1) ) * 0.5_wp
261            vslp(ji,jj,jk) = vslp(ji,jj,jk) * ( vmask(ji+1,jj,jk) + vmask(ji-1,jj,jk  ) ) * 0.5_wp   &
262               &                            * ( vmask(ji  ,jj,jk) + vmask(ji  ,jj,jk+1) ) * 0.5_wp
263         END_2D
264      END DO
265
266
267      ! II.  slopes at w point           | wslpi = mij( d/di( prd ) / d/dz( prd )
268      ! ===========================      | wslpj = mij( d/dj( prd ) / d/dz( prd )
269      !
270      DO_3D_00_00( 2, jpkm1 )
271         !                                  !* Local vertical density gradient evaluated from N^2
272         zbw = zm1_2g * pn2 (ji,jj,jk) * ( prd (ji,jj,jk) + prd (ji,jj,jk-1) + 2. )
273         !                                  !* Slopes at w point
274         !                                        ! i- & j-gradient of density at w-points
275         zci = MAX(  umask(ji-1,jj,jk  ) + umask(ji,jj,jk  )           &
276            &      + umask(ji-1,jj,jk-1) + umask(ji,jj,jk-1) , zeps  ) * e1t(ji,jj)
277         zcj = MAX(  vmask(ji,jj-1,jk  ) + vmask(ji,jj,jk-1)           &
278            &      + vmask(ji,jj-1,jk-1) + vmask(ji,jj,jk  ) , zeps  ) * e2t(ji,jj)
279         zai =    (  zgru (ji-1,jj,jk  ) + zgru (ji,jj,jk-1)           &
280            &      + zgru (ji-1,jj,jk-1) + zgru (ji,jj,jk  )   ) / zci * wmask (ji,jj,jk)
281         zaj =    (  zgrv (ji,jj-1,jk  ) + zgrv (ji,jj,jk-1)           &
282            &      + zgrv (ji,jj-1,jk-1) + zgrv (ji,jj,jk  )   ) / zcj * wmask (ji,jj,jk)
283         !                                        ! bound the slopes: abs(zw.)<= 1/100 and zb..<0.
284         !                                        ! + kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
285         zbi = MIN( zbw ,- 100._wp* ABS( zai ) , -7.e+3_wp/e3w(ji,jj,jk,Kmm)* ABS( zai )  )
286         zbj = MIN( zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/e3w(ji,jj,jk,Kmm)* ABS( zaj )  )
287         !                                        ! wslpi and wslpj with ML flattening (output in zwz and zww, resp.)
288         zfk = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) )   ! zfk=1 in the ML otherwise zfk=0
289         zck = ( gdepw(ji,jj,jk,Kmm) - gdepw(ji,jj,mikt(ji,jj),Kmm) ) / MAX( hmlp(ji,jj) - gdepw(ji,jj,mikt(ji,jj),Kmm), 10._wp )
290         zwz(ji,jj,jk) = (  zai / ( zbi - zeps ) * ( 1._wp - zfk ) + zck * wslpiml(ji,jj) * zfk  ) * wmask(ji,jj,jk)
291         zww(ji,jj,jk) = (  zaj / ( zbj - zeps ) * ( 1._wp - zfk ) + zck * wslpjml(ji,jj) * zfk  ) * wmask(ji,jj,jk)
292
293!!gm  modif to suppress omlmask....  (as in Griffies operator)
294!               !                                         ! jk must be >= ML level for zfk=1. otherwise  zfk=0.
295!               zfk = REAL( 1 - 1/(1 + jk / nmln(ji+1,jj)), wp )
296!               zck = gdepw(ji,jj,jk)    / MAX( hmlp(ji,jj), 10. )
297!               zwz(ji,jj,jk) = ( zfk * zai / ( zbi - zeps ) + ( 1._wp - zfk ) * wslpiml(ji,jj) * zck ) * tmask(ji,jj,jk)
298!               zww(ji,jj,jk) = ( zfk * zaj / ( zbj - zeps ) + ( 1._wp - zfk ) * wslpjml(ji,jj) * zck ) * tmask(ji,jj,jk)
299!!gm end modif
300      END_3D
301      CALL lbc_lnk_multi( 'ldfslp', zwz, 'T', -1.,  zww, 'T', -1. )      ! lateral boundary conditions
302      !
303      !                                           !* horizontal Shapiro filter
304      DO jk = 2, jpkm1
305         DO_2D_00_00
306            zcofw = wmask(ji,jj,jk) * z1_16
307            wslpi(ji,jj,jk) = (         zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     &
308                 &               +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     &
309                 &               + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)     &
310                 &               +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )   &
311                 &               + 4.*  zwz(ji  ,jj  ,jk)                         ) * zcofw
312
313            wslpj(ji,jj,jk) = (         zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)     &
314                 &               +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)     &
315                 &               + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)     &
316                 &               +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )   &
317                 &               + 4.*  zww(ji  ,jj  ,jk)                         ) * zcofw
318         END_2D
319         DO jj = 3, jpj-2                               ! other rows
320            DO ji = fs_2, fs_jpim1   ! vector opt.
321               zcofw = wmask(ji,jj,jk) * z1_16
322               wslpi(ji,jj,jk) = (         zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     &
323                    &               +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     &
324                    &               + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)     &
325                    &               +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )   &
326                    &               + 4.*  zwz(ji  ,jj  ,jk)                         ) * zcofw
327
328               wslpj(ji,jj,jk) = (         zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)     &
329                    &               +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)     &
330                    &               + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)     &
331                    &               +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )   &
332                    &               + 4.*  zww(ji  ,jj  ,jk)                         ) * zcofw
333            END DO
334         END DO
335         !                                        !* decrease in vicinity of topography
336         DO_2D_00_00
337            zck =   ( umask(ji,jj,jk) + umask(ji-1,jj,jk) )   &
338               &  * ( vmask(ji,jj,jk) + vmask(ji,jj-1,jk) ) * 0.25
339            wslpi(ji,jj,jk) = wslpi(ji,jj,jk) * zck
340            wslpj(ji,jj,jk) = wslpj(ji,jj,jk) * zck
341         END_2D
342      END DO
343
344      ! IV. Lateral boundary conditions
345      ! ===============================
346      CALL lbc_lnk_multi( 'ldfslp', uslp , 'U', -1. , vslp , 'V', -1. , wslpi, 'W', -1., wslpj, 'W', -1. )
347
348      IF(sn_cfctl%l_prtctl) THEN
349         CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp  - u : ', tab3d_2=vslp,  clinfo2=' v : ', kdim=jpk)
350         CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk)
351      ENDIF
352      !
353      IF( ln_timing )   CALL timing_stop('ldf_slp')
354      !
355   END SUBROUTINE ldf_slp
356
357
358   SUBROUTINE ldf_slp_triad ( kt, Kbb, Kmm )
359      !!----------------------------------------------------------------------
360      !!                 ***  ROUTINE ldf_slp_triad  ***
361      !!
362      !! ** Purpose :   Compute the squared slopes of neutral surfaces (slope
363      !!      of iso-pycnal surfaces referenced locally) (ln_traldf_triad=T)
364      !!      at W-points using the Griffies quarter-cells.
365      !!
366      !! ** Method  :   calculates alpha and beta at T-points
367      !!
368      !! ** Action : - triadi_g, triadj_g   T-pts i- and j-slope triads relative to geopot. (used for eiv)
369      !!             - triadi , triadj    T-pts i- and j-slope triads relative to model-coordinate
370      !!             - wslp2              squared slope of neutral surfaces at w-points.
371      !!----------------------------------------------------------------------
372      INTEGER, INTENT( in ) ::   kt             ! ocean time-step index
373      INTEGER , INTENT(in)  ::   Kbb, Kmm       ! ocean time level indices
374      !!
375      INTEGER  ::   ji, jj, jk, jl, ip, jp, kp  ! dummy loop indices
376      INTEGER  ::   iku, ikv                    ! local integer
377      REAL(wp) ::   zfacti, zfactj              ! local scalars
378      REAL(wp) ::   znot_thru_surface           ! local scalars
379      REAL(wp) ::   zdit, zdis, zdkt, zbu, zbti, zisw
380      REAL(wp) ::   zdjt, zdjs, zdks, zbv, zbtj, zjsw
381      REAL(wp) ::   zdxrho_raw, zti_coord, zti_raw, zti_lim, zti_g_raw, zti_g_lim
382      REAL(wp) ::   zdyrho_raw, ztj_coord, ztj_raw, ztj_lim, ztj_g_raw, ztj_g_lim
383      REAL(wp) ::   zdzrho_raw
384      REAL(wp) ::   zbeta0, ze3_e1, ze3_e2
385      REAL(wp), DIMENSION(jpi,jpj)     ::   z1_mlbw
386      REAL(wp), DIMENSION(jpi,jpj,jpk,0:1) ::   zdxrho , zdyrho, zdzrho     ! Horizontal and vertical density gradients
387      REAL(wp), DIMENSION(jpi,jpj,0:1,0:1) ::   zti_mlb, ztj_mlb            ! for Griffies operator only
388      !!----------------------------------------------------------------------
389      !
390      IF( ln_timing )   CALL timing_start('ldf_slp_triad')
391      !
392      !--------------------------------!
393      !  Some preliminary calculation  !
394      !--------------------------------!
395      !
396      DO jl = 0, 1                            !==  unmasked before density i- j-, k-gradients  ==!
397         !
398         ip = jl   ;   jp = jl                ! guaranteed nonzero gradients ( absolute value larger than repsln)
399         DO_3D_10_10( 1, jpkm1 )
400            zdit = ( ts(ji+1,jj,jk,jp_tem,Kbb) - ts(ji,jj,jk,jp_tem,Kbb) )    ! i-gradient of T & S at u-point
401            zdis = ( ts(ji+1,jj,jk,jp_sal,Kbb) - ts(ji,jj,jk,jp_sal,Kbb) )
402            zdjt = ( ts(ji,jj+1,jk,jp_tem,Kbb) - ts(ji,jj,jk,jp_tem,Kbb) )    ! j-gradient of T & S at v-point
403            zdjs = ( ts(ji,jj+1,jk,jp_sal,Kbb) - ts(ji,jj,jk,jp_sal,Kbb) )
404            zdxrho_raw = ( - rab_b(ji+ip,jj   ,jk,jp_tem) * zdit + rab_b(ji+ip,jj   ,jk,jp_sal) * zdis ) * r1_e1u(ji,jj)
405            zdyrho_raw = ( - rab_b(ji   ,jj+jp,jk,jp_tem) * zdjt + rab_b(ji   ,jj+jp,jk,jp_sal) * zdjs ) * r1_e2v(ji,jj)
406            zdxrho(ji+ip,jj   ,jk,1-ip) = SIGN(  MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw  )   ! keep the sign
407            zdyrho(ji   ,jj+jp,jk,1-jp) = SIGN(  MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw  )
408         END_3D
409         !
410         IF( ln_zps .AND. l_grad_zps ) THEN     ! partial steps: correction of i- & j-grad on bottom
411            DO_2D_10_10
412               iku  = mbku(ji,jj)          ;   ikv  = mbkv(ji,jj)             ! last ocean level (u- & v-points)
413               zdit = gtsu(ji,jj,jp_tem)   ;   zdjt = gtsv(ji,jj,jp_tem)      ! i- & j-gradient of Temperature
414               zdis = gtsu(ji,jj,jp_sal)   ;   zdjs = gtsv(ji,jj,jp_sal)      ! i- & j-gradient of Salinity
415               zdxrho_raw = ( - rab_b(ji+ip,jj   ,iku,jp_tem) * zdit + rab_b(ji+ip,jj   ,iku,jp_sal) * zdis ) * r1_e1u(ji,jj)
416               zdyrho_raw = ( - rab_b(ji   ,jj+jp,ikv,jp_tem) * zdjt + rab_b(ji   ,jj+jp,ikv,jp_sal) * zdjs ) * r1_e2v(ji,jj)
417               zdxrho(ji+ip,jj   ,iku,1-ip) = SIGN( MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw )   ! keep the sign
418               zdyrho(ji   ,jj+jp,ikv,1-jp) = SIGN( MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw )
419            END_2D
420         ENDIF
421         !
422      END DO
423
424      DO kp = 0, 1                            !==  unmasked before density i- j-, k-gradients  ==!
425         DO_3D_11_11( 1, jpkm1 )
426            IF( jk+kp > 1 ) THEN        ! k-gradient of T & S a jk+kp
427               zdkt = ( ts(ji,jj,jk+kp-1,jp_tem,Kbb) - ts(ji,jj,jk+kp,jp_tem,Kbb) )
428               zdks = ( ts(ji,jj,jk+kp-1,jp_sal,Kbb) - ts(ji,jj,jk+kp,jp_sal,Kbb) )
429            ELSE
430               zdkt = 0._wp                                             ! 1st level gradient set to zero
431               zdks = 0._wp
432            ENDIF
433            zdzrho_raw = ( - rab_b(ji,jj,jk   ,jp_tem) * zdkt & 
434                       &   + rab_b(ji,jj,jk   ,jp_sal) * zdks &
435                       & ) / e3w(ji,jj,jk+kp,Kmm) 
436            zdzrho(ji,jj,jk,kp) = - MIN( - repsln , zdzrho_raw )    ! force zdzrho >= repsln
437         END_3D
438      END DO
439      !
440      DO_2D_11_11
441         jk = MIN( nmln(ji,jj), mbkt(ji,jj) ) + 1     ! MIN in case ML depth is the ocean depth
442         z1_mlbw(ji,jj) = 1._wp / gdepw(ji,jj,jk,Kmm)
443      END_2D
444      !
445      !                                       !==  intialisations to zero  ==!
446      !
447      wslp2  (:,:,:)     = 0._wp              ! wslp2 will be cumulated 3D field set to zero
448      triadi_g(:,:,1,:,:) = 0._wp   ;   triadi_g(:,:,jpk,:,:) = 0._wp   ! set surface and bottom slope to zero
449      triadj_g(:,:,1,:,:) = 0._wp   ;   triadj_g(:,:,jpk,:,:) = 0._wp
450      !!gm _iso set to zero missing
451      triadi  (:,:,1,:,:) = 0._wp   ;   triadj  (:,:,jpk,:,:) = 0._wp   ! set surface and bottom slope to zero
452      triadj  (:,:,1,:,:) = 0._wp   ;   triadj  (:,:,jpk,:,:) = 0._wp
453
454      !-------------------------------------!
455      !  Triads just below the Mixed Layer  !
456      !-------------------------------------!
457      !
458      DO jl = 0, 1                            ! calculate slope of the 4 triads immediately ONE level below mixed-layer base
459         DO kp = 0, 1                         ! with only the slope-max limit   and   MASKED
460            DO_2D_10_10
461               ip = jl   ;   jp = jl
462               !
463               jk = nmln(ji+ip,jj) + 1
464               IF( jk > mbkt(ji+ip,jj) ) THEN   ! ML reaches bottom
465                  zti_mlb(ji+ip,jj   ,1-ip,kp) = 0.0_wp
466               ELSE                             
467                  ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth)
468                  zti_g_raw = (  zdxrho(ji+ip,jj,jk-kp,1-ip) / zdzrho(ji+ip,jj,jk-kp,kp)      &
469                     &          - ( gdept(ji+1,jj,jk-kp,Kmm) - gdept(ji,jj,jk-kp,Kmm) ) * r1_e1u(ji,jj)  ) * umask(ji,jj,jk)
470                  ze3_e1    =  e3w(ji+ip,jj,jk-kp,Kmm) * r1_e1u(ji,jj) 
471                  zti_mlb(ji+ip,jj   ,1-ip,kp) = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e1  , ABS( zti_g_raw ) ), zti_g_raw )
472               ENDIF
473               !
474               jk = nmln(ji,jj+jp) + 1
475               IF( jk >  mbkt(ji,jj+jp) ) THEN  !ML reaches bottom
476                  ztj_mlb(ji   ,jj+jp,1-jp,kp) = 0.0_wp
477               ELSE
478                  ztj_g_raw = (  zdyrho(ji,jj+jp,jk-kp,1-jp) / zdzrho(ji,jj+jp,jk-kp,kp)      &
479                     &      - ( gdept(ji,jj+1,jk-kp,Kmm) - gdept(ji,jj,jk-kp,Kmm) ) / e2v(ji,jj)  ) * vmask(ji,jj,jk)
480                  ze3_e2    =  e3w(ji,jj+jp,jk-kp,Kmm) / e2v(ji,jj)
481                  ztj_mlb(ji   ,jj+jp,1-jp,kp) = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e2  , ABS( ztj_g_raw ) ), ztj_g_raw )
482               ENDIF
483            END_2D
484         END DO
485      END DO
486
487      !-------------------------------------!
488      !  Triads with surface limits         !
489      !-------------------------------------!
490      !
491      DO kp = 0, 1                            ! k-index of triads
492         DO jl = 0, 1
493            ip = jl   ;   jp = jl             ! i- and j-indices of triads (i-k and j-k planes)
494            DO jk = 1, jpkm1
495               ! Must mask contribution to slope from dz/dx at constant s for triads jk=1,kp=0 that poke up though ocean surface
496               znot_thru_surface = REAL( 1-1/(jk+kp), wp )  !jk+kp=1,=0.; otherwise=1.0
497               DO_2D_10_10
498                  !
499                  ! Calculate slope relative to geopotentials used for GM skew fluxes
500                  ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth)
501                  ! Limit by slope *relative to geopotentials* by rn_slpmax, and mask by psi-point
502                  ! masked by umask taken at the level of dz(rho)
503                  !
504                  ! raw slopes: unmasked unbounded slopes (relative to geopotential (zti_g) and model surface (zti)
505                  !
506                  zti_raw   = zdxrho(ji+ip,jj   ,jk,1-ip) / zdzrho(ji+ip,jj   ,jk,kp)                   ! unmasked
507                  ztj_raw   = zdyrho(ji   ,jj+jp,jk,1-jp) / zdzrho(ji   ,jj+jp,jk,kp)
508                  !
509                  ! Must mask contribution to slope for triad jk=1,kp=0 that poke up though ocean surface
510                  zti_coord = znot_thru_surface * ( gdept(ji+1,jj  ,jk,Kmm) - gdept(ji,jj,jk,Kmm) ) * r1_e1u(ji,jj)
511                  ztj_coord = znot_thru_surface * ( gdept(ji  ,jj+1,jk,Kmm) - gdept(ji,jj,jk,Kmm) ) * r1_e2v(ji,jj)     ! unmasked
512                  zti_g_raw = zti_raw - zti_coord      ! ref to geopot surfaces
513                  ztj_g_raw = ztj_raw - ztj_coord
514                  ! additional limit required in bilaplacian case
515                  ze3_e1    = e3w(ji+ip,jj   ,jk+kp,Kmm) * r1_e1u(ji,jj)
516                  ze3_e2    = e3w(ji   ,jj+jp,jk+kp,Kmm) * r1_e2v(ji,jj)
517                  ! NB: hard coded factor 5 (can be a namelist parameter...)
518                  zti_g_lim = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e1, ABS( zti_g_raw ) ), zti_g_raw )
519                  ztj_g_lim = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e2, ABS( ztj_g_raw ) ), ztj_g_raw )
520                  !
521                  ! Below  ML use limited zti_g as is & mask
522                  ! Inside ML replace by linearly reducing sx_mlb towards surface & mask
523                  !
524                  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
525                  zfactj = REAL( 1 - 1/(1 + (jk+kp-1)/nmln(ji,jj+jp)), wp )  ! must be .ge. nmln(ji,jj) for zfact=1
526                  !                                                          !                   otherwise  zfact=0
527                  zti_g_lim =          ( zfacti   * zti_g_lim                       &
528                     &      + ( 1._wp - zfacti ) * zti_mlb(ji+ip,jj,1-ip,kp)   &
529                     &                           * gdepw(ji+ip,jj,jk+kp,Kmm) * z1_mlbw(ji+ip,jj) ) * umask(ji,jj,jk+kp)
530                  ztj_g_lim =          ( zfactj   * ztj_g_lim                       &
531                     &      + ( 1._wp - zfactj ) * ztj_mlb(ji,jj+jp,1-jp,kp)   &
532                     &                           * gdepw(ji,jj+jp,jk+kp,Kmm) * z1_mlbw(ji,jj+jp) ) * vmask(ji,jj,jk+kp)
533                  !
534                  triadi_g(ji+ip,jj   ,jk,1-ip,kp) = zti_g_lim
535                  triadj_g(ji   ,jj+jp,jk,1-jp,kp) = ztj_g_lim
536                  !
537                  ! Get coefficients of isoneutral diffusion tensor
538                  ! 1. Utilise gradients *relative* to s-coordinate, so add t-point slopes (*subtract* depth gradients)
539                  ! 2. We require that isoneutral diffusion  gives no vertical buoyancy flux
540                  !     i.e. 33 term = (real slope* 31, 13 terms)
541                  ! To do this, retain limited sx**2  in vertical flux, but divide by real slope for 13/31 terms
542                  ! Equivalent to tapering A_iso = sx_limited**2/(real slope)**2
543                  !
544                  zti_lim  = ( zti_g_lim + zti_coord ) * umask(ji,jj,jk+kp)    ! remove coordinate slope => relative to coordinate surfaces
545                  ztj_lim  = ( ztj_g_lim + ztj_coord ) * vmask(ji,jj,jk+kp)
546                  !
547                  IF( ln_triad_iso ) THEN
548                     zti_raw = zti_lim*zti_lim / zti_raw
549                     ztj_raw = ztj_lim*ztj_lim / ztj_raw
550                     zti_raw = SIGN( MIN( ABS(zti_lim), ABS( zti_raw ) ), zti_raw )
551                     ztj_raw = SIGN( MIN( ABS(ztj_lim), ABS( ztj_raw ) ), ztj_raw )
552                     zti_lim = zfacti * zti_lim + ( 1._wp - zfacti ) * zti_raw
553                     ztj_lim = zfactj * ztj_lim + ( 1._wp - zfactj ) * ztj_raw
554                  ENDIF
555                  !                                      ! switching triad scheme
556                  zisw = (1._wp - rn_sw_triad ) + rn_sw_triad    &
557                     &            * 2._wp * ABS( 0.5_wp - kp - ( 0.5_wp - ip ) * SIGN( 1._wp , zdxrho(ji+ip,jj,jk,1-ip) )  )
558                  zjsw = (1._wp - rn_sw_triad ) + rn_sw_triad    &
559                     &            * 2._wp * ABS( 0.5_wp - kp - ( 0.5_wp - jp ) * SIGN( 1._wp , zdyrho(ji,jj+jp,jk,1-jp) )  )
560                  !
561                  triadi(ji+ip,jj   ,jk,1-ip,kp) = zti_lim * zisw
562                  triadj(ji   ,jj+jp,jk,1-jp,kp) = ztj_lim * zjsw
563                  !
564                  zbu  = e1e2u(ji   ,jj   ) * e3u(ji   ,jj   ,jk   ,Kmm)
565                  zbv  = e1e2v(ji   ,jj   ) * e3v(ji   ,jj   ,jk   ,Kmm)
566                  zbti = e1e2t(ji+ip,jj   ) * e3w(ji+ip,jj   ,jk+kp,Kmm)
567                  zbtj = e1e2t(ji   ,jj+jp) * e3w(ji   ,jj+jp,jk+kp,Kmm)
568                  !
569                  wslp2(ji+ip,jj,jk+kp) = wslp2(ji+ip,jj,jk+kp) + 0.25_wp * zbu / zbti * zti_g_lim*zti_g_lim      ! masked
570                  wslp2(ji,jj+jp,jk+kp) = wslp2(ji,jj+jp,jk+kp) + 0.25_wp * zbv / zbtj * ztj_g_lim*ztj_g_lim
571               END_2D
572            END DO
573         END DO
574      END DO
575      !
576      wslp2(:,:,1) = 0._wp                ! force the surface wslp to zero
577
578      CALL lbc_lnk( 'ldfslp', wslp2, 'W', 1. )      ! lateral boundary confition on wslp2 only   ==>>> gm : necessary ? to be checked
579      !
580      IF( ln_timing )   CALL timing_stop('ldf_slp_triad')
581      !
582   END SUBROUTINE ldf_slp_triad
583
584
585   SUBROUTINE ldf_slp_mxl( prd, pn2, p_gru, p_grv, p_dzr, Kmm )
586      !!----------------------------------------------------------------------
587      !!                  ***  ROUTINE ldf_slp_mxl  ***
588      !!
589      !! ** Purpose :   Compute the slopes of iso-neutral surface just below
590      !!              the mixed layer.
591      !!
592      !! ** Method  :   The slope in the i-direction is computed at u- & w-points
593      !!              (uslpml, wslpiml) and the slope in the j-direction is computed
594      !!              at v- and w-points (vslpml, wslpjml) with the same bounds as
595      !!              in ldf_slp.
596      !!
597      !! ** Action  :   uslpml, wslpiml :  i- &  j-slopes of neutral surfaces
598      !!                vslpml, wslpjml    just below the mixed layer
599      !!                omlmask         :  mixed layer mask
600      !!----------------------------------------------------------------------
601      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   prd            ! in situ density
602      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   pn2            ! Brunt-Vaisala frequency (locally ref.)
603      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_gru, p_grv   ! i- & j-gradient of density (u- & v-pts)
604      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_dzr          ! z-gradient of density      (T-point)
605      INTEGER , INTENT(in)                   ::   Kmm            ! ocean time level indices
606      !!
607      INTEGER  ::   ji , jj , jk                   ! dummy loop indices
608      INTEGER  ::   iku, ikv, ik, ikm1             ! local integers
609      REAL(wp) ::   zeps, zm1_g, zm1_2g, z1_slpmax ! local scalars
610      REAL(wp) ::   zci, zfi, zau, zbu, zai, zbi   !   -      -
611      REAL(wp) ::   zcj, zfj, zav, zbv, zaj, zbj   !   -      -
612      REAL(wp) ::   zck, zfk,      zbw             !   -      -
613      !!----------------------------------------------------------------------
614      !
615      zeps   =  1.e-20_wp        !==   Local constant initialization   ==!
616      zm1_g  = -1.0_wp / grav
617      zm1_2g = -0.5_wp / grav
618      z1_slpmax = 1._wp / rn_slpmax
619      !
620      uslpml (1,:) = 0._wp      ;      uslpml (jpi,:) = 0._wp
621      vslpml (1,:) = 0._wp      ;      vslpml (jpi,:) = 0._wp
622      wslpiml(1,:) = 0._wp      ;      wslpiml(jpi,:) = 0._wp
623      wslpjml(1,:) = 0._wp      ;      wslpjml(jpi,:) = 0._wp
624      !
625      !                                            !==   surface mixed layer mask   !
626      DO_3D_11_11( 1, jpk )
627         ik = nmln(ji,jj) - 1
628         IF( jk <= ik ) THEN   ;   omlmask(ji,jj,jk) = 1._wp
629         ELSE                  ;   omlmask(ji,jj,jk) = 0._wp
630         ENDIF
631      END_3D
632
633
634      ! Slopes of isopycnal surfaces just before bottom of mixed layer
635      ! --------------------------------------------------------------
636      ! The slope are computed as in the 3D case.
637      ! A key point here is the definition of the mixed layer at u- and v-points.
638      ! It is assumed to be the maximum of the two neighbouring T-point mixed layer depth.
639      ! Otherwise, a n2 value inside the mixed layer can be involved in the computation
640      ! of the slope, resulting in a too steep diagnosed slope and thus a spurious eddy
641      ! induce velocity field near the base of the mixed layer.
642      !-----------------------------------------------------------------------
643      !
644      DO_2D_00_00
645         !                        !==   Slope at u- & v-points just below the Mixed Layer   ==!
646         !
647         !                        !- vertical density gradient for u- and v-slopes (from dzr at T-point)
648         iku = MIN(  MAX( 1, nmln(ji,jj) , nmln(ji+1,jj) ) , jpkm1  )   ! ML (MAX of T-pts, bound by jpkm1)
649         ikv = MIN(  MAX( 1, nmln(ji,jj) , nmln(ji,jj+1) ) , jpkm1  )   !
650         zbu = 0.5_wp * ( p_dzr(ji,jj,iku) + p_dzr(ji+1,jj  ,iku) )
651         zbv = 0.5_wp * ( p_dzr(ji,jj,ikv) + p_dzr(ji  ,jj+1,ikv) )
652         !                        !- horizontal density gradient at u- & v-points
653         zau = p_gru(ji,jj,iku) * r1_e1u(ji,jj)
654         zav = p_grv(ji,jj,ikv) * r1_e2v(ji,jj)
655         !                        !- bound the slopes: abs(zw.)<= 1/100 and zb..<0
656         !                           kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
657         zbu = MIN(  zbu , - z1_slpmax * ABS( zau ) , -7.e+3_wp/e3u(ji,jj,iku,Kmm)* ABS( zau )  )
658         zbv = MIN(  zbv , - z1_slpmax * ABS( zav ) , -7.e+3_wp/e3v(ji,jj,ikv,Kmm)* ABS( zav )  )
659         !                        !- Slope at u- & v-points (uslpml, vslpml)
660         uslpml(ji,jj) = zau / ( zbu - zeps ) * umask(ji,jj,iku)
661         vslpml(ji,jj) = zav / ( zbv - zeps ) * vmask(ji,jj,ikv)
662         !
663         !                        !==   i- & j-slopes at w-points just below the Mixed Layer   ==!
664         !
665         ik   = MIN( nmln(ji,jj) + 1, jpk )
666         ikm1 = MAX( 1, ik-1 )
667         !                        !- vertical density gradient for w-slope (from N^2)
668         zbw = zm1_2g * pn2 (ji,jj,ik) * ( prd (ji,jj,ik) + prd (ji,jj,ikm1) + 2. )
669         !                        !- horizontal density i- & j-gradient at w-points
670         zci = MAX(   umask(ji-1,jj,ik  ) + umask(ji,jj,ik  )           &
671            &       + umask(ji-1,jj,ikm1) + umask(ji,jj,ikm1) , zeps  ) * e1t(ji,jj)
672         zcj = MAX(   vmask(ji,jj-1,ik  ) + vmask(ji,jj,ik  )           &
673            &       + vmask(ji,jj-1,ikm1) + vmask(ji,jj,ikm1) , zeps  ) * e2t(ji,jj)
674         zai =    (   p_gru(ji-1,jj,ik  ) + p_gru(ji,jj,ik)             &
675            &       + p_gru(ji-1,jj,ikm1) + p_gru(ji,jj,ikm1  )  ) / zci  * tmask(ji,jj,ik)
676         zaj =    (   p_grv(ji,jj-1,ik  ) + p_grv(ji,jj,ik  )           &
677            &       + p_grv(ji,jj-1,ikm1) + p_grv(ji,jj,ikm1)  ) / zcj  * tmask(ji,jj,ik)
678         !                        !- bound the slopes: abs(zw.)<= 1/100 and zb..<0.
679         !                           kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
680         zbi = MIN(  zbw , -100._wp* ABS( zai ) , -7.e+3_wp/e3w(ji,jj,ik,Kmm)* ABS( zai )  )
681         zbj = MIN(  zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/e3w(ji,jj,ik,Kmm)* ABS( zaj )  )
682         !                        !- i- & j-slope at w-points (wslpiml, wslpjml)
683         wslpiml(ji,jj) = zai / ( zbi - zeps ) * tmask (ji,jj,ik)
684         wslpjml(ji,jj) = zaj / ( zbj - zeps ) * tmask (ji,jj,ik)
685      END_2D
686      !!gm this lbc_lnk should be useless....
687      CALL lbc_lnk_multi( 'ldfslp', uslpml , 'U', -1. , vslpml , 'V', -1. , wslpiml, 'W', -1. , wslpjml, 'W', -1. ) 
688      !
689   END SUBROUTINE ldf_slp_mxl
690
691
692   SUBROUTINE ldf_slp_init
693      !!----------------------------------------------------------------------
694      !!                  ***  ROUTINE ldf_slp_init  ***
695      !!
696      !! ** Purpose :   Initialization for the isopycnal slopes computation
697      !!
698      !! ** Method  :   
699      !!----------------------------------------------------------------------
700      INTEGER ::   ji, jj, jk   ! dummy loop indices
701      INTEGER ::   ierr         ! local integer
702      !!----------------------------------------------------------------------
703      !
704      IF(lwp) THEN
705         WRITE(numout,*)
706         WRITE(numout,*) 'ldf_slp_init : direction of lateral mixing'
707         WRITE(numout,*) '~~~~~~~~~~~~'
708      ENDIF
709      !
710      ALLOCATE( ah_wslp2(jpi,jpj,jpk) , akz(jpi,jpj,jpk) , STAT=ierr )
711      IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate ah_slp2 or akz' )
712      !
713      IF( ln_traldf_triad ) THEN        ! Griffies operator : triad of slopes
714         IF(lwp) WRITE(numout,*) '   ==>>>   triad) operator (Griffies)'
715         ALLOCATE( triadi_g(jpi,jpj,jpk,0:1,0:1) , triadj_g(jpi,jpj,jpk,0:1,0:1) ,     &
716            &      triadi  (jpi,jpj,jpk,0:1,0:1) , triadj  (jpi,jpj,jpk,0:1,0:1) ,     &
717            &      wslp2   (jpi,jpj,jpk)                                         , STAT=ierr )
718         IF( ierr > 0      )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Griffies operator slope' )
719         IF( ln_dynldf_iso )   CALL ctl_stop( 'ldf_slp_init: Griffies operator on momentum not supported' )
720         !
721      ELSE                             ! Madec operator : slopes at u-, v-, and w-points
722         IF(lwp) WRITE(numout,*) '   ==>>>   iso operator (Madec)'
723         ALLOCATE( omlmask(jpi,jpj,jpk) ,                                                                        &
724            &      uslp(jpi,jpj,jpk) , uslpml(jpi,jpj) , wslpi(jpi,jpj,jpk) , wslpiml(jpi,jpj) ,     &
725            &      vslp(jpi,jpj,jpk) , vslpml(jpi,jpj) , wslpj(jpi,jpj,jpk) , wslpjml(jpi,jpj) , STAT=ierr )
726         IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Madec operator slope ' )
727
728         ! Direction of lateral diffusion (tracers and/or momentum)
729         ! ------------------------------
730         uslp (:,:,:) = 0._wp   ;   uslpml (:,:) = 0._wp      ! set the slope to zero (even in s-coordinates)
731         vslp (:,:,:) = 0._wp   ;   vslpml (:,:) = 0._wp
732         wslpi(:,:,:) = 0._wp   ;   wslpiml(:,:) = 0._wp
733         wslpj(:,:,:) = 0._wp   ;   wslpjml(:,:) = 0._wp
734
735         !!gm I no longer understand this.....
736!!gm         IF( (ln_traldf_hor .OR. ln_dynldf_hor) .AND. .NOT. (.NOT.ln_linssh .AND. ln_rstart) ) THEN
737!            IF(lwp)   WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces'
738!
739!            ! geopotential diffusion in s-coordinates on tracers and/or momentum
740!            ! The slopes of s-surfaces are computed once (no call to ldfslp in step)
741!            ! The slopes for momentum diffusion are i- or j- averaged of those on tracers
742!
743!            ! set the slope of diffusion to the slope of s-surfaces
744!            !      ( c a u t i o n : minus sign as dep has positive value )
745!            DO jk = 1, jpk
746!               DO jj = 2, jpjm1
747!                  DO ji = fs_2, fs_jpim1   ! vector opt.
748!                     uslp (ji,jj,jk) = - ( gdept(ji+1,jj,jk,Kmm) - gdept(ji ,jj ,jk,Kmm) ) * r1_e1u(ji,jj) * umask(ji,jj,jk)
749!                     vslp (ji,jj,jk) = - ( gdept(ji,jj+1,jk,Kmm) - gdept(ji ,jj ,jk,Kmm) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk)
750!                     wslpi(ji,jj,jk) = - ( gdepw(ji+1,jj,jk,Kmm) - gdepw(ji-1,jj,jk,Kmm) ) * r1_e1t(ji,jj) * wmask(ji,jj,jk) * 0.5
751!                     wslpj(ji,jj,jk) = - ( gdepw(ji,jj+1,jk,Kmm) - gdepw(ji,jj-1,jk,Kmm) ) * r1_e2t(ji,jj) * wmask(ji,jj,jk) * 0.5
752!                  END DO
753!               END DO
754!            END DO
755!            CALL lbc_lnk_multi( 'ldfslp', uslp , 'U', -1. ; CALL lbc_lnk( 'ldfslp', vslp , 'V', -1.,  wslpi, 'W', -1.,  wslpj, 'W', -1. )
756!!gm         ENDIF
757      ENDIF
758      !
759   END SUBROUTINE ldf_slp_init
760
761   !!======================================================================
762END MODULE ldfslp
Note: See TracBrowser for help on using the repository browser.