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/2019/dev_r11943_MERGE_2019/src/OCE/LDF – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/LDF/ldfdyn.F90 @ 12252

Last change on this file since 12252 was 12252, checked in by smasson, 4 years ago

rev12240_dev_r11943_MERGE_2019: same as [12251], merge trunk 12072:12248, all sette tests ok, GYRE_PISCES, AMM12, ISOMIP, VORTEX intentical to 12236

  • Property svn:keywords set to Id
File size: 31.2 KB
RevLine 
[3]1MODULE ldfdyn
2   !!======================================================================
3   !!                       ***  MODULE  ldfdyn  ***
4   !! Ocean physics:  lateral viscosity coefficient
5   !!=====================================================================
[1601]6   !! History :  OPA  ! 1997-07  (G. Madec)  multi dimensional coefficients
7   !!   NEMO     1.0  ! 2002-09  (G. Madec)  F90: Free form and module
[5836]8   !!            3.7  ! 2014-01  (F. Lemarie, G. Madec)  restructuration/simplification of ahm specification,
9   !!                 !                                  add velocity dependent coefficient and optional read in file
[1601]10   !!----------------------------------------------------------------------
[3]11
12   !!----------------------------------------------------------------------
[5836]13   !!   ldf_dyn_init  : initialization, namelist read, and parameters control
14   !!   ldf_dyn       : update lateral eddy viscosity coefficients at each time step
[3]15   !!----------------------------------------------------------------------
16   USE oce             ! ocean dynamics and tracers   
17   USE dom_oce         ! ocean space and time domain
18   USE phycst          ! physical constants
[9490]19   USE ldfslp          ! lateral diffusion: slopes of mixing orientation
[5836]20   USE ldfc1d_c2d      ! lateral diffusion: 1D and 2D cases
21   !
[3]22   USE in_out_manager  ! I/O manager
[5836]23   USE iom             ! I/O module for ehanced bottom friction file
24   USE timing          ! Timing
[3]25   USE lib_mpp         ! distribued memory computing library
26   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
27
28   IMPLICIT NONE
29   PRIVATE
30
[5836]31   PUBLIC   ldf_dyn_init   ! called by nemogcm.F90
32   PUBLIC   ldf_dyn        ! called by step.F90
[3]33
[9490]34   !                                    !!* Namelist namdyn_ldf : lateral mixing on momentum *
[9526]35   LOGICAL , PUBLIC ::   ln_dynldf_OFF   !: No operator (i.e. no explicit diffusion)
[5836]36   LOGICAL , PUBLIC ::   ln_dynldf_lap   !: laplacian operator
37   LOGICAL , PUBLIC ::   ln_dynldf_blp   !: bilaplacian operator
38   LOGICAL , PUBLIC ::   ln_dynldf_lev   !: iso-level direction
39   LOGICAL , PUBLIC ::   ln_dynldf_hor   !: horizontal (geopotential) direction
[9490]40!  LOGICAL , PUBLIC ::   ln_dynldf_iso   !: iso-neutral direction                        (see ldfslp)
[5836]41   INTEGER , PUBLIC ::   nn_ahm_ijk_t    !: choice of time & space variations of the lateral eddy viscosity coef.
[9490]42   !                                        !  time invariant coefficients:  aht = 1/2  Ud*Ld   (lap case)
43   !                                           !                             bht = 1/12 Ud*Ld^3 (blp case)
44   REAL(wp), PUBLIC ::   rn_Uv                 !: lateral viscous velocity  [m/s]
45   REAL(wp), PUBLIC ::   rn_Lv                 !: lateral viscous length    [m]
46   !                                        ! Smagorinsky viscosity  (nn_ahm_ijk_t = 32)
47   REAL(wp), PUBLIC ::   rn_csmc               !: Smagorinsky constant of proportionality
48   REAL(wp), PUBLIC ::   rn_minfac             !: Multiplicative factor of theorectical minimum Smagorinsky viscosity
49   REAL(wp), PUBLIC ::   rn_maxfac             !: Multiplicative factor of theorectical maximum Smagorinsky viscosity
50   !                                        ! iso-neutral laplacian (ln_dynldf_lap=ln_dynldf_iso=T)
51   REAL(wp), PUBLIC ::   rn_ahm_b              !: lateral laplacian background eddy viscosity  [m2/s]
[3]52
[9490]53   !                                    !!* Parameter to control the type of lateral viscous operator
54   INTEGER, PARAMETER, PUBLIC ::   np_ERROR  =-10                       !: error in setting the operator
55   INTEGER, PARAMETER, PUBLIC ::   np_no_ldf = 00                       !: without operator (i.e. no lateral viscous trend)
56   !                          !!      laplacian     !    bilaplacian    !
57   INTEGER, PARAMETER, PUBLIC ::   np_lap    = 10   ,   np_blp    = 20  !: iso-level operator
58   INTEGER, PARAMETER, PUBLIC ::   np_lap_i  = 11                       !: iso-neutral or geopotential operator
59   !
60   INTEGER           , PUBLIC ::   nldf_dyn         !: type of lateral diffusion used defined from ln_dynldf_... (namlist logicals)
61   LOGICAL           , PUBLIC ::   l_ldfdyn_time    !: flag for time variation of the lateral eddy viscosity coef.
[5836]62
[9490]63   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ahmt, ahmf   !: eddy viscosity coef. at T- and F-points [m2/s or m4/s]
[10784]64   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dtensq       !: horizontal tension squared         (Smagorinsky only)
65   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dshesq       !: horizontal shearing strain squared (Smagorinsky only)
[7646]66   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   esqt, esqf   !: Square of the local gridscale (e1e2/(e1+e2))**2           
[5836]67
[9490]68   REAL(wp) ::   r1_2    = 0.5_wp            ! =1/2
[5836]69   REAL(wp) ::   r1_4    = 0.25_wp           ! =1/4
[7646]70   REAL(wp) ::   r1_8    = 0.125_wp          ! =1/8
[9490]71   REAL(wp) ::   r1_12   = 1._wp / 12._wp    ! =1/12
[5836]72   REAL(wp) ::   r1_288  = 1._wp / 288._wp   ! =1/( 12^2 * 2 )
73
[3]74   !! * Substitutions
[5836]75#  include "vectopt_loop_substitute.h90"
[3]76   !!----------------------------------------------------------------------
[9598]77   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[1152]78   !! $Id$
[10068]79   !! Software governed by the CeCILL license (see ./LICENSE)
[3]80   !!----------------------------------------------------------------------
81CONTAINS
82
83   SUBROUTINE ldf_dyn_init
84      !!----------------------------------------------------------------------
85      !!                  ***  ROUTINE ldf_dyn_init  ***
86      !!                   
87      !! ** Purpose :   set the horizontal ocean dynamics physics
88      !!
[5836]89      !! ** Method  :   the eddy viscosity coef. specification depends on:
90      !!              - the operator:
91      !!             ln_dynldf_lap = T     laplacian operator
92      !!             ln_dynldf_blp = T   bilaplacian operator
93      !!              - the parameter nn_ahm_ijk_t:
94      !!    nn_ahm_ijk_t  =  0 => = constant
95      !!                  = 10 => = F(z) :     = constant with a reduction of 1/4 with depth
96      !!                  =-20 => = F(i,j)     = shape read in 'eddy_viscosity.nc' file
97      !!                  = 20    = F(i,j)     = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
98      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_viscosity.nc'  file
99      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
100      !!                  = 31    = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
101      !!                                                           or |u|e^3/12 bilaplacian operator )
[7646]102      !!                  = 32    = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
103      !!                                                           (   L^2|D|      laplacian operator
104      !!                                                           or  L^4|D|/8  bilaplacian operator )
[3]105      !!----------------------------------------------------------------------
[9490]106      INTEGER  ::   ji, jj, jk                     ! dummy loop indices
107      INTEGER  ::   ioptio, ierr, inum, ios, inn   ! local integer
108      REAL(wp) ::   zah0, zah_max, zUfac           ! local scalar
109      CHARACTER(len=5) ::   cl_Units               ! units (m2/s or m4/s)
110      !!
[9526]111      NAMELIST/namdyn_ldf/ ln_dynldf_OFF, ln_dynldf_lap, ln_dynldf_blp,   &   ! type of operator
112         &                 ln_dynldf_lev, ln_dynldf_hor, ln_dynldf_iso,   &   ! acting direction of the operator
113         &                 nn_ahm_ijk_t , rn_Uv    , rn_Lv,   rn_ahm_b,   &   ! lateral eddy coefficient
114         &                 rn_csmc      , rn_minfac    , rn_maxfac            ! Smagorinsky settings
[5836]115      !!----------------------------------------------------------------------
116      !
[4147]117      READ  ( numnam_ref, namdyn_ldf, IOSTAT = ios, ERR = 901)
[11536]118901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namdyn_ldf in reference namelist' )
[3]119
[4147]120      READ  ( numnam_cfg, namdyn_ldf, IOSTAT = ios, ERR = 902 )
[11536]121902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namdyn_ldf in configuration namelist' )
[4624]122      IF(lwm) WRITE ( numond, namdyn_ldf )
[4147]123
[1601]124      IF(lwp) THEN                      ! Parameter print
[3]125         WRITE(numout,*)
126         WRITE(numout,*) 'ldf_dyn : lateral momentum physics'
127         WRITE(numout,*) '~~~~~~~'
[4147]128         WRITE(numout,*) '   Namelist namdyn_ldf : set lateral mixing parameters'
[5836]129         !
130         WRITE(numout,*) '      type :'
[9526]131         WRITE(numout,*) '         no explicit diffusion                ln_dynldf_OFF = ', ln_dynldf_OFF
[5836]132         WRITE(numout,*) '         laplacian operator                   ln_dynldf_lap = ', ln_dynldf_lap
133         WRITE(numout,*) '         bilaplacian operator                 ln_dynldf_blp = ', ln_dynldf_blp
134         !
135         WRITE(numout,*) '      direction of action :'
136         WRITE(numout,*) '         iso-level                            ln_dynldf_lev = ', ln_dynldf_lev
137         WRITE(numout,*) '         horizontal (geopotential)            ln_dynldf_hor = ', ln_dynldf_hor
138         WRITE(numout,*) '         iso-neutral                          ln_dynldf_iso = ', ln_dynldf_iso
139         !
140         WRITE(numout,*) '      coefficients :'
141         WRITE(numout,*) '         type of time-space variation         nn_ahm_ijk_t  = ', nn_ahm_ijk_t
[9490]142         WRITE(numout,*) '         lateral viscous velocity  (if cst)      rn_Uv      = ', rn_Uv, ' m/s'
143         WRITE(numout,*) '         lateral viscous length    (if cst)      rn_Lv      = ', rn_Lv, ' m'
144         WRITE(numout,*) '         background viscosity (iso-lap case)     rn_ahm_b   = ', rn_ahm_b, ' m2/s'
145         !
[9019]146         WRITE(numout,*) '      Smagorinsky settings (nn_ahm_ijk_t  = 32) :'
[7646]147         WRITE(numout,*) '         Smagorinsky coefficient              rn_csmc       = ', rn_csmc
[9490]148         WRITE(numout,*) '         factor multiplier for eddy visc.'
149         WRITE(numout,*) '            lower limit (default 1.0)         rn_minfac    = ', rn_minfac
150         WRITE(numout,*) '            upper limit (default 1.0)         rn_maxfac    = ', rn_maxfac
[3]151      ENDIF
152
[9490]153      !
154      !           !==  type of lateral operator used  ==!   (set nldf_dyn)
155      !           !=====================================!
156      !
157      nldf_dyn = np_ERROR
158      ioptio = 0
[9526]159      IF( ln_dynldf_OFF ) THEN   ;   nldf_dyn = np_no_ldf   ;   ioptio = ioptio + 1   ;   ENDIF
160      IF( ln_dynldf_lap ) THEN   ;                              ioptio = ioptio + 1   ;   ENDIF
161      IF( ln_dynldf_blp ) THEN   ;                              ioptio = ioptio + 1   ;   ENDIF
162      IF( ioptio /= 1   )   CALL ctl_stop( 'dyn_ldf_init: use ONE of the 3 operator options (NONE/lap/blp)' )
[9490]163      !
[9526]164      IF(.NOT.ln_dynldf_OFF ) THEN     !==  direction ==>> type of operator  ==!
[9490]165         ioptio = 0
166         IF( ln_dynldf_lev )   ioptio = ioptio + 1
167         IF( ln_dynldf_hor )   ioptio = ioptio + 1
168         IF( ln_dynldf_iso )   ioptio = ioptio + 1
169         IF( ioptio /= 1   )   CALL ctl_stop( 'dyn_ldf_init: use ONE of the 3 direction options (level/hor/iso)' )
170         !
171         !                             ! Set nldf_dyn, the type of lateral diffusion, from ln_dynldf_... logicals
172         ierr = 0
173         IF( ln_dynldf_lap ) THEN         ! laplacian operator
174            IF( ln_zco ) THEN                   ! z-coordinate
175               IF ( ln_dynldf_lev )   nldf_dyn = np_lap     ! iso-level = horizontal (no rotation)
176               IF ( ln_dynldf_hor )   nldf_dyn = np_lap     ! iso-level = horizontal (no rotation)
177               IF ( ln_dynldf_iso )   nldf_dyn = np_lap_i   ! iso-neutral            (   rotation)
178            ENDIF
179            IF( ln_zps ) THEN                   ! z-coordinate with partial step
180               IF ( ln_dynldf_lev )   nldf_dyn = np_lap     ! iso-level              (no rotation)
181               IF ( ln_dynldf_hor )   nldf_dyn = np_lap     ! iso-level              (no rotation)
182               IF ( ln_dynldf_iso )   nldf_dyn = np_lap_i   ! iso-neutral            (   rotation)
183            ENDIF
184            IF( ln_sco ) THEN                   ! s-coordinate
185               IF ( ln_dynldf_lev )   nldf_dyn = np_lap     ! iso-level = horizontal (no rotation)
186               IF ( ln_dynldf_hor )   nldf_dyn = np_lap_i   ! horizontal             (   rotation)
187               IF ( ln_dynldf_iso )   nldf_dyn = np_lap_i   ! iso-neutral            (   rotation)
188            ENDIF
189         ENDIF
190         !
191         IF( ln_dynldf_blp ) THEN         ! bilaplacian operator
192            IF( ln_zco ) THEN                   ! z-coordinate
193               IF( ln_dynldf_lev )   nldf_dyn = np_blp   ! iso-level = horizontal (no rotation)
194               IF( ln_dynldf_hor )   nldf_dyn = np_blp   ! iso-level = horizontal (no rotation)
195               IF( ln_dynldf_iso )   ierr = 2            ! iso-neutral            (   rotation)
196            ENDIF
197            IF( ln_zps ) THEN                   ! z-coordinate with partial step
198               IF( ln_dynldf_lev )   nldf_dyn = np_blp   ! iso-level              (no rotation)
199               IF( ln_dynldf_hor )   nldf_dyn = np_blp   ! iso-level              (no rotation)
200               IF( ln_dynldf_iso )   ierr = 2            ! iso-neutral            (   rotation)
201            ENDIF
202            IF( ln_sco ) THEN                   ! s-coordinate
203               IF( ln_dynldf_lev )   nldf_dyn = np_blp   ! iso-level              (no rotation)
204               IF( ln_dynldf_hor )   ierr = 2            ! horizontal             (   rotation)
205               IF( ln_dynldf_iso )   ierr = 2            ! iso-neutral            (   rotation)
206            ENDIF
207         ENDIF
208         !
209         IF( ierr == 2 )   CALL ctl_stop( 'rotated bi-laplacian operator does not exist' )
210         !
211         IF( nldf_dyn == np_lap_i )   l_ldfslp = .TRUE.  ! rotation require the computation of the slopes
212         !
[3]213      ENDIF
[5836]214      !
[9490]215      IF(lwp) THEN
216         WRITE(numout,*)
217         SELECT CASE( nldf_dyn )
218         CASE( np_no_ldf )   ;   WRITE(numout,*) '   ==>>>   NO lateral viscosity'
219         CASE( np_lap    )   ;   WRITE(numout,*) '   ==>>>   iso-level laplacian operator'
220         CASE( np_lap_i  )   ;   WRITE(numout,*) '   ==>>>   rotated laplacian operator with iso-level background'
221         CASE( np_blp    )   ;   WRITE(numout,*) '   ==>>>   iso-level bi-laplacian operator'
222         END SELECT
223         WRITE(numout,*)
[3]224      ENDIF
[9490]225     
[1601]226      !
[9490]227      !           !==  Space/time variation of eddy coefficients  ==!
228      !           !=================================================!
[5836]229      !
[9490]230      l_ldfdyn_time = .FALSE.                ! no time variation except in case defined below
[1601]231      !
[9526]232      IF( ln_dynldf_OFF ) THEN
[9490]233         IF(lwp) WRITE(numout,*) '   ==>>>   No viscous operator selected. ahmt and ahmf are not allocated'
234         RETURN
[5836]235         !
[9490]236      ELSE                                   !==  a lateral diffusion operator is used  ==!
[5836]237         !
[9490]238         !                                         ! allocate the ahm arrays
239         ALLOCATE( ahmt(jpi,jpj,jpk) , ahmf(jpi,jpj,jpk) , STAT=ierr )
240         IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate arrays')
241         !
[10784]242         ahmt(:,:,:) = 0._wp                       ! init to 0 needed
243         ahmf(:,:,:) = 0._wp
[9490]244         !
245         !                                         ! value of lap/blp eddy mixing coef.
246         IF(     ln_dynldf_lap ) THEN   ;   zUfac = r1_2 *rn_Uv   ;   inn = 1   ;   cl_Units = ' m2/s'   !   laplacian
247         ELSEIF( ln_dynldf_blp ) THEN   ;   zUfac = r1_12*rn_Uv   ;   inn = 3   ;   cl_Units = ' m4/s'   ! bilaplacian
248         ENDIF
249         zah0    = zUfac *    rn_Lv**inn              ! mixing coefficient
250         zah_max = zUfac * (ra*rad)**inn              ! maximum reachable coefficient (value at the Equator)
251         !
252         SELECT CASE(  nn_ahm_ijk_t  )             !* Specification of space-time variations of ahmt, ahmf
253         !
[5836]254         CASE(   0  )      !==  constant  ==!
[9490]255            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity. = constant = ', zah0, cl_Units
256            ahmt(:,:,1:jpkm1) = zah0
257            ahmf(:,:,1:jpkm1) = zah0
[5836]258            !
259         CASE(  10  )      !==  fixed profile  ==!
[9490]260            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( depth )'
261            IF(lwp) WRITE(numout,*) '           surface viscous coef. = constant = ', zah0, cl_Units
262            ahmt(:,:,1) = zah0                        ! constant surface value
263            ahmf(:,:,1) = zah0
264            CALL ldf_c1d( 'DYN', ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )
[5836]265            !
266         CASE ( -20 )      !== fixed horizontal shape read in file  ==!
[9490]267            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F(i,j) read in eddy_viscosity.nc file'
[5836]268            CALL iom_open( 'eddy_viscosity_2D.nc', inum )
269            CALL iom_get ( inum, jpdom_data, 'ahmt_2d', ahmt(:,:,1) )
270            CALL iom_get ( inum, jpdom_data, 'ahmf_2d', ahmf(:,:,1) )
271            CALL iom_close( inum )
272            DO jk = 2, jpkm1
[9490]273               ahmt(:,:,jk) = ahmt(:,:,1)
274               ahmf(:,:,jk) = ahmf(:,:,1)
[5836]275            END DO
276            !
277         CASE(  20  )      !== fixed horizontal shape  ==!
[9490]278            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( e1, e2 ) or F( e1^3, e2^3 ) (lap. or blp. case)'
279            IF(lwp) WRITE(numout,*) '           using a fixed viscous velocity = ', rn_Uv  ,' m/s   and   Lv = Max(e1,e2)'
280            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
281            CALL ldf_c2d( 'DYN', zUfac      , inn        , ahmt, ahmf )         ! surface value proportional to scale factor^inn
[5836]282            !
283         CASE( -30  )      !== fixed 3D shape read in file  ==!
[9490]284            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F(i,j,k) read in eddy_viscosity_3D.nc file'
[5836]285            CALL iom_open( 'eddy_viscosity_3D.nc', inum )
286            CALL iom_get ( inum, jpdom_data, 'ahmt_3d', ahmt )
287            CALL iom_get ( inum, jpdom_data, 'ahmf_3d', ahmf )
288            CALL iom_close( inum )
289            !
290         CASE(  30  )       !==  fixed 3D shape  ==!
[9490]291            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth )'
292            IF(lwp) WRITE(numout,*) '           using a fixed viscous velocity = ', rn_Uv  ,' m/s   and   Ld = Max(e1,e2)'
293            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
294            CALL ldf_c2d( 'DYN', zUfac      , inn        , ahmt, ahmf )         ! surface value proportional to scale factor^inn
295            CALL ldf_c1d( 'DYN', ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )  ! reduction with depth
[5836]296            !
297         CASE(  31  )       !==  time varying 3D field  ==!
[9490]298            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth , time )'
299            IF(lwp) WRITE(numout,*) '           proportional to the local velocity : 1/2 |u|e (lap) or 1/12 |u|e^3 (blp)'
[5836]300            !
301            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
302            !
[7646]303         CASE(  32  )       !==  time varying 3D field  ==!
[9490]304            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth , time )'
305            IF(lwp) WRITE(numout,*) '           proportional to the local deformation rate and gridscale (Smagorinsky)'
[7646]306            !
307            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
308            !
[9490]309            !                          ! allocate arrays used in ldf_dyn.
[10784]310            ALLOCATE( dtensq(jpi,jpj,jpk) , dshesq(jpi,jpj,jpk) , esqt(jpi,jpj) , esqf(jpi,jpj) , STAT=ierr )
[7646]311            IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate Smagorinsky arrays')
312            !
[10784]313            DO jj = 1, jpj             ! Set local gridscale values
314               DO ji = 1, jpi
[11653]315                  esqt(ji,jj) = ( 2._wp * e1e2t(ji,jj) / ( e1t(ji,jj) + e2t(ji,jj) ) )**2 
316                  esqf(ji,jj) = ( 2._wp * e1e2f(ji,jj) / ( e1f(ji,jj) + e2f(ji,jj) ) )**2 
[7646]317               END DO
318            END DO
319            !
[5836]320         CASE DEFAULT
321            CALL ctl_stop('ldf_dyn_init: wrong choice for nn_ahm_ijk_t, the type of space-time variation of ahm')
322         END SELECT
323         !
[9490]324         IF( .NOT.l_ldfdyn_time ) THEN             !* No time variation
325            IF(     ln_dynldf_lap ) THEN                 !   laplacian operator (mask only)
326               ahmt(:,:,1:jpkm1) =       ahmt(:,:,1:jpkm1)   * tmask(:,:,1:jpkm1)
327               ahmf(:,:,1:jpkm1) =       ahmf(:,:,1:jpkm1)   * fmask(:,:,1:jpkm1)
328            ELSEIF( ln_dynldf_blp ) THEN                 ! bilaplacian operator (square root + mask)
329               ahmt(:,:,1:jpkm1) = SQRT( ahmt(:,:,1:jpkm1) ) * tmask(:,:,1:jpkm1)
330               ahmf(:,:,1:jpkm1) = SQRT( ahmf(:,:,1:jpkm1) ) * fmask(:,:,1:jpkm1)
331            ENDIF
[5836]332         ENDIF
333         !
[71]334      ENDIF
[1601]335      !
[5836]336   END SUBROUTINE ldf_dyn_init
[71]337
338
[11949]339   SUBROUTINE ldf_dyn( kt, Kbb )
[3]340      !!----------------------------------------------------------------------
[5836]341      !!                  ***  ROUTINE ldf_dyn  ***
342      !!
343      !! ** Purpose :   update at kt the momentum lateral mixing coeff. (ahmt and ahmf)
[3]344      !!
[5836]345      !! ** Method  :   time varying eddy viscosity coefficients:
[3]346      !!
[5836]347      !!    nn_ahm_ijk_t = 31    ahmt, ahmf = F(i,j,k,t) = F(local velocity)
348      !!                         ( |u|e /12  or  |u|e^3/12 for laplacian or bilaplacian operator )
[3]349      !!
[7646]350      !!    nn_ahm_ijk_t = 32    ahmt, ahmf = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
351      !!                         ( L^2|D|    or  L^4|D|/8  for laplacian or bilaplacian operator )
352      !!
353      !! ** note    :    in BLP cases the sqrt of the eddy coef is returned, since bilaplacian is en re-entrant laplacian
354      !! ** action  :    ahmt, ahmf   updated at each time step
[3]355      !!----------------------------------------------------------------------
[5836]356      INTEGER, INTENT(in) ::   kt   ! time step index
[11949]357      INTEGER, INTENT(in) ::   Kbb  ! ocean time level indices
[1601]358      !
[5836]359      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[10784]360      REAL(wp) ::   zu2pv2_ij_p1, zu2pv2_ij, zu2pv2_ij_m1, zemax   ! local scalar (option 31)
361      REAL(wp) ::   zcmsmag, zstabf_lo, zstabf_up, zdelta, zdb     ! local scalar (option 32)
[5836]362      !!----------------------------------------------------------------------
363      !
[9019]364      IF( ln_timing )   CALL timing_start('ldf_dyn')
[5836]365      !
366      SELECT CASE(  nn_ahm_ijk_t  )       !== Eddy vicosity coefficients ==!
367      !
368      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
369         !
[7646]370         IF( ln_dynldf_lap   ) THEN        ! laplacian operator : |u| e /12 = |u/144| e
[5836]371            DO jk = 1, jpkm1
372               DO jj = 2, jpjm1
373                  DO ji = fs_2, fs_jpim1
[11949]374                     zu2pv2_ij    = uu(ji  ,jj  ,jk,Kbb) * uu(ji  ,jj  ,jk,Kbb) + vv(ji  ,jj  ,jk,Kbb) * vv(ji  ,jj  ,jk,Kbb)
375                     zu2pv2_ij_m1 = uu(ji-1,jj  ,jk,Kbb) * uu(ji-1,jj  ,jk,Kbb) + vv(ji  ,jj-1,jk,Kbb) * vv(ji  ,jj-1,jk,Kbb)
[10784]376                     zemax = MAX( e1t(ji,jj) , e2t(ji,jj) )
377                     ahmt(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zemax * tmask(ji,jj,jk)      ! 288= 12*12 * 2
[5836]378                  END DO
379               END DO
[10784]380               DO jj = 1, jpjm1
381                  DO ji = 1, fs_jpim1
[11949]382                     zu2pv2_ij_p1 = uu(ji  ,jj+1,jk, Kbb) * uu(ji  ,jj+1,jk, Kbb) + vv(ji+1,jj  ,jk, Kbb) * vv(ji+1,jj  ,jk, Kbb)
383                     zu2pv2_ij    = uu(ji  ,jj  ,jk, Kbb) * uu(ji  ,jj  ,jk, Kbb) + vv(ji  ,jj  ,jk, Kbb) * vv(ji  ,jj  ,jk, Kbb)
[10784]384                     zemax = MAX( e1f(ji,jj) , e2f(ji,jj) )
385                     ahmf(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zemax * fmask(ji,jj,jk)      ! 288= 12*12 * 2
386                  END DO
387               END DO
[5836]388            END DO
389         ELSEIF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( |u| e^3 /12 ) = sqrt( |u/144| e ) * e
390            DO jk = 1, jpkm1
391               DO jj = 2, jpjm1
392                  DO ji = fs_2, fs_jpim1
[11949]393                     zu2pv2_ij    = uu(ji  ,jj  ,jk,Kbb) * uu(ji  ,jj  ,jk,Kbb) + vv(ji  ,jj  ,jk,Kbb) * vv(ji  ,jj  ,jk,Kbb)
394                     zu2pv2_ij_m1 = uu(ji-1,jj  ,jk,Kbb) * uu(ji-1,jj  ,jk,Kbb) + vv(ji  ,jj-1,jk,Kbb) * vv(ji  ,jj-1,jk,Kbb)
[10784]395                     zemax = MAX( e1t(ji,jj) , e2t(ji,jj) )
396                     ahmt(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zemax  ) * zemax * tmask(ji,jj,jk)
[5836]397                  END DO
398               END DO
[10784]399               DO jj = 1, jpjm1
400                  DO ji = 1, fs_jpim1
[11949]401                     zu2pv2_ij_p1 = uu(ji  ,jj+1,jk, Kbb) * uu(ji  ,jj+1,jk, Kbb) + vv(ji+1,jj  ,jk, Kbb) * vv(ji+1,jj  ,jk, Kbb)
402                     zu2pv2_ij    = uu(ji  ,jj  ,jk, Kbb) * uu(ji  ,jj  ,jk, Kbb) + vv(ji  ,jj  ,jk, Kbb) * vv(ji  ,jj  ,jk, Kbb)
[10784]403                     zemax = MAX( e1f(ji,jj) , e2f(ji,jj) )
404                     ahmf(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zemax  ) * zemax * fmask(ji,jj,jk)
405                  END DO
406               END DO
[5836]407            END DO
408         ENDIF
409         !
[10425]410         CALL lbc_lnk_multi( 'ldfdyn', ahmt, 'T', 1.,  ahmf, 'F', 1. )
[5836]411         !
[7646]412         !
413      CASE(  32  )       !==  time varying 3D field  ==!   = F( local deformation rate and gridscale ) (Smagorinsky)
414         !
415         IF( ln_dynldf_lap .OR. ln_dynldf_blp  ) THEN        ! laplacian operator : (C_smag/pi)^2 L^2 |D|
416            !
[10784]417            zcmsmag   = (rn_csmc/rpi)**2                                            ! (C_smag/pi)^2
[12252]418            zstabf_lo = rn_minfac * rn_minfac / ( 2._wp * 12._wp * 12._wp * zcmsmag ) ! lower limit stability factor scaling
[10784]419            zstabf_up = rn_maxfac / ( 4._wp * zcmsmag * 2._wp * rdt )               ! upper limit stability factor scaling
[7646]420            IF( ln_dynldf_blp ) zstabf_lo = ( 16._wp / 9._wp ) * zstabf_lo          ! provide |U|L^3/12 lower limit instead
421            !                                                                       ! of |U|L^3/16 in blp case
422            DO jk = 1, jpkm1
423               !
[10784]424               DO jj = 2, jpjm1
425                  DO ji = 2, jpim1
[11949]426                     zdb =    ( uu(ji,jj,jk,Kbb) * r1_e2u(ji,jj) -  uu(ji-1,jj,jk,Kbb) * r1_e2u(ji-1,jj) )  &
427                          &                      * r1_e1t(ji,jj) * e2t(ji,jj)                           &
428                          & - ( vv(ji,jj,jk,Kbb) * r1_e1v(ji,jj) -  vv(ji,jj-1,jk,Kbb) * r1_e1v(ji,jj-1) )  &
429                          &                      * r1_e2t(ji,jj) * e1t(ji,jj)
[10784]430                     dtensq(ji,jj,jk) = zdb * zdb * tmask(ji,jj,jk)
[7646]431                  END DO
432               END DO
433               !
434               DO jj = 1, jpjm1
435                  DO ji = 1, jpim1
[11949]436                     zdb =   (  uu(ji,jj+1,jk,Kbb) * r1_e1u(ji,jj+1) -  uu(ji,jj,jk,Kbb) * r1_e1u(ji,jj) )  &
437                          &                        * r1_e2f(ji,jj)   * e1f(ji,jj)                       &
438                          & + ( vv(ji+1,jj,jk,Kbb) * r1_e2v(ji+1,jj) -  vv(ji,jj,jk,Kbb) * r1_e2v(ji,jj) )  &
439                          &                        * r1_e1f(ji,jj)   * e2f(ji,jj)
[10784]440                     dshesq(ji,jj,jk) = zdb * zdb * fmask(ji,jj,jk)
[7646]441                  END DO
442               END DO
443               !
[10784]444            END DO
445            !
446            CALL lbc_lnk_multi( 'ldfdyn', dtensq, 'T', 1. )  ! lbc_lnk on dshesq not needed
447            !
448            DO jk = 1, jpkm1
449              !
450               DO jj = 2, jpjm1                                ! T-point value
[7646]451                  DO ji = fs_2, fs_jpim1
452                     !
[11949]453                     zu2pv2_ij    = uu(ji  ,jj  ,jk,Kbb) * uu(ji  ,jj  ,jk,Kbb) + vv(ji  ,jj  ,jk,Kbb) * vv(ji  ,jj  ,jk,Kbb)
454                     zu2pv2_ij_m1 = uu(ji-1,jj  ,jk,Kbb) * uu(ji-1,jj  ,jk,Kbb) + vv(ji  ,jj-1,jk,Kbb) * vv(ji  ,jj-1,jk,Kbb)
[10784]455                     !
[7646]456                     zdelta         = zcmsmag * esqt(ji,jj)                                        ! L^2 * (C_smag/pi)^2
[10784]457                     ahmt(ji,jj,jk) = zdelta * SQRT(          dtensq(ji  ,jj,jk) +                         &
458                        &                            r1_4 * ( dshesq(ji  ,jj,jk) + dshesq(ji  ,jj-1,jk) +  &
459                        &                                     dshesq(ji-1,jj,jk) + dshesq(ji-1,jj-1,jk) ) )
460                     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
461                     ahmt(ji,jj,jk) = MIN( ahmt(ji,jj,jk),                                    zdelta * zstabf_up )   ! Impose upper limit == maxfac  * L^2/(4*2dt)
462                     !
463                  END DO
464               END DO
465               !
466               DO jj = 1, jpjm1                                ! F-point value
467                  DO ji = 1, fs_jpim1
468                     !
[11949]469                     zu2pv2_ij_p1 = uu(ji  ,jj+1,jk, kbb) * uu(ji  ,jj+1,jk, kbb) + vv(ji+1,jj  ,jk, kbb) * vv(ji+1,jj  ,jk, kbb)
470                     zu2pv2_ij    = uu(ji  ,jj  ,jk, kbb) * uu(ji  ,jj  ,jk, kbb) + vv(ji  ,jj  ,jk, kbb) * vv(ji  ,jj  ,jk, kbb)
[10784]471                     !
[7646]472                     zdelta         = zcmsmag * esqf(ji,jj)                                        ! L^2 * (C_smag/pi)^2
[10784]473                     ahmf(ji,jj,jk) = zdelta * SQRT(          dshesq(ji  ,jj,jk) +                         &
474                        &                            r1_4 * ( dtensq(ji  ,jj,jk) + dtensq(ji  ,jj+1,jk) +  &
475                        &                                     dtensq(ji+1,jj,jk) + dtensq(ji+1,jj+1,jk) ) )
476                     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
477                     ahmf(ji,jj,jk) = MIN( ahmf(ji,jj,jk),                                    zdelta * zstabf_up )   ! Impose upper limit == maxfac  * L^2/(4*2dt)
[7646]478                     !
479                  END DO
480               END DO
[10784]481               !
[7646]482            END DO
483            !
484         ENDIF
485         !
486         IF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( (C_smag/pi)^2 L^4 |D|/8)
[9490]487            !                          !                      = sqrt( A_lap_smag L^2/8 )
488            !                          ! stability limits already applied to laplacian values
489            !                          ! effective default limits are 1/12 |U|L^3 < B_hm < 1//(32*2dt) L^4
[7646]490            DO jk = 1, jpkm1
491               DO jj = 2, jpjm1
492                  DO ji = fs_2, fs_jpim1
[9490]493                     ahmt(ji,jj,jk) = SQRT( r1_8 * esqt(ji,jj) * ahmt(ji,jj,jk) )
[10784]494                  END DO
495               END DO
496               DO jj = 1, jpjm1
497                  DO ji = 1, fs_jpim1
[9490]498                     ahmf(ji,jj,jk) = SQRT( r1_8 * esqf(ji,jj) * ahmf(ji,jj,jk) )
[7646]499                  END DO
500               END DO
501            END DO
502            !
503         ENDIF
504         !
[10425]505         CALL lbc_lnk_multi( 'ldfdyn', ahmt, 'T', 1. , ahmf, 'F', 1. )
[7646]506         !
[5836]507      END SELECT
508      !
509      CALL iom_put( "ahmt_2d", ahmt(:,:,1) )   ! surface u-eddy diffusivity coeff.
510      CALL iom_put( "ahmf_2d", ahmf(:,:,1) )   ! surface v-eddy diffusivity coeff.
511      CALL iom_put( "ahmt_3d", ahmt(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
512      CALL iom_put( "ahmf_3d", ahmf(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
513      !
[9019]514      IF( ln_timing )   CALL timing_stop('ldf_dyn')
[5836]515      !
516   END SUBROUTINE ldf_dyn
[3]517
518   !!======================================================================
519END MODULE ldfdyn
Note: See TracBrowser for help on using the repository browser.