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 @ 2435

Last change on this file since 2435 was 2435, checked in by cetlod, 13 years ago

Improve the 1D vertical configuration in v3.3beta

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