New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
ldfslp.F90 in branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/RK3_SRC/LDF – NEMO

source: branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/RK3_SRC/LDF/ldfslp.F90 @ 8568

Last change on this file since 8568 was 8568, checked in by gm, 7 years ago

#1911 (ENHANCE-09): PART I.2 - _NONE option + remove zts + see associated wiki page

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