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

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 2344

Last change on this file since 2344 was 2344, checked in by gm, 14 years ago

v3.3beta: Suppress obsolete key_mpp_omp

  • Property svn:keywords set to Id
File size: 50.6 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     0.5  ! 2002-10  (G. Madec)  Free form, F90
10   !!            1.0  ! 2005-10  (A. Beckmann)  correction for s-coordinates
11   !!            3.3  ! 2006-10  (C. Harris, G. Nurser)  add ldf_slp_grif (Griffies operator)
12   !!----------------------------------------------------------------------
13#if   defined key_ldfslp   ||   defined key_esopa
14   !!----------------------------------------------------------------------
15   !!   'key_ldfslp'                      Rotation of lateral mixing tensor
16   !!----------------------------------------------------------------------
17   !!   ldf_slp      : compute the slopes of neutral surface
18   !!   ldf_slp_mxl  : compute the slopes of iso-neutral surface
19   !!   ldf_slp_init : initialization of the slopes computation
20   !!----------------------------------------------------------------------
21   USE oce             ! ocean dynamics and tracers
22   USE dom_oce         ! ocean space and time domain
23   USE ldftra_oce
24   USE ldfdyn_oce
25   USE phycst          ! physical constants
26   USE zdfmxl          ! mixed layer depth
27   USE eosbn2
28   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
29   USE in_out_manager  ! I/O manager
30   USE prtctl          ! Print control
31
32   IMPLICIT NONE
33   PRIVATE
34
35   PUBLIC   ldf_slp         ! routine called by step.F90
36   PUBLIC   ldf_slp_init    ! routine called by opa.F90
37   PUBLIC   ldf_slp_grif    ! routine called by step.F90
38
39   LOGICAL , PUBLIC, PARAMETER              ::   lk_ldfslp = .TRUE.   !: slopes flag
40   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   uslp, wslpi          !: i_slope at U- and W-points
41   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   vslp, wslpj          !: j-slope at V- and W-points
42   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   wslp2                !: wslp**2 from Griffies quarter cells
43   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   alpha, beta          !: alpha,beta at T points
44   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   tfw,sfw,ftu,fsu,ftv,fsv,ftud,fsud,ftvd,fsvd
45   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   psix_eiv
46   REAL(wp), PUBLIC, DIMENSION(jpi,jpj,jpk) ::   psiy_eiv
47   
48   REAL(wp), DIMENSION(jpi,jpj,jpk) ::   omlmask           ! mask of the surface mixed layer at T-pt
49   REAL(wp), DIMENSION(jpi,jpj)     ::   uslpml, wslpiml   ! i_slope at U- and W-points just below the mixed layer
50   REAL(wp), DIMENSION(jpi,jpj)     ::   vslpml, wslpjml   ! j_slope at V- and W-points just below the mixed layer
51
52   !! * Substitutions
53#  include "domzgr_substitute.h90"
54#  include "ldftra_substitute.h90"
55#  include "ldfeiv_substitute.h90"
56#  include "vectopt_loop_substitute.h90"
57   !!----------------------------------------------------------------------
58   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
59   !! $Id$
60   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
61   !!----------------------------------------------------------------------
62CONTAINS
63
64   SUBROUTINE ldf_slp( kt, prd, pn2 )
65      !!----------------------------------------------------------------------
66      !!                 ***  ROUTINE ldf_slp  ***
67      !!
68      !! ** Purpose :   Compute the slopes of neutral surface (slope of isopycnal
69      !!              surfaces referenced locally) ('key_traldfiso').
70      !!
71      !! ** Method  :   The slope in the i-direction is computed at U- and
72      !!      W-points (uslp, wslpi) and the slope in the j-direction is
73      !!      computed at V- and W-points (vslp, wslpj).
74      !!      They are bounded by 1/100 over the whole ocean, and within the
75      !!      surface layer they are bounded by the distance to the surface
76      !!      ( slope<= depth/l  where l is the length scale of horizontal
77      !!      diffusion (here, aht=2000m2/s ==> l=20km with a typical velocity
78      !!      of 10cm/s)
79      !!        A horizontal shapiro filter is applied to the slopes
80      !!        ln_sco=T, s-coordinate, add to the previously computed slopes
81      !!      the slope of the model level surface.
82      !!        macro-tasked on horizontal slab (jk-loop)  (2, jpk-1)
83      !!      [slopes already set to zero at level 1, and to zero or the ocean
84      !!      bottom slope (ln_sco=T) at level jpk in inildf]
85      !!
86      !! ** Action : - uslp, wslpi, and vslp, wslpj, the i- and  j-slopes
87      !!               of now neutral surfaces at u-, w- and v- w-points, resp.
88      !!----------------------------------------------------------------------
89      USE oce , zgru  => ua   ! use ua as workspace
90      USE oce , zgrv  => va   ! use va as workspace
91      USE oce , zwy   => ta   ! use ta as workspace
92      USE oce , zwz   => sa   ! use sa as workspace
93      !!
94      INTEGER , INTENT(in)                         ::   kt    ! ocean time-step index
95      REAL(wp), INTENT(in), DIMENSION(jpi,jpj,jpk) ::   prd   ! in situ density
96      REAL(wp), INTENT(in), DIMENSION(jpi,jpj,jpk) ::   pn2   ! Brunt-Vaisala frequency (locally ref.)
97      !!
98      INTEGER  ::   ji , jj , jk    ! dummy loop indices
99      INTEGER  ::   ii0, ii1, iku   ! temporary integer
100      INTEGER  ::   ij0, ij1, ikv   ! temporary integer
101      REAL(wp) ::   zeps, zmg, zm05g, zalpha        ! temporary scalars
102      REAL(wp) ::   zcoef1, zcoef2, zcoef3          !    -         -
103      REAL(wp) ::   zcofu , zcofv , zcofw           !    -         -
104      REAL(wp) ::   zau, zbu, zai, zbi, z1u, z1wu   !    -         -
105      REAL(wp) ::   zav, zbv, zaj, zbj, z1v, z1wv   !
106      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zww     ! 3D workspace
107      !!----------------------------------------------------------------------
108     
109      zeps  =  1.e-20                             ! Local constant initialization
110      zmg   = -1.0 / grav
111      zm05g = -0.5 / grav
112      !
113      zww(:,:,:) = 0.e0
114      zwz(:,:,:) = 0.e0
115      !                                           ! horizontal density gradient computation
116      DO jk = 1, jpk
117         DO jj = 1, jpjm1
118            DO ji = 1, fs_jpim1   ! vector opt.
119               zgru(ji,jj,jk) = umask(ji,jj,jk) * ( prd(ji+1,jj  ,jk) - prd(ji,jj,jk) ) 
120               zgrv(ji,jj,jk) = vmask(ji,jj,jk) * ( prd(ji  ,jj+1,jk) - prd(ji,jj,jk) ) 
121            END DO
122         END DO
123      END DO
124      IF( ln_zps ) THEN                           ! partial steps correction at the bottom ocean level
125# if defined key_vectopt_loop 
126         DO jj = 1, 1
127            DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
128# else
129         DO jj = 1, jpjm1
130            DO ji = 1, jpim1
131# endif
132               iku = MIN ( mbathy(ji,jj), mbathy(ji+1,jj) ) - 1               ! last ocean level
133               ikv = MIN ( mbathy(ji,jj), mbathy(ji,jj+1) ) - 1
134               zgru(ji,jj,iku) = gru(ji,jj) 
135               zgrv(ji,jj,ikv) = grv(ji,jj)               
136            END DO
137         END DO
138      ENDIF
139     
140      CALL ldf_slp_mxl( prd, pn2 )                ! Slopes of isopycnal surfaces just below the mixed layer
141
142     
143      ! I.  slopes at u and v point      | uslp = d/di( prd ) / d/dz( prd )
144      ! ===========================      | vslp = d/dj( prd ) / d/dz( prd )
145      !               
146      !                                           !* Local vertical density gradient evaluated from N^2
147      DO jk = 2, jpkm1                            !  zwy = d/dz(prd)= - ( prd ) / grav * mk(pn2) -- at t point
148         DO jj = 1, jpj
149            DO ji = 1, jpi
150               zwy(ji,jj,jk) = zmg * ( prd(ji,jj,jk) + 1. ) * (    pn2  (ji,jj,jk) + pn2  (ji,jj,jk+1)     )   &
151                  &                                         / MAX( tmask(ji,jj,jk) + tmask(ji,jj,jk+1), 1. )
152            END DO
153         END DO
154      END DO
155      DO jk = 2, jpkm1                            !* Slopes at u and v points
156         DO jj = 2, jpjm1
157            DO ji = fs_2, fs_jpim1   ! vector opt.
158               ! horizontal and vertical density gradient at u- and v-points
159               zau = 1. / e1u(ji,jj) * zgru(ji,jj,jk)
160               zav = 1. / e2v(ji,jj) * zgrv(ji,jj,jk)
161               zbu = 0.5 * ( zwy(ji,jj,jk) + zwy(ji+1,jj  ,jk) )
162               zbv = 0.5 * ( zwy(ji,jj,jk) + zwy(ji  ,jj+1,jk) )
163               ! bound the slopes: abs(zw.)<= 1/100 and zb..<0
164               !                   kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
165               zbu = MIN( zbu, -100.*ABS( zau ), -7.e+3/fse3u(ji,jj,jk)*ABS( zau ) )
166               zbv = MIN( zbv, -100.*ABS( zav ), -7.e+3/fse3v(ji,jj,jk)*ABS( zav ) )
167               ! uslp and vslp output in zwz and zww, resp.
168               zalpha = MAX( omlmask(ji,jj,jk), omlmask(ji+1,jj,jk) )
169               zwz (ji,jj,jk) = ( ( 1. - zalpha) * zau / ( zbu - zeps )                                           &
170                  &                    + zalpha  * uslpml(ji,jj)                                                  &
171                  &                              * 0.5 * ( fsdept(ji+1,jj,jk)+fsdept(ji,jj,jk)-fse3u(ji,jj,1) )   &
172                  &                              / MAX( hmlpt(ji,jj), hmlpt(ji+1,jj), 5. ) ) * umask(ji,jj,jk)
173               zalpha = MAX( omlmask(ji,jj,jk), omlmask(ji,jj+1,jk) )
174               zww (ji,jj,jk) = ( ( 1. - zalpha) * zav / ( zbv - zeps )                                           &
175                  &                    + zalpha  * vslpml(ji,jj)                                                  &
176                  &                              * 0.5 * ( fsdept(ji,jj+1,jk)+fsdept(ji,jj,jk)-fse3v(ji,jj,1) )   &
177                  &                              / MAX( hmlpt(ji,jj), hmlpt(ji,jj+1), 5. ) ) * vmask(ji,jj,jk)
178            END DO
179         END DO
180      END DO
181      CALL lbc_lnk( zwz, 'U', -1. )   ;   CALL lbc_lnk( zww, 'V', -1. )      ! lateral boundary conditions
182      !
183      zcofu = 1. / 16.                            !* horizontal Shapiro filter
184      zcofv = 1. / 16.
185      DO jk = 2, jpkm1
186         DO jj = 2, jpjm1, jpj-3                        ! rows jj=2 and =jpjm1 only
187            DO ji = 2, jpim1 
188               uslp(ji,jj,jk) = zcofu * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      &
189                  &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      &
190                  &                       + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)      &
191                  &                       +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )    &
192                  &                       + 4.*  zwz(ji  ,jj  ,jk)                       )
193               vslp(ji,jj,jk) = zcofv * (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)      &
194                  &                       +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)      &
195                  &                       + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)      &
196                  &                       +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )    &
197                  &                       + 4.*  zww(ji,jj    ,jk)                       )
198            END DO
199         END DO
200         DO jj = 3, jpj-2                               ! other rows
201            DO ji = fs_2, fs_jpim1   ! vector opt.
202               uslp(ji,jj,jk) = zcofu * (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)      &
203                  &                       +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)      &
204                  &                       + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)      &
205                  &                       +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )    &
206                  &                       + 4.*  zwz(ji  ,jj  ,jk)                       )
207               vslp(ji,jj,jk) = zcofv * (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)      &
208                  &                       +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)      &
209                  &                       + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)      &
210                  &                       +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )    &
211                  &                       + 4.*  zww(ji,jj    ,jk)                       )
212            END DO
213         END DO
214         !                                        !* decrease along coastal boundaries
215         DO jj = 2, jpjm1
216            DO ji = fs_2, fs_jpim1   ! vector opt.
217               z1u  = ( umask(ji,jj+1,jk) + umask(ji,jj-1,jk) )*.5
218               z1v  = ( vmask(ji+1,jj,jk) + vmask(ji-1,jj,jk) )*.5
219               z1wu = ( umask(ji,jj,jk)   + umask(ji,jj,jk+1) )*.5
220               z1wv = ( vmask(ji,jj,jk)   + vmask(ji,jj,jk+1) )*.5
221               uslp(ji,jj,jk) = uslp(ji,jj,jk) * z1u * z1wu
222               vslp(ji,jj,jk) = vslp(ji,jj,jk) * z1v * z1wv
223            END DO
224         END DO
225      END DO
226
227
228      ! II.  slopes at w point           | wslpi = mij( d/di( prd ) / d/dz( prd )
229      ! ===========================      | wslpj = mij( d/dj( prd ) / d/dz( prd )
230      !               
231      !                                           !* Local vertical density gradient evaluated from N^2
232      DO jk = 2, jpkm1                            !  zwy = d/dz(prd)= - mk ( prd ) / grav * pn2 -- at w point
233         DO jj = 1, jpj
234            DO ji = 1, jpi
235               zwy(ji,jj,jk) = zm05g * pn2 (ji,jj,jk) * ( prd (ji,jj,jk) + prd (ji,jj,jk-1) + 2. )
236            END DO
237         END DO
238      END DO
239      DO jk = 2, jpkm1                            !* Slopes at w point
240         DO jj = 2, jpjm1
241            DO ji = fs_2, fs_jpim1   ! vector opt.
242               !                                        ! horizontal density i-gradient at w-points
243               zcoef1 = MAX( zeps, umask(ji-1,jj,jk  )+umask(ji,jj,jk  )    &
244                  &               +umask(ji-1,jj,jk-1)+umask(ji,jj,jk-1) )
245               zcoef1 = 1. / ( zcoef1 * e1t (ji,jj) )
246               zai = zcoef1 * (  zgru(ji  ,jj,jk  ) + zgru(ji  ,jj,jk-1)   &
247                  &            + zgru(ji-1,jj,jk-1) + zgru(ji-1,jj,jk  ) ) * tmask (ji,jj,jk)
248               !                                        ! horizontal density j-gradient at w-points
249               zcoef2 = MAX( zeps, vmask(ji,jj-1,jk  )+vmask(ji,jj,jk-1)   &
250                  &               +vmask(ji,jj-1,jk-1)+vmask(ji,jj,jk  ) )
251               zcoef2 = 1.0 / ( zcoef2 *  e2t (ji,jj) )
252               zaj = zcoef2 * (  zgrv(ji,jj  ,jk  ) + zgrv(ji,jj  ,jk-1)   &
253                  &            + zgrv(ji,jj-1,jk-1) + zgrv(ji,jj-1,jk  ) ) * tmask (ji,jj,jk)
254               !                                        ! bound the slopes: abs(zw.)<= 1/100 and zb..<0.
255               !                                        ! static instability: kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
256               zbi = MIN( zwy (ji,jj,jk),- 100.*ABS(zai), -7.e+3/fse3w(ji,jj,jk)*ABS(zai) )
257               zbj = MIN( zwy (ji,jj,jk), -100.*ABS(zaj), -7.e+3/fse3w(ji,jj,jk)*ABS(zaj) )
258               !                                        ! wslpi and wslpj output in zwz and zww, resp.
259               zalpha = MAX( omlmask(ji,jj,jk), omlmask(ji,jj,jk-1) )
260               zcoef3 = fsdepw(ji,jj,jk) / MAX( hmlp(ji,jj), 10. )
261               zwz(ji,jj,jk) = (     zai / ( zbi - zeps)  * ( 1. - zalpha )   &
262                  &             + zcoef3 * wslpiml(ji,jj) *        zalpha   ) * tmask (ji,jj,jk)
263               zww(ji,jj,jk) = (     zaj / ( zbj - zeps)  * ( 1. - zalpha )   &
264                  &             + zcoef3 * wslpjml(ji,jj) *        zalpha   ) * tmask (ji,jj,jk)
265            END DO
266         END DO
267      END DO
268      CALL lbc_lnk( zwz, 'T', -1. )   ;    CALL lbc_lnk( zww, 'T', -1. )      ! lateral boundary conditions on zwz and zww
269      !
270      !                                           !* horizontal Shapiro filter
271      DO jk = 2, jpkm1
272         DO jj = 2, jpjm1, jpj-3                        ! rows jj=2 and =jpjm1
273            DO ji = 2, jpim1
274               zcofw = tmask(ji,jj,jk) / 16.
275               wslpi(ji,jj,jk) = (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     &
276                  &                +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     &
277                  &                + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)     &
278                  &                +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )   &
279                  &                + 4.*  zwz(ji  ,jj  ,jk)                        ) * zcofw
280
281               wslpj(ji,jj,jk) = (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)     &
282                  &                +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)     &
283                  &                + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)     &
284                  &                +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )   &
285                  &                + 4.*  zww(ji  ,jj  ,jk)                        ) * zcofw
286            END DO
287         END DO 
288         DO jj = 3, jpj-2                               ! other rows
289            DO ji = fs_2, fs_jpim1   ! vector opt.
290               zcofw = tmask(ji,jj,jk) / 16.
291               wslpi(ji,jj,jk) = (        zwz(ji-1,jj-1,jk) + zwz(ji+1,jj-1,jk)     &
292                  &                +      zwz(ji-1,jj+1,jk) + zwz(ji+1,jj+1,jk)     &
293                  &                + 2.*( zwz(ji  ,jj-1,jk) + zwz(ji-1,jj  ,jk)     &
294                  &                +      zwz(ji+1,jj  ,jk) + zwz(ji  ,jj+1,jk) )   &
295                  &                + 4.*  zwz(ji  ,jj  ,jk)                        ) * zcofw
296
297               wslpj(ji,jj,jk) = (        zww(ji-1,jj-1,jk) + zww(ji+1,jj-1,jk)     &
298                  &                +      zww(ji-1,jj+1,jk) + zww(ji+1,jj+1,jk)     &
299                  &                + 2.*( zww(ji  ,jj-1,jk) + zww(ji-1,jj  ,jk)     &
300                  &                +      zww(ji+1,jj  ,jk) + zww(ji  ,jj+1,jk) )   &
301                  &                + 4.*  zww(ji  ,jj  ,jk)                        ) * zcofw
302            END DO
303         END DO
304         !                                        !* decrease along coastal boundaries
305         DO jj = 2, jpjm1
306            DO ji = fs_2, fs_jpim1   ! vector opt.
307               z1u = ( umask(ji,jj,jk) + umask(ji-1,jj,jk) ) *.5
308               z1v = ( vmask(ji,jj,jk) + vmask(ji,jj-1,jk) ) *.5
309               wslpi(ji,jj,jk) = wslpi(ji,jj,jk) * z1u * z1v
310               wslpj(ji,jj,jk) = wslpj(ji,jj,jk) * z1u * z1v
311            END DO
312         END DO
313      END DO
314         
315      ! III.  Specific grid points     
316      ! ===========================
317      !               
318      IF( cp_cfg == "orca" .AND. jp_cfg == 4 ) THEN     !  ORCA_R4 configuration: horizontal diffusion in specific area
319         !                                                    ! Gibraltar Strait
320         ij0 =  50   ;   ij1 =  53
321         ii0 =  69   ;   ii1 =  71   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
322         ij0 =  51   ;   ij1 =  53
323         ii0 =  68   ;   ii1 =  71   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
324         ii0 =  69   ;   ii1 =  71   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
325         ii0 =  69   ;   ii1 =  71   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
326         !
327         !                                                    ! Mediterrannean Sea
328         ij0 =  49   ;   ij1 =  56
329         ii0 =  71   ;   ii1 =  90   ;   uslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
330         ij0 =  50   ;   ij1 =  56
331         ii0 =  70   ;   ii1 =  90   ;   vslp ( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
332         ii0 =  71   ;   ii1 =  90   ;   wslpi( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
333         ii0 =  71   ;   ii1 =  90   ;   wslpj( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) , : ) = 0.e0
334      ENDIF
335
336     
337      ! IV. Lateral boundary conditions
338      ! ===============================
339      CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. )
340      CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. )
341
342
343      IF(ln_ctl) THEN
344         CALL prt_ctl(tab3d_1=uslp , clinfo1=' slp  - u : ', tab3d_2=vslp,  clinfo2=' v : ', kdim=jpk)
345         CALL prt_ctl(tab3d_1=wslpi, clinfo1=' slp  - wi: ', tab3d_2=wslpj, clinfo2=' wj: ', kdim=jpk)
346      ENDIF
347      !
348   END SUBROUTINE ldf_slp
349
350
351   SUBROUTINE ldf_slp_grif ( kt )
352     !!----------------------------------------------------------------------
353     !!                 ***  ROUTINE ldf_slp_grif  ***
354     !!
355     !! ** Purpose :   Compute the squared slopes of neutral surfaces (slope
356     !!      of iso-pycnal surfaces referenced locally) ('key_traldfiso')
357     !!      at W-points using the Griffies quarter-cells.  Also calculates
358     !!      alpha and beta at T-points for use in the Griffies isopycnal
359     !!      scheme.
360     !!
361     !! ** Method  :   The slope in the i-direction is computed at U- and
362     !!      W-points (uslp, wslpi) and the slope in the j-direction is
363     !!      computed at V- and W-points (vslp, wslpj).
364     !!
365     !! ** Action : - alpha, beta
366     !!               wslp2 squared slope of neutral surfaces at w-points.
367     !!----------------------------------------------------------------------
368     USE oce, zdit  => ua   ! use ua as workspace
369     USE oce, zdis  => va   ! use va as workspace
370     USE oce, zdjt  => ta   ! use ta as workspace
371     USE oce, zdjs  => sa   ! use sa as workspace
372     !!
373     INTEGER, INTENT( in ) ::   kt         ! ocean time-step index
374     !!
375     INTEGER  ::   ji, jj, jk, ip, jp, kp  ! dummy loop indices
376     INTEGER  ::   iku, ikv                ! local integer
377     REAL(wp) ::   zt, zs, zh, zt2, zsp5, zp1t1           ! local scalars
378     REAL(wp) ::     zdenr, zrhotmp, zdndt, zdddt         !   -      -
379     REAL(wp) ::     zdnds, zddds, znum, zden             !   -      -
380     REAL(wp) ::     zslope, za_sxe, zslopec, zdsloper    !   -      -
381     REAL(wp) ::     zfact, zepsln, zatempw,zatempu,zatempv    !   -      -
382     REAL(wp) ::     ze1ur, ze2vr, ze3wr, zdxt, zdxs, zdyt, zdys, zdzt, zdzs, zvolf 
383     REAL(wp) ::     zr_slpmax, zdxrho, zdyrho, zabs_dzrho
384     REAL(wp), DIMENSION(jpi,jpj,jpk,0:1,0:1) ::   zsx,zsy
385     REAL(wp), DIMENSION(jpi,jpj    ,0:1,0:1) ::   zsx_ml_base, zsy_ml_base
386     REAL(wp), DIMENSION(jpi,jpj,jpk)         ::   zdkt, zdks
387     REAL(wp), DIMENSION(jpi,jpj) ::   zr_ml_basew
388     !!----------------------------------------------------------------------
389
390     ! 0. Local constant initialization
391     ! --------------------------------
392     zr_slpmax = 1.0_wp/slpmax
393
394     ! zslopec=0.004
395     ! zdsloper=1000.0
396     zepsln=1e-25
397
398     SELECT CASE ( nn_eos )
399
400     CASE ( 0 )               ! Jackett and McDougall (1994) formulation
401
402        !                                                ! ===============
403        DO jk = 1, jpk                                   ! Horizontal slab
404           !                                             ! ===============
405           DO jj = 1, jpjm1
406              DO ji = 1, fs_jpim1
407                 zt = tb(ji,jj,jk)          ! potential temperature
408                 zs = sb(ji,jj,jk) - 35.0   ! salinity anomaly (s-35)
409                 zh = fsdept(ji,jj,jk)      ! depth in meters
410
411                 beta(ji,jj,jk) = ( ( -0.415613e-09 * zt + 0.555579e-07 ) * zt      &
412                      &                            - 0.301985e-05 ) * zt      &
413                      &   + 0.785567e-03                                      &
414                      &   + (     0.515032e-08 * zs                           &
415                      &         + 0.788212e-08 * zt - 0.356603e-06 ) * zs     &
416                      &   +(  (   0.121551e-17 * zh                           &
417                      &         - 0.602281e-15 * zs                           &
418                      &         - 0.175379e-14 * zt + 0.176621e-12 ) * zh     &
419                      &                             + 0.408195e-10   * zs     &
420                      &     + ( - 0.213127e-11 * zt + 0.192867e-09 ) * zt     &
421                      &                             - 0.121555e-07 ) * zh
422
423                 alpha(ji,jj,jk) = - beta(ji,jj,jk) *                       &
424                      &     (((( - 0.255019e-07 * zt + 0.298357e-05 ) * zt   &
425                      &                               - 0.203814e-03 ) * zt   &
426                      &                               + 0.170907e-01 ) * zt   &
427                      &   + 0.665157e-01                                      &
428                      &   +     ( - 0.678662e-05 * zs                         &
429                      &           - 0.846960e-04 * zt + 0.378110e-02 ) * zs   &
430                      &   +   ( ( - 0.302285e-13 * zh                         &
431                      &           - 0.251520e-11 * zs                         &
432                      &           + 0.512857e-12 * zt * zt           ) * zh   &
433                      &           - 0.164759e-06 * zs                         &
434                      &        +(   0.791325e-08 * zt - 0.933746e-06 ) * zt   &
435                      &                               + 0.380374e-04 ) * zh)
436
437              END DO
438           END DO
439        END DO
440
441     CASE ( 1 )
442
443        alpha(:,:,:)=-rn_alpha
444        beta(:,:,:)=0.0
445
446     CASE ( 2 )
447
448        alpha(:,:,:)=-rn_alpha
449        beta (:,:,:)=rn_beta
450
451     CASE ( 3 )
452
453        DO jk =1, jpk
454           DO jj = 1, jpjm1
455              DO ji = 1, fs_jpim1
456                 zt = tb(ji,jj,jk)
457                 zs = sb(ji,jj,jk)
458                 zh = fsdept(ji,jj,jk)
459                 zt2 = zt**2
460                 zsp5 = sqrt(ABS(zs))
461                 zp1t1=zh*zt
462                 znum=9.99843699e+02+zt*(7.35212840e+00+zt*(-5.45928211e-02+3.98476704e-04*zt)) &
463                      +zs*(2.96938239e+00-7.23268813e-03*zt+2.12382341e-03*zs)  &
464                      +zh*(1.04004591e-02+1.03970529e-07*zt2+5.18761880e-06*zs+ &
465                      zh*(-3.24041825e-08-1.23869360e-11*zt2))
466                 zden=1.00000000e+00+zt*(7.28606739e-03+zt*(-4.60835542e-05+zt*(3.68390573e-07+zt*1.80809186e-10))) &
467                      +zs*(2.14691708e-03+zt*(-9.27062484e-06-1.78343643e-10*zt2)+zsp5*(4.76534122e-06+1.63410736e-09*zt2)) &
468                      + zh*(5.30848875e-06+zh*zt*(-3.03175128e-16*zt2-1.27934137e-17*zh))
469                 zdenr=1.0/zden
470                 zrhotmp=znum*zdenr
471                 zdndt=7.35212840e+00+zt*(-1.091856422e-01+1.195430112e-03*zt)-7.23268813e-03*zs &
472                      +zp1t1*(2.07941058e-07-2.4773872e-11*zh)
473                 zdddt=7.28606739e-03+zt*(-9.21671084e-05+zt*(1.105171719e-06+7.23236744e-10*zt))  &
474                      +zs*(-9.27062484e-06+zt*(-5.35030929e-10*zt+3.26821472e-09*zsp5))  &
475                      +zh*zh*(-9.09525384e-16*zt2-1.27934137e-17*zh)
476                 zdnds=2.96938239e+00-7.23268813e-03*zt+2*2.12382341e-03*zs+5.18761880e-06*zh
477                 zddds=2.14691708e-03+zt*(-9.27062484e-06-1.78343643e-10*zt2)+zsp5*(7.14801183e-06+2.45116104e-09*zt2)
478                 alpha(ji,jj,jk)=(zdndt-zrhotmp*zdddt)*zdenr
479                 beta(ji,jj,jk)=zdenr*(zdnds-zrhotmp*zddds)
480
481              END DO
482           END DO
483        END DO
484
485     CASE DEFAULT
486
487        IF(lwp) WRITE(numout,cform_err)
488        IF(lwp) WRITE(numout,*) '          bad flag value for nn_eos = ', nn_eos
489        nstop = nstop + 1
490
491     END SELECT
492
493     CALL lbc_lnk( alpha, 'T', 1. )
494     CALL lbc_lnk( beta, 'T', 1. )
495
496
497     ! ---------------------------------------
498     ! 1. Horizontal tracer gradients at T-level jk
499     ! ---------------------------------------
500     DO jk = 1, jpkm1
501        DO jj = 1, jpjm1
502           DO ji = 1, fs_jpim1   ! vector opt.
503              ! i-gradient of T and S at jj
504              zdit (ji,jj,jk) = ( tb(ji+1,jj,jk)-tb(ji,jj,jk) ) * umask(ji,jj,jk)
505              zdis (ji,jj,jk) = ( sb(ji+1,jj,jk)-sb(ji,jj,jk) ) * umask(ji,jj,jk)
506              ! j-gradient of T and S at jj
507              zdjt (ji,jj,jk) = ( tb(ji,jj+1,jk)-tb(ji,jj,jk) ) * vmask(ji,jj,jk)
508              zdjs (ji,jj,jk) = ( sb(ji,jj+1,jk)-sb(ji,jj,jk) ) * vmask(ji,jj,jk)
509           END DO
510        END DO
511     END DO
512
513     IF( ln_zps ) THEN      ! partial steps correction at the last level
514# if defined key_vectopt_loop
515     DO jj = 1, 1
516        DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
517# else
518           DO jj = 1, jpjm1
519              DO ji = 1, jpim1
520# endif
521                 ! last ocean level
522                 iku  = MIN( mbathy(ji,jj), mbathy(ji+1,jj  ) ) - 1
523                 ikv  = MIN( mbathy(ji,jj), mbathy(ji  ,jj+1) ) - 1
524                 ! i-gradient of T and S
525                 zdit (ji,jj,iku) = gtsu(ji,jj,jp_tem)
526                 zdis (ji,jj,iku) = gtsu(ji,jj,jp_sal)
527                 ! j-gradient of T and S
528                 zdjt (ji,jj,ikv) = gtsv(ji,jj,jp_tem)
529                 zdjs (ji,jj,ikv) = gtsv(ji,jj,jp_sal)
530              END DO
531           END DO
532        ENDIF
533
534        ! ---------------------------------------
535        ! 1. Vertical tracer gradient at w-level jk
536        ! ---------------------------------------
537        DO jk = 2, jpk
538           zdkt(:,:,jk) = ( tb(:,:,jk-1) - tb(:,:,jk) ) * tmask(:,:,jk)
539           zdks(:,:,jk) = ( sb(:,:,jk-1) - sb(:,:,jk) ) * tmask(:,:,jk)
540        END DO
541
542        zdkt(:,:,1) = 0.0
543        zdks(:,:,1) = 0.0
544
545        ! ---------------------------------------
546        ! Depth of the w-point below ML base
547        ! ---------------------------------------
548        DO jj = 1, jpj
549           DO ji = 1, jpi
550              jk = nmln(ji,jj)
551              zr_ml_basew(ji,jj)=1.0/fsdepw(ji,jj,jk+1)
552           END DO
553        END DO
554
555
556        wslp2(:,:,:)=0.0
557        tfw(:,:,:) = 0.0
558        sfw(:,:,:) = 0.0
559        ftu(:,:,:) = 0.0
560        fsu(:,:,:) = 0.0
561        ftv(:,:,:) = 0.0
562        fsv(:,:,:) = 0.0
563        ftud(:,:,:) = 0.0
564        fsud(:,:,:) = 0.0
565        ftvd(:,:,:) = 0.0
566        fsvd(:,:,:) = 0.0
567        psix_eiv(:,:,:) = 0.0
568        psiy_eiv(:,:,:) = 0.0
569
570        ! ----------------------------------------------------------------------
571        ! x-z plane
572        ! ----------------------------------------------------------------------
573
574        ! calculate limited triad x-slopes zsx in interior (1=<jk=<jpk-1)
575        DO ip=0,1
576           DO kp=0,1
577
578              DO jk = 1, jpkm1
579                 DO jj = 1, jpjm1
580                    DO ji = 1, fs_jpim1   ! vector opt.
581
582                       ze1ur=1.0/e1u(ji,jj)
583                       zdxt=zdit(ji,jj,jk)*ze1ur
584                       zdxs=zdis(ji,jj,jk)*ze1ur
585
586                       ze3wr=1.0/fse3w(ji+ip,jj,jk+kp)
587                       zdzt=zdkt(ji+ip,jj,jk+kp)*ze3wr
588                       zdzs=zdks(ji+ip,jj,jk+kp)*ze3wr
589                       ! Calculate the horizontal and vertical density gradient
590                       zdxrho = alpha(ji+ip,jj,jk)*zdxt+beta(ji+ip,jj,jk)*zdxs
591                       zabs_dzrho = ABS(alpha(ji+ip,jj,jk)*zdzt+beta(ji+ip,jj,jk)*zdzs)+zepsln
592                       ! Limit by slpmax, and mask by psi-point
593                       zsx(ji+ip,jj,jk,1-ip,kp) = umask(ji,jj,jk+kp) &
594                            & *zdxrho/MAX(zabs_dzrho,zr_slpmax*ABS(zdxrho))
595                    END DO
596                 END DO
597              END DO
598
599           END DO
600        END DO
601
602        ! calculate slope of triad immediately below mixed-layer base
603        DO ip=0,1
604           DO kp=0,1
605              DO jj = 1, jpjm1
606                 DO ji = 1, fs_jpim1
607                    jk = nmln(ji+ip,jj)
608                    zsx_ml_base(ji+ip,jj,1-ip,kp)=zsx(ji+ip,jj,jk+1-kp,1-ip,kp)
609                 END DO
610              END DO
611           END DO
612        END DO
613
614        ! Below ML use limited zsx as is
615        ! Inside ML replace by linearly reducing zsx_ml_base towards surface
616        DO ip=0,1
617           DO kp=0,1
618
619              DO jk = 1, jpkm1
620
621                 DO jj = 1, jpjm1
622
623                    DO ji = 1, fs_jpim1   ! vector opt.
624                       ! k index of uppermost point(s) of triad is jk+kp-1
625                       ! must be .ge. nmln(ji,jj) for zfact=1.
626                       !                   otherwise  zfact=0.
627                       zfact = 1 - 1/(1 + (jk+kp-1)/nmln(ji+ip,jj))
628                       zsx(ji+ip,jj,jk,1-ip,kp) = zfact*zsx(ji+ip,jj,jk,1-ip,kp) + &
629                            & (1.0_wp-zfact)*(fsdepw(ji+ip,jj,jk+kp)*zr_ml_basew(ji+ip,jj))*zsx_ml_base(ji+ip,jj,1-ip,kp) 
630                    END DO
631
632                 END DO
633
634              END DO
635           END DO
636        END DO
637
638        ! Use zsx to calculate fluxes and save averaged slopes psix_eiv at psi-points
639        DO ip=0,1
640           DO kp=0,1
641
642              DO jk = 1, jpkm1
643
644                 DO jj = 1, jpjm1
645
646                    DO ji = 1, fs_jpim1   ! vector opt.
647
648                       ze1ur=1.0/e1u(ji,jj)
649                       zdxt=zdit(ji,jj,jk)*ze1ur
650                       zdxs=zdis(ji,jj,jk)*ze1ur
651
652                       ze3wr=1.0/fse3w(ji+ip,jj,jk+kp)
653                       zdzt=zdkt(ji+ip,jj,jk+kp)*ze3wr
654                       zdzs=zdks(ji+ip,jj,jk+kp)*ze3wr
655                       zslope=zsx(ji+ip,jj,jk,1-ip,kp)
656
657                       zvolf = 0.25_wp*e1u(ji,jj)*e2u(ji,jj)*fse3u(ji,jj,jk)
658
659                       ftu(ji,jj,jk)= ftu(ji,jj,jk)+zslope*zdzt*zvolf*ze1ur
660                       fsu(ji,jj,jk)= fsu(ji,jj,jk)+zslope*zdzs*zvolf*ze1ur
661                       ftud(ji,jj,jk)=ftud(ji,jj,jk)+fsahtu(ji,jj,jk)*zdxt*zvolf*ze1ur
662                       fsud(ji,jj,jk)=fsud(ji,jj,jk)+fsahtu(ji,jj,jk)*zdxs*zvolf*ze1ur
663                       tfw(ji+ip,jj,jk+kp)=tfw(ji+ip,jj,jk+kp)+(zvolf*ze3wr)*zslope*zdxt
664                       sfw(ji+ip,jj,jk+kp)=sfw(ji+ip,jj,jk+kp)+(zvolf*ze3wr)*zslope*zdxs
665                       wslp2(ji+ip,jj,jk+kp)=wslp2(ji+ip,jj,jk+kp)+ &
666                            & ((zvolf*ze3wr)/(e1t(ji+ip,jj)*e2t(ji+ip,jj)))*(zslope)**2
667                       psix_eiv(ji,jj,jk+kp) = psix_eiv(ji,jj,jk+kp) + 0.25_wp*zslope
668
669                    END DO
670                 END DO
671
672              END DO
673           END DO
674        END DO
675
676        ! ----------------------------------------------------------------------
677        ! y-z plane
678        ! ----------------------------------------------------------------------
679
680        ! calculate limited triad y-slopes zsy in interior (1=<jk=<jpk-1)
681        DO jp=0,1
682           DO kp=0,1
683
684              DO jk = 1, jpkm1
685                 DO jj = 1, jpjm1
686                    DO ji = 1, fs_jpim1   ! vector opt.
687
688                       ze2vr=1.0/e2v(ji,jj)
689                       zdyt=zdjt(ji,jj,jk)*ze2vr
690                       zdys=zdjs(ji,jj,jk)*ze2vr
691
692                       ze3wr=1.0/fse3w(ji,jj+jp,jk+kp)
693                       zdzt=zdkt(ji,jj+jp,jk+kp)*ze3wr
694                       zdzs=zdks(ji,jj+jp,jk+kp)*ze3wr
695                       ! Calculate the horizontal and vertical density gradient
696                       zdyrho = alpha(ji,jj+jp,jk)*zdyt+beta(ji,jj+jp,jk)*zdys
697                       zabs_dzrho = ABS(alpha(ji,jj+jp,jk)*zdzt+beta(ji,jj+jp,jk)*zdzs)+zepsln
698                       ! Limit by slpmax, and mask by psi-point
699                       zsy(ji,jj+jp,jk,1-jp,kp) = vmask(ji,jj,jk+kp) &
700                            & *zdyrho/MAX(zabs_dzrho,zr_slpmax*ABS(zdyrho))
701                    END DO
702                 END DO
703              END DO
704
705           END DO
706        END DO
707
708        ! calculate slope of triad immediately below mixed-layer base
709        DO jp=0,1
710           DO kp=0,1
711              DO jj = 1, jpjm1
712                 DO ji = 1, fs_jpim1
713                    jk = nmln(ji,jj+jp)
714                    zsy_ml_base(ji,jj+jp,1-jp,kp)=zsy(ji,jj+jp,jk+1-kp,1-jp,kp)
715                 END DO
716              END DO
717           END DO
718        END DO
719
720        ! Below ML use limited zsy as is
721        ! Inside ML replace by linearly reducing zsy_ml_base towards surface
722        DO jp=0,1
723           DO kp=0,1
724              DO jk = 1, jpkm1
725                 DO jj = 1, jpjm1
726                    DO ji = 1, fs_jpim1   ! vector opt.
727                       ! k index of uppermost point(s) of triad is jk+kp-1
728                       ! must be .ge. nmln(ji,jj) for zfact=1.
729                       !                   otherwise  zfact=0.
730                       zfact = 1 - 1/(1 + (jk+kp-1)/nmln(ji,jj+jp))
731                       zsy(ji,jj+jp,jk,1-jp,kp) = zfact*zsy(ji,jj+jp,jk,1-jp,kp) + &
732                           & (1.0_wp-zfact)*(fsdepw(ji,jj+jp,jk+kp)*zr_ml_basew(ji,jj+jp))*zsy_ml_base(ji,jj+jp,1-jp,kp) 
733                    END DO
734
735                 END DO
736
737              END DO
738           END DO
739        END DO
740
741        ! Use zsy to calculate fluxes and save averaged slopes psiy_eiv at psi-points
742        DO jp=0,1
743           DO kp=0,1
744              DO jk = 1, jpkm1
745                 DO jj = 1, jpjm1
746                    DO ji = 1, fs_jpim1   ! vector opt.
747
748                       ze2vr=1.0/e2v(ji,jj)
749                       zdyt=zdjt(ji,jj,jk)*ze2vr
750                       zdys=zdjs(ji,jj,jk)*ze2vr
751
752                       ze3wr=1.0/fse3w(ji,jj+jp,jk+kp)
753                       zdzt=zdkt(ji,jj+jp,jk+kp)*ze3wr
754                       zdzs=zdks(ji,jj+jp,jk+kp)*ze3wr
755                       zslope=zsy(ji,jj+jp,jk,1-jp,kp)
756
757                       zvolf = 0.25_wp*e1v(ji,jj)*e2v(ji,jj)*fse3v(ji,jj,jk)
758
759                       ftv(ji,jj,jk)= ftv(ji,jj,jk)+zslope*zdzt*zvolf*ze2vr
760                       fsv(ji,jj,jk)= fsv(ji,jj,jk)+zslope*zdzs*zvolf*ze2vr
761                       ftvd(ji,jj,jk)=ftvd(ji,jj,jk)+fsahtv(ji,jj,jk)*zdyt*zvolf*ze2vr
762                       fsvd(ji,jj,jk)=fsvd(ji,jj,jk)+fsahtv(ji,jj,jk)*zdys*zvolf*ze2vr
763                       tfw(ji,jj+jp,jk+kp)=tfw(ji,jj+jp,jk+kp)+(zvolf*ze3wr)*zslope*zdyt
764                       sfw(ji,jj+jp,jk+kp)=sfw(ji,jj+jp,jk+kp)+(zvolf*ze3wr)*zslope*zdys
765                       wslp2(ji,jj+jp,jk+kp)=wslp2(ji,jj+jp,jk+kp)+ &
766                            & ((zvolf*ze3wr)/(e1t(ji,jj+jp)*e2t(ji,jj+jp)))*(zslope)**2
767                       psiy_eiv(ji,jj,jk+kp) = psiy_eiv(ji,jj,jk+kp) + 0.25_wp*zslope
768
769                    END DO
770                 END DO
771              END DO
772           END DO
773        END DO
774
775        tfw(:,:,1)=0.e0
776        sfw(:,:,1)=0.e0
777        wslp2(:,:,1)=0.e0
778
779        CALL lbc_lnk( wslp2, 'W', 1. )
780        CALL lbc_lnk( tfw, 'W', 1. )
781        CALL lbc_lnk( sfw, 'W', 1. )
782        CALL lbc_lnk( ftu, 'U', -1. )
783        CALL lbc_lnk( fsu, 'U', -1. )
784        CALL lbc_lnk( ftv, 'V', -1. )
785        CALL lbc_lnk( fsv, 'V', -1. )
786        CALL lbc_lnk( ftud, 'U', -1. )
787        CALL lbc_lnk( fsud, 'U', -1. )
788        CALL lbc_lnk( ftvd, 'V', -1. )
789        CALL lbc_lnk( fsvd, 'V', -1. )
790        CALL lbc_lnk( psix_eiv, 'U', -1. )
791        CALL lbc_lnk( psiy_eiv, 'V', -1. )
792      !
793   END SUBROUTINE ldf_slp_grif
794
795
796   SUBROUTINE ldf_slp_mxl( prd, pn2 )
797      !!----------------------------------------------------------------------
798      !!                  ***  ROUTINE ldf_slp_mxl  ***
799      !!
800      !! ** Purpose :   Compute the slopes of iso-neutral surface just below
801      !!              the mixed layer.
802      !!
803      !! ** Method :
804      !!      The slope in the i-direction is computed at u- and w-points
805      !!   (uslp, wslpi) and the slope in the j-direction is computed at
806      !!   v- and w-points (vslp, wslpj).
807      !!   They are bounded by 1/100 over the whole ocean, and within the
808      !!   surface layer they are bounded by the distance to the surface
809      !!   ( slope<= depth/l  where l is the length scale of horizontal
810      !!   diffusion (here, aht=2000m2/s ==> l=20km with a typical velocity
811      !!   of 10cm/s)
812      !!
813      !! ** Action :   Compute uslp, wslpi, and vslp, wslpj, the i- and  j-slopes
814      !!   of now neutral surfaces at u-, w- and v- w-points, resp.
815      !!----------------------------------------------------------------------
816      USE oce , zgru  => ua   ! ua, va used as workspace and set to hor.
817      USE oce , zgrv  => va   ! density gradient in ldf_slp
818      !!
819      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) ::   prd   ! in situ density
820      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) ::   pn2   ! Brunt-Vaisala frequency (locally ref.)
821      !!
822      INTEGER  ::   ji, jj, jk           ! dummy loop indices
823      INTEGER  ::   ik, ikm1             ! temporary integers
824      REAL(wp) ::   zeps, zmg, zm05g     ! temporary scalars
825      REAL(wp) ::   zcoef1, zcoef2       !    -         -
826      REAL(wp) ::   zau, zbu, zai, zbi   !    -         -
827      REAL(wp) ::   zav, zbv, zaj, zbj   !    -         -
828      REAL(wp), DIMENSION(jpi,jpj) ::   zwy   ! 2D workspace
829      !!----------------------------------------------------------------------
830
831      zeps  =  1.e-20                    ! Local constant initialization
832      zmg   = -1.0 / grav
833      zm05g = -0.5 / grav
834      !
835      uslpml (1,:) = 0.e0      ;      uslpml (jpi,:) = 0.e0
836      vslpml (1,:) = 0.e0      ;      vslpml (jpi,:) = 0.e0
837      wslpiml(1,:) = 0.e0      ;      wslpiml(jpi,:) = 0.e0
838      wslpjml(1,:) = 0.e0      ;      wslpjml(jpi,:) = 0.e0
839
840      !                                  ! surface mixed layer mask
841      DO jk = 1, jpk                          ! =1 inside the mixed layer, =0 otherwise
842# if defined key_vectopt_loop
843         DO jj = 1, 1
844            DO ji = 1, jpij   ! vector opt. (forced unrolling)
845# else
846         DO jj = 1, jpj
847            DO ji = 1, jpi
848# endif
849               ik = nmln(ji,jj) - 1
850               IF( jk <= ik ) THEN   ;   omlmask(ji,jj,jk) = 1.e0
851               ELSE                  ;   omlmask(ji,jj,jk) = 0.e0
852               ENDIF
853            END DO
854         END DO
855      END DO
856
857
858      ! Slopes of isopycnal surfaces just before bottom of mixed layer
859      ! --------------------------------------------------------------
860      ! The slope are computed as in the 3D case.
861      ! A key point here is the definition of the mixed layer at u- and v-points.
862      ! It is assumed to be the maximum of the two neighbouring T-point mixed layer depth.
863      ! Otherwise, a n2 value inside the mixed layer can be involved in the computation
864      ! of the slope, resulting in a too steep diagnosed slope and thus a spurious eddy
865      ! induce velocity field near the base of the mixed layer.
866      !-----------------------------------------------------------------------
867      !
868      zwy(:,jpj) = 0.e0            !* vertical density gradient for u-slope (from N^2)
869      zwy(jpi,:) = 0.e0
870# if defined key_vectopt_loop
871      DO jj = 1, 1
872         DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
873# else
874      DO jj = 1, jpjm1
875         DO ji = 1, jpim1
876# endif
877            ik = MAX( 1, nmln(ji,jj) , nmln(ji+1,jj) )       ! avoid spurious recirculation
878            ik = MIN( ik, jpkm1 )                            ! if ik = jpk take jpkm1 values
879            zwy(ji,jj) = zmg * ( prd(ji,jj,ik) + 1. ) * (    pn2  (ji,jj,ik) + pn2  (ji,jj,ik+1)     )   &
880               &                                      / MAX( tmask(ji,jj,ik) + tmask(ji,jj,ik+1), 1. )
881         END DO
882      END DO
883      CALL lbc_lnk( zwy, 'U', 1. )      ! lateral boundary conditions NO sign change
884
885      !                            !* Slope at u points
886# if defined key_vectopt_loop
887      DO jj = 1, 1
888         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
889# else
890      DO jj = 2, jpjm1
891         DO ji = 2, jpim1
892# endif
893            ! horizontal and vertical density gradient at u-points
894            ik = MAX( 1, nmln(ji,jj) , nmln(ji+1,jj) )
895            ik = MIN( ik, jpkm1 )
896            zau = 1./ e1u(ji,jj) * zgru(ji,jj,ik)
897            zbu = 0.5*( zwy(ji,jj) + zwy(ji+1,jj) )
898            ! bound the slopes: abs(zw.)<= 1/100 and zb..<0
899            !                         kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
900            zbu = MIN( zbu, -100.*ABS(zau), -7.e+3/fse3u(ji,jj,ik)*ABS(zau) )
901            ! uslpml
902            uslpml (ji,jj) = zau / ( zbu - zeps ) * umask (ji,jj,ik)
903         END DO
904      END DO
905      CALL lbc_lnk( uslpml, 'U', -1. )      ! lateral boundary conditions (i-gradient => sign change)
906
907      !                            !* vertical density gradient for v-slope (from N^2)
908# if defined key_vectopt_loop
909      DO jj = 1, 1
910         DO ji = 1, jpij-jpi   ! vector opt. (forced unrolling)
911# else
912      DO jj = 1, jpjm1
913         DO ji = 1, jpim1
914# endif
915            ik = MAX( 1, nmln(ji,jj) , nmln(ji,jj+1) )
916            ik = MIN( ik, jpkm1 )
917            zwy(ji,jj) = zmg * ( prd(ji,jj,ik) + 1. ) * (    pn2  (ji,jj,ik) + pn2  (ji,jj,ik+1)     )   &
918               &                                      / MAX( tmask(ji,jj,ik) + tmask(ji,jj,ik+1), 1. )
919         END DO
920      END DO
921      CALL lbc_lnk( zwy, 'V', 1. )      ! lateral boundary conditions NO sign change
922
923      !                            !* Slope at v points
924# if defined key_vectopt_loop
925      DO jj = 1, 1
926         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
927# else
928      DO jj = 2, jpjm1
929         DO ji = 2, jpim1
930# endif
931            ! horizontal and vertical density gradient at v-points
932            ik = MAX( 1, nmln(ji,jj) , nmln(ji,jj+1) )
933            ik = MIN( ik,jpkm1 )
934            zav = 1./ e2v(ji,jj) * zgrv(ji,jj,ik)
935            zbv = 0.5*( zwy(ji,jj) + zwy(ji,jj+1) )
936            ! bound the slopes: abs(zw.)<= 1/100 and zb..<0
937            !                         kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
938            zbv = MIN( zbv, -100.*ABS(zav), -7.e+3/fse3v(ji,jj,ik)*ABS( zav ) )
939            ! vslpml
940            vslpml (ji,jj) = zav / ( zbv - zeps ) * vmask (ji,jj,ik)
941         END DO
942      END DO
943      CALL lbc_lnk( vslpml, 'V', -1. )      ! lateral boundary conditions (j-gradient => sign change)
944
945
946      !                            !* vertical density gradient for w-slope (from N^2)
947# if defined key_vectopt_loop
948      DO jj = 1, 1
949         DO ji = 1, jpij   ! vector opt. (forced unrolling)
950# else
951      DO jj = 1, jpj
952         DO ji = 1, jpi
953# endif
954            ik = nmln(ji,jj) + 1
955            ik = MIN( ik, jpk )
956            ikm1 = MAX ( 1, ik-1)
957            zwy (ji,jj) = zm05g * pn2 (ji,jj,ik) *     &
958               &             ( prd (ji,jj,ik) + prd (ji,jj,ikm1) + 2. )
959         END DO
960      END DO
961
962      !                            !* Slopes at w points
963# if defined key_vectopt_loop
964      DO jj = 1, 1
965         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
966# else
967      DO jj = 2, jpjm1
968         DO ji = 2, jpim1
969# endif
970            ik = nmln(ji,jj) + 1
971            ik = MIN( ik, jpk )
972            ikm1 = MAX ( 1, ik-1 )
973            ! horizontal density i-gradient at w-points
974            zcoef1 = MAX( zeps, umask(ji-1,jj,ik  )+umask(ji,jj,ik  )   &
975               &               +umask(ji-1,jj,ikm1)+umask(ji,jj,ikm1) )
976            zcoef1 = 1. / ( zcoef1 * e1t (ji,jj) )
977            zai = zcoef1 * (  zgru(ji  ,jj,ik  ) + zgru(ji  ,jj,ikm1)   &
978               &            + zgru(ji-1,jj,ikm1) + zgru(ji-1,jj,ik  ) ) * tmask (ji,jj,ik)
979            ! horizontal density j-gradient at w-points
980            zcoef2 = MAX( zeps, vmask(ji,jj-1,ik  )+vmask(ji,jj,ikm1)    &
981               &               +vmask(ji,jj-1,ikm1)+vmask(ji,jj,ik  ) )
982            zcoef2 = 1.0 / ( zcoef2 *  e2t (ji,jj) )
983            zaj = zcoef2 * (  zgrv(ji,jj  ,ik  ) + zgrv(ji,jj  ,ikm1)   &
984               &            + zgrv(ji,jj-1,ikm1) + zgrv(ji,jj-1,ik  ) ) * tmask (ji,jj,ik)
985            ! bound the slopes: abs(zw.)<= 1/100 and zb..<0.
986            !                   static instability:
987            !                   kxz max= ah slope max =< e1 e3 /(pi**2 2 dt)
988            zbi = MIN ( zwy (ji,jj),- 100.*ABS(zai), -7.e+3/fse3w(ji,jj,ik)*ABS(zai) )
989            zbj = MIN ( zwy (ji,jj), -100.*ABS(zaj), -7.e+3/fse3w(ji,jj,ik)*ABS(zaj) )
990            ! wslpiml and wslpjml
991            wslpiml (ji,jj) = zai / ( zbi - zeps) * tmask (ji,jj,ik)
992            wslpjml (ji,jj) = zaj / ( zbj - zeps) * tmask (ji,jj,ik)
993         END DO
994      END DO
995      CALL lbc_lnk( wslpiml, 'W', -1. )   ;   CALL lbc_lnk( wslpjml, 'W', -1. )      ! lateral boundary conditions
996      !
997   END SUBROUTINE ldf_slp_mxl
998
999
1000   SUBROUTINE ldf_slp_init
1001      !!----------------------------------------------------------------------
1002      !!                  ***  ROUTINE ldf_slp_init  ***
1003      !!
1004      !! ** Purpose :   Initialization for the isopycnal slopes computation
1005      !!
1006      !! ** Method  :   read the nammbf namelist and check the parameter
1007      !!      values called by tra_dmp at the first timestep (nit000)
1008      !!----------------------------------------------------------------------
1009      INTEGER ::   ji, jj, jk   ! dummy loop indices
1010      !!----------------------------------------------------------------------
1011     
1012      IF(lwp) THEN   
1013         WRITE(numout,*)
1014         WRITE(numout,*) 'ldf_slp : direction of lateral mixing'
1015         WRITE(numout,*) '~~~~~~~'
1016      ENDIF
1017
1018      ! Direction of lateral diffusion (tracers and/or momentum)
1019      ! ------------------------------
1020      ! set the slope to zero (even in s-coordinates)
1021
1022      uslp (:,:,:) = 0.e0
1023      vslp (:,:,:) = 0.e0
1024      wslpi(:,:,:) = 0.e0
1025      wslpj(:,:,:) = 0.e0
1026
1027      uslpml (:,:) = 0.e0
1028      vslpml (:,:) = 0.e0
1029      wslpiml(:,:) = 0.e0
1030      wslpjml(:,:) = 0.e0
1031
1032      IF( (ln_traldf_hor .OR. ln_dynldf_hor) .AND. .NOT. (lk_vvl .AND. ln_rstart) ) THEN
1033         IF(lwp) THEN
1034            WRITE(numout,*) '          Horizontal mixing in s-coordinate: slope = slope of s-surfaces'
1035         ENDIF
1036
1037         ! geopotential diffusion in s-coordinates on tracers and/or momentum
1038         ! The slopes of s-surfaces are computed once (no call to ldfslp in step)
1039         ! The slopes for momentum diffusion are i- or j- averaged of those on tracers
1040
1041         ! set the slope of diffusion to the slope of s-surfaces
1042         !      ( c a u t i o n : minus sign as fsdep has positive value )
1043         DO jk = 1, jpk
1044            DO jj = 2, jpjm1
1045               DO ji = fs_2, fs_jpim1   ! vector opt.
1046                  uslp (ji,jj,jk) = -1./e1u(ji,jj) * ( fsdept(ji+1,jj,jk) - fsdept(ji ,jj ,jk) ) * umask(ji,jj,jk)
1047                  vslp (ji,jj,jk) = -1./e2v(ji,jj) * ( fsdept(ji,jj+1,jk) - fsdept(ji ,jj ,jk) ) * vmask(ji,jj,jk)
1048                  wslpi(ji,jj,jk) = -1./e1t(ji,jj) * ( fsdepw(ji+1,jj,jk) - fsdepw(ji-1,jj,jk) ) * tmask(ji,jj,jk) * 0.5
1049                  wslpj(ji,jj,jk) = -1./e2t(ji,jj) * ( fsdepw(ji,jj+1,jk) - fsdepw(ji,jj-1,jk) ) * tmask(ji,jj,jk) * 0.5
1050               END DO
1051            END DO
1052         END DO
1053         ! Lateral boundary conditions on the slopes
1054         CALL lbc_lnk( uslp , 'U', -1. )      ;      CALL lbc_lnk( vslp , 'V', -1. )
1055         CALL lbc_lnk( wslpi, 'W', -1. )      ;      CALL lbc_lnk( wslpj, 'W', -1. )
1056      ENDIF
1057      !
1058   END SUBROUTINE ldf_slp_init
1059
1060#else
1061   !!------------------------------------------------------------------------
1062   !!   Dummy module :                 NO Rotation of lateral mixing tensor
1063   !!------------------------------------------------------------------------
1064   LOGICAL, PUBLIC, PARAMETER ::   lk_ldfslp = .FALSE.    !: slopes flag
1065CONTAINS
1066   SUBROUTINE ldf_slp( kt, prd, pn2 )        ! Dummy routine
1067      INTEGER, INTENT(in) :: kt 
1068      REAL, DIMENSION(:,:,:), INTENT(in) :: prd, pn2
1069      WRITE(*,*) 'ldf_slp: You should not have seen this print! error?', kt, prd(1,1,1), pn2(1,1,1)
1070   END SUBROUTINE ldf_slp
1071   SUBROUTINE ldf_slp_init       ! Dummy routine
1072   END SUBROUTINE ldf_slp_init
1073#endif
1074
1075   !!======================================================================
1076END MODULE ldfslp
Note: See TracBrowser for help on using the repository browser.