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/2020/dev_r12527_Gurvan_ShallowWater/src/SWE – NEMO

source: NEMO/branches/2020/dev_r12527_Gurvan_ShallowWater/src/SWE/ldfdyn.F90 @ 14071

Last change on this file since 14071 was 13416, checked in by gm, 4 years ago

First dev of BVP, cond 3.4, gridF fixes

  • Property svn:keywords set to Id
File size: 32.1 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   USE usrdef_nam , ONLY : ln_dynldf_lap_PM
29   !
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC   ldf_dyn_init   ! called by nemogcm.F90
34   PUBLIC   ldf_dyn        ! called by step.F90
35
36   !                                    !!* Namelist namdyn_ldf : lateral mixing on momentum *
37   LOGICAL , PUBLIC ::   ln_dynldf_OFF   !: No operator (i.e. no explicit diffusion)
38   LOGICAL , PUBLIC ::   ln_dynldf_lap   !: laplacian operator
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   REAL(wp), PUBLIC ::   rn_Uv                 !: lateral viscous velocity  [m/s]
47   REAL(wp), PUBLIC ::   rn_Lv                 !: lateral viscous length    [m]
48   !                                        ! Smagorinsky viscosity  (nn_ahm_ijk_t = 32)
49   REAL(wp), PUBLIC ::   rn_csmc               !: Smagorinsky constant of proportionality
50   REAL(wp), PUBLIC ::   rn_minfac             !: Multiplicative factor of theorectical minimum Smagorinsky viscosity
51   REAL(wp), PUBLIC ::   rn_maxfac             !: Multiplicative factor of theorectical maximum Smagorinsky viscosity
52   !                                        ! iso-neutral laplacian (ln_dynldf_lap=ln_dynldf_iso=T)
53   REAL(wp), PUBLIC ::   rn_ahm_b              !: lateral laplacian background eddy viscosity  [m2/s]
54
55   !                                    !!* Parameter to control the type of lateral viscous operator
56   INTEGER, PARAMETER, PUBLIC ::   np_ERROR  =-10                       !: error in setting the operator
57   INTEGER, PARAMETER, PUBLIC ::   np_no_ldf = 00                       !: without operator (i.e. no lateral viscous trend)
58   !                          !!      laplacian     !    bilaplacian    !
59   INTEGER, PARAMETER, PUBLIC ::   np_lap    = 10   ,   np_blp    = 20  !: iso-level operator
60   INTEGER, PARAMETER, PUBLIC ::   np_lap_i  = 11                       !: iso-neutral or geopotential operator
61   !
62   INTEGER           , PUBLIC ::   nldf_dyn         !: type of lateral diffusion used defined from ln_dynldf_... (namlist logicals)
63   LOGICAL           , PUBLIC ::   l_ldfdyn_time    !: flag for time variation of the lateral eddy viscosity coef.
64
65   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ahmt, ahmf   !: eddy viscosity coef. at T- and F-points [m2/s or m4/s]
66   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dtensq       !: horizontal tension squared         (Smagorinsky only)
67   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   dshesq       !: horizontal shearing strain squared (Smagorinsky only)
68   REAL(wp),         ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   esqt, esqf   !: Square of the local gridscale (e1e2/(e1+e2))**2           
69
70   REAL(wp) ::   r1_2    = 0.5_wp            ! =1/2
71   REAL(wp) ::   r1_4    = 0.25_wp           ! =1/4
72   REAL(wp) ::   r1_8    = 0.125_wp          ! =1/8
73   REAL(wp) ::   r1_12   = 1._wp / 12._wp    ! =1/12
74   REAL(wp) ::   r1_288  = 1._wp / 288._wp   ! =1/( 12^2 * 2 )
75
76   !! * Substitutions
77#  include "do_loop_substitute.h90"
78   !!----------------------------------------------------------------------
79   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
80   !! $Id$
81   !! Software governed by the CeCILL license (see ./LICENSE)
82   !!----------------------------------------------------------------------
83CONTAINS
84
85   SUBROUTINE ldf_dyn_init
86      !!----------------------------------------------------------------------
87      !!                  ***  ROUTINE ldf_dyn_init  ***
88      !!                   
89      !! ** Purpose :   set the horizontal ocean dynamics physics
90      !!
91      !! ** Method  :   the eddy viscosity coef. specification depends on:
92      !!              - the operator:
93      !!             ln_dynldf_lap = T     laplacian operator
94      !!             ln_dynldf_blp = T   bilaplacian operator
95      !!              - the parameter nn_ahm_ijk_t:
96      !!    nn_ahm_ijk_t  =  0 => = constant
97      !!                  = 10 => = F(z) :     = constant with a reduction of 1/4 with depth
98      !!                  =-20 => = F(i,j)     = shape read in 'eddy_viscosity.nc' file
99      !!                  = 20    = F(i,j)     = F(e1,e2) or F(e1^3,e2^3) (lap or bilap case)
100      !!                  =-30 => = F(i,j,k)   = shape read in 'eddy_viscosity.nc'  file
101      !!                  = 30    = F(i,j,k)   = 2D (case 20) + decrease with depth (case 10)
102      !!                  = 31    = F(i,j,k,t) = F(local velocity) (  |u|e  /12   laplacian operator
103      !!                                                           or |u|e^3/12 bilaplacian operator )
104      !!                  = 32    = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
105      !!                                                           (   L^2|D|      laplacian operator
106      !!                                                           or  L^4|D|/8  bilaplacian operator )
107      !!----------------------------------------------------------------------
108      INTEGER  ::   ji, jj, jk                     ! dummy loop indices
109      INTEGER  ::   ioptio, ierr, inum, ios, inn   ! local integer
110      REAL(wp) ::   zah0, zah_max, zUfac           ! local scalar
111      REAL(wp) ::   zsum                           ! local scalar
112      CHARACTER(len=5) ::   cl_Units               ! units (m2/s or m4/s)
113      !!
114      NAMELIST/namdyn_ldf/ ln_dynldf_OFF, ln_dynldf_lap, ln_dynldf_blp,   &   ! type of operator
115         &                 ln_dynldf_lev, ln_dynldf_hor, ln_dynldf_iso,   &   ! acting direction of the operator
116         &                 nn_ahm_ijk_t , rn_Uv    , rn_Lv,   rn_ahm_b,   &   ! lateral eddy coefficient
117         &                 rn_csmc      , rn_minfac    , rn_maxfac            ! Smagorinsky settings
118      !!----------------------------------------------------------------------
119      !
120      READ  ( numnam_ref, namdyn_ldf, IOSTAT = ios, ERR = 901)
121901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namdyn_ldf in reference namelist' )
122
123      READ  ( numnam_cfg, namdyn_ldf, IOSTAT = ios, ERR = 902 )
124902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namdyn_ldf in configuration namelist' )
125      IF(lwm) WRITE ( numond, namdyn_ldf )
126
127      IF(lwp) THEN                      ! Parameter print
128         WRITE(numout,*)
129         WRITE(numout,*) 'ldf_dyn : lateral momentum physics'
130         WRITE(numout,*) '~~~~~~~'
131         WRITE(numout,*) '   Namelist namdyn_ldf : set lateral mixing parameters'
132         !
133         WRITE(numout,*) '      type :'
134         WRITE(numout,*) '         no explicit diffusion                ln_dynldf_OFF = ', ln_dynldf_OFF
135         WRITE(numout,*) '         laplacian operator                   ln_dynldf_lap = ', ln_dynldf_lap
136         WRITE(numout,*) '         bilaplacian operator                 ln_dynldf_blp = ', ln_dynldf_blp
137         !
138         WRITE(numout,*) '      direction of action :'
139         WRITE(numout,*) '         iso-level                            ln_dynldf_lev = ', ln_dynldf_lev
140         WRITE(numout,*) '         horizontal (geopotential)            ln_dynldf_hor = ', ln_dynldf_hor
141         WRITE(numout,*) '         iso-neutral                          ln_dynldf_iso = ', ln_dynldf_iso
142         !
143         WRITE(numout,*) '      coefficients :'
144         WRITE(numout,*) '         type of time-space variation         nn_ahm_ijk_t  = ', nn_ahm_ijk_t
145         WRITE(numout,*) '         lateral viscous velocity  (if cst)      rn_Uv      = ', rn_Uv, ' m/s'
146         WRITE(numout,*) '         lateral viscous length    (if cst)      rn_Lv      = ', rn_Lv, ' m'
147         WRITE(numout,*) '         background viscosity (iso-lap case)     rn_ahm_b   = ', rn_ahm_b, ' m2/s'
148         !
149         WRITE(numout,*) '      Smagorinsky settings (nn_ahm_ijk_t  = 32) :'
150         WRITE(numout,*) '         Smagorinsky coefficient              rn_csmc       = ', rn_csmc
151         WRITE(numout,*) '         factor multiplier for eddy visc.'
152         WRITE(numout,*) '            lower limit (default 1.0)         rn_minfac    = ', rn_minfac
153         WRITE(numout,*) '            upper limit (default 1.0)         rn_maxfac    = ', rn_maxfac
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      ENDIF
228     
229      !
230      !           !==  Space/time variation of eddy coefficients  ==!
231      !           !=================================================!
232      !
233      l_ldfdyn_time = .FALSE.                ! no time variation except in case defined below
234      !
235      IF( ln_dynldf_OFF ) THEN
236         IF(lwp) WRITE(numout,*) '   ==>>>   No viscous operator selected. ahmt and ahmf are not allocated'
237         RETURN
238         !
239      ELSE                                   !==  a lateral diffusion operator is used  ==!
240         !
241         !                                         ! allocate the ahm arrays
242         ALLOCATE( ahmt(jpi,jpj,jpk) , ahmf(jpi,jpj,jpk) , STAT=ierr )
243         IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate arrays')
244         !
245         ahmt(:,:,:) = 0._wp                       ! init to 0 needed
246         ahmf(:,:,:) = 0._wp
247         !
248         !                                         ! value of lap/blp eddy mixing coef.
249         IF(     ln_dynldf_lap ) THEN   ;   zUfac = r1_2 *rn_Uv   ;   inn = 1   ;   cl_Units = ' m2/s'   !   laplacian
250         ELSEIF( ln_dynldf_blp ) THEN   ;   zUfac = r1_12*rn_Uv   ;   inn = 3   ;   cl_Units = ' m4/s'   ! bilaplacian
251         ENDIF
252         zah0    = zUfac *    rn_Lv**inn              ! mixing coefficient
253         zah_max = zUfac * (ra*rad)**inn              ! maximum reachable coefficient (value at the Equator)
254         !
255         SELECT CASE(  nn_ahm_ijk_t  )             !* Specification of space-time variations of ahmt, ahmf
256         !
257         CASE(   0  )      !==  constant  ==!
258            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity. = constant = ', zah0, cl_Units
259            ahmt(:,:,1:jpkm1) = zah0
260            ahmf(:,:,1:jpkm1) = zah0
261            !
262         CASE(  10  )      !==  fixed profile  ==!
263            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( depth )'
264            IF(lwp) WRITE(numout,*) '           surface viscous coef. = constant = ', zah0, cl_Units
265            ahmt(:,:,1) = zah0                        ! constant surface value
266            ahmf(:,:,1) = zah0
267            CALL ldf_c1d( 'DYN', ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )
268            !
269         CASE ( -20 )      !== fixed horizontal shape read in file  ==!
270            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F(i,j) read in eddy_viscosity.nc file'
271            CALL iom_open( 'eddy_viscosity_2D.nc', inum )
272            CALL iom_get ( inum, jpdom_data, 'ahmt_2d', ahmt(:,:,1) )
273            CALL iom_get ( inum, jpdom_data, 'ahmf_2d', ahmf(:,:,1) )
274            CALL iom_close( inum )
275            DO jk = 2, jpkm1
276               ahmt(:,:,jk) = ahmt(:,:,1)
277               ahmf(:,:,jk) = ahmf(:,:,1)
278            END DO
279            !
280         CASE(  20  )      !== fixed horizontal shape  ==!
281            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( e1, e2 ) or F( e1^3, e2^3 ) (lap. or blp. case)'
282            IF(lwp) WRITE(numout,*) '           using a fixed viscous velocity = ', rn_Uv  ,' m/s   and   Lv = Max(e1,e2)'
283            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
284            CALL ldf_c2d( 'DYN', zUfac      , inn        , ahmt, ahmf )         ! surface value proportional to scale factor^inn
285            !
286         CASE( -30  )      !== fixed 3D shape read in file  ==!
287            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F(i,j,k) read in eddy_viscosity_3D.nc file'
288            CALL iom_open( 'eddy_viscosity_3D.nc', inum )
289            CALL iom_get ( inum, jpdom_data, 'ahmt_3d', ahmt )
290            CALL iom_get ( inum, jpdom_data, 'ahmf_3d', ahmf )
291            CALL iom_close( inum )
292            !
293         CASE(  30  )       !==  fixed 3D shape  ==!
294            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth )'
295            IF(lwp) WRITE(numout,*) '           using a fixed viscous velocity = ', rn_Uv  ,' m/s   and   Ld = Max(e1,e2)'
296            IF(lwp) WRITE(numout,*) '           maximum reachable coefficient (at the Equator) = ', zah_max, cl_Units, '  for e1=1°)'
297            CALL ldf_c2d( 'DYN', zUfac      , inn        , ahmt, ahmf )         ! surface value proportional to scale factor^inn
298            CALL ldf_c1d( 'DYN', ahmt(:,:,1), ahmf(:,:,1), ahmt, ahmf )  ! reduction with depth
299            !
300         CASE(  31  )       !==  time varying 3D field  ==!
301            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth , time )'
302            IF(lwp) WRITE(numout,*) '           proportional to the local velocity : 1/2 |u|e (lap) or 1/12 |u|e^3 (blp)'
303            !
304            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
305            !
306         CASE(  32  )       !==  time varying 3D field  ==!
307            IF(lwp) WRITE(numout,*) '   ==>>>   eddy viscosity = F( latitude, longitude, depth , time )'
308            IF(lwp) WRITE(numout,*) '           proportional to the local deformation rate and gridscale (Smagorinsky)'
309            !
310            l_ldfdyn_time = .TRUE.     ! will be calculated by call to ldf_dyn routine in step.F90
311            !
312            !                          ! allocate arrays used in ldf_dyn.
313            ALLOCATE( dtensq(jpi,jpj,jpk) , dshesq(jpi,jpj,jpk) , esqt(jpi,jpj) , esqf(jpi,jpj) , STAT=ierr )
314            IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'ldf_dyn_init: failed to allocate Smagorinsky arrays')
315            !
316            DO_2D_11_11
317               esqt(ji,jj) = ( 2._wp * e1e2t(ji,jj) / ( e1t(ji,jj) + e2t(ji,jj) ) )**2 
318               esqf(ji,jj) = ( 2._wp * e1e2f(ji,jj) / ( e1f(ji,jj) + e2f(ji,jj) ) )**2 
319            END_2D
320            !
321         CASE DEFAULT
322            CALL ctl_stop('ldf_dyn_init: wrong choice for nn_ahm_ijk_t, the type of space-time variation of ahm')
323         END SELECT
324         !
325         IF( .NOT.l_ldfdyn_time ) THEN             !* No time variation
326            IF(     ln_dynldf_lap ) THEN                 !   laplacian operator (mask only)
327!!an          !
328            WRITE(numout,*) '   ln_dynldf_lap_PM = ',ln_dynldf_lap_PM 
329               IF(     ln_dynldf_lap_PM ) THEN                 !   laplacian operator (mask only)
330!! mask tension at the coast (equivalent of ghostpoints at T)
331                  DO jk = 1, jpk
332                     DO jj = 1, jpjm1
333                        DO ji = 1, jpim1      ! NO vector opt.
334                           ! si sum(fmask)==3 = mouillé (on touche pas)
335                           ! si sum = 2 alors on met a 0 zsum = fmask + fmask + fmask,.. et si zsum < 2 -> 0 sinon = 1
336                           zsum =   fmask(ji,jj  ,jk) + fmask(ji+1,jj  ,jk)   &
337                              &   + fmask(ji,jj+1,jk) + fmask(ji+1,jj+1,jk)
338                           IF ( zsum < 2._wp )   ahmt(ji,jj,jk) = 0
339                           !
340                           !ahmt(ji,jj,jk) = ahmt(ji,jj,jk) * fmask(ji,jj  ,jk) * fmask(ji+1,jj  ,jk)   &
341                           !   &                            * fmask(ji,jj+1,jk) * fmask(ji+1,jj+1,jk)
342                        END DO
343                     END DO
344                  END DO
345                  ahmt(jpi,:,1:jpkm1) =  0._wp
346                  ahmt(:,jpj,1:jpkm1) =  0._wp
347                  WRITE(numout,*) '  ahmt x0'
348!! apply no slip at the coast (ssfmask = 1 within the domain and at the coast contrary to fmask in free slip)
349                   DO jk = 1, jpkm1
350                      ahmf(:,:,jk) =    ahmf(:,:,jk) * ( 2._wp * ssfmask(:,:) - fmask(:,:,jk) )
351                   END DO
352                   WRITE(numout,*) '  ahmf x2'
353               ELSE
354               ! classic boundary condition on the viscosity coefficient
355                  ahmt(:,:,1:jpkm1) =       ahmt(:,:,1:jpkm1)   * tmask(:,:,1:jpkm1)
356                  WRITE(numout,*) ' ahmt tmasked '
357                  ahmf(:,:,1:jpkm1) =       ahmf(:,:,1:jpkm1)   * fmask(:,:,1:jpkm1)
358                  WRITE(numout,*) ' ahmf fmasked '
359               ENDIF
360!!an         !                 
361            ELSEIF( ln_dynldf_blp ) THEN                 ! bilaplacian operator (square root + mask)
362               ahmt(:,:,1:jpkm1) = SQRT( ahmt(:,:,1:jpkm1) ) * tmask(:,:,1:jpkm1)
363               ahmf(:,:,1:jpkm1) = SQRT( ahmf(:,:,1:jpkm1) ) * fmask(:,:,1:jpkm1)
364            ENDIF
365         ENDIF
366         !
367      ENDIF
368      !
369   END SUBROUTINE ldf_dyn_init
370
371
372   SUBROUTINE ldf_dyn( kt, Kbb )
373      !!----------------------------------------------------------------------
374      !!                  ***  ROUTINE ldf_dyn  ***
375      !!
376      !! ** Purpose :   update at kt the momentum lateral mixing coeff. (ahmt and ahmf)
377      !!
378      !! ** Method  :   time varying eddy viscosity coefficients:
379      !!
380      !!    nn_ahm_ijk_t = 31    ahmt, ahmf = F(i,j,k,t) = F(local velocity)
381      !!                         ( |u|e /12  or  |u|e^3/12 for laplacian or bilaplacian operator )
382      !!
383      !!    nn_ahm_ijk_t = 32    ahmt, ahmf = F(i,j,k,t) = F(local deformation rate and gridscale) (D and L) (Smagorinsky) 
384      !!                         ( L^2|D|    or  L^4|D|/8  for laplacian or bilaplacian operator )
385      !!
386      !! ** note    :    in BLP cases the sqrt of the eddy coef is returned, since bilaplacian is en re-entrant laplacian
387      !! ** action  :    ahmt, ahmf   updated at each time step
388      !!----------------------------------------------------------------------
389      INTEGER, INTENT(in) ::   kt   ! time step index
390      INTEGER, INTENT(in) ::   Kbb  ! ocean time level indices
391      !
392      INTEGER  ::   ji, jj, jk   ! dummy loop indices
393      REAL(wp) ::   zu2pv2_ij_p1, zu2pv2_ij, zu2pv2_ij_m1, zemax   ! local scalar (option 31)
394      REAL(wp) ::   zcmsmag, zstabf_lo, zstabf_up, zdelta, zdb     ! local scalar (option 32)
395      !!----------------------------------------------------------------------
396      !
397      IF( ln_timing )   CALL timing_start('ldf_dyn')
398      !
399      SELECT CASE(  nn_ahm_ijk_t  )       !== Eddy vicosity coefficients ==!
400      !
401      CASE(  31  )       !==  time varying 3D field  ==!   = F( local velocity )
402         !
403         IF( ln_dynldf_lap   ) THEN        ! laplacian operator : |u| e /12 = |u/144| e
404            DO jk = 1, jpkm1
405               DO_2D_00_00
406                  zu2pv2_ij    = uu(ji  ,jj  ,jk,Kbb) * uu(ji  ,jj  ,jk,Kbb) + vv(ji  ,jj  ,jk,Kbb) * vv(ji  ,jj  ,jk,Kbb)
407                  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)
408                  zemax = MAX( e1t(ji,jj) , e2t(ji,jj) )
409                  ahmt(ji,jj,jk) = SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zemax * tmask(ji,jj,jk)      ! 288= 12*12 * 2
410               END_2D
411               DO_2D_10_10
412                  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)
413                  zu2pv2_ij    = uu(ji  ,jj  ,jk, Kbb) * uu(ji  ,jj  ,jk, Kbb) + vv(ji  ,jj  ,jk, Kbb) * vv(ji  ,jj  ,jk, Kbb)
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_2D
417            END DO
418         ELSEIF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( |u| e^3 /12 ) = sqrt( |u/144| e ) * e
419            DO jk = 1, jpkm1
420               DO_2D_00_00
421                  zu2pv2_ij    = uu(ji  ,jj  ,jk,Kbb) * uu(ji  ,jj  ,jk,Kbb) + vv(ji  ,jj  ,jk,Kbb) * vv(ji  ,jj  ,jk,Kbb)
422                  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)
423                  zemax = MAX( e1t(ji,jj) , e2t(ji,jj) )
424                  ahmt(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_m1) * r1_288 ) * zemax  ) * zemax * tmask(ji,jj,jk)
425               END_2D
426               DO_2D_10_10
427                  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)
428                  zu2pv2_ij    = uu(ji  ,jj  ,jk, Kbb) * uu(ji  ,jj  ,jk, Kbb) + vv(ji  ,jj  ,jk, Kbb) * vv(ji  ,jj  ,jk, Kbb)
429                  zemax = MAX( e1f(ji,jj) , e2f(ji,jj) )
430                  ahmf(ji,jj,jk) = SQRT(  SQRT( (zu2pv2_ij + zu2pv2_ij_p1) * r1_288 ) * zemax  ) * zemax * fmask(ji,jj,jk)
431               END_2D
432            END DO
433         ENDIF
434         !
435         CALL lbc_lnk_multi( 'ldfdyn', ahmt, 'T', 1.,  ahmf, 'F', 1. )
436         !
437         !
438      CASE(  32  )       !==  time varying 3D field  ==!   = F( local deformation rate and gridscale ) (Smagorinsky)
439         !
440         IF( ln_dynldf_lap .OR. ln_dynldf_blp  ) THEN        ! laplacian operator : (C_smag/pi)^2 L^2 |D|
441            !
442            zcmsmag   = (rn_csmc/rpi)**2                                            ! (C_smag/pi)^2
443            zstabf_lo = rn_minfac * rn_minfac / ( 2._wp * 12._wp * 12._wp * zcmsmag ) ! lower limit stability factor scaling
444            zstabf_up = rn_maxfac / ( 4._wp * zcmsmag * 2._wp * rn_Dt )               ! upper limit stability factor scaling
445            IF( ln_dynldf_blp ) zstabf_lo = ( 16._wp / 9._wp ) * zstabf_lo          ! provide |U|L^3/12 lower limit instead
446            !                                                                       ! of |U|L^3/16 in blp case
447            DO jk = 1, jpkm1
448               !
449               DO_2D_00_00
450                  zdb =    ( uu(ji,jj,jk,Kbb) * r1_e2u(ji,jj) -  uu(ji-1,jj,jk,Kbb) * r1_e2u(ji-1,jj) )  &
451                       &                      * r1_e1t(ji,jj) * e2t(ji,jj)                           &
452                       & - ( vv(ji,jj,jk,Kbb) * r1_e1v(ji,jj) -  vv(ji,jj-1,jk,Kbb) * r1_e1v(ji,jj-1) )  &
453                       &                      * r1_e2t(ji,jj) * e1t(ji,jj)
454                  dtensq(ji,jj,jk) = zdb * zdb * tmask(ji,jj,jk)
455               END_2D
456               !
457               DO_2D_10_10
458                  zdb =   (  uu(ji,jj+1,jk,Kbb) * r1_e1u(ji,jj+1) -  uu(ji,jj,jk,Kbb) * r1_e1u(ji,jj) )  &
459                       &                        * r1_e2f(ji,jj)   * e1f(ji,jj)                       &
460                       & + ( vv(ji+1,jj,jk,Kbb) * r1_e2v(ji+1,jj) -  vv(ji,jj,jk,Kbb) * r1_e2v(ji,jj) )  &
461                       &                        * r1_e1f(ji,jj)   * e2f(ji,jj)
462                  dshesq(ji,jj,jk) = zdb * zdb * fmask(ji,jj,jk)
463               END_2D
464               !
465            END DO
466            !
467            CALL lbc_lnk_multi( 'ldfdyn', dtensq, 'T', 1. )  ! lbc_lnk on dshesq not needed
468            !
469            DO jk = 1, jpkm1
470              !
471               DO_2D_00_00
472                  !
473                  zu2pv2_ij    = uu(ji  ,jj  ,jk,Kbb) * uu(ji  ,jj  ,jk,Kbb) + vv(ji  ,jj  ,jk,Kbb) * vv(ji  ,jj  ,jk,Kbb)
474                  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)
475                  !
476                  zdelta         = zcmsmag * esqt(ji,jj)                                        ! L^2 * (C_smag/pi)^2
477                  ahmt(ji,jj,jk) = zdelta * SQRT(          dtensq(ji  ,jj,jk) +                         &
478                     &                            r1_4 * ( dshesq(ji  ,jj,jk) + dshesq(ji  ,jj-1,jk) +  &
479                     &                                     dshesq(ji-1,jj,jk) + dshesq(ji-1,jj-1,jk) ) )
480                  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
481                  ahmt(ji,jj,jk) = MIN( ahmt(ji,jj,jk),                                    zdelta * zstabf_up )   ! Impose upper limit == maxfac  * L^2/(4*2dt)
482                  !
483               END_2D
484               !
485               DO_2D_10_10
486                  !
487                  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)
488                  zu2pv2_ij    = uu(ji  ,jj  ,jk, kbb) * uu(ji  ,jj  ,jk, kbb) + vv(ji  ,jj  ,jk, kbb) * vv(ji  ,jj  ,jk, kbb)
489                  !
490                  zdelta         = zcmsmag * esqf(ji,jj)                                        ! L^2 * (C_smag/pi)^2
491                  ahmf(ji,jj,jk) = zdelta * SQRT(          dshesq(ji  ,jj,jk) +                         &
492                     &                            r1_4 * ( dtensq(ji  ,jj,jk) + dtensq(ji  ,jj+1,jk) +  &
493                     &                                     dtensq(ji+1,jj,jk) + dtensq(ji+1,jj+1,jk) ) )
494                  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
495                  ahmf(ji,jj,jk) = MIN( ahmf(ji,jj,jk),                                    zdelta * zstabf_up )   ! Impose upper limit == maxfac  * L^2/(4*2dt)
496                  !
497               END_2D
498               !
499            END DO
500            !
501         ENDIF
502         !
503         IF( ln_dynldf_blp ) THEN      ! bilaplacian operator : sqrt( (C_smag/pi)^2 L^4 |D|/8)
504            !                          !                      = sqrt( A_lap_smag L^2/8 )
505            !                          ! stability limits already applied to laplacian values
506            !                          ! effective default limits are 1/12 |U|L^3 < B_hm < 1//(32*2dt) L^4
507            DO jk = 1, jpkm1
508               DO_2D_00_00
509                  ahmt(ji,jj,jk) = SQRT( r1_8 * esqt(ji,jj) * ahmt(ji,jj,jk) )
510               END_2D
511               DO_2D_10_10
512                  ahmf(ji,jj,jk) = SQRT( r1_8 * esqf(ji,jj) * ahmf(ji,jj,jk) )
513               END_2D
514            END DO
515            !
516         ENDIF
517         !
518         CALL lbc_lnk_multi( 'ldfdyn', ahmt, 'T', 1. , ahmf, 'F', 1. )
519         !
520      END SELECT
521      !
522      CALL iom_put( "ahmt_2d", ahmt(:,:,1) )   ! surface u-eddy diffusivity coeff.
523      CALL iom_put( "ahmf_2d", ahmf(:,:,1) )   ! surface v-eddy diffusivity coeff.
524      CALL iom_put( "ahmt_3d", ahmt(:,:,:) )   ! 3D      u-eddy diffusivity coeff.
525      CALL iom_put( "ahmf_3d", ahmf(:,:,:) )   ! 3D      v-eddy diffusivity coeff.
526      !
527      IF( ln_timing )   CALL timing_stop('ldf_dyn')
528      !
529   END SUBROUTINE ldf_dyn
530
531   !!======================================================================
532END MODULE ldfdyn
Note: See TracBrowser for help on using the repository browser.