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

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/LDF/ldfslp.F90 @ 2678

Last change on this file since 2678 was 2678, checked in by rblod, 13 years ago

Phasing branch dev_r2586_dynamic_mem with revision 2675 off the trunk

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