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 @ 12252

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

rev12240_dev_r11943_MERGE_2019: same as [12251], merge trunk 12072:12248, all sette tests ok, GYRE_PISCES, AMM12, ISOMIP, VORTEX intentical to 12236

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