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.
ldfdyn.F90 in branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/OPA_SRC/LDF – NEMO

source: branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/OPA_SRC/LDF/ldfdyn.F90 @ 8016

Last change on this file since 8016 was 8016, checked in by timgraham, 7 years ago

Delete some remaining "USE wrk_array" lines

  • Property svn:keywords set to Id
File size: 24.9 KB
Line 
1MODULE ldfdyn
2   !!======================================================================
3   !!                       ***  MODULE  ldfdyn  ***
4   !! Ocean physics:  lateral viscosity coefficient
5   !!=====================================================================
6   !! History :  OPA  ! 1997-07  (G. Madec)  multi dimensional coefficients
7   !!   NEMO     1.0  ! 2002-09  (G. Madec)  F90: Free form and module
8   !!            3.7  ! 2014-01  (F. Lemarie, G. Madec)  restructuration/simplification of ahm specification,
9   !!                 !                                  add velocity dependent coefficient and optional read in file
10   !!----------------------------------------------------------------------
11
12   !!----------------------------------------------------------------------
13   !!   ldf_dyn_init  : initialization, namelist read, and parameters control
14   !!   ldf_dyn       : update lateral eddy viscosity coefficients at each time step
15   !!----------------------------------------------------------------------
16   USE oce             ! ocean dynamics and tracers   
17   USE dom_oce         ! ocean space and time domain
18   USE phycst          ! physical constants
19   USE ldfc1d_c2d      ! lateral diffusion: 1D and 2D cases
20   !
21   USE in_out_manager  ! I/O manager
22   USE iom             ! I/O module for ehanced bottom friction file
23   USE timing          ! Timing
24   USE lib_mpp         ! distribued memory computing library
25   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
26
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC   ldf_dyn_init   ! called by nemogcm.F90
31   PUBLIC   ldf_dyn        ! called by step.F90
32
33   !                                                !!* Namelist namdyn_ldf : lateral mixing on momentum *
34   LOGICAL , PUBLIC ::   ln_dynldf_lap   !: laplacian operator
35   LOGICAL , PUBLIC ::   ln_dynldf_blp   !: bilaplacian operator
36   LOGICAL , PUBLIC ::   ln_dynldf_lev   !: iso-level direction
37   LOGICAL , PUBLIC ::   ln_dynldf_hor   !: horizontal (geopotential) direction
38   LOGICAL , PUBLIC ::   ln_dynldf_iso   !: iso-neutral direction
39   INTEGER , PUBLIC ::   nn_ahm_ijk_t    !: choice of time & space variations of the lateral eddy viscosity coef.
40   REAL(wp), PUBLIC ::   rn_ahm_0        !: lateral laplacian eddy viscosity            [m2/s]
41   REAL(wp), PUBLIC ::   rn_ahm_b        !: lateral laplacian background eddy viscosity [m2/s]
42   REAL(wp), PUBLIC ::   rn_bhm_0        !: lateral bilaplacian eddy viscosity          [m4/s]
43                                         !! If nn_ahm_ijk_t = 32 a time and space varying Smagorinsky viscosity
44                                         !! will be computed.
45   REAL(wp), PUBLIC ::   rn_csmc         !: Smagorinsky constant of proportionality
46   REAL(wp), PUBLIC ::   rn_minfac       !: Multiplicative factor of theorectical minimum Smagorinsky viscosity
47   REAL(wp), PUBLIC ::   rn_maxfac       !: Multiplicative factor of theorectical maximum Smagorinsky viscosity
48
49   LOGICAL , PUBLIC ::   l_ldfdyn_time   !: flag for time variation of the lateral eddy viscosity coef.
50
51   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ahmt, ahmf   !: eddy diffusivity coef. at U- and V-points   [m2/s or m4/s]
52   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dtensq       !: horizontal tension squared         (Smagorinsky only)
53   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dshesq       !: horizontal shearing strain squared (Smagorinsky only)
54   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   esqt, esqf   !: Square of the local gridscale (e1e2/(e1+e2))**2           
55
56   REAL(wp) ::   r1_12   = 1._wp / 12._wp    ! =1/12
57   REAL(wp) ::   r1_4    = 0.25_wp           ! =1/4
58   REAL(wp) ::   r1_8    = 0.125_wp          ! =1/8
59   REAL(wp) ::   r1_288  = 1._wp / 288._wp   ! =1/( 12^2 * 2 )
60
61   !! * Substitutions
62#  include "vectopt_loop_substitute.h90"
63   !!----------------------------------------------------------------------
64   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
65   !! $Id$
66   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
67   !!----------------------------------------------------------------------
68CONTAINS
69
70   SUBROUTINE ldf_dyn_init
71      !!----------------------------------------------------------------------
72      !!                  ***  ROUTINE ldf_dyn_init  ***
73      !!                   
74      !! ** Purpose :   set the horizontal ocean dynamics physics
75      !!
76      !! ** Method  :   the eddy viscosity coef. specification depends on:
77      !!              - the operator:
78      !!             ln_dynldf_lap = T     laplacian operator
79      !!             ln_dynldf_blp = T   bilaplacian operator
80      !!              - the parameter nn_ahm_ijk_t:
81      !!    nn_ahm_ijk_t  =  0 => = constant
82      !!                  = 10 => = F(z) :     = constant with a reduction of 1/4 with depth
83      !!                  =-20 => = F(i,j)     = shape read in 'eddy_viscosity.nc' file
84      !!                  = 20    = F(i,j)     = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
85      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_viscosity.nc'  file
86      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
87      !!                  = 31    = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
88      !!                                                           or |u|e^3/12 bilaplacian operator )
89      !!                  = 32    = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
90      !!                                                           (   L^2|D|      laplacian operator
91      !!                                                           or  L^4|D|/8  bilaplacian operator )
92      !!----------------------------------------------------------------------
93      INTEGER  ::   ji, jj, jk        ! dummy loop indices
94      INTEGER  ::   ierr, inum, ios   ! local integer
95      REAL(wp) ::   zah0              ! local scalar
96      !
97      NAMELIST/namdyn_ldf/ ln_dynldf_lap, ln_dynldf_blp,                  &
98         &                 ln_dynldf_lev, ln_dynldf_hor, ln_dynldf_iso,   &
99         &                 nn_ahm_ijk_t , rn_ahm_0, rn_ahm_b, rn_bhm_0,   &
100         &                 rn_csmc      , rn_minfac, rn_maxfac
101      !!----------------------------------------------------------------------
102      !
103      REWIND( numnam_ref )              ! Namelist namdyn_ldf in reference namelist : Lateral physics
104      READ  ( numnam_ref, namdyn_ldf, IOSTAT = ios, ERR = 901)
105901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn_ldf in reference namelist', lwp )
106
107      REWIND( numnam_cfg )              ! Namelist namdyn_ldf in configuration namelist : Lateral physics
108      READ  ( numnam_cfg, namdyn_ldf, IOSTAT = ios, ERR = 902 )
109902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn_ldf in configuration namelist', lwp )
110      IF(lwm) WRITE ( numond, namdyn_ldf )
111
112      IF(lwp) THEN                      ! Parameter print
113         WRITE(numout,*)
114         WRITE(numout,*) 'ldf_dyn : lateral momentum physics'
115         WRITE(numout,*) '~~~~~~~'
116         WRITE(numout,*) '   Namelist namdyn_ldf : set lateral mixing parameters'
117         !
118         WRITE(numout,*) '      type :'
119         WRITE(numout,*) '         laplacian operator                   ln_dynldf_lap = ', ln_dynldf_lap
120         WRITE(numout,*) '         bilaplacian operator                 ln_dynldf_blp = ', ln_dynldf_blp
121         !
122         WRITE(numout,*) '      direction of action :'
123         WRITE(numout,*) '         iso-level                            ln_dynldf_lev = ', ln_dynldf_lev
124         WRITE(numout,*) '         horizontal (geopotential)            ln_dynldf_hor = ', ln_dynldf_hor
125         WRITE(numout,*) '         iso-neutral                          ln_dynldf_iso = ', ln_dynldf_iso
126         !
127         WRITE(numout,*) '      coefficients :'
128         WRITE(numout,*) '         type of time-space variation         nn_ahm_ijk_t  = ', nn_ahm_ijk_t
129         WRITE(numout,*) '         lateral laplacian eddy viscosity     rn_ahm_0      = ', rn_ahm_0, ' m2/s'
130         WRITE(numout,*) '         background viscosity (iso case)      rn_ahm_b      = ', rn_ahm_b, ' m2/s'
131         WRITE(numout,*) '         lateral bilaplacian eddy viscosity   rn_bhm_0      = ', rn_bhm_0, ' m4/s'
132         WRITE(numout,*) '      smagorinsky settings (nn_ahm_ijk_t  = 32) :'
133         WRITE(numout,*) '         Smagorinsky coefficient              rn_csmc       = ', rn_csmc
134         WRITE(numout,*) '         factor multiplier for theorectical lower limit for '
135         WRITE(numout,*) '           Smagorinsky eddy visc (def. 1.0)   rn_minfac    = ', rn_minfac
136         WRITE(numout,*) '         factor multiplier for theorectical lower upper for '
137         WRITE(numout,*) '           Smagorinsky eddy visc (def. 1.0)   rn_maxfac    = ', rn_maxfac
138      ENDIF
139
140      !                                ! Parameter control
141      IF( .NOT.ln_dynldf_lap .AND. .NOT.ln_dynldf_blp ) THEN
142         IF(lwp) WRITE(numout,*) '   No viscous operator selected. ahmt and ahmf are not allocated'
143         l_ldfdyn_time = .FALSE.
144         RETURN
145      ENDIF
146      !
147      IF( ln_dynldf_blp .AND. ln_dynldf_iso ) THEN     ! iso-neutral bilaplacian not implemented
148         CALL ctl_stop( 'dyn_ldf_init: iso-neutral bilaplacian not coded yet' ) 
149      ENDIF
150
151      ! ... Space/Time variation of eddy coefficients
152      !                                               ! allocate the ahm arrays
153      ALLOCATE( ahmt(jpi,jpj,jpk) , ahmf(jpi,jpj,jpk) , STAT=ierr )
154      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate arrays')
155      !
156      ahmt(:,:,jpk) = 0._wp                           ! last level always 0 
157      ahmf(:,:,jpk) = 0._wp
158      !
159      !                                               ! value of eddy mixing coef.
160      IF    ( ln_dynldf_lap ) THEN   ;   zah0 =      rn_ahm_0         !   laplacian operator
161      ELSEIF( ln_dynldf_blp ) THEN   ;   zah0 = ABS( rn_bhm_0 )       ! bilaplacian operator
162      ELSE                                                                  ! NO viscous  operator
163         CALL ctl_warn( 'ldf_dyn_init: No lateral viscous operator used ' )
164      ENDIF
165      !
166      l_ldfdyn_time = .FALSE.                         ! no time variation except in case defined below
167      !
168      IF( ln_dynldf_lap .OR. ln_dynldf_blp ) THEN     ! only if a lateral diffusion operator is used
169         !
170         SELECT CASE(  nn_ahm_ijk_t  )                ! Specification of space time variations of ahmt, ahmf
171         !
172         CASE(   0  )      !==  constant  ==!
173            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = constant '
174            ahmt(:,:,:) = zah0 * tmask(:,:,:)
175            ahmf(:,:,:) = zah0 * fmask(:,:,:)
176            !
177         CASE(  10  )      !==  fixed profile  ==!
178            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F( depth )'
179            ahmt(:,:,1) = zah0 * tmask(:,:,1)                      ! constant surface value
180            ahmf(:,:,1) = zah0 * fmask(:,:,1)
181            CALL ldf_c1d( 'DYN', r1_4, ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )
182            !
183         CASE ( -20 )      !== fixed horizontal shape read in file  ==!
184            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F(i,j) read in eddy_viscosity.nc file'
185            CALL iom_open( 'eddy_viscosity_2D.nc', inum )
186            CALL iom_get ( inum, jpdom_data, 'ahmt_2d', ahmt(:,:,1) )
187            CALL iom_get ( inum, jpdom_data, 'ahmf_2d', ahmf(:,:,1) )
188            CALL iom_close( inum )
189!!gm Question : info for LAP or BLP case  to take into account the SQRT in the bilaplacian case ???
190!!              do we introduce a scaling by the max value of the array, and then multiply by zah0 ????
191!!              better:  check that the max is <=1  i.e. it is a shape from 0 to 1, not a coef that has physical dimension
192            DO jk = 2, jpkm1
193               ahmt(:,:,jk) = ahmt(:,:,1) * tmask(:,:,jk)
194               ahmf(:,:,jk) = ahmf(:,:,1) * fmask(:,:,jk)
195            END DO
196            !
197         CASE(  20  )      !== fixed horizontal shape  ==!
198            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F( e1, e2 ) or F( e1^3, e2^3 ) (lap. or blp. case)'
199            IF( ln_dynldf_lap )   CALL ldf_c2d( 'DYN', 'LAP', zah0, ahmt, ahmf )    ! surface value proportional to scale factor
200            IF( ln_dynldf_blp )   CALL ldf_c2d( 'DYN', 'BLP', zah0, ahmt, ahmf )    ! surface value proportional to scale factor^3
201            !
202         CASE( -30  )      !== fixed 3D shape read in file  ==!
203            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F(i,j,k) read in eddy_diffusivity_3D.nc file'
204            CALL iom_open( 'eddy_viscosity_3D.nc', inum )
205            CALL iom_get ( inum, jpdom_data, 'ahmt_3d', ahmt )
206            CALL iom_get ( inum, jpdom_data, 'ahmf_3d', ahmf )
207            CALL iom_close( inum )
208!!gm Question : info for LAP or BLP case  to take into account the SQRT in the bilaplacian case ????
209!!              do we introduce a scaling by the max value of the array, and then multiply by zah0 ????
210            DO jk = 1, jpkm1
211               ahmt(:,:,jk) = ahmt(:,:,jk) * tmask(:,:,jk)
212               ahmf(:,:,jk) = ahmf(:,:,jk) * fmask(:,:,jk)
213            END DO
214            !
215         CASE(  30  )       !==  fixed 3D shape  ==!
216            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F( latitude, longitude, depth )'
217            IF( ln_dynldf_lap )   CALL ldf_c2d( 'DYN', 'LAP', zah0, ahmt, ahmf )    ! surface value proportional to scale factor
218            IF( ln_dynldf_blp )   CALL ldf_c2d( 'DYN', 'BLP', zah0, ahmt, ahmf )    ! surface value proportional to scale factor
219            !                                                    ! reduction with depth
220            CALL ldf_c1d( 'DYN', r1_4, ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )
221            !
222         CASE(  31  )       !==  time varying 3D field  ==!
223            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F( latitude, longitude, depth , time )'
224            IF(lwp) WRITE(numout,*) '                                proportional to the velocity : |u|e/12 or |u|e^3/12'
225            !
226            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
227            !
228         CASE(  32  )       !==  time varying 3D field  ==!
229            IF(lwp) WRITE(numout,*) '          momentum mixing coef. = F( latitude, longitude, depth , time )'
230            IF(lwp) WRITE(numout,*) '             proportional to the local deformation rate and gridscale (Smagorinsky)'
231            IF(lwp) WRITE(numout,*) '                                                             : L^2|D| or L^4|D|/8'
232            !
233            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
234            !
235            ! allocate arrays used in ldf_dyn.
236            ALLOCATE( dtensq(jpi,jpj) , dshesq(jpi,jpj) , esqt(jpi,jpj) ,  esqf(jpi,jpj) , STAT=ierr )
237            IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate Smagorinsky arrays')
238            !
239            ! Set local gridscale values
240            DO jj = 2, jpjm1
241               DO ji = fs_2, fs_jpim1
242                  esqt(ji,jj) = ( e1e2t(ji,jj) /( e1t(ji,jj) + e2t(ji,jj) ) )**2 
243                  esqf(ji,jj) = ( e1e2f(ji,jj) /( e1f(ji,jj) + e2f(ji,jj) ) )**2 
244               END DO
245            END DO
246            !
247         CASE DEFAULT
248            CALL ctl_stop('ldf_dyn_init: wrong choice for nn_ahm_ijk_t, the type of space-time variation of ahm')
249         END SELECT
250         !
251         IF( ln_dynldf_blp .AND. .NOT. l_ldfdyn_time ) THEN       ! bilapcian and no time variation:
252            ahmt(:,:,:) = SQRT( ahmt(:,:,:) )                     ! take the square root of the coefficient
253            ahmf(:,:,:) = SQRT( ahmf(:,:,:) )
254         ENDIF
255         !
256      ENDIF
257      !
258   END SUBROUTINE ldf_dyn_init
259
260
261   SUBROUTINE ldf_dyn( kt )
262      !!----------------------------------------------------------------------
263      !!                  ***  ROUTINE ldf_dyn  ***
264      !!
265      !! ** Purpose :   update at kt the momentum lateral mixing coeff. (ahmt and ahmf)
266      !!
267      !! ** Method  :   time varying eddy viscosity coefficients:
268      !!
269      !!    nn_ahm_ijk_t = 31    ahmt, ahmf = F(i,j,k,t) = F(local velocity)
270      !!                         ( |u|e /12  or  |u|e^3/12 for laplacian or bilaplacian operator )
271      !!
272      !!    nn_ahm_ijk_t = 32    ahmt, ahmf = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
273      !!                         ( L^2|D|    or  L^4|D|/8  for laplacian or bilaplacian operator )
274      !!
275      !! ** note    :    in BLP cases the sqrt of the eddy coef is returned, since bilaplacian is en re-entrant laplacian
276      !! ** action  :    ahmt, ahmf   updated at each time step
277      !!----------------------------------------------------------------------
278      INTEGER, INTENT(in) ::   kt   ! time step index
279      !
280      INTEGER  ::   ji, jj, jk   ! dummy loop indices
281      REAL(wp) ::   zu2pv2_ij_p1, zu2pv2_ij, zu2pv2_ij_m1, zetmax, zefmax   ! local scalar
282      REAL(wp) ::   zcmsmag, zstabf_lo, zstabf_up, zdelta, zdb              ! local scalar
283      !!----------------------------------------------------------------------
284      !
285      IF( nn_timing == 1 )  CALL timing_start('ldf_dyn')
286      !
287      SELECT CASE(  nn_ahm_ijk_t  )       !== Eddy vicosity coefficients ==!
288      !
289      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
290         !
291         IF( ln_dynldf_lap   ) THEN        ! laplacian operator : |u| e /12 = |u/144| e
292            DO jk = 1, jpkm1
293               DO jj = 2, jpjm1
294                  DO ji = fs_2, fs_jpim1
295                     zu2pv2_ij_p1 = ub(ji  ,jj+1,jk) * ub(ji  ,jj+1,jk) + vb(ji+1,jj  ,jk) * vb(ji+1,jj  ,jk)
296                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
297                     zu2pv2_ij_m1 = ub(ji-1,jj  ,jk) * ub(ji-1,jj  ,jk) + vb(ji  ,jj-1,jk) * vb(ji  ,jj-1,jk)
298                     zetmax = MAX( e1t(ji,jj) , e2t(ji,jj) )
299                     zefmax = MAX( e1f(ji,jj) , e2f(ji,jj) )
300                     ahmt(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zetmax * tmask(ji,jj,jk)      ! 288= 12*12 * 2
301                     ahmf(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zefmax * fmask(ji,jj,jk)
302                  END DO
303               END DO
304            END DO
305         ELSEIF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( |u| e^3 /12 ) = sqrt( |u/144| e ) * e
306            DO jk = 1, jpkm1
307               DO jj = 2, jpjm1
308                  DO ji = fs_2, fs_jpim1
309                     zu2pv2_ij_p1 = ub(ji  ,jj+1,jk) * ub(ji  ,jj+1,jk) + vb(ji+1,jj  ,jk) * vb(ji+1,jj  ,jk)
310                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
311                     zu2pv2_ij_m1 = ub(ji-1,jj  ,jk) * ub(ji-1,jj  ,jk) + vb(ji  ,jj-1,jk) * vb(ji  ,jj-1,jk)
312                     zetmax = MAX( e1t(ji,jj) , e2t(ji,jj) )
313                     zefmax = MAX( e1f(ji,jj) , e2f(ji,jj) )
314                     ahmt(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zetmax  ) * zetmax * tmask(ji,jj,jk)
315                     ahmf(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zefmax  ) * zefmax * fmask(ji,jj,jk)
316                  END DO
317               END DO
318            END DO
319         ENDIF
320         !
321         CALL lbc_lnk( ahmt, 'T', 1. )   ;   CALL lbc_lnk( ahmf, 'F', 1. )
322         !
323         !
324      CASE(  32  )       !==  time varying 3D field  ==!   = F( local deformation rate and gridscale ) (Smagorinsky)
325         !
326         IF( ln_dynldf_lap .OR. ln_dynldf_blp  ) THEN        ! laplacian operator : (C_smag/pi)^2 L^2 |D|
327            !
328            zcmsmag = (rn_csmc/rpi)**2                                              ! (C_smag/pi)^2
329            zstabf_lo  = rn_minfac * rn_minfac / ( 2._wp * 4._wp * zcmsmag )        ! lower limit stability factor scaling
330            zstabf_up  = rn_maxfac / ( 4._wp * zcmsmag * 2._wp * rdt )              ! upper limit stability factor scaling
331            IF( ln_dynldf_blp ) zstabf_lo = ( 16._wp / 9._wp ) * zstabf_lo          ! provide |U|L^3/12 lower limit instead
332            !                                                                       ! of |U|L^3/16 in blp case
333            DO jk = 1, jpkm1
334               !
335               DO jj = 2, jpj
336                  DO ji = 2, jpi
337                     zdb = ( (  ub(ji,jj,jk) * r1_e2u(ji,jj) -  ub(ji-1,jj,jk) * r1_e2u(ji-1,jj) )  &
338                          &                  * r1_e1t(ji,jj) * e2t(ji,jj)                           &
339                          & - ( vb(ji,jj,jk) * r1_e1v(ji,jj) -  vb(ji,jj-1,jk) * r1_e1v(ji,jj-1) )  &
340                          &                  * r1_e2t(ji,jj) * e1t(ji,jj)    ) * tmask(ji,jj,jk)
341                     dtensq(ji,jj) = zdb*zdb
342                  END DO
343               END DO
344               !
345               DO jj = 1, jpjm1
346                  DO ji = 1, jpim1
347                     zdb = ( (  ub(ji,jj+1,jk) * r1_e1u(ji,jj+1) -  ub(ji,jj,jk) * r1_e1u(ji,jj) )  &
348                          &                    * r1_e2f(ji,jj)   * e1f(ji,jj)                       &
349                          & + ( vb(ji+1,jj,jk) * r1_e2v(ji+1,jj) -  vb(ji,jj,jk) * r1_e2v(ji,jj) )  &
350                          &                    * r1_e1f(ji,jj)   * e2f(ji,jj)  ) * fmask(ji,jj,jk)
351                     dshesq(ji,jj) = zdb*zdb
352                  END DO
353               END DO
354               !
355               DO jj = 2, jpjm1
356                  DO ji = fs_2, fs_jpim1
357                     !
358                     zu2pv2_ij_p1 = ub(ji  ,jj+1,jk) * ub(ji  ,jj+1,jk) + vb(ji+1,jj  ,jk) * vb(ji+1,jj  ,jk)
359                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
360                     zu2pv2_ij_m1 = ub(ji-1,jj  ,jk) * ub(ji-1,jj  ,jk) + vb(ji  ,jj-1,jk) * vb(ji  ,jj-1,jk)
361                                                     ! T-point value
362                     zdelta         = zcmsmag * esqt(ji,jj)                                        ! L^2 * (C_smag/pi)^2
363                     ahmt(ji,jj,jk) = zdelta * sqrt(          dtensq(ji,jj)   +                        &
364                                     &               r1_4 * ( dshesq(ji,jj)   + dshesq(ji,jj-1)   +    &
365                                     &                        dshesq(ji-1,jj) + dshesq(ji-1,jj-1) ) )
366                     ahmt(ji,jj,jk) =   MAX( ahmt(ji,jj,jk),   &
367                                     &   SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * zdelta * zstabf_lo ) ) ! Impose lower limit == minfac  * |U|L/2
368                     ahmt(ji,jj,jk) = MIN( ahmt(ji,jj,jk), zdelta * zstabf_up )                    ! Impose upper limit == maxfac  * L^2/(4*2dt)
369                                                     ! F-point value
370                     zdelta         = zcmsmag * esqf(ji,jj)                                        ! L^2 * (C_smag/pi)^2
371                     ahmf(ji,jj,jk) = zdelta * sqrt(          dshesq(ji,jj)   +                        &
372                                     &               r1_4 * ( dtensq(ji,jj)   + dtensq(ji,jj+1)   +    &
373                                     &                        dtensq(ji+1,jj) + dtensq(ji+1,jj+1) ) )
374                     ahmf(ji,jj,jk) =   MAX( ahmf(ji,jj,jk),   &
375                                     &   SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * zdelta * zstabf_lo ) ) ! Impose lower limit == minfac  * |U|L/2
376                     ahmf(ji,jj,jk) = MIN( ahmf(ji,jj,jk), zdelta * zstabf_up )                    ! Impose upper limit == maxfac  * L^2/(4*2dt)
377                     !
378                  END DO
379               END DO
380            END DO
381            !
382         ENDIF
383         !
384         IF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( (C_smag/pi)^2 L^4 |D|/8)
385                                       !                      = sqrt( A_lap_smag L^2/8 )
386                                       ! stability limits already applied to laplacian values
387                                       ! effective default limits are |U|L^3/12 < B_hm < L^4/(32*2dt)
388            !
389            DO jk = 1, jpkm1
390               !
391               DO jj = 2, jpjm1
392                  DO ji = fs_2, fs_jpim1
393                     !
394                     ahmt(ji,jj,jk) = sqrt( r1_8 * esqt(ji,jj) * ahmt(ji,jj,jk) )
395                     ahmf(ji,jj,jk) = sqrt( r1_8 * esqf(ji,jj) * ahmf(ji,jj,jk) )
396                     !
397                  END DO
398               END DO
399            END DO
400            !
401         ENDIF
402         !
403         CALL lbc_lnk( ahmt, 'T', 1. )   ;   CALL lbc_lnk( ahmf, 'F', 1. )
404         !
405      END SELECT
406      !
407      CALL iom_put( "ahmt_2d", ahmt(:,:,1) )   ! surface u-eddy diffusivity coeff.
408      CALL iom_put( "ahmf_2d", ahmf(:,:,1) )   ! surface v-eddy diffusivity coeff.
409      CALL iom_put( "ahmt_3d", ahmt(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
410      CALL iom_put( "ahmf_3d", ahmf(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
411      !
412      IF( nn_timing == 1 )  CALL timing_stop('ldf_dyn')
413      !
414   END SUBROUTINE ldf_dyn
415
416   !!======================================================================
417END MODULE ldfdyn
Note: See TracBrowser for help on using the repository browser.