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

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
ldfslp.F90 in branches/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 9616

Last change on this file since 9616 was 9176, checked in by andmirek, 6 years ago

#2001: OMP directives

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