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 NEMO/trunk/src/OCE/LDF – NEMO

source: NEMO/trunk/src/OCE/LDF/ldfslp.F90 @ 14433

Last change on this file since 14433 was 14433, checked in by smasson, 3 years ago

trunk: merge dev_r14312_MPI_Interface into the trunk, #2598

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