New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
ldfslp.F90 in branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 4596

Last change on this file since 4596 was 4596, checked in by gm, 10 years ago

#1260: LDF simplification + bilap iso-neutral for TRA and GYRE

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