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 @ 12788

Last change on this file since 12788 was 12788, checked in by davestorkey, 4 years ago

UKMO/NEMO_4.0.1_biharmonic_GM : New calculation for zmu with new rn_bgm_msc namelist parameter.

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