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

source: trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 4605

Last change on this file since 4605 was 4488, checked in by rfurner, 10 years ago

fixes to enable proper calulation of slopes for geopotential diffusion in scoords

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