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

source: branches/DEV_r2191_3partymerge2010/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 2208

Last change on this file since 2208 was 2208, checked in by rblod, 14 years ago

Put FCM NEMO code changes in DEV_r2191_3partymerge2010 branch

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