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

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

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

Last change on this file since 3116 was 3116, checked in by cetlod, 13 years ago

dev_NEMO_MERGE_2011: add in changes dev_NOC_UKMO_MERGE developments

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