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_merge_2017/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 9124

Last change on this file since 9124 was 9124, checked in by gm, 6 years ago

dev_merge_2017: ln_timing instead of nn_timing + restricted timing to nemo_init and routine called by step in OPA_SRC

  • Property svn:keywords set to Id
File size: 48.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 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_multi( zwz, 'U', -1.,  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_multi( zwz, 'T', -1.,  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_multi( uslp , 'U', -1. , vslp , 'V', -1. , wslpi, 'W', -1., wslpj, 'W', -1. )
365
366      IF(ln_ctl) THEN
367         CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp  - u : ', tab3d_2=vslp,  clinfo2=' v : ', kdim=jpk)
368         CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk)
369      ENDIF
370      !
371      IF( ln_timing )   CALL timing_stop('ldf_slp')
372      !
373   END SUBROUTINE ldf_slp
374
375
376   SUBROUTINE ldf_slp_triad ( kt )
377      !!----------------------------------------------------------------------
378      !!                 ***  ROUTINE ldf_slp_triad  ***
379      !!
380      !! ** Purpose :   Compute the squared slopes of neutral surfaces (slope
381      !!      of iso-pycnal surfaces referenced locally) (ln_traldf_triad=T)
382      !!      at W-points using the Griffies quarter-cells.
383      !!
384      !! ** Method  :   calculates alpha and beta at T-points
385      !!
386      !! ** Action : - triadi_g, triadj_g   T-pts i- and j-slope triads relative to geopot. (used for eiv)
387      !!             - triadi , triadj    T-pts i- and j-slope triads relative to model-coordinate
388      !!             - wslp2              squared slope of neutral surfaces at w-points.
389      !!----------------------------------------------------------------------
390      INTEGER, INTENT( in ) ::   kt             ! ocean time-step index
391      !!
392      INTEGER  ::   ji, jj, jk, jl, ip, jp, kp  ! dummy loop indices
393      INTEGER  ::   iku, ikv                    ! local integer
394      REAL(wp) ::   zfacti, zfactj              ! local scalars
395      REAL(wp) ::   znot_thru_surface           ! local scalars
396      REAL(wp) ::   zdit, zdis, zdkt, zbu, zbti, zisw
397      REAL(wp) ::   zdjt, zdjs, zdks, zbv, zbtj, zjsw
398      REAL(wp) ::   zdxrho_raw, zti_coord, zti_raw, zti_lim, zti_g_raw, zti_g_lim
399      REAL(wp) ::   zdyrho_raw, ztj_coord, ztj_raw, ztj_lim, ztj_g_raw, ztj_g_lim
400      REAL(wp) ::   zdzrho_raw
401      REAL(wp) ::   zbeta0, ze3_e1, ze3_e2
402      REAL(wp), DIMENSION(jpi,jpj)     ::   z1_mlbw
403      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zalbet
404      REAL(wp), DIMENSION(jpi,jpj,jpk,0:1) ::   zdxrho , zdyrho, zdzrho     ! Horizontal and vertical density gradients
405      REAL(wp), DIMENSION(jpi,jpj,0:1,0:1) ::   zti_mlb, ztj_mlb            ! for Griffies operator only
406      !!----------------------------------------------------------------------
407      !
408      IF( ln_timing )   CALL timing_start('ldf_slp_triad')
409      !
410      !--------------------------------!
411      !  Some preliminary calculation  !
412      !--------------------------------!
413      !
414      DO jl = 0, 1                            !==  unmasked before density i- j-, k-gradients  ==!
415         !
416         ip = jl   ;   jp = jl                ! guaranteed nonzero gradients ( absolute value larger than repsln)
417         DO jk = 1, jpkm1                     ! done each pair of triad
418            DO jj = 1, jpjm1                  ! NB: not masked ==>  a minimum value is set
419               DO ji = 1, fs_jpim1            ! vector opt.
420                  zdit = ( tsb(ji+1,jj,jk,jp_tem) - tsb(ji,jj,jk,jp_tem) )    ! i-gradient of T & S at u-point
421                  zdis = ( tsb(ji+1,jj,jk,jp_sal) - tsb(ji,jj,jk,jp_sal) )
422                  zdjt = ( tsb(ji,jj+1,jk,jp_tem) - tsb(ji,jj,jk,jp_tem) )    ! j-gradient of T & S at v-point
423                  zdjs = ( tsb(ji,jj+1,jk,jp_sal) - tsb(ji,jj,jk,jp_sal) )
424                  zdxrho_raw = ( - rab_b(ji+ip,jj   ,jk,jp_tem) * zdit + rab_b(ji+ip,jj   ,jk,jp_sal) * zdis ) * r1_e1u(ji,jj)
425                  zdyrho_raw = ( - rab_b(ji   ,jj+jp,jk,jp_tem) * zdjt + rab_b(ji   ,jj+jp,jk,jp_sal) * zdjs ) * r1_e2v(ji,jj)
426                  zdxrho(ji+ip,jj   ,jk,1-ip) = SIGN(  MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw  )   ! keep the sign
427                  zdyrho(ji   ,jj+jp,jk,1-jp) = SIGN(  MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw  )
428               END DO
429            END DO
430         END DO
431         !
432         IF( ln_zps .AND. l_grad_zps ) THEN     ! partial steps: correction of i- & j-grad on bottom
433            DO jj = 1, jpjm1
434               DO ji = 1, jpim1
435                  iku  = mbku(ji,jj)          ;   ikv  = mbkv(ji,jj)             ! last ocean level (u- & v-points)
436                  zdit = gtsu(ji,jj,jp_tem)   ;   zdjt = gtsv(ji,jj,jp_tem)      ! i- & j-gradient of Temperature
437                  zdis = gtsu(ji,jj,jp_sal)   ;   zdjs = gtsv(ji,jj,jp_sal)      ! i- & j-gradient of Salinity
438                  zdxrho_raw = ( - rab_b(ji+ip,jj   ,iku,jp_tem) * zdit + rab_b(ji+ip,jj   ,iku,jp_sal) * zdis ) * r1_e1u(ji,jj)
439                  zdyrho_raw = ( - rab_b(ji   ,jj+jp,ikv,jp_tem) * zdjt + rab_b(ji   ,jj+jp,ikv,jp_sal) * zdjs ) * r1_e2v(ji,jj)
440                  zdxrho(ji+ip,jj   ,iku,1-ip) = SIGN( MAX( repsln, ABS( zdxrho_raw ) ), zdxrho_raw )   ! keep the sign
441                  zdyrho(ji   ,jj+jp,ikv,1-jp) = SIGN( MAX( repsln, ABS( zdyrho_raw ) ), zdyrho_raw )
442               END DO
443            END DO
444         ENDIF
445         !
446      END DO
447
448      DO kp = 0, 1                            !==  unmasked before density i- j-, k-gradients  ==!
449         DO jk = 1, jpkm1                     ! done each pair of triad
450            DO jj = 1, jpj                    ! NB: not masked ==>  a minimum value is set
451               DO ji = 1, jpi                 ! vector opt.
452                  IF( jk+kp > 1 ) THEN        ! k-gradient of T & S a jk+kp
453                     zdkt = ( tsb(ji,jj,jk+kp-1,jp_tem) - tsb(ji,jj,jk+kp,jp_tem) )
454                     zdks = ( tsb(ji,jj,jk+kp-1,jp_sal) - tsb(ji,jj,jk+kp,jp_sal) )
455                  ELSE
456                     zdkt = 0._wp                                             ! 1st level gradient set to zero
457                     zdks = 0._wp
458                  ENDIF
459                  zdzrho_raw = ( - zalbet(ji,jj,jk) * zdkt + zbeta0*zdks ) / e3w_n(ji,jj,jk+kp)
460                  zdzrho(ji,jj,jk,kp) = - MIN( - repsln , zdzrho_raw )    ! force zdzrho >= repsln
461                 END DO
462            END DO
463         END DO
464      END DO
465      !
466      DO jj = 1, jpj                          !==  Reciprocal depth of the w-point below ML base  ==!
467         DO ji = 1, jpi
468            jk = MIN( nmln(ji,jj), mbkt(ji,jj) ) + 1     ! MIN in case ML depth is the ocean depth
469            z1_mlbw(ji,jj) = 1._wp / gdepw_n(ji,jj,jk)
470         END DO
471      END DO
472      !
473      !                                       !==  intialisations to zero  ==!
474      !
475      wslp2  (:,:,:)     = 0._wp              ! wslp2 will be cumulated 3D field set to zero
476      triadi_g(:,:,1,:,:) = 0._wp   ;   triadi_g(:,:,jpk,:,:) = 0._wp   ! set surface and bottom slope to zero
477      triadj_g(:,:,1,:,:) = 0._wp   ;   triadj_g(:,:,jpk,:,:) = 0._wp
478      !!gm _iso set to zero missing
479      triadi  (:,:,1,:,:) = 0._wp   ;   triadj  (:,:,jpk,:,:) = 0._wp   ! set surface and bottom slope to zero
480      triadj  (:,:,1,:,:) = 0._wp   ;   triadj  (:,:,jpk,:,:) = 0._wp
481
482      !-------------------------------------!
483      !  Triads just below the Mixed Layer  !
484      !-------------------------------------!
485      !
486      DO jl = 0, 1                            ! calculate slope of the 4 triads immediately ONE level below mixed-layer base
487         DO kp = 0, 1                         ! with only the slope-max limit   and   MASKED
488            DO jj = 1, jpjm1
489               DO ji = 1, fs_jpim1
490                  ip = jl   ;   jp = jl
491                  !
492                  jk = nmln(ji+ip,jj) + 1
493                  IF( jk > mbkt(ji+ip,jj) ) THEN   ! ML reaches bottom
494                     zti_mlb(ji+ip,jj   ,1-ip,kp) = 0.0_wp
495                  ELSE                             
496                     ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth)
497                     zti_g_raw = (  zdxrho(ji+ip,jj,jk-kp,1-ip) / zdzrho(ji+ip,jj,jk-kp,kp)      &
498                        &          - ( gdept_n(ji+1,jj,jk-kp) - gdept_n(ji,jj,jk-kp) ) * r1_e1u(ji,jj)  ) * umask(ji,jj,jk)
499                     ze3_e1    =  e3w_n(ji+ip,jj,jk-kp) * r1_e1u(ji,jj) 
500                     zti_mlb(ji+ip,jj   ,1-ip,kp) = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e1  , ABS( zti_g_raw ) ), zti_g_raw )
501                  ENDIF
502                  !
503                  jk = nmln(ji,jj+jp) + 1
504                  IF( jk >  mbkt(ji,jj+jp) ) THEN  !ML reaches bottom
505                     ztj_mlb(ji   ,jj+jp,1-jp,kp) = 0.0_wp
506                  ELSE
507                     ztj_g_raw = (  zdyrho(ji,jj+jp,jk-kp,1-jp) / zdzrho(ji,jj+jp,jk-kp,kp)      &
508                        &      - ( gdept_n(ji,jj+1,jk-kp) - gdept_n(ji,jj,jk-kp) ) / e2v(ji,jj)  ) * vmask(ji,jj,jk)
509                     ze3_e2    =  e3w_n(ji,jj+jp,jk-kp) / e2v(ji,jj)
510                     ztj_mlb(ji   ,jj+jp,1-jp,kp) = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e2  , ABS( ztj_g_raw ) ), ztj_g_raw )
511                  ENDIF
512               END DO
513            END DO
514         END DO
515      END DO
516
517      !-------------------------------------!
518      !  Triads with surface limits         !
519      !-------------------------------------!
520      !
521      DO kp = 0, 1                            ! k-index of triads
522         DO jl = 0, 1
523            ip = jl   ;   jp = jl             ! i- and j-indices of triads (i-k and j-k planes)
524            DO jk = 1, jpkm1
525               ! Must mask contribution to slope from dz/dx at constant s for triads jk=1,kp=0 that poke up though ocean surface
526               znot_thru_surface = REAL( 1-1/(jk+kp), wp )  !jk+kp=1,=0.; otherwise=1.0
527               DO jj = 1, jpjm1
528                  DO ji = 1, fs_jpim1         ! vector opt.
529                     !
530                     ! Calculate slope relative to geopotentials used for GM skew fluxes
531                     ! Add s-coordinate slope at t-points (do this by *subtracting* gradient of depth)
532                     ! Limit by slope *relative to geopotentials* by rn_slpmax, and mask by psi-point
533                     ! masked by umask taken at the level of dz(rho)
534                     !
535                     ! raw slopes: unmasked unbounded slopes (relative to geopotential (zti_g) and model surface (zti)
536                     !
537                     zti_raw   = zdxrho(ji+ip,jj   ,jk,1-ip) / zdzrho(ji+ip,jj   ,jk,kp)                   ! unmasked
538                     ztj_raw   = zdyrho(ji   ,jj+jp,jk,1-jp) / zdzrho(ji   ,jj+jp,jk,kp)
539                     !
540                     ! Must mask contribution to slope for triad jk=1,kp=0 that poke up though ocean surface
541                     zti_coord = znot_thru_surface * ( gdept_n(ji+1,jj  ,jk) - gdept_n(ji,jj,jk) ) * r1_e1u(ji,jj)
542                     ztj_coord = znot_thru_surface * ( gdept_n(ji  ,jj+1,jk) - gdept_n(ji,jj,jk) ) * r1_e2v(ji,jj)     ! unmasked
543                     zti_g_raw = zti_raw - zti_coord      ! ref to geopot surfaces
544                     ztj_g_raw = ztj_raw - ztj_coord
545                     ! additional limit required in bilaplacian case
546                     ze3_e1    = e3w_n(ji+ip,jj   ,jk+kp) * r1_e1u(ji,jj)
547                     ze3_e2    = e3w_n(ji   ,jj+jp,jk+kp) * r1_e2v(ji,jj)
548                     ! NB: hard coded factor 5 (can be a namelist parameter...)
549                     zti_g_lim = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e1, ABS( zti_g_raw ) ), zti_g_raw )
550                     ztj_g_lim = SIGN( MIN( rn_slpmax, 5.0_wp * ze3_e2, ABS( ztj_g_raw ) ), ztj_g_raw )
551                     !
552                     ! Below  ML use limited zti_g as is & mask
553                     ! Inside ML replace by linearly reducing sx_mlb towards surface & mask
554                     !
555                     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
556                     zfactj = REAL( 1 - 1/(1 + (jk+kp-1)/nmln(ji,jj+jp)), wp )  ! must be .ge. nmln(ji,jj) for zfact=1
557                     !                                                          !                   otherwise  zfact=0
558                     zti_g_lim =          ( zfacti   * zti_g_lim                       &
559                        &      + ( 1._wp - zfacti ) * zti_mlb(ji+ip,jj,1-ip,kp)   &
560                        &                           * gdepw_n(ji+ip,jj,jk+kp) * z1_mlbw(ji+ip,jj) ) * umask(ji,jj,jk+kp)
561                     ztj_g_lim =          ( zfactj   * ztj_g_lim                       &
562                        &      + ( 1._wp - zfactj ) * ztj_mlb(ji,jj+jp,1-jp,kp)   &
563                        &                           * gdepw_n(ji,jj+jp,jk+kp) * z1_mlbw(ji,jj+jp) ) * vmask(ji,jj,jk+kp)
564                     !
565                     triadi_g(ji+ip,jj   ,jk,1-ip,kp) = zti_g_lim
566                     triadj_g(ji   ,jj+jp,jk,1-jp,kp) = ztj_g_lim
567                     !
568                     ! Get coefficients of isoneutral diffusion tensor
569                     ! 1. Utilise gradients *relative* to s-coordinate, so add t-point slopes (*subtract* depth gradients)
570                     ! 2. We require that isoneutral diffusion  gives no vertical buoyancy flux
571                     !     i.e. 33 term = (real slope* 31, 13 terms)
572                     ! To do this, retain limited sx**2  in vertical flux, but divide by real slope for 13/31 terms
573                     ! Equivalent to tapering A_iso = sx_limited**2/(real slope)**2
574                     !
575                     zti_lim  = ( zti_g_lim + zti_coord ) * umask(ji,jj,jk+kp)    ! remove coordinate slope => relative to coordinate surfaces
576                     ztj_lim  = ( ztj_g_lim + ztj_coord ) * vmask(ji,jj,jk+kp)
577                     !
578                     IF( ln_triad_iso ) THEN
579                        zti_raw = zti_lim*zti_lim / zti_raw
580                        ztj_raw = ztj_lim*ztj_lim / ztj_raw
581                        zti_raw = SIGN( MIN( ABS(zti_lim), ABS( zti_raw ) ), zti_raw )
582                        ztj_raw = SIGN( MIN( ABS(ztj_lim), ABS( ztj_raw ) ), ztj_raw )
583                        zti_lim = zfacti * zti_lim + ( 1._wp - zfacti ) * zti_raw
584                        ztj_lim = zfactj * ztj_lim + ( 1._wp - zfactj ) * ztj_raw
585                     ENDIF
586                     !                                      ! switching triad scheme
587                     zisw = (rn_sw_triad - 1._wp ) + rn_sw_triad    &
588                        &            * 2._wp * ABS( 0.5_wp - kp - ( 0.5_wp - ip ) * SIGN( 1._wp , zdxrho(ji+ip,jj,jk,1-ip) )  )
589                     zjsw = (rn_sw_triad - 1._wp ) + rn_sw_triad    &
590                        &            * 2._wp * ABS( 0.5_wp - kp - ( 0.5_wp - jp ) * SIGN( 1._wp , zdyrho(ji,jj+jp,jk,1-jp) )  )
591                     !
592                     triadi(ji+ip,jj   ,jk,1-ip,kp) = zti_lim * zisw
593                     triadj(ji   ,jj+jp,jk,1-jp,kp) = ztj_lim * zjsw
594                     !
595                     zbu  = e1e2u(ji   ,jj   ) * e3u_n(ji   ,jj   ,jk   )
596                     zbv  = e1e2v(ji   ,jj   ) * e3v_n(ji   ,jj   ,jk   )
597                     zbti = e1e2t(ji+ip,jj   ) * e3w_n(ji+ip,jj   ,jk+kp)
598                     zbtj = e1e2t(ji   ,jj+jp) * e3w_n(ji   ,jj+jp,jk+kp)
599                     !
600                     wslp2(ji+ip,jj,jk+kp) = wslp2(ji+ip,jj,jk+kp) + 0.25_wp * zbu / zbti * zti_g_lim*zti_g_lim      ! masked
601                     wslp2(ji,jj+jp,jk+kp) = wslp2(ji,jj+jp,jk+kp) + 0.25_wp * zbv / zbtj * ztj_g_lim*ztj_g_lim
602                  END DO
603               END DO
604            END DO
605         END DO
606      END DO
607      !
608      wslp2(:,:,1) = 0._wp                ! force the surface wslp to zero
609
610      CALL lbc_lnk( wslp2, 'W', 1. )      ! lateral boundary confition on wslp2 only   ==>>> gm : necessary ? to be checked
611      !
612      IF( ln_timing )   CALL timing_stop('ldf_slp_triad')
613      !
614   END SUBROUTINE ldf_slp_triad
615
616
617   SUBROUTINE ldf_slp_mxl( prd, pn2, p_gru, p_grv, p_dzr )
618      !!----------------------------------------------------------------------
619      !!                  ***  ROUTINE ldf_slp_mxl  ***
620      !!
621      !! ** Purpose :   Compute the slopes of iso-neutral surface just below
622      !!              the mixed layer.
623      !!
624      !! ** Method  :   The slope in the i-direction is computed at u- & w-points
625      !!              (uslpml, wslpiml) and the slope in the j-direction is computed
626      !!              at v- and w-points (vslpml, wslpjml) with the same bounds as
627      !!              in ldf_slp.
628      !!
629      !! ** Action  :   uslpml, wslpiml :  i- &  j-slopes of neutral surfaces
630      !!                vslpml, wslpjml    just below the mixed layer
631      !!                omlmask         :  mixed layer mask
632      !!----------------------------------------------------------------------
633      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   prd            ! in situ density
634      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   pn2            ! Brunt-Vaisala frequency (locally ref.)
635      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_gru, p_grv   ! i- & j-gradient of density (u- & v-pts)
636      REAL(wp), DIMENSION(:,:,:), INTENT(in) ::   p_dzr          ! z-gradient of density      (T-point)
637      !!
638      INTEGER  ::   ji , jj , jk                   ! dummy loop indices
639      INTEGER  ::   iku, ikv, ik, ikm1             ! local integers
640      REAL(wp) ::   zeps, zm1_g, zm1_2g, z1_slpmax ! local scalars
641      REAL(wp) ::   zci, zfi, zau, zbu, zai, zbi   !   -      -
642      REAL(wp) ::   zcj, zfj, zav, zbv, zaj, zbj   !   -      -
643      REAL(wp) ::   zck, zfk,      zbw             !   -      -
644      !!----------------------------------------------------------------------
645      !
646      zeps   =  1.e-20_wp        !==   Local constant initialization   ==!
647      zm1_g  = -1.0_wp / grav
648      zm1_2g = -0.5_wp / grav
649      z1_slpmax = 1._wp / rn_slpmax
650      !
651      uslpml (1,:) = 0._wp      ;      uslpml (jpi,:) = 0._wp
652      vslpml (1,:) = 0._wp      ;      vslpml (jpi,:) = 0._wp
653      wslpiml(1,:) = 0._wp      ;      wslpiml(jpi,:) = 0._wp
654      wslpjml(1,:) = 0._wp      ;      wslpjml(jpi,:) = 0._wp
655      !
656      !                                            !==   surface mixed layer mask   !
657      DO jk = 1, jpk                               ! =1 inside the mixed layer, =0 otherwise
658         DO jj = 1, jpj
659            DO ji = 1, jpi
660               ik = nmln(ji,jj) - 1
661               IF( jk <= ik ) THEN   ;   omlmask(ji,jj,jk) = 1._wp
662               ELSE                  ;   omlmask(ji,jj,jk) = 0._wp
663               ENDIF
664            END DO
665         END DO
666      END DO
667
668
669      ! Slopes of isopycnal surfaces just before bottom of mixed layer
670      ! --------------------------------------------------------------
671      ! The slope are computed as in the 3D case.
672      ! A key point here is the definition of the mixed layer at u- and v-points.
673      ! It is assumed to be the maximum of the two neighbouring T-point mixed layer depth.
674      ! Otherwise, a n2 value inside the mixed layer can be involved in the computation
675      ! of the slope, resulting in a too steep diagnosed slope and thus a spurious eddy
676      ! induce velocity field near the base of the mixed layer.
677      !-----------------------------------------------------------------------
678      !
679      DO jj = 2, jpjm1
680         DO ji = 2, jpim1
681            !                        !==   Slope at u- & v-points just below the Mixed Layer   ==!
682            !
683            !                        !- vertical density gradient for u- and v-slopes (from dzr at T-point)
684            iku = MIN(  MAX( 1, nmln(ji,jj) , nmln(ji+1,jj) ) , jpkm1  )   ! ML (MAX of T-pts, bound by jpkm1)
685            ikv = MIN(  MAX( 1, nmln(ji,jj) , nmln(ji,jj+1) ) , jpkm1  )   !
686            zbu = 0.5_wp * ( p_dzr(ji,jj,iku) + p_dzr(ji+1,jj  ,iku) )
687            zbv = 0.5_wp * ( p_dzr(ji,jj,ikv) + p_dzr(ji  ,jj+1,ikv) )
688            !                        !- horizontal density gradient at u- & v-points
689            zau = p_gru(ji,jj,iku) * r1_e1u(ji,jj)
690            zav = p_grv(ji,jj,ikv) * r1_e2v(ji,jj)
691            !                        !- bound the slopes: abs(zw.)<= 1/100 and zb..<0
692            !                           kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
693            zbu = MIN(  zbu , - z1_slpmax * ABS( zau ) , -7.e+3_wp/e3u_n(ji,jj,iku)* ABS( zau )  )
694            zbv = MIN(  zbv , - z1_slpmax * ABS( zav ) , -7.e+3_wp/e3v_n(ji,jj,ikv)* ABS( zav )  )
695            !                        !- Slope at u- & v-points (uslpml, vslpml)
696            uslpml(ji,jj) = zau / ( zbu - zeps ) * umask(ji,jj,iku)
697            vslpml(ji,jj) = zav / ( zbv - zeps ) * vmask(ji,jj,ikv)
698            !
699            !                        !==   i- & j-slopes at w-points just below the Mixed Layer   ==!
700            !
701            ik   = MIN( nmln(ji,jj) + 1, jpk )
702            ikm1 = MAX( 1, ik-1 )
703            !                        !- vertical density gradient for w-slope (from N^2)
704            zbw = zm1_2g * pn2 (ji,jj,ik) * ( prd (ji,jj,ik) + prd (ji,jj,ikm1) + 2. )
705            !                        !- horizontal density i- & j-gradient at w-points
706            zci = MAX(   umask(ji-1,jj,ik  ) + umask(ji,jj,ik  )           &
707               &       + umask(ji-1,jj,ikm1) + umask(ji,jj,ikm1) , zeps  ) * e1t(ji,jj)
708            zcj = MAX(   vmask(ji,jj-1,ik  ) + vmask(ji,jj,ik  )           &
709               &       + vmask(ji,jj-1,ikm1) + vmask(ji,jj,ikm1) , zeps  ) * e2t(ji,jj)
710            zai =    (   p_gru(ji-1,jj,ik  ) + p_gru(ji,jj,ik)             &
711               &       + p_gru(ji-1,jj,ikm1) + p_gru(ji,jj,ikm1  )  ) / zci  * tmask(ji,jj,ik)
712            zaj =    (   p_grv(ji,jj-1,ik  ) + p_grv(ji,jj,ik  )           &
713               &       + p_grv(ji,jj-1,ikm1) + p_grv(ji,jj,ikm1)  ) / zcj  * tmask(ji,jj,ik)
714            !                        !- bound the slopes: abs(zw.)<= 1/100 and zb..<0.
715            !                           kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
716            zbi = MIN(  zbw , -100._wp* ABS( zai ) , -7.e+3_wp/e3w_n(ji,jj,ik)* ABS( zai )  )
717            zbj = MIN(  zbw , -100._wp* ABS( zaj ) , -7.e+3_wp/e3w_n(ji,jj,ik)* ABS( zaj )  )
718            !                        !- i- & j-slope at w-points (wslpiml, wslpjml)
719            wslpiml(ji,jj) = zai / ( zbi - zeps ) * tmask (ji,jj,ik)
720            wslpjml(ji,jj) = zaj / ( zbj - zeps ) * tmask (ji,jj,ik)
721         END DO
722      END DO
723      !!gm this lbc_lnk should be useless....
724      CALL lbc_lnk_multi( uslpml , 'U', -1. , vslpml , 'V', -1. , wslpiml, 'W', -1. , wslpjml, 'W', -1. ) 
725      !
726   END SUBROUTINE ldf_slp_mxl
727
728
729   SUBROUTINE ldf_slp_init
730      !!----------------------------------------------------------------------
731      !!                  ***  ROUTINE ldf_slp_init  ***
732      !!
733      !! ** Purpose :   Initialization for the isopycnal slopes computation
734      !!
735      !! ** Method  :   
736      !!----------------------------------------------------------------------
737      INTEGER ::   ji, jj, jk   ! dummy loop indices
738      INTEGER ::   ierr         ! local integer
739      !!----------------------------------------------------------------------
740      !
741      IF(lwp) THEN
742         WRITE(numout,*)
743         WRITE(numout,*) 'ldf_slp_init : direction of lateral mixing'
744         WRITE(numout,*) '~~~~~~~~~~~~'
745      ENDIF
746      !
747      ALLOCATE( ah_wslp2(jpi,jpj,jpk) , akz(jpi,jpj,jpk) , STAT=ierr )
748      IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate ah_slp2 or akz' )
749      !
750      IF( ln_traldf_triad ) THEN        ! Griffies operator : triad of slopes
751         IF(lwp) WRITE(numout,*) '              Griffies (triad) operator initialisation'
752         ALLOCATE( triadi_g(jpi,jpj,jpk,0:1,0:1) , triadj_g(jpi,jpj,jpk,0:1,0:1) ,     &
753            &      triadi  (jpi,jpj,jpk,0:1,0:1) , triadj  (jpi,jpj,jpk,0:1,0:1) ,     &
754            &      wslp2   (jpi,jpj,jpk)                                         , STAT=ierr )
755         IF( ierr > 0      )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Griffies operator slope' )
756         IF( ln_dynldf_iso )   CALL ctl_stop( 'ldf_slp_init: Griffies operator on momentum not supported' )
757         !
758      ELSE                             ! Madec operator : slopes at u-, v-, and w-points
759         IF(lwp) WRITE(numout,*) '              Madec operator initialisation'
760         ALLOCATE( omlmask(jpi,jpj,jpk) ,                                                                        &
761            &      uslp(jpi,jpj,jpk) , uslpml(jpi,jpj) , wslpi(jpi,jpj,jpk) , wslpiml(jpi,jpj) ,     &
762            &      vslp(jpi,jpj,jpk) , vslpml(jpi,jpj) , wslpj(jpi,jpj,jpk) , wslpjml(jpi,jpj) , STAT=ierr )
763         IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'ldf_slp_init : unable to allocate Madec operator slope ' )
764
765         ! Direction of lateral diffusion (tracers and/or momentum)
766         ! ------------------------------
767         uslp (:,:,:) = 0._wp   ;   uslpml (:,:) = 0._wp      ! set the slope to zero (even in s-coordinates)
768         vslp (:,:,:) = 0._wp   ;   vslpml (:,:) = 0._wp
769         wslpi(:,:,:) = 0._wp   ;   wslpiml(:,:) = 0._wp
770         wslpj(:,:,:) = 0._wp   ;   wslpjml(:,:) = 0._wp
771
772         !!gm I no longer understand this.....
773!!gm         IF( (ln_traldf_hor .OR. ln_dynldf_hor) .AND. .NOT. (.NOT.ln_linssh .AND. ln_rstart) ) THEN
774!            IF(lwp)   WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces'
775!
776!            ! geopotential diffusion in s-coordinates on tracers and/or momentum
777!            ! The slopes of s-surfaces are computed once (no call to ldfslp in step)
778!            ! The slopes for momentum diffusion are i- or j- averaged of those on tracers
779!
780!            ! set the slope of diffusion to the slope of s-surfaces
781!            !      ( c a u t i o n : minus sign as dep has positive value )
782!            DO jk = 1, jpk
783!               DO jj = 2, jpjm1
784!                  DO ji = fs_2, fs_jpim1   ! vector opt.
785!                     uslp (ji,jj,jk) = - ( gdept_n(ji+1,jj,jk) - gdept_n(ji ,jj ,jk) ) * r1_e1u(ji,jj) * umask(ji,jj,jk)
786!                     vslp (ji,jj,jk) = - ( gdept_n(ji,jj+1,jk) - gdept_n(ji ,jj ,jk) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk)
787!                     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
788!                     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
789!                  END DO
790!               END DO
791!            END DO
792!            CALL lbc_lnk_multi( uslp , 'U', -1. ; CALL lbc_lnk( vslp , 'V', -1.,  wslpi, 'W', -1.,  wslpj, 'W', -1. )
793!!gm         ENDIF
794      ENDIF
795      !
796   END SUBROUTINE ldf_slp_init
797
798   !!======================================================================
799END MODULE ldfslp
Note: See TracBrowser for help on using the repository browser.