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

Last change on this file since 3229 was 3229, checked in by charris, 10 years ago

Added timing calls to most significant routines in LDF, SBC and ZDF.

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