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

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

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

Last change on this file since 3393 was 3294, checked in by rblod, 12 years ago

Merge of 3.4beta into the trunk

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