Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
sx_cfl.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34submodule(opr_sx) sx_cfl
35 implicit none
36
37contains
38
39 module function opr_sx_cfl(dt, u, v, w, xh, coef, nelv) result(cfl)
40 type(space_t), intent(in) :: Xh
41 type(coef_t), intent(in) :: coef
42 integer, intent(in) :: nelv
43 real(kind=rp), intent(in) :: dt
44 real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz, nelv) :: u, v, w
45 real(kind=rp) :: cfl
46
47 select case (xh%lx)
48 case (14)
49 cfl = sx_cfl_lx14(dt, u, v, w, &
50 coef%drdx, coef%dsdx, coef%dtdx, &
51 coef%drdy, coef%dsdy, coef%dtdy, &
52 coef%drdz, coef%dsdz, coef%dtdz, &
53 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
54 coef%jacinv, nelv)
55 case (13)
56 cfl = sx_cfl_lx13(dt, u, v, w, &
57 coef%drdx, coef%dsdx, coef%dtdx, &
58 coef%drdy, coef%dsdy, coef%dtdy, &
59 coef%drdz, coef%dsdz, coef%dtdz, &
60 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
61 coef%jacinv, nelv)
62 case (12)
63 cfl = sx_cfl_lx12(dt, u, v, w, &
64 coef%drdx, coef%dsdx, coef%dtdx, &
65 coef%drdy, coef%dsdy, coef%dtdy, &
66 coef%drdz, coef%dsdz, coef%dtdz, &
67 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
68 coef%jacinv, nelv)
69 case (11)
70 cfl = sx_cfl_lx11(dt, u, v, w, &
71 coef%drdx, coef%dsdx, coef%dtdx, &
72 coef%drdy, coef%dsdy, coef%dtdy, &
73 coef%drdz, coef%dsdz, coef%dtdz, &
74 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
75 coef%jacinv, nelv)
76 case (10)
77 cfl = sx_cfl_lx10(dt, u, v, w, &
78 coef%drdx, coef%dsdx, coef%dtdx, &
79 coef%drdy, coef%dsdy, coef%dtdy, &
80 coef%drdz, coef%dsdz, coef%dtdz, &
81 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
82 coef%jacinv, nelv)
83 case (9)
84 cfl = sx_cfl_lx9(dt, u, v, w, &
85 coef%drdx, coef%dsdx, coef%dtdx, &
86 coef%drdy, coef%dsdy, coef%dtdy, &
87 coef%drdz, coef%dsdz, coef%dtdz, &
88 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
89 coef%jacinv, nelv)
90 case (8)
91 cfl = sx_cfl_lx8(dt, u, v, w, &
92 coef%drdx, coef%dsdx, coef%dtdx, &
93 coef%drdy, coef%dsdy, coef%dtdy, &
94 coef%drdz, coef%dsdz, coef%dtdz, &
95 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
96 coef%jacinv, nelv)
97 case (7)
98 cfl = sx_cfl_lx7(dt, u, v, w, &
99 coef%drdx, coef%dsdx, coef%dtdx, &
100 coef%drdy, coef%dsdy, coef%dtdy, &
101 coef%drdz, coef%dsdz, coef%dtdz, &
102 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
103 coef%jacinv, nelv)
104 case (6)
105 cfl = sx_cfl_lx6(dt, u, v, w, &
106 coef%drdx, coef%dsdx, coef%dtdx, &
107 coef%drdy, coef%dsdy, coef%dtdy, &
108 coef%drdz, coef%dsdz, coef%dtdz, &
109 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
110 coef%jacinv, nelv)
111 case (5)
112 cfl = sx_cfl_lx5(dt, u, v, w, &
113 coef%drdx, coef%dsdx, coef%dtdx, &
114 coef%drdy, coef%dsdy, coef%dtdy, &
115 coef%drdz, coef%dsdz, coef%dtdz, &
116 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
117 coef%jacinv, nelv)
118 case (4)
119 cfl = sx_cfl_lx4(dt, u, v, w, &
120 coef%drdx, coef%dsdx, coef%dtdx, &
121 coef%drdy, coef%dsdy, coef%dtdy, &
122 coef%drdz, coef%dsdz, coef%dtdz, &
123 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
124 coef%jacinv, nelv)
125 case (3)
126 cfl = sx_cfl_lx3(dt, u, v, w, &
127 coef%drdx, coef%dsdx, coef%dtdx, &
128 coef%drdy, coef%dsdy, coef%dtdy, &
129 coef%drdz, coef%dsdz, coef%dtdz, &
130 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
131 coef%jacinv, nelv)
132 case (2)
133 cfl = sx_cfl_lx2(dt, u, v, w, &
134 coef%drdx, coef%dsdx, coef%dtdx, &
135 coef%drdy, coef%dsdy, coef%dtdy, &
136 coef%drdz, coef%dsdz, coef%dtdz, &
137 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
138 coef%jacinv, nelv)
139 case default
140 cfl = sx_cfl_lx(dt, u, v, w, &
141 coef%drdx, coef%dsdx, coef%dtdx, &
142 coef%drdy, coef%dsdy, coef%dtdy, &
143 coef%drdz, coef%dsdz, coef%dtdz, &
144 xh%dr_inv, xh%ds_inv, xh%dt_inv, &
145 coef%jacinv, nelv, xh%lx)
146 end select
147
148 end function opr_sx_cfl
149
150 function sx_cfl_lx(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
151 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
152 jacinv, nelv, lx) result(cfl)
153 integer, intent(in) :: nelv, lx
154 real(kind=rp), intent(in) :: dt
155 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
156 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
157 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
158 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
159 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
160 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
161 real(kind=rp) :: cflr, cfls, cflt, cflm
162 real(kind=rp) :: ur, us, ut
163 real(kind=rp) :: cfl
164 integer :: i, j, k, e
165 cfl = 0d0
166
167 do k = 1, lx
168 do j = 1, lx
169 do i = 1, lx
170 do e = 1, nelv
171 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
172 + v(i,j,k,e)*drdy(i,j,k,e) &
173 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
174 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
175 + v(i,j,k,e)*dsdy(i,j,k,e) &
176 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
177 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
178 + v(i,j,k,e)*dtdy(i,j,k,e) &
179 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
180
181 cflr = abs(dt*ur*dr_inv(i))
182 cfls = abs(dt*us*ds_inv(j))
183 cflt = abs(dt*ut*dt_inv(k))
184
185 cflm = cflr + cfls + cflt
186 cfl = max(cfl, cflm)
187 end do
188 end do
189 end do
190 end do
191
192 end function sx_cfl_lx
193
194 function sx_cfl_lx14(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
195 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
196 jacinv, nelv) result(cfl)
197 integer, parameter :: lx = 14
198 integer, intent(in) :: nelv
199 real(kind=rp), intent(in) :: dt
200 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
201 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
202 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
203 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
204 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
205 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
206 real(kind=rp) :: cflr, cfls, cflt, cflm
207 real(kind=rp) :: ur, us, ut
208 real(kind=rp) :: cfl
209 integer :: i, j, k, e
210 cfl = 0d0
211
212 do k = 1, lx
213 do j = 1, lx
214 do i = 1, lx
215 do e = 1, nelv
216 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
217 + v(i,j,k,e)*drdy(i,j,k,e) &
218 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
219 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
220 + v(i,j,k,e)*dsdy(i,j,k,e) &
221 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
222 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
223 + v(i,j,k,e)*dtdy(i,j,k,e) &
224 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
225
226 cflr = abs(dt*ur*dr_inv(i))
227 cfls = abs(dt*us*ds_inv(j))
228 cflt = abs(dt*ut*dt_inv(k))
229
230 cflm = cflr + cfls + cflt
231 cfl = max(cfl, cflm)
232 end do
233 end do
234 end do
235 end do
236
237 end function sx_cfl_lx14
238
239 function sx_cfl_lx13(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
240 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
241 jacinv, nelv) result(cfl)
242 integer, parameter :: lx = 13
243 integer, intent(in) :: nelv
244 real(kind=rp), intent(in) :: dt
245 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
246 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
247 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
248 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
249 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
250 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
251 real(kind=rp) :: cflr, cfls, cflt, cflm
252 real(kind=rp) :: ur, us, ut
253 real(kind=rp) :: cfl
254 integer :: i, j, k, e
255 cfl = 0d0
256
257 do k = 1, lx
258 do j = 1, lx
259 do i = 1, lx
260 do e = 1, nelv
261 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
262 + v(i,j,k,e)*drdy(i,j,k,e) &
263 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
264 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
265 + v(i,j,k,e)*dsdy(i,j,k,e) &
266 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
267 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
268 + v(i,j,k,e)*dtdy(i,j,k,e) &
269 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
270
271 cflr = abs(dt*ur*dr_inv(i))
272 cfls = abs(dt*us*ds_inv(j))
273 cflt = abs(dt*ut*dt_inv(k))
274
275 cflm = cflr + cfls + cflt
276 cfl = max(cfl, cflm)
277 end do
278 end do
279 end do
280 end do
281
282 end function sx_cfl_lx13
283
284 function sx_cfl_lx12(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
285 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
286 jacinv, nelv) result(cfl)
287 integer, parameter :: lx = 12
288 integer, intent(in) :: nelv
289 real(kind=rp), intent(in) :: dt
290 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
291 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
292 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
293 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
294 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
295 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
296 real(kind=rp) :: cflr, cfls, cflt, cflm
297 real(kind=rp) :: ur, us, ut
298 real(kind=rp) :: cfl
299 integer :: i, j, k, e
300 cfl = 0d0
301
302 do k = 1, lx
303 do j = 1, lx
304 do i = 1, lx
305 do e = 1, nelv
306 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
307 + v(i,j,k,e)*drdy(i,j,k,e) &
308 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
309 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
310 + v(i,j,k,e)*dsdy(i,j,k,e) &
311 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
312 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
313 + v(i,j,k,e)*dtdy(i,j,k,e) &
314 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
315
316 cflr = abs(dt*ur*dr_inv(i))
317 cfls = abs(dt*us*ds_inv(j))
318 cflt = abs(dt*ut*dt_inv(k))
319
320 cflm = cflr + cfls + cflt
321 cfl = max(cfl, cflm)
322 end do
323 end do
324 end do
325 end do
326
327 end function sx_cfl_lx12
328
329 function sx_cfl_lx11(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
330 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
331 jacinv, nelv) result(cfl)
332 integer, parameter :: lx = 11
333 integer, intent(in) :: nelv
334 real(kind=rp), intent(in) :: dt
335 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
336 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
337 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
338 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
339 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
340 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
341 real(kind=rp) :: cflr, cfls, cflt, cflm
342 real(kind=rp) :: ur, us, ut
343 real(kind=rp) :: cfl
344 integer :: i, j, k, e
345 cfl = 0d0
346
347 do k = 1, lx
348 do j = 1, lx
349 do i = 1, lx
350 do e = 1, nelv
351 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
352 + v(i,j,k,e)*drdy(i,j,k,e) &
353 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
354 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
355 + v(i,j,k,e)*dsdy(i,j,k,e) &
356 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
357 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
358 + v(i,j,k,e)*dtdy(i,j,k,e) &
359 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
360
361 cflr = abs(dt*ur*dr_inv(i))
362 cfls = abs(dt*us*ds_inv(j))
363 cflt = abs(dt*ut*dt_inv(k))
364
365 cflm = cflr + cfls + cflt
366 cfl = max(cfl, cflm)
367 end do
368 end do
369 end do
370 end do
371
372 end function sx_cfl_lx11
373
374 function sx_cfl_lx10(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
375 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
376 jacinv, nelv) result(cfl)
377 integer, parameter :: lx = 10
378 integer, intent(in) :: nelv
379 real(kind=rp), intent(in) :: dt
380 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
381 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
382 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
383 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
384 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
385 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
386 real(kind=rp) :: cflr, cfls, cflt, cflm
387 real(kind=rp) :: ur, us, ut
388 real(kind=rp) :: cfl
389 integer :: i, j, k, e
390 cfl = 0d0
391
392 do k = 1, lx
393 do j = 1, lx
394 do i = 1, lx
395 do e = 1, nelv
396 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
397 + v(i,j,k,e)*drdy(i,j,k,e) &
398 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
399 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
400 + v(i,j,k,e)*dsdy(i,j,k,e) &
401 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
402 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
403 + v(i,j,k,e)*dtdy(i,j,k,e) &
404 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
405
406 cflr = abs(dt*ur*dr_inv(i))
407 cfls = abs(dt*us*ds_inv(j))
408 cflt = abs(dt*ut*dt_inv(k))
409
410 cflm = cflr + cfls + cflt
411 cfl = max(cfl, cflm)
412 end do
413 end do
414 end do
415 end do
416
417 end function sx_cfl_lx10
418
419 function sx_cfl_lx9(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
420 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
421 jacinv, nelv) result(cfl)
422 integer, parameter :: lx = 9
423 integer, intent(in) :: nelv
424 real(kind=rp), intent(in) :: dt
425 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
426 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
427 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
428 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
429 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
430 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
431 real(kind=rp) :: cflr, cfls, cflt, cflm
432 real(kind=rp) :: ur, us, ut
433 real(kind=rp) :: cfl
434 integer :: i, j, k, e
435 cfl = 0d0
436
437 do k = 1, lx
438 do j = 1, lx
439 do i = 1, lx
440 do e = 1, nelv
441 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
442 + v(i,j,k,e)*drdy(i,j,k,e) &
443 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
444 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
445 + v(i,j,k,e)*dsdy(i,j,k,e) &
446 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
447 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
448 + v(i,j,k,e)*dtdy(i,j,k,e) &
449 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
450
451 cflr = abs(dt*ur*dr_inv(i))
452 cfls = abs(dt*us*ds_inv(j))
453 cflt = abs(dt*ut*dt_inv(k))
454
455 cflm = cflr + cfls + cflt
456 cfl = max(cfl, cflm)
457 end do
458 end do
459 end do
460 end do
461
462 end function sx_cfl_lx9
463
464 function sx_cfl_lx8(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
465 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
466 jacinv, nelv) result(cfl)
467 integer, parameter :: lx = 8
468 integer, intent(in) :: nelv
469 real(kind=rp), intent(in) :: dt
470 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
471 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
472 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
473 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
474 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
475 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
476 real(kind=rp) :: cflr, cfls, cflt, cflm
477 real(kind=rp) :: ur, us, ut
478 real(kind=rp) :: cfl
479 integer :: i, j, k, e
480 cfl = 0d0
481
482 do k = 1, lx
483 do j = 1, lx
484 do i = 1, lx
485 do e = 1, nelv
486 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
487 + v(i,j,k,e)*drdy(i,j,k,e) &
488 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
489 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
490 + v(i,j,k,e)*dsdy(i,j,k,e) &
491 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
492 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
493 + v(i,j,k,e)*dtdy(i,j,k,e) &
494 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
495
496 cflr = abs(dt*ur*dr_inv(i))
497 cfls = abs(dt*us*ds_inv(j))
498 cflt = abs(dt*ut*dt_inv(k))
499
500 cflm = cflr + cfls + cflt
501 cfl = max(cfl, cflm)
502 end do
503 end do
504 end do
505 end do
506
507 end function sx_cfl_lx8
508
509 function sx_cfl_lx7(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
510 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
511 jacinv, nelv) result(cfl)
512 integer, parameter :: lx = 7
513 integer, intent(in) :: nelv
514 real(kind=rp), intent(in) :: dt
515 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
516 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
517 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
518 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
519 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
520 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
521 real(kind=rp) :: cflr, cfls, cflt, cflm
522 real(kind=rp) :: ur, us, ut
523 real(kind=rp) :: cfl
524 integer :: i, j, k, e
525 cfl = 0d0
526
527 do k = 1, lx
528 do j = 1, lx
529 do i = 1, lx
530 do e = 1, nelv
531 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
532 + v(i,j,k,e)*drdy(i,j,k,e) &
533 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
534 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
535 + v(i,j,k,e)*dsdy(i,j,k,e) &
536 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
537 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
538 + v(i,j,k,e)*dtdy(i,j,k,e) &
539 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
540
541 cflr = abs(dt*ur*dr_inv(i))
542 cfls = abs(dt*us*ds_inv(j))
543 cflt = abs(dt*ut*dt_inv(k))
544
545 cflm = cflr + cfls + cflt
546 cfl = max(cfl, cflm)
547 end do
548 end do
549 end do
550 end do
551
552 end function sx_cfl_lx7
553
554 function sx_cfl_lx6(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
555 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
556 jacinv, nelv) result(cfl)
557 integer, parameter :: lx = 6
558 integer, intent(in) :: nelv
559 real(kind=rp), intent(in) :: dt
560 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
561 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
562 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
563 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
564 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
565 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
566 real(kind=rp) :: cflr, cfls, cflt, cflm
567 real(kind=rp) :: ur, us, ut
568 real(kind=rp) :: cfl
569 integer :: i, j, k, e
570 cfl = 0d0
571
572 do k = 1, lx
573 do j = 1, lx
574 do i = 1, lx
575 do e = 1, nelv
576 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
577 + v(i,j,k,e)*drdy(i,j,k,e) &
578 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
579 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
580 + v(i,j,k,e)*dsdy(i,j,k,e) &
581 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
582 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
583 + v(i,j,k,e)*dtdy(i,j,k,e) &
584 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
585
586 cflr = abs(dt*ur*dr_inv(i))
587 cfls = abs(dt*us*ds_inv(j))
588 cflt = abs(dt*ut*dt_inv(k))
589
590 cflm = cflr + cfls + cflt
591 cfl = max(cfl, cflm)
592 end do
593 end do
594 end do
595 end do
596
597 end function sx_cfl_lx6
598
599 function sx_cfl_lx5(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
600 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
601 jacinv, nelv) result(cfl)
602 integer, parameter :: lx = 5
603 integer, intent(in) :: nelv
604 real(kind=rp), intent(in) :: dt
605 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
606 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
607 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
608 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
609 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
610 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
611 real(kind=rp) :: cflr, cfls, cflt, cflm
612 real(kind=rp) :: ur, us, ut
613 real(kind=rp) :: cfl
614 integer :: i, j, k, e
615 cfl = 0d0
616
617 do k = 1, lx
618 do j = 1, lx
619 do i = 1, lx
620 do e = 1, nelv
621 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
622 + v(i,j,k,e)*drdy(i,j,k,e) &
623 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
624 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
625 + v(i,j,k,e)*dsdy(i,j,k,e) &
626 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
627 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
628 + v(i,j,k,e)*dtdy(i,j,k,e) &
629 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
630
631 cflr = abs(dt*ur*dr_inv(i))
632 cfls = abs(dt*us*ds_inv(j))
633 cflt = abs(dt*ut*dt_inv(k))
634
635 cflm = cflr + cfls + cflt
636 cfl = max(cfl, cflm)
637 end do
638 end do
639 end do
640 end do
641
642 end function sx_cfl_lx5
643
644 function sx_cfl_lx4(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
645 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
646 jacinv, nelv) result(cfl)
647 integer, parameter :: lx = 4
648 integer, intent(in) :: nelv
649 real(kind=rp), intent(in) :: dt
650 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
651 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
652 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
653 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
654 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
655 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
656 real(kind=rp) :: cflr, cfls, cflt, cflm
657 real(kind=rp) :: ur, us, ut
658 real(kind=rp) :: cfl
659 integer :: i, j, k, e
660 cfl = 0d0
661
662 do k = 1, lx
663 do j = 1, lx
664 do i = 1, lx
665 do e = 1, nelv
666 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
667 + v(i,j,k,e)*drdy(i,j,k,e) &
668 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
669 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
670 + v(i,j,k,e)*dsdy(i,j,k,e) &
671 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
672 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
673 + v(i,j,k,e)*dtdy(i,j,k,e) &
674 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
675
676 cflr = abs(dt*ur*dr_inv(i))
677 cfls = abs(dt*us*ds_inv(j))
678 cflt = abs(dt*ut*dt_inv(k))
679
680 cflm = cflr + cfls + cflt
681 cfl = max(cfl, cflm)
682 end do
683 end do
684 end do
685 end do
686
687 end function sx_cfl_lx4
688
689 function sx_cfl_lx3(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
690 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
691 jacinv, nelv) result(cfl)
692 integer, parameter :: lx = 3
693 integer, intent(in) :: nelv
694 real(kind=rp), intent(in) :: dt
695 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
696 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
697 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
698 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
699 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
700 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
701 real(kind=rp) :: cflr, cfls, cflt, cflm
702 real(kind=rp) :: ur, us, ut
703 real(kind=rp) :: cfl
704 integer :: i, j, k, e
705 cfl = 0d0
706
707 do k = 1, lx
708 do j = 1, lx
709 do i = 1, lx
710 do e = 1, nelv
711 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
712 + v(i,j,k,e)*drdy(i,j,k,e) &
713 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
714 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
715 + v(i,j,k,e)*dsdy(i,j,k,e) &
716 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
717 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
718 + v(i,j,k,e)*dtdy(i,j,k,e) &
719 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
720
721 cflr = abs(dt*ur*dr_inv(i))
722 cfls = abs(dt*us*ds_inv(j))
723 cflt = abs(dt*ut*dt_inv(k))
724
725 cflm = cflr + cfls + cflt
726 cfl = max(cfl, cflm)
727 end do
728 end do
729 end do
730 end do
731
732 end function sx_cfl_lx3
733
734 function sx_cfl_lx2(dt, u, v, w, drdx, dsdx, dtdx, drdy, dsdy, dtdy, &
735 drdz, dsdz, dtdz, dr_inv, ds_inv, dt_inv, &
736 jacinv, nelv) result(cfl)
737 integer, parameter :: lx = 2
738 integer, intent(in) :: nelv
739 real(kind=rp), intent(in) :: dt
740 real(kind=rp), dimension(lx, lx, lx, nelv) :: u, v, w
741 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdx, dsdx, dtdx
742 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdy, dsdy, dtdy
743 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: drdz, dsdz, dtdz
744 real(kind=rp), dimension(lx), intent(in) :: dr_inv, ds_inv, dt_inv
745 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: jacinv
746 real(kind=rp) :: cflr, cfls, cflt, cflm
747 real(kind=rp) :: ur, us, ut
748 real(kind=rp) :: cfl
749 integer :: i, j, k, e
750 cfl = 0d0
751
752 do k = 1, lx
753 do j = 1, lx
754 do i = 1, lx
755 do e = 1, nelv
756 ur = ( u(i,j,k,e)*drdx(i,j,k,e) &
757 + v(i,j,k,e)*drdy(i,j,k,e) &
758 + w(i,j,k,e)*drdz(i,j,k,e) ) * jacinv(i,j,k,e)
759 us = ( u(i,j,k,e)*dsdx(i,j,k,e) &
760 + v(i,j,k,e)*dsdy(i,j,k,e) &
761 + w(i,j,k,e)*dsdz(i,j,k,e) ) * jacinv(i,j,k,e)
762 ut = ( u(i,j,k,e)*dtdx(i,j,k,e) &
763 + v(i,j,k,e)*dtdy(i,j,k,e) &
764 + w(i,j,k,e)*dtdz(i,j,k,e) ) * jacinv(i,j,k,e)
765
766 cflr = abs(dt*ur*dr_inv(i))
767 cfls = abs(dt*us*ds_inv(j))
768 cflt = abs(dt*ut*dt_inv(k))
769
770 cflm = cflr + cfls + cflt
771 cfl = max(cfl, cflm)
772 end do
773 end do
774 end do
775 end do
776
777 end function sx_cfl_lx2
778
779end submodule sx_cfl
Operators SX-Aurora backend.
Definition opr_sx.f90:2
#define max(a, b)
Definition tensor.cu:40