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

source: branches/2015/dev_r5151_UKMO_ISF/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 5956

Last change on this file since 5956 was 5956, checked in by mathiot, 8 years ago

ISF : merged trunk (5936) into branch

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