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

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 5866

Last change on this file since 5866 was 5866, checked in by gm, 8 years ago

#1613: vvl by default: add ln_linssh and remove key_vvl

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