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

source: branches/2014/dev_r4879_UKMO_NOC_MERGE/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 4925

Last change on this file since 4925 was 4812, checked in by mathiot, 10 years ago

UKMO02 ISF: correction of some TOP/ and OFF/ files (compilation and SETTE issues) + correction of bdyvol + slope of iso. beneath ice shelf + top friction in dynzdf_imp + minor change in zdfddm

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