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 NEMO/branches/UKMO/NEMO_4.0.1_biharmonic_GM/src/OCE/LDF – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.1_biharmonic_GM/src/OCE/LDF/ldfdyn.F90 @ 12776

Last change on this file since 12776 was 12776, checked in by davestorkey, 5 years ago

UKMO/NEMO_4.0.1_biharmonic_GM : Allow option of running with standard viscosity and biharmonic GM. This version bit compares with previous version.

File size: 39.7 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   !!            4.0  ! 2020-02  (C. Wilson, ...) add bhm coefficient for bi-Laplacian GM implementation via momentum     
11   !!----------------------------------------------------------------------
12
13   !!----------------------------------------------------------------------
14   !!   ldf_dyn_init  : initialization, namelist read, and parameters control
15   !!   ldf_dyn       : update lateral eddy viscosity coefficients at each time step
16   !!----------------------------------------------------------------------
17   USE oce             ! ocean dynamics and tracers   
18   USE dom_oce         ! ocean space and time domain
19   USE phycst          ! physical constants
20   USE ldfslp          ! lateral diffusion: slopes of mixing orientation
21   USE ldfc1d_c2d      ! lateral diffusion: 1D and 2D cases
22   !
23   USE in_out_manager  ! I/O manager
24   USE iom             ! I/O module for ehanced bottom friction file
25   USE timing          ! Timing
26   USE lib_mpp         ! distribued memory computing library
27   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   ldf_dyn_init   ! called by nemogcm.F90
33   PUBLIC   ldf_dyn        ! called by step.F90
34
35   !                                    !!* Namelist namdyn_ldf : lateral mixing on momentum *
36   LOGICAL , PUBLIC ::   ln_dynldf_OFF   !: No operator (i.e. no explicit diffusion)
37   LOGICAL , PUBLIC ::   ln_dynldf_lap   !: laplacian operator
38   LOGICAL , PUBLIC ::   ln_dynldf_bgm   !: bi-Laplacian GM implementation via momentum
39   LOGICAL , PUBLIC ::   ln_dynldf_blp   !: bilaplacian operator
40   LOGICAL , PUBLIC ::   ln_dynldf_lev   !: iso-level direction
41   LOGICAL , PUBLIC ::   ln_dynldf_hor   !: horizontal (geopotential) direction
42!  LOGICAL , PUBLIC ::   ln_dynldf_iso   !: iso-neutral direction                        (see ldfslp)
43   INTEGER , PUBLIC ::   nn_ahm_ijk_t    !: choice of time & space variations of the lateral eddy viscosity coef.
44   !                                        !  time invariant coefficients:  aht = 1/2  Ud*Ld   (lap case)
45   !                                           !                             bht = 1/12 Ud*Ld^3 (blp case)
46   INTEGER , PUBLIC ::   nn_bhm_ijk_t    !: choice of time & space variations of the lateral bi-Laplacian GM coef.
47   REAL(wp), PUBLIC ::   rn_Uv                 !: lateral viscous velocity  [m/s]
48   REAL(wp), PUBLIC ::   rn_Lv                 !: lateral viscous length    [m]
49   !                                        ! Smagorinsky viscosity  (nn_ahm_ijk_t = 32)
50   REAL(wp), PUBLIC ::   rn_csmc               !: Smagorinsky constant of proportionality
51   REAL(wp), PUBLIC ::   rn_minfac             !: Multiplicative factor of theorectical minimum Smagorinsky viscosity
52   REAL(wp), PUBLIC ::   rn_maxfac             !: Multiplicative factor of theorectical maximum Smagorinsky viscosity
53   !                                        ! iso-neutral laplacian (ln_dynldf_lap=ln_dynldf_iso=T)
54   REAL(wp), PUBLIC ::   rn_ahm_b              !: lateral laplacian background eddy viscosity  [m2/s]
55   REAL(wp), PUBLIC ::   rn_bhm_b              !: lateral bi-Laplacian GM background eddy viscosity  [m4/s]
56   REAL(wp), PUBLIC ::   rn_bgmzcrit           !: critical depth (>=0) for linear tapering of bi-Lap. GM coef. to zero at surface [m]
57   !                                    !!* Parameter to control the type of lateral viscous operator
58   INTEGER, PARAMETER, PUBLIC ::   np_ERROR  =-10                       !: error in setting the operator
59   INTEGER, PARAMETER, PUBLIC ::   np_no_ldf = 00                       !: without operator (i.e. no lateral viscous trend)
60   !                          !!      laplacian     !    bilaplacian    !
61   INTEGER, PARAMETER, PUBLIC ::   np_lap    = 10   ,   np_blp    = 20  !: iso-level operator
62   INTEGER, PARAMETER, PUBLIC ::   np_lap_i  = 11                       !: iso-neutral or geopotential operator
63   INTEGER           , PUBLIC ::   nldf_dyn         !: type of lateral diffusion used defined from ln_dynldf_... (namlist logicals)
64   LOGICAL           , PUBLIC ::   l_ldfdyn_time    !: flag for time variation of the lateral eddy viscosity coef.
65
66   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ahmt, ahmf   !: eddy viscosity coef. at T- and F-points [m2/s or m4/s]
67   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   bhm          !: eddy viscosity coef. for bi-Laplacian GM at W-points (cf. avm) [m4/s]
68   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dtensq       !: horizontal tension squared         (Smagorinsky only)
69   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dshesq       !: horizontal shearing strain squared (Smagorinsky only)
70   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   esqt, esqf   !: Square of the local gridscale (e1e2/(e1+e2))**2           
71
72   REAL(wp) ::   r1_2    = 0.5_wp            ! =1/2
73   REAL(wp) ::   r1_4    = 0.25_wp           ! =1/4
74   REAL(wp) ::   r1_8    = 0.125_wp          ! =1/8
75   REAL(wp) ::   r1_12   = 1._wp / 12._wp    ! =1/12
76   REAL(wp) ::   r1_288  = 1._wp / 288._wp   ! =1/( 12^2 * 2 )
77
78   !! * Substitutions
79#  include "vectopt_loop_substitute.h90"
80   !!----------------------------------------------------------------------
81   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
82   !! $Id: ldfdyn.F90 11653 2019-10-04 12:34:02Z davestorkey $
83   !! Software governed by the CeCILL license (see ./LICENSE)
84   !!----------------------------------------------------------------------
85CONTAINS
86
87   SUBROUTINE ldf_dyn_init
88      !!----------------------------------------------------------------------
89      !!                  ***  ROUTINE ldf_dyn_init  ***
90      !!                   
91      !! ** Purpose :   set the horizontal ocean dynamics physics
92      !!
93      !! ** Method  :   the eddy viscosity coef. specification depends on:
94      !!              - the operator:
95      !!             ln_dynldf_lap = T     laplacian operator
96      !!             ln_dynldf_blp = T   bilaplacian operator
97      !!             ln_dynldf_bgm = T   bi-Laplacian GM operator
98      !!              - the parameter nn_ahm_ijk_t:
99      !!    nn_ahm_ijk_t  =  0 => = constant
100      !!                  = 10 => = F(z) :     = constant with a reduction of 1/4 with depth
101      !!                  =-20 => = F(i,j)     = shape read in 'eddy_viscosity.nc' file
102      !!                  = 20    = F(i,j)     = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
103      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_viscosity.nc'  file
104      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
105      !!                  = 31    = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
106      !!                                                           or |u|e^3/12 bilaplacian operator )
107      !!                  = 32    = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
108      !!                                                           (   L^2|D|      laplacian operator
109      !!                                                           or  L^4|D|/8  bilaplacian operator )
110      !!                - the parameter nn_bhm_ijk_t:
111      !!    nn_bhm_ijk_t  = 11 => = F(z) :     = constant, with BCs set to zero at top and bottom
112      !!    nn_bhm_ijk_t  = 12 => = F(z) :     = constant in interior, linear taper to zero at surface,
113      !!                                         for depths < rn_bgmzcrit, with BCs set to zero at top and bottom
114      !!    nn_bhm_ijk_t  = 13 => = F(z) :     = bhm0 X dx^2 X dz^2, with BCs set to zero at top and bottom
115      !!----------------------------------------------------------------------
116      INTEGER  ::   ji, jj, jk                     ! dummy loop indices
117      INTEGER  ::   ioptio, ierr, inum, ios, inn   ! local integer
118      INTEGER  ::   ikw                            ! local integer
119      REAL(wp) ::   zah0, zah_max, zUfac           ! local scalar
120      REAL(wp) ::   bhm0, lhscrit, lhscritmax      ! local scalar
121      CHARACTER(len=5) ::   cl_Units               ! units (m2/s or m4/s)
122      !!
123      NAMELIST/namdyn_ldf/ ln_dynldf_OFF, ln_dynldf_lap, ln_dynldf_blp,   &   ! type of operator
124         &                 ln_dynldf_bgm,                                 &   ! type of operator
125         &                 ln_dynldf_lev, ln_dynldf_hor, ln_dynldf_iso,   &   ! acting direction of the operator
126         &                 nn_ahm_ijk_t , rn_Uv    , rn_Lv,   rn_ahm_b,   &   ! lateral eddy coefficient
127         &                 nn_bhm_ijk_t, rn_bhm_b, rn_bgmzcrit,           &   ! lateral bi-Lap. GM coefficient
128         &                 rn_csmc      , rn_minfac    , rn_maxfac            ! Smagorinsky settings
129      !!----------------------------------------------------------------------
130      !
131      REWIND( numnam_ref )              ! Namelist namdyn_ldf in reference namelist : Lateral physics
132      READ  ( numnam_ref, namdyn_ldf, IOSTAT = ios, ERR = 901)
133901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namdyn_ldf in reference namelist' )
134
135      REWIND( numnam_cfg )              ! Namelist namdyn_ldf in configuration namelist : Lateral physics
136      READ  ( numnam_cfg, namdyn_ldf, IOSTAT = ios, ERR = 902 )
137902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namdyn_ldf in configuration namelist' )
138      IF(lwm) WRITE ( numond, namdyn_ldf )
139
140      IF(lwp) THEN                      ! Parameter print
141         WRITE(numout,*)
142         WRITE(numout,*) 'ldf_dyn : lateral momentum physics'
143         WRITE(numout,*) '~~~~~~~'
144         WRITE(numout,*) '   Namelist namdyn_ldf : set lateral mixing parameters'
145         !
146         WRITE(numout,*) '      type :'
147         WRITE(numout,*) '         no explicit diffusion                ln_dynldf_OFF = ', ln_dynldf_OFF
148         WRITE(numout,*) '         laplacian operator                   ln_dynldf_lap = ', ln_dynldf_lap
149         WRITE(numout,*) '         bilaplacian operator                 ln_dynldf_blp = ', ln_dynldf_blp
150         WRITE(numout,*) '         bi-Laplacian GM operator             ln_dynldf_bgm = ', ln_dynldf_bgm
151         WRITE(numout,*) '         critical depth for taper (if used)   rn_bgmzcrit   = ', rn_bgmzcrit
152         !
153         WRITE(numout,*) '      direction of action :'
154         WRITE(numout,*) '         iso-level                            ln_dynldf_lev = ', ln_dynldf_lev
155         WRITE(numout,*) '         horizontal (geopotential)            ln_dynldf_hor = ', ln_dynldf_hor
156         WRITE(numout,*) '         iso-neutral                          ln_dynldf_iso = ', ln_dynldf_iso
157         !
158         WRITE(numout,*) '      coefficients :'
159         WRITE(numout,*) '         type of time-space variation         nn_ahm_ijk_t  = ', nn_ahm_ijk_t
160         WRITE(numout,*) '         lateral viscous velocity  (if cst)      rn_Uv      = ', rn_Uv, ' m/s'
161         WRITE(numout,*) '         lateral viscous length    (if cst)      rn_Lv      = ', rn_Lv, ' m'
162         WRITE(numout,*) '         background viscosity (iso-lap case)     rn_ahm_b   = ', rn_ahm_b, ' m2/s'
163         WRITE(numout,*) '         bi-Laplacian GM background viscosity    rn_bhm_b   = ', rn_bhm_b, ' m4/s'
164         WRITE(numout,*) '      Smagorinsky settings (nn_ahm_ijk_t  = 32) :'
165         WRITE(numout,*) '         Smagorinsky coefficient              rn_csmc       = ', rn_csmc
166         WRITE(numout,*) '         factor multiplier for eddy visc.'
167         WRITE(numout,*) '            lower limit (default 1.0)         rn_minfac    = ', rn_minfac
168         WRITE(numout,*) '            upper limit (default 1.0)         rn_maxfac    = ', rn_maxfac
169      ENDIF
170
171      !
172      !           !==  type of lateral operator used  ==!   (set nldf_dyn)
173      !           !=====================================!
174      !
175      nldf_dyn = np_ERROR
176      ioptio = 0
177      !CW exclude combination of lap+blp (as in existing code), but allow other combinations
178      IF( ln_dynldf_OFF ) THEN   ;   nldf_dyn = np_no_ldf   ;   ioptio = ioptio + 1   ;   ENDIF
179      IF( ln_dynldf_lap ) THEN   ;                              ioptio = ioptio + 5   ;   ENDIF
180      IF( ln_dynldf_blp ) THEN   ;                              ioptio = ioptio + 5   ;   ENDIF
181      !
182      IF( (ioptio < 1).OR.(ioptio > 6) )   CALL ctl_stop( 'dyn_ldf_init: use at least ONE of the 4 operator options (NONE/lap/blp) but not (lap+blp)' )
183      !
184
185      IF(.NOT.ln_dynldf_OFF ) THEN     !==  direction ==>> type of operator  ==!
186         ioptio = 0
187         IF( ln_dynldf_lev )   ioptio = ioptio + 1
188         IF( ln_dynldf_hor )   ioptio = ioptio + 1
189         IF( ln_dynldf_iso )   ioptio = ioptio + 1
190         IF( ioptio /= 1   )   CALL ctl_stop( 'dyn_ldf_init: use ONE of the 3 direction options (level/hor/iso)' )
191         !
192         !                             ! Set nldf_dyn, the type of lateral diffusion, from ln_dynldf_... logicals
193         ierr = 0
194         IF( ln_dynldf_lap ) THEN         ! laplacian operator
195            IF( ln_zco ) THEN                   ! z-coordinate
196               IF ( ln_dynldf_lev )   nldf_dyn = np_lap     ! iso-level = horizontal (no rotation)
197               IF ( ln_dynldf_hor )   nldf_dyn = np_lap     ! iso-level = horizontal (no rotation)
198               IF ( ln_dynldf_iso )   nldf_dyn = np_lap_i   ! iso-neutral            (   rotation)
199            ENDIF
200            IF( ln_zps ) THEN                   ! z-coordinate with partial step
201               IF ( ln_dynldf_lev )   nldf_dyn = np_lap     ! iso-level              (no rotation)
202               IF ( ln_dynldf_hor )   nldf_dyn = np_lap     ! iso-level              (no rotation)
203               IF ( ln_dynldf_iso )   nldf_dyn = np_lap_i   ! iso-neutral            (   rotation)
204            ENDIF
205            IF( ln_sco ) THEN                   ! s-coordinate
206               IF ( ln_dynldf_lev )   nldf_dyn = np_lap     ! iso-level = horizontal (no rotation)
207               IF ( ln_dynldf_hor )   nldf_dyn = np_lap_i   ! horizontal             (   rotation)
208               IF ( ln_dynldf_iso )   nldf_dyn = np_lap_i   ! iso-neutral            (   rotation)
209            ENDIF
210         ENDIF
211         !
212         IF( ln_dynldf_blp ) THEN         ! bilaplacian operator
213            IF( ln_zco ) THEN                   ! z-coordinate
214               IF( ln_dynldf_lev )   nldf_dyn = np_blp   ! iso-level = horizontal (no rotation)
215               IF( ln_dynldf_hor )   nldf_dyn = np_blp   ! iso-level = horizontal (no rotation)
216               IF( ln_dynldf_iso )   ierr = 2            ! iso-neutral            (   rotation)
217            ENDIF
218            IF( ln_zps ) THEN                   ! z-coordinate with partial step
219               IF( ln_dynldf_lev )   nldf_dyn = np_blp   ! iso-level              (no rotation)
220               IF( ln_dynldf_hor )   nldf_dyn = np_blp   ! iso-level              (no rotation)
221               IF( ln_dynldf_iso )   ierr = 2            ! iso-neutral            (   rotation)
222            ENDIF
223            IF( ln_sco ) THEN                   ! s-coordinate
224               IF( ln_dynldf_lev )   nldf_dyn = np_blp   ! iso-level              (no rotation)
225               IF( ln_dynldf_hor )   ierr = 2            ! horizontal             (   rotation)
226               IF( ln_dynldf_iso )   ierr = 2            ! iso-neutral            (   rotation)
227            ENDIF
228         ENDIF
229         !
230
231         IF( ierr == 2 )   CALL ctl_stop( 'rotated bi-laplacian operator does not exist' )
232         !
233
234         IF( nldf_dyn == np_lap_i )   l_ldfslp = .TRUE.  ! rotation require the computation of the slopes
235         !
236      ENDIF
237      !
238      IF(lwp) THEN
239         WRITE(numout,*)
240         SELECT CASE( nldf_dyn )
241         CASE( np_no_ldf )   ;   WRITE(numout,*) '   ==>>>   NO lateral viscosity'
242         CASE( np_lap    )   ;   WRITE(numout,*) '   ==>>>   iso-level laplacian operator'
243         CASE( np_lap_i  )   ;   WRITE(numout,*) '   ==>>>   rotated laplacian operator with iso-level background'
244         CASE( np_blp    )   ;   WRITE(numout,*) '   ==>>>   iso-level bi-laplacian operator'
245         END SELECT
246         WRITE(numout,*)
247      ENDIF
248     
249      !
250      !           !==  Space/time variation of eddy coefficients  ==!
251      !           !=================================================!
252      !
253      l_ldfdyn_time = .FALSE.                ! no time variation except in case defined below
254      !
255      IF( ln_dynldf_OFF ) THEN
256         !
257         IF(lwp) WRITE(numout,*) '   ==>>>   No viscous operator selected. ahmt, ahmf are not allocated'
258         !
259      ELSE                                   !==  a lateral diffusion operator is used  ==!
260         !
261         !                                         ! allocate the ahm arrays
262         ALLOCATE( ahmt(jpi,jpj,jpk) , ahmf(jpi,jpj,jpk) , bhm(jpi,jpj,jpk) , STAT=ierr )
263         !
264         IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate arrays')
265         !
266         ahmt(:,:,:) = 0._wp                       ! init to 0 needed
267         ahmf(:,:,:) = 0._wp
268         !                                        ! value of lap/blp eddy mixing coef.
269         IF(     ln_dynldf_lap ) THEN   ;   zUfac = r1_2 *rn_Uv   ;   inn = 1   ;   cl_Units = ' m2/s'   !   laplacian
270         ELSEIF( ln_dynldf_blp ) THEN   ;   zUfac = r1_12*rn_Uv   ;   inn = 3   ;   cl_Units = ' m4/s'   ! bilaplacian
271         ENDIF
272         zah0    = zUfac *    rn_Lv**inn              ! mixing coefficient
273         zah_max = zUfac * (ra*rad)**inn              ! maximum reachable coefficient (value at the Equator)
274         !
275         SELECT CASE(  nn_ahm_ijk_t  )             !* Specification of space-time variations of ahmt, ahmf
276         !
277         CASE(   0  )      !==  constant  ==!
278            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity. = constant = ', zah0, cl_Units
279            ahmt(:,:,1:jpkm1) = zah0
280            ahmf(:,:,1:jpkm1) = zah0
281            !
282         CASE(  10  )      !==  fixed profile  ==!
283            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( depth )'
284            IF(lwp) WRITE(numout,*) '           surface viscous coef. = constant = ', zah0, cl_Units
285            ahmt(:,:,1) = zah0                        ! constant surface value
286            ahmf(:,:,1) = zah0
287            CALL ldf_c1d( 'DYN', ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )
288            !       
289         CASE ( -20 )      !== fixed horizontal shape read in file  ==!
290            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F(i,j) read in eddy_viscosity.nc file'
291            CALL iom_open( 'eddy_viscosity_2D.nc', inum )
292            CALL iom_get ( inum, jpdom_data, 'ahmt_2d', ahmt(:,:,1) )
293            CALL iom_get ( inum, jpdom_data, 'ahmf_2d', ahmf(:,:,1) )
294            CALL iom_close( inum )
295            DO jk = 2, jpkm1
296               ahmt(:,:,jk) = ahmt(:,:,1)
297               ahmf(:,:,jk) = ahmf(:,:,1)
298            END DO
299            !
300         CASE(  20  )      !== fixed horizontal shape  ==!
301            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( e1, e2 ) or F( e1^3, e2^3 ) (lap. or blp. case)'
302            IF(lwp) WRITE(numout,*) '           using a fixed viscous velocity = ', rn_Uv  ,' m/s   and   Lv = Max(e1,e2)'
303            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
304            CALL ldf_c2d( 'DYN', zUfac      , inn        , ahmt, ahmf )         ! surface value proportional to scale factor^inn
305            !
306         CASE( -30  )      !== fixed 3D shape read in file  ==!
307            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F(i,j,k) read in eddy_viscosity_3D.nc file'
308            CALL iom_open( 'eddy_viscosity_3D.nc', inum )
309            CALL iom_get ( inum, jpdom_data, 'ahmt_3d', ahmt )
310            CALL iom_get ( inum, jpdom_data, 'ahmf_3d', ahmf )
311            CALL iom_close( inum )
312            !
313         CASE(  30  )       !==  fixed 3D shape  ==!
314            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth )'
315            IF(lwp) WRITE(numout,*) '           using a fixed viscous velocity = ', rn_Uv  ,' m/s   and   Ld = Max(e1,e2)'
316            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
317            CALL ldf_c2d( 'DYN', zUfac      , inn        , ahmt, ahmf )         ! surface value proportional to scale factor^inn
318            CALL ldf_c1d( 'DYN', ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )  ! reduction with depth
319            !
320         CASE(  31  )       !==  time varying 3D field  ==!
321            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth , time )'
322            IF(lwp) WRITE(numout,*) '           proportional to the local velocity : 1/2 |u|e (lap) or 1/12 |u|e^3 (blp)'
323            !
324            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
325            !
326         CASE(  32  )       !==  time varying 3D field  ==!
327            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth , time )'
328            IF(lwp) WRITE(numout,*) '           proportional to the local deformation rate and gridscale (Smagorinsky)'
329            !
330            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
331            !
332            !                          ! allocate arrays used in ldf_dyn.
333            ALLOCATE( dtensq(jpi,jpj,jpk) , dshesq(jpi,jpj,jpk) , esqt(jpi,jpj) , esqf(jpi,jpj) , STAT=ierr )
334            IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate Smagorinsky arrays')
335            !
336            DO jj = 1, jpj             ! Set local gridscale values
337               DO ji = 1, jpi
338                  esqt(ji,jj) = ( 2._wp * e1e2t(ji,jj) / ( e1t(ji,jj) + e2t(ji,jj) ) )**2 
339                  esqf(ji,jj) = ( 2._wp * e1e2f(ji,jj) / ( e1f(ji,jj) + e2f(ji,jj) ) )**2 
340               END DO
341            END DO
342            !
343         CASE DEFAULT
344            CALL ctl_stop('ldf_dyn_init: wrong choice for nn_ahm_ijk_t, the type of space-time variation of ahm')
345         END SELECT
346
347         !
348         IF( .NOT.l_ldfdyn_time ) THEN             !* No time variation
349            IF(     ln_dynldf_lap ) THEN                 !   laplacian operator (mask only)
350               ahmt(:,:,1:jpkm1) =       ahmt(:,:,1:jpkm1)   * tmask(:,:,1:jpkm1)
351               ahmf(:,:,1:jpkm1) =       ahmf(:,:,1:jpkm1)   * fmask(:,:,1:jpkm1)
352            ELSEIF( ln_dynldf_blp ) THEN                 ! bilaplacian operator (square root + mask)
353               ahmt(:,:,1:jpkm1) = SQRT( ahmt(:,:,1:jpkm1) ) * tmask(:,:,1:jpkm1)
354               ahmf(:,:,1:jpkm1) = SQRT( ahmf(:,:,1:jpkm1) ) * fmask(:,:,1:jpkm1)
355            ENDIF
356         ENDIF
357         !
358      ENDIF
359      !
360!CW   Add separate structure options for bi-Laplacian GM, to allow it to be used in combination with other types of diffusion, above.
361      IF( ln_dynldf_bgm ) THEN
362         IF(lwp) WRITE(numout,*) '  ==>>> Initialising bi-Laplacian GM eddy viscosity coefficient.'
363         ALLOCATE( bhm(jpi,jpj,jpk) , STAT=ierr )
364         bhm(:,:,:)  = 0._wp                       ! init to 0 needed
365         bhm0 = rn_bhm_b       ;             ;   cl_Units = ' m4/s'   ! bi-Laplacian GM
366         SELECT CASE (nn_bhm_ijk_t)
367            !CW Define bi-Laplacian GM diffusivity and test the related computational stability criterion   
368
369         CASE(  11  )      !==  fixed profile: constant except for zero at top and bottom, so that eddy-induced velocity, w*=0  ==!
370            IF(lwp) WRITE(numout,*) '   ==>>>   bi-Laplacian GM eddy viscosity = constant in interior, zero at top and bottom'
371            IF(lwp) WRITE(numout,*) '           interior viscous coef. = constant = ', bhm0, cl_Units
372
373            !     First set to constant in interior, all levels
374            DO jj = 2, jpjm1
375               DO ji = 2, jpim1
376                  bhm(ji,jj,:) = bhm0
377               ENDDO
378            ENDDO
379
380            !     Surface BC : set to zero
381            bhm(:,:,1)   = 0._wp
382            !     Flat bottom BC : set to zero
383            bhm(:,:,jpk) = 0._wp 
384            !     Variable bathymetry case: diffusive fluxes masked in dyn_ldf_bgm
385
386            !CW Test criterion for two grid-length mode for BGM scheme with leapfrog timestepping
387            ! \frac{{\kappa \Delta t}{(\Delta z)^2 (\Delta x)^2}< 1/32
388            ! test at T-point - may need to revise this choice!!!!
389
390            ! Find domain maximum value of LHS, then test inequality.
391
392            ! initialise
393            lhscritmax=0._wp
394
395            DO jj = 2, jpjm1
396               DO ji = 2, jpim1
397                  ikw = mbkt(ji,jj) !bottom last wet T-level
398                  DO jk = 2, ikw
399                     lhscrit=bhm0*rdt/(e3t_n(ji,jj,jk)**2*e1t(ji,jj)**2)
400                     IF(lhscrit .gt. lhscritmax) lhscritmax=lhscrit
401                  ENDDO
402               ENDDO
403            ENDDO
404
405            IF ( lhscritmax .ge. 1._wp/32._wp) THEN
406               IF(lwp) THEN
407                  WRITE(numout,*)
408                  WRITE(numout,*) ' WARNING : Bi-Laplacian GM eddy viscosity is not',              &
409                       &            ' consistent with the model time step and grid setup: ', &
410                       &            ' LHS=',lhscrit, & 
411                       &            'Should be < 0.03125 in order to avoid computational instability (Mike Bell, Dave Storkey)'
412                  WRITE(numout,*) ' ========='
413                  WRITE(numout,*)
414               ENDIF
415               ! CW: warn, don't stop, as we may wish to violate this theoretical limit.
416               ! CALL ctl_stop( 'ldf_dyn_init', 'Inconsistent Bi-Lap. GM diffusivity')
417            ELSE
418               IF(lwp) THEN
419                  WRITE(numout,*)
420                  WRITE(numout,*) ' INFORMATION : Bi-Laplacian GM eddy viscosity is ',              &
421                       &            ' consistent with the model time step and grid setup: ', &
422                       &            ' LHS=',lhscrit, & 
423                       &            'Should be < 0.03125 in order to avoid computational instability (Mike Bell, Dave Storkey)'
424                  WRITE(numout,*) ' ========='
425                  WRITE(numout,*)
426               ENDIF
427            ENDIF
428            !-----------------------------
429         CASE(  12  )      !==  fixed profile: constant in interior, zero at bottom and linearly-tapering to zero at top, over depth shallower than rn_bgmzcrit  ==!
430            IF(lwp) WRITE(numout,*) '   ==>>>   bi-Laplacian GM eddy viscosity = constant in interior, zero at bottom and'
431            IF(lwp) WRITE(numout,*) '           linearly-tapering to zero at top, over depth shallower than ', rn_bgmzcrit, ' metres.'
432            IF(lwp) WRITE(numout,*) '           Interior viscous coef. = constant = ', bhm0, cl_Units
433
434            !     Impose linear taper from interior value to zero at zero depth, for depths < rn_bgmzcrit
435            !     N.B. Test criterion for two grid-length mode for BGM scheme with leapfrog timestepping not implemented for this case yet.
436
437            DO jk = 1, jpk
438
439               IF (gdept_1d(jk) .lt. rn_bgmzcrit) THEN
440
441                  DO jj = 2, jpjm1
442                     DO ji = 2, jpim1
443                        bhm(ji,jj,jk) = bhm0 * gdepw_1d(jk)/rn_bgmzcrit
444                     ENDDO
445                  ENDDO
446
447               ELSE
448
449                  ! constant at depth rn_bgmzcrit and larger
450                  DO jj = 2, jpjm1
451                     DO ji = 2, jpim1
452                        bhm(ji,jj,jk) = bhm0
453                     ENDDO
454                  ENDDO
455
456               ENDIF
457
458            ENDDO
459
460            !     Surface BC : set to zero
461            bhm(:,:,1)   = 0._wp 
462            !     Flat bottom BC : set to zero
463            bhm(:,:,jpk) = 0._wp 
464            !     Variable bathymetry case: diffusive fluxes masked in dyn_ldf_bgm
465
466            !--------------------------------
467         CASE(  13  )      !==  fixed profile: steady profile of the form bhm0 X dx^2 X dz^2 in interior, zero at bottom and top  ==!
468
469            cl_Units = ' 1/s'
470
471            IF(lwp) WRITE(numout,*) '   ==>>>   bi-Laplacian GM eddy viscosity : steady profile of the form'
472            IF(lwp) WRITE(numout,*) '           bhm0 X dx^2 X dz^2 in interior, zero boundary conds. at bottom and top.'
473            IF(lwp) WRITE(numout,*) '           In this case, bhm0 is the constant of proportionality,'
474            IF(lwp) WRITE(numout,*) '           dimensioned as [diffusivity]/L^4, i.e. bhm0=', bhm0, cl_Units
475
476            cl_Units = ' m4/s'
477
478            !     N.B. Test criterion for two grid-length mode for BGM scheme with leapfrog timestepping not implemented for this case yet.
479
480            DO jk = 1, jpk
481               ! constant at depth rn_bgmzcrit and larger
482               DO jj = 2, jpjm1
483                  DO ji = 2, jpim1
484                     bhm(ji,jj,jk) = bhm0*e1t(ji,jj)**2*e3w_n(ji,jj,jk)**2
485                  ENDDO
486               ENDDO
487            ENDDO
488
489            !     Surface BC : set to zero
490            bhm(:,:,1)   = 0._wp 
491            !     Flat bottom BC : set to zero
492            bhm(:,:,jpk) = 0._wp 
493            !     Variable bathymetry case: diffusive fluxes masked in dyn_ldf_bgm
494
495            !--------------------------------
496
497         CASE DEFAULT
498            CALL ctl_stop('ldf_dyn_init: wrong choice for nn_bhm_ijk_t, the type of space-time variation of bhm')
499         END SELECT
500         !
501         ! mask and lbc_lnk
502         bhm(:,:,1:jpkm1) =       bhm(:,:,1:jpkm1)   * wmask(:,:,1:jpkm1)
503         CALL lbc_lnk( 'ldf_dyn_init', bhm , 'W', 1. ) 
504         !
505      ENDIF ! ln_dynldf_bgm
506      !
507    END SUBROUTINE ldf_dyn_init
508
509
510    SUBROUTINE ldf_dyn( kt )
511      !!----------------------------------------------------------------------
512      !!                  ***  ROUTINE ldf_dyn  ***
513      !!
514      !! ** Purpose :   update at kt the momentum lateral mixing coeff. (ahmt and ahmf)
515      !!
516      !! ** Method  :   time varying eddy viscosity coefficients:
517      !!
518      !!    nn_ahm_ijk_t = 31    ahmt, ahmf = F(i,j,k,t) = F(local velocity)
519      !!                         ( |u|e /12  or  |u|e^3/12 for laplacian or bilaplacian operator )
520      !!
521      !!    nn_ahm_ijk_t = 32    ahmt, ahmf = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
522      !!                         ( L^2|D|    or  L^4|D|/8  for laplacian or bilaplacian operator )
523      !!
524      !! ** note    :    in BLP cases the sqrt of the eddy coef is returned, since bilaplacian is en re-entrant laplacian
525      !! ** action  :    ahmt, ahmf   updated at each time step
526      !!----------------------------------------------------------------------
527      INTEGER, INTENT(in) ::   kt   ! time step index
528      !
529      INTEGER  ::   ji, jj, jk   ! dummy loop indices
530      REAL(wp) ::   zu2pv2_ij_p1, zu2pv2_ij, zu2pv2_ij_m1, zemax   ! local scalar (option 31)
531      REAL(wp) ::   zcmsmag, zstabf_lo, zstabf_up, zdelta, zdb     ! local scalar (option 32)
532      !!----------------------------------------------------------------------
533      !
534      IF( ln_timing )   CALL timing_start('ldf_dyn')
535      !
536      SELECT CASE(  nn_ahm_ijk_t  )       !== Eddy vicosity coefficients ==!
537         !
538      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
539         !
540         IF( ln_dynldf_lap   ) THEN        ! laplacian operator : |u| e /12 = |u/144| e
541            DO jk = 1, jpkm1
542               DO jj = 2, jpjm1
543                  DO ji = fs_2, fs_jpim1
544                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
545                     zu2pv2_ij_m1 = ub(ji-1,jj  ,jk) * ub(ji-1,jj  ,jk) + vb(ji  ,jj-1,jk) * vb(ji  ,jj-1,jk)
546                     zemax = MAX( e1t(ji,jj) , e2t(ji,jj) )
547                     ahmt(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zemax * tmask(ji,jj,jk)      ! 288= 12*12 * 2
548                  END DO
549               END DO
550               DO jj = 1, jpjm1
551                  DO ji = 1, fs_jpim1
552                     zu2pv2_ij_p1 = ub(ji  ,jj+1,jk) * ub(ji  ,jj+1,jk) + vb(ji+1,jj  ,jk) * vb(ji+1,jj  ,jk)
553                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
554                     zemax = MAX( e1f(ji,jj) , e2f(ji,jj) )
555                     ahmf(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zemax * fmask(ji,jj,jk)      ! 288= 12*12 * 2
556                  END DO
557               END DO
558            END DO
559         ELSEIF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( |u| e^3 /12 ) = sqrt( |u/144| e ) * e
560            DO jk = 1, jpkm1
561               DO jj = 2, jpjm1
562                  DO ji = fs_2, fs_jpim1
563                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
564                     zu2pv2_ij_m1 = ub(ji-1,jj  ,jk) * ub(ji-1,jj  ,jk) + vb(ji  ,jj-1,jk) * vb(ji  ,jj-1,jk)
565                     zemax = MAX( e1t(ji,jj) , e2t(ji,jj) )
566                     ahmt(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zemax  ) * zemax * tmask(ji,jj,jk)
567                  END DO
568               END DO
569               DO jj = 1, jpjm1
570                  DO ji = 1, fs_jpim1
571                     zu2pv2_ij_p1 = ub(ji  ,jj+1,jk) * ub(ji  ,jj+1,jk) + vb(ji+1,jj  ,jk) * vb(ji+1,jj  ,jk)
572                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
573                     zemax = MAX( e1f(ji,jj) , e2f(ji,jj) )
574                     ahmf(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zemax  ) * zemax * fmask(ji,jj,jk)
575                  END DO
576               END DO
577            END DO
578         ENDIF
579         !
580         CALL lbc_lnk_multi( 'ldfdyn', ahmt, 'T', 1.,  ahmf, 'F', 1. )
581         !
582         !
583      CASE(  32  )       !==  time varying 3D field  ==!   = F( local deformation rate and gridscale ) (Smagorinsky)
584         !
585         IF( ln_dynldf_lap .OR. ln_dynldf_blp  ) THEN        ! laplacian operator : (C_smag/pi)^2 L^2 |D|
586            !
587            zcmsmag   = (rn_csmc/rpi)**2                                            ! (C_smag/pi)^2
588            zstabf_lo = rn_minfac * rn_minfac / ( 2._wp * 4._wp * zcmsmag )         ! lower limit stability factor scaling
589            zstabf_up = rn_maxfac / ( 4._wp * zcmsmag * 2._wp * rdt )               ! upper limit stability factor scaling
590            IF( ln_dynldf_blp ) zstabf_lo = ( 16._wp / 9._wp ) * zstabf_lo          ! provide |U|L^3/12 lower limit instead
591            !                                                                       ! of |U|L^3/16 in blp case
592            DO jk = 1, jpkm1
593               !
594               DO jj = 2, jpjm1
595                  DO ji = 2, jpim1
596                     zdb =   ( ub(ji,jj,jk) * r1_e2u(ji,jj) - ub(ji-1,jj,jk) * r1_e2u(ji-1,jj) ) * r1_e1t(ji,jj) * e2t(ji,jj)  &
597                          &  - ( vb(ji,jj,jk) * r1_e1v(ji,jj) - vb(ji,jj-1,jk) * r1_e1v(ji,jj-1) ) * r1_e2t(ji,jj) * e1t(ji,jj)
598                     dtensq(ji,jj,jk) = zdb * zdb * tmask(ji,jj,jk)
599                  END DO
600               END DO
601               !
602               DO jj = 1, jpjm1
603                  DO ji = 1, jpim1
604                     zdb =   ( ub(ji,jj+1,jk) * r1_e1u(ji,jj+1) - ub(ji,jj,jk) * r1_e1u(ji,jj) ) * r1_e2f(ji,jj) * e1f(ji,jj)  &
605                          &  + ( vb(ji+1,jj,jk) * r1_e2v(ji+1,jj) - vb(ji,jj,jk) * r1_e2v(ji,jj) ) * r1_e1f(ji,jj) * e2f(ji,jj)
606                     dshesq(ji,jj,jk) = zdb * zdb * fmask(ji,jj,jk)
607                  END DO
608               END DO
609               !
610            END DO
611            !
612            CALL lbc_lnk_multi( 'ldfdyn', dtensq, 'T', 1. )  ! lbc_lnk on dshesq not needed
613            !
614            DO jk = 1, jpkm1
615               !
616               DO jj = 2, jpjm1                                ! T-point value
617                  DO ji = fs_2, fs_jpim1
618                     !
619                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
620                     zu2pv2_ij_m1 = ub(ji-1,jj  ,jk) * ub(ji-1,jj  ,jk) + vb(ji  ,jj-1,jk) * vb(ji  ,jj-1,jk)
621                     !
622                     zdelta         = zcmsmag * esqt(ji,jj)                                        ! L^2 * (C_smag/pi)^2
623                     ahmt(ji,jj,jk) = zdelta * SQRT(          dtensq(ji  ,jj,jk) +                         &
624                          &                            r1_4 * ( dshesq(ji  ,jj,jk) + dshesq(ji  ,jj-1,jk) +  &
625                          &                                     dshesq(ji-1,jj,jk) + dshesq(ji-1,jj-1,jk) ) )
626                     ahmt(ji,jj,jk) = MAX( ahmt(ji,jj,jk), SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * zdelta * zstabf_lo ) ) ! Impose lower limit == minfac  * |U|L/2
627                     ahmt(ji,jj,jk) = MIN( ahmt(ji,jj,jk),                                    zdelta * zstabf_up )   ! Impose upper limit == maxfac  * L^2/(4*2dt)
628                     !
629                  END DO
630               END DO
631               !
632               DO jj = 1, jpjm1                                ! F-point value
633                  DO ji = 1, fs_jpim1
634                     !
635                     zu2pv2_ij_p1 = ub(ji  ,jj+1,jk) * ub(ji  ,jj+1,jk) + vb(ji+1,jj  ,jk) * vb(ji+1,jj  ,jk)
636                     zu2pv2_ij    = ub(ji  ,jj  ,jk) * ub(ji  ,jj  ,jk) + vb(ji  ,jj  ,jk) * vb(ji  ,jj  ,jk)
637                     !
638                     zdelta         = zcmsmag * esqf(ji,jj)                                        ! L^2 * (C_smag/pi)^2
639                     ahmf(ji,jj,jk) = zdelta * SQRT(          dshesq(ji  ,jj,jk) +                         &
640                          &                            r1_4 * ( dtensq(ji  ,jj,jk) + dtensq(ji  ,jj+1,jk) +  &
641                          &                                     dtensq(ji+1,jj,jk) + dtensq(ji+1,jj+1,jk) ) )
642                     ahmf(ji,jj,jk) = MAX( ahmf(ji,jj,jk), SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * zdelta * zstabf_lo ) ) ! Impose lower limit == minfac  * |U|L/2
643                     ahmf(ji,jj,jk) = MIN( ahmf(ji,jj,jk),                                    zdelta * zstabf_up )   ! Impose upper limit == maxfac  * L^2/(4*2dt)
644                     !
645                  END DO
646               END DO
647               !
648            END DO
649            !
650         ENDIF
651         !
652         IF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( (C_smag/pi)^2 L^4 |D|/8)
653            !                          !                      = sqrt( A_lap_smag L^2/8 )
654            !                          ! stability limits already applied to laplacian values
655            !                          ! effective default limits are 1/12 |U|L^3 < B_hm < 1//(32*2dt) L^4
656            DO jk = 1, jpkm1
657               DO jj = 2, jpjm1
658                  DO ji = fs_2, fs_jpim1
659                     ahmt(ji,jj,jk) = SQRT( r1_8 * esqt(ji,jj) * ahmt(ji,jj,jk) )
660                  END DO
661               END DO
662               DO jj = 1, jpjm1
663                  DO ji = 1, fs_jpim1
664                     ahmf(ji,jj,jk) = SQRT( r1_8 * esqf(ji,jj) * ahmf(ji,jj,jk) )
665                  END DO
666               END DO
667            END DO
668            !
669         ENDIF
670         !
671         CALL lbc_lnk_multi( 'ldfdyn', ahmt, 'T', 1. , ahmf, 'F', 1. )
672         !
673      END SELECT
674      !
675      CALL iom_put( "ahmt_2d", ahmt(:,:,1) )   ! surface u-eddy diffusivity coeff.
676      CALL iom_put( "ahmf_2d", ahmf(:,:,1) )   ! surface v-eddy diffusivity coeff.
677      CALL iom_put( "ahmt_3d", ahmt(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
678      CALL iom_put( "ahmf_3d", ahmf(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
679      !
680      IF( ln_timing )   CALL timing_stop('ldf_dyn')
681      !
682   END SUBROUTINE ldf_dyn
683
684   !!======================================================================
685END MODULE ldfdyn
Note: See TracBrowser for help on using the repository browser.