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