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

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 3221

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