New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
ldfslp.F90 in NEMO/branches/2019/ENHANCE-02_ISF_nemo/src/OCE/LDF – NEMO

source: NEMO/branches/2019/ENHANCE-02_ISF_nemo/src/OCE/LDF/ldfslp.F90 @ 11395

Last change on this file since 11395 was 11395, checked in by mathiot, 5 years ago

ENHANCE-02_ISF_nemo : Initial commit isf simplification (add ISF directory, moved isf routine in and split isf cavity and isf parametrisation, ...) (ticket #2142)

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