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

source: branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 5038

Last change on this file since 5038 was 5038, checked in by jamesharle, 9 years ago

Merging branch with HEAD of the trunk

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