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/2021/dev_r14312_MOI-01_MOI-GRIFFIESTRIADSMOOTHING/src/OCE/LDF – NEMO

source: NEMO/branches/2021/dev_r14312_MOI-01_MOI-GRIFFIESTRIADSMOOTHING/src/OCE/LDF/ldfslp.F90 @ 14357

Last change on this file since 14357 was 14357, checked in by cbricaud, 3 years ago

add in dev_r14312_MOI-01_MOI-GRIFFIESTRIADSMOOTHING branch ISF top limit condition for triads (ticket #2601)

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