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

source: NEMO/branches/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/LDF/ldfdyn.F90 @ 10986

Last change on this file since 10986 was 10986, checked in by andmirek, 5 years ago

GMED 462 add flush

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