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_r12072_MERGE_OPTION2_2019/src/OCE/LDF – NEMO

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/OCE/LDF/ldfslp.F90 @ 12251

Last change on this file since 12251 was 12251, checked in by smasson, 4 years ago

rev12232_dev_r12072_MERGE_OPTION2_2019: merge trunk 12072:12248, all sette tests ok, GYRE_PISCES, AMM12, ISOMIP, VORTEX intentical to 12210

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