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

source: trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 3558

Last change on this file since 3558 was 3558, checked in by rblod, 11 years ago

Fix issues when using key_nosignedzeo, see ticket #996

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