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

source: branches/2014/dev_CNRS_2014/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 4896

Last change on this file since 4896 was 4896, checked in by cetlod, 10 years ago

2014/dev_CNRS_2014 : merge the 1st branch onto dev_CNRS_2014, see ticket #1415

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