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/dev_r11233_AGRIF-05_jchanut_vert_coord_interp/src/OCE/LDF – NEMO

source: NEMO/branches/2019/dev_r11233_AGRIF-05_jchanut_vert_coord_interp/src/OCE/LDF/ldfslp.F90 @ 11574

Last change on this file since 11574 was 10425, checked in by smasson, 5 years ago

trunk: merge back dev_r10164_HPC09_ESIWACE_PREP_MERGE@10424 into the trunk

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