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

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90

Last change on this file was 11101, checked in by frrh, 5 years ago

Merge changes from Met Office GMED ticket 450 to reduce unnecessary
text output from NEMO.
This output, which is typically not switchable, is rarely of interest
in normal (non-debugging) runs and simply redunantley consumes extra
file space.
Further, the presence of this text output has been shown to
significantly degrade performance of models which are run during
Met Office HPC RAID (disk) checks.
The new code introduces switches which are configurable via the
changes made in the associated Met Office MOCI ticket 399.

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