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

source: branches/2015/dev_r5803_NOC_WAD/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 5870

Last change on this file since 5870 was 5870, checked in by acc, 8 years ago

Branch 2015/dev_r5803_NOC_WAD. Merge in trunk changes from 5803 to 5869 in preparation for merge. Also tidied and reorganised some wetting and drying code. Renamed wadlmt.F90 to wetdry.F90. Wetting drying code changes restricted to domzgr.F90, domvvl.F90 nemogcm.F90 sshwzv.F90, dynspg_ts.F90, wetdry.F90 and dynhpg.F90. Code passes full SETTE tests with ln_wd=.false.. Still awaiting test case for checking with ln_wd=.false.

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