Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
convect_scalar.f90
Go to the documentation of this file.
1! Copyright (c) 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_cpu) cpu_convect_scalar
35 use num_types, only : rp
36 use math, only : col2
37 implicit none
38
39contains
40
41 module subroutine opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, coef_gll, &
42 coef_gl, gll_to_gl)
43 type(space_t), intent(in) :: Xh_GL
44 type(space_t), intent(in) :: Xh_GLL
45 type(coef_t), intent(in) :: coef_GLL
46 type(coef_t), intent(in) :: coef_GL
47 type(interpolator_t), intent(inout) :: GLL_to_GL
48 real(kind=rp), intent(inout) :: &
49 du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
50 real(kind=rp), intent(inout) :: &
51 u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
52 real(kind=rp), intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
53 associate(dx => xh_gl%dx, dy => xh_gl%dy, dz => xh_gl%dz, &
54 lx => xh_gl%lx, nelv => coef_gl%msh%nelv)
55
56 select case (lx)
57 case (18)
58 call cpu_convect_scalar_lx18(du, u, c, dx, dy, dz, &
59 xh_gll, coef_gll, gll_to_gl, nelv)
60 case (17)
61 call cpu_convect_scalar_lx17(du, u, c, dx, dy, dz, &
62 xh_gll, coef_gll, gll_to_gl, nelv)
63 case (16)
64 call cpu_convect_scalar_lx16(du, u, c, dx, dy, dz, &
65 xh_gll, coef_gll, gll_to_gl, nelv)
66 case (15)
67 call cpu_convect_scalar_lx15(du, u, c, dx, dy, dz, &
68 xh_gll, coef_gll, gll_to_gl, nelv)
69 case (14)
70 call cpu_convect_scalar_lx14(du, u, c, dx, dy, dz, &
71 xh_gll, coef_gll, gll_to_gl, nelv)
72 case (13)
73 call cpu_convect_scalar_lx13(du, u, c, dx, dy, dz, &
74 xh_gll, coef_gll, gll_to_gl, nelv)
75 case (12)
76 call cpu_convect_scalar_lx12(du, u, c, dx, dy, dz, &
77 xh_gll, coef_gll, gll_to_gl, nelv)
78 case (11)
79 call cpu_convect_scalar_lx11(du, u, c, dx, dy, dz, &
80 xh_gll, coef_gll, gll_to_gl, nelv)
81 case (10)
82 call cpu_convect_scalar_lx10(du, u, c, dx, dy, dz, &
83 xh_gll, coef_gll, gll_to_gl, nelv)
84 case (9)
85 call cpu_convect_scalar_lx9(du, u, c, dx, dy, dz, &
86 xh_gll, coef_gll, gll_to_gl, nelv)
87 case (8)
88 call cpu_convect_scalar_lx8(du, u, c, dx, dy, dz, &
89 xh_gll, coef_gll, gll_to_gl, nelv)
90 case (7)
91 call cpu_convect_scalar_lx7(du, u, c, dx, dy, dz, &
92 xh_gll, coef_gll, gll_to_gl, nelv)
93 case (6)
94 call cpu_convect_scalar_lx6(du, u, c, dx, dy, dz, &
95 xh_gll, coef_gll, gll_to_gl, nelv)
96 case (5)
97 call cpu_convect_scalar_lx5(du, u, c, dx, dy, dz, &
98 xh_gll, coef_gll, gll_to_gl, nelv)
99 case (4)
100 call cpu_convect_scalar_lx4(du, u, c, dx, dy, dz, &
101 xh_gll, coef_gll, gll_to_gl, nelv)
102 case (3)
103 call cpu_convect_scalar_lx3(du, u, c, dx, dy, dz, &
104 xh_gll, coef_gll, gll_to_gl, nelv)
105 case (2)
106 call cpu_convect_scalar_lx2(du, u, c, dx, dy, dz, &
107 xh_gll, coef_gll, gll_to_gl, nelv)
108 case default
109 call cpu_convect_scalar_lx(du, u, c, dx, dy, dz, &
110 xh_gll, coef_gll, gll_to_gl, nelv, lx)
111 end select
112 end associate
113
114 end subroutine opr_cpu_convect_scalar
115
116 subroutine cpu_convect_scalar_lx(du, u, c, dx, dy, dz, Xh_GLL, &
117 coef_GLL, GLL_to_GL, nelv, lx)
118 integer, intent(in) :: nelv, lx
119 type(space_t), intent(in) :: Xh_GLL
120 type(coef_t), intent(in) :: coef_GLL
121 type(interpolator_t), intent(inout) :: GLL_to_GL
122 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
123 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
124 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
125 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
126 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
127 real(kind=rp) :: ud(lx*lx*lx)
128 real(kind=rp) :: tmp
129 integer :: e, i, j, k, l, idx, n_GLL
130
131 n_gll = nelv * xh_gll%lxyz
132
133 do e = 1, nelv
134 do j = 1, lx * lx
135 do i = 1, lx
136 tmp = 0.0_rp
137 do k = 1, lx
138 tmp = tmp + dx(i,k) * u(k,j,1,e)
139 end do
140 ur(i,j,1) = tmp
141 end do
142 end do
143
144 do k = 1, lx
145 do j = 1, lx
146 do i = 1, lx
147 tmp = 0.0_rp
148 do l = 1, lx
149 tmp = tmp + dy(j,l) * u(i,l,k,e)
150 end do
151 us(i,j,k) = tmp
152 end do
153 end do
154 end do
155
156 do k = 1, lx
157 do i = 1, lx*lx
158 tmp = 0.0_rp
159 do l = 1, lx
160 tmp = tmp + dz(k,l) * u(i,1,l,e)
161 end do
162 ut(i,1,k) = tmp
163 end do
164 end do
165
166 do i = 1, lx * lx * lx
167 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
168 + c(i,e,3) * ut(i,1,1)
169 end do
170 idx = (e-1) * xh_gll%lxyz+1
171 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
172 end do
173 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
174 call col2(du, coef_gll%Binv, n_gll)
175
176 end subroutine cpu_convect_scalar_lx
177
178 subroutine cpu_convect_scalar_lx18(du, u, c, dx, dy, dz, Xh_GLL, &
179 coef_GLL, GLL_to_GL, nelv)
180 integer, parameter :: lx = 18
181 integer, intent(in) :: nelv
182 type(space_t), intent(in) :: Xh_GLL
183 type(coef_t), intent(in) :: coef_GLL
184 type(interpolator_t), intent(inout) :: GLL_to_GL
185 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
186 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
187 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
188 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
189 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
190 real(kind=rp) :: ud(lx*lx*lx)
191 real(kind=rp) :: tmp
192 integer :: e, i, j, k, l, idx, n_GLL
193
194 n_gll = nelv * xh_gll%lxyz
195
196 do e = 1, nelv
197 do j = 1, lx * lx
198 do i = 1, lx
199 tmp = 0.0_rp
200 do k = 1, lx
201 tmp = tmp + dx(i,k) * u(k,j,1,e)
202 end do
203 ur(i,j,1) = tmp
204 end do
205 end do
206
207 do k = 1, lx
208 do j = 1, lx
209 do i = 1, lx
210 tmp = 0.0_rp
211 do l = 1, lx
212 tmp = tmp + dy(j,l) * u(i,l,k,e)
213 end do
214 us(i,j,k) = tmp
215 end do
216 end do
217 end do
218
219 do k = 1, lx
220 do i = 1, lx*lx
221 tmp = 0.0_rp
222 do l = 1, lx
223 tmp = tmp + dz(k,l) * u(i,1,l,e)
224 end do
225 ut(i,1,k) = tmp
226 end do
227 end do
228
229 do i = 1, lx * lx * lx
230 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
231 + c(i,e,3) * ut(i,1,1)
232 end do
233 idx = (e-1) * xh_gll%lxyz+1
234 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
235 end do
236 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
237 call col2(du, coef_gll%Binv, n_gll)
238
239 end subroutine cpu_convect_scalar_lx18
240
241 subroutine cpu_convect_scalar_lx17(du, u, c, dx, dy, dz, Xh_GLL, &
242 coef_GLL, GLL_to_GL, nelv)
243 integer, parameter :: lx = 17
244 integer, intent(in) :: nelv
245 type(space_t), intent(in) :: Xh_GLL
246 type(coef_t), intent(in) :: coef_GLL
247 type(interpolator_t), intent(inout) :: GLL_to_GL
248 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
249 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
250 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
251 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
252 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
253 real(kind=rp) :: ud(lx*lx*lx)
254 real(kind=rp) :: tmp
255 integer :: e, i, j, k, l, idx, n_GLL
256
257 n_gll = nelv * xh_gll%lxyz
258
259 do e = 1, nelv
260 do j = 1, lx * lx
261 do i = 1, lx
262 tmp = 0.0_rp
263 do k = 1, lx
264 tmp = tmp + dx(i,k) * u(k,j,1,e)
265 end do
266 ur(i,j,1) = tmp
267 end do
268 end do
269
270 do k = 1, lx
271 do j = 1, lx
272 do i = 1, lx
273 tmp = 0.0_rp
274 do l = 1, lx
275 tmp = tmp + dy(j,l) * u(i,l,k,e)
276 end do
277 us(i,j,k) = tmp
278 end do
279 end do
280 end do
281
282 do k = 1, lx
283 do i = 1, lx*lx
284 tmp = 0.0_rp
285 do l = 1, lx
286 tmp = tmp + dz(k,l) * u(i,1,l,e)
287 end do
288 ut(i,1,k) = tmp
289 end do
290 end do
291
292 do i = 1, lx * lx * lx
293 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
294 + c(i,e,3) * ut(i,1,1)
295 end do
296 idx = (e-1) * xh_gll%lxyz+1
297 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
298 end do
299 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
300 call col2(du, coef_gll%Binv, n_gll)
301
302 end subroutine cpu_convect_scalar_lx17
303
304 subroutine cpu_convect_scalar_lx16(du, u, c, dx, dy, dz, Xh_GLL, &
305 coef_GLL, GLL_to_GL, nelv)
306 integer, parameter :: lx = 16
307 integer, intent(in) :: nelv
308 type(space_t), intent(in) :: Xh_GLL
309 type(coef_t), intent(in) :: coef_GLL
310 type(interpolator_t), intent(inout) :: GLL_to_GL
311 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
312 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
313 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
314 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
315 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
316 real(kind=rp) :: ud(lx*lx*lx)
317 real(kind=rp) :: tmp
318 integer :: e, i, j, k, l, idx, n_GLL
319
320 n_gll = nelv * xh_gll%lxyz
321
322 do e = 1, nelv
323 do j = 1, lx * lx
324 do i = 1, lx
325 tmp = 0.0_rp
326 do k = 1, lx
327 tmp = tmp + dx(i,k) * u(k,j,1,e)
328 end do
329 ur(i,j,1) = tmp
330 end do
331 end do
332
333 do k = 1, lx
334 do j = 1, lx
335 do i = 1, lx
336 tmp = 0.0_rp
337 do l = 1, lx
338 tmp = tmp + dy(j,l) * u(i,l,k,e)
339 end do
340 us(i,j,k) = tmp
341 end do
342 end do
343 end do
344
345 do k = 1, lx
346 do i = 1, lx*lx
347 tmp = 0.0_rp
348 do l = 1, lx
349 tmp = tmp + dz(k,l) * u(i,1,l,e)
350 end do
351 ut(i,1,k) = tmp
352 end do
353 end do
354
355 do i = 1, lx * lx * lx
356 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
357 + c(i,e,3) * ut(i,1,1)
358 end do
359 idx = (e-1) * xh_gll%lxyz+1
360 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
361 end do
362 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
363 call col2(du, coef_gll%Binv, n_gll)
364
365 end subroutine cpu_convect_scalar_lx16
366
367 subroutine cpu_convect_scalar_lx15(du, u, c, dx, dy, dz, Xh_GLL, &
368 coef_GLL, GLL_to_GL, nelv)
369 integer, parameter :: lx = 15
370 integer, intent(in) :: nelv
371 type(space_t), intent(in) :: Xh_GLL
372 type(coef_t), intent(in) :: coef_GLL
373 type(interpolator_t), intent(inout) :: GLL_to_GL
374 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
375 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
376 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
377 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
378 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
379 real(kind=rp) :: ud(lx*lx*lx)
380 real(kind=rp) :: tmp
381 integer :: e, i, j, k, l, idx, n_GLL
382
383 n_gll = nelv * xh_gll%lxyz
384
385 do e = 1, nelv
386 do j = 1, lx * lx
387 do i = 1, lx
388 tmp = 0.0_rp
389 do k = 1, lx
390 tmp = tmp + dx(i,k) * u(k,j,1,e)
391 end do
392 ur(i,j,1) = tmp
393 end do
394 end do
395
396 do k = 1, lx
397 do j = 1, lx
398 do i = 1, lx
399 tmp = 0.0_rp
400 do l = 1, lx
401 tmp = tmp + dy(j,l) * u(i,l,k,e)
402 end do
403 us(i,j,k) = tmp
404 end do
405 end do
406 end do
407
408 do k = 1, lx
409 do i = 1, lx*lx
410 tmp = 0.0_rp
411 do l = 1, lx
412 tmp = tmp + dz(k,l) * u(i,1,l,e)
413 end do
414 ut(i,1,k) = tmp
415 end do
416 end do
417
418 do i = 1, lx * lx * lx
419 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
420 + c(i,e,3) * ut(i,1,1)
421 end do
422 idx = (e-1) * xh_gll%lxyz+1
423 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
424 end do
425 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
426 call col2(du, coef_gll%Binv, n_gll)
427
428 end subroutine cpu_convect_scalar_lx15
429
430 subroutine cpu_convect_scalar_lx14(du, u, c, dx, dy, dz, Xh_GLL, &
431 coef_GLL, GLL_to_GL, nelv)
432 integer, parameter :: lx = 14
433 integer, intent(in) :: nelv
434 type(space_t), intent(in) :: Xh_GLL
435 type(coef_t), intent(in) :: coef_GLL
436 type(interpolator_t), intent(inout) :: GLL_to_GL
437 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
438 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
439 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
440 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
441 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
442 real(kind=rp) :: ud(lx*lx*lx)
443 real(kind=rp) :: tmp
444 integer :: e, i, j, k, l, idx, n_GLL
445
446 n_gll = nelv * xh_gll%lxyz
447
448 do e = 1, nelv
449 do j = 1, lx * lx
450 do i = 1, lx
451 tmp = 0.0_rp
452 do k = 1, lx
453 tmp = tmp + dx(i,k) * u(k,j,1,e)
454 end do
455 ur(i,j,1) = tmp
456 end do
457 end do
458
459 do k = 1, lx
460 do j = 1, lx
461 do i = 1, lx
462 tmp = 0.0_rp
463 do l = 1, lx
464 tmp = tmp + dy(j,l) * u(i,l,k,e)
465 end do
466 us(i,j,k) = tmp
467 end do
468 end do
469 end do
470
471 do k = 1, lx
472 do i = 1, lx*lx
473 tmp = 0.0_rp
474 do l = 1, lx
475 tmp = tmp + dz(k,l) * u(i,1,l,e)
476 end do
477 ut(i,1,k) = tmp
478 end do
479 end do
480
481 do i = 1, lx * lx * lx
482 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
483 + c(i,e,3) * ut(i,1,1)
484 end do
485 idx = (e-1) * xh_gll%lxyz+1
486 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
487 end do
488 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
489 call col2(du, coef_gll%Binv, n_gll)
490
491 end subroutine cpu_convect_scalar_lx14
492
493 subroutine cpu_convect_scalar_lx13(du, u, c, dx, dy, dz, Xh_GLL, &
494 coef_GLL, GLL_to_GL, nelv)
495 integer, parameter :: lx = 13
496 integer, intent(in) :: nelv
497 type(space_t), intent(in) :: Xh_GLL
498 type(coef_t), intent(in) :: coef_GLL
499 type(interpolator_t), intent(inout) :: GLL_to_GL
500 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
501 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
502 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
503 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
504 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
505 real(kind=rp) :: ud(lx*lx*lx)
506 real(kind=rp) :: tmp
507 integer :: e, i, j, k, l, idx, n_GLL
508
509 n_gll = nelv * xh_gll%lxyz
510
511 do e = 1, nelv
512 do j = 1, lx * lx
513 do i = 1, lx
514 tmp = 0.0_rp
515 do k = 1, lx
516 tmp = tmp + dx(i,k) * u(k,j,1,e)
517 end do
518 ur(i,j,1) = tmp
519 end do
520 end do
521
522 do k = 1, lx
523 do j = 1, lx
524 do i = 1, lx
525 tmp = 0.0_rp
526 do l = 1, lx
527 tmp = tmp + dy(j,l) * u(i,l,k,e)
528 end do
529 us(i,j,k) = tmp
530 end do
531 end do
532 end do
533
534 do k = 1, lx
535 do i = 1, lx*lx
536 tmp = 0.0_rp
537 do l = 1, lx
538 tmp = tmp + dz(k,l) * u(i,1,l,e)
539 end do
540 ut(i,1,k) = tmp
541 end do
542 end do
543
544 do i = 1, lx * lx * lx
545 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
546 + c(i,e,3) * ut(i,1,1)
547 end do
548 idx = (e-1) * xh_gll%lxyz+1
549 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
550 end do
551 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
552 call col2(du, coef_gll%Binv, n_gll)
553
554 end subroutine cpu_convect_scalar_lx13
555
556 subroutine cpu_convect_scalar_lx12(du, u, c, dx, dy, dz, Xh_GLL, &
557 coef_GLL, GLL_to_GL, nelv)
558 integer, parameter :: lx = 12
559 integer, intent(in) :: nelv
560 type(space_t), intent(in) :: Xh_GLL
561 type(coef_t), intent(in) :: coef_GLL
562 type(interpolator_t), intent(inout) :: GLL_to_GL
563 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
564 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
565 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
566 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
567 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
568 real(kind=rp) :: ud(lx*lx*lx)
569 real(kind=rp) :: tmp
570 integer :: e, i, j, k, l, idx, n_GLL
571
572 n_gll = nelv * xh_gll%lxyz
573
574 do e = 1, nelv
575 do j = 1, lx * lx
576 do i = 1, lx
577 tmp = 0.0_rp
578 do k = 1, lx
579 tmp = tmp + dx(i,k) * u(k,j,1,e)
580 end do
581 ur(i,j,1) = tmp
582 end do
583 end do
584
585 do k = 1, lx
586 do j = 1, lx
587 do i = 1, lx
588 tmp = 0.0_rp
589 do l = 1, lx
590 tmp = tmp + dy(j,l) * u(i,l,k,e)
591 end do
592 us(i,j,k) = tmp
593 end do
594 end do
595 end do
596
597 do k = 1, lx
598 do i = 1, lx*lx
599 tmp = 0.0_rp
600 do l = 1, lx
601 tmp = tmp + dz(k,l) * u(i,1,l,e)
602 end do
603 ut(i,1,k) = tmp
604 end do
605 end do
606
607 do i = 1, lx * lx * lx
608 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
609 + c(i,e,3) * ut(i,1,1)
610 end do
611 idx = (e-1) * xh_gll%lxyz+1
612 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
613 end do
614 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
615 call col2(du, coef_gll%Binv, n_gll)
616
617 end subroutine cpu_convect_scalar_lx12
618
619 subroutine cpu_convect_scalar_lx11(du, u, c, dx, dy, dz, Xh_GLL, &
620 coef_GLL, GLL_to_GL, nelv)
621 integer, parameter :: lx = 11
622 integer, intent(in) :: nelv
623 type(space_t), intent(in) :: Xh_GLL
624 type(coef_t), intent(in) :: coef_GLL
625 type(interpolator_t), intent(inout) :: GLL_to_GL
626 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
627 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
628 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
629 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
630 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
631 real(kind=rp) :: ud(lx*lx*lx)
632 real(kind=rp) :: tmp
633 integer :: e, i, j, k, l, idx, n_GLL
634
635 n_gll = nelv * xh_gll%lxyz
636
637 do e = 1, nelv
638 do j = 1, lx * lx
639 do i = 1, lx
640 tmp = 0.0_rp
641 do k = 1, lx
642 tmp = tmp + dx(i,k) * u(k,j,1,e)
643 end do
644 ur(i,j,1) = tmp
645 end do
646 end do
647
648 do k = 1, lx
649 do j = 1, lx
650 do i = 1, lx
651 tmp = 0.0_rp
652 do l = 1, lx
653 tmp = tmp + dy(j,l) * u(i,l,k,e)
654 end do
655 us(i,j,k) = tmp
656 end do
657 end do
658 end do
659
660 do k = 1, lx
661 do i = 1, lx*lx
662 tmp = 0.0_rp
663 do l = 1, lx
664 tmp = tmp + dz(k,l) * u(i,1,l,e)
665 end do
666 ut(i,1,k) = tmp
667 end do
668 end do
669
670 do i = 1, lx * lx * lx
671 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
672 + c(i,e,3) * ut(i,1,1)
673 end do
674 idx = (e-1) * xh_gll%lxyz+1
675 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
676 end do
677 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
678 call col2(du, coef_gll%Binv, n_gll)
679
680 end subroutine cpu_convect_scalar_lx11
681
682 subroutine cpu_convect_scalar_lx10(du, u, c, dx, dy, dz, Xh_GLL, &
683 coef_GLL, GLL_to_GL, nelv)
684 integer, parameter :: lx = 10
685 integer, intent(in) :: nelv
686 type(space_t), intent(in) :: Xh_GLL
687 type(coef_t), intent(in) :: coef_GLL
688 type(interpolator_t), intent(inout) :: GLL_to_GL
689 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
690 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
691 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
692 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
693 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
694 real(kind=rp) :: ud(lx*lx*lx)
695 real(kind=rp) :: tmp
696 integer :: e, i, j, k, l, idx, n_GLL
697
698 n_gll = nelv * xh_gll%lxyz
699
700 do e = 1, nelv
701 do j = 1, lx * lx
702 do i = 1, lx
703 tmp = 0.0_rp
704 do k = 1, lx
705 tmp = tmp + dx(i,k) * u(k,j,1,e)
706 end do
707 ur(i,j,1) = tmp
708 end do
709 end do
710
711 do k = 1, lx
712 do j = 1, lx
713 do i = 1, lx
714 tmp = 0.0_rp
715 do l = 1, lx
716 tmp = tmp + dy(j,l) * u(i,l,k,e)
717 end do
718 us(i,j,k) = tmp
719 end do
720 end do
721 end do
722
723 do k = 1, lx
724 do i = 1, lx*lx
725 tmp = 0.0_rp
726 do l = 1, lx
727 tmp = tmp + dz(k,l) * u(i,1,l,e)
728 end do
729 ut(i,1,k) = tmp
730 end do
731 end do
732
733 do i = 1, lx * lx * lx
734 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
735 + c(i,e,3) * ut(i,1,1)
736 end do
737 idx = (e-1) * xh_gll%lxyz+1
738 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
739 end do
740 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
741 call col2(du, coef_gll%Binv, n_gll)
742
743 end subroutine cpu_convect_scalar_lx10
744
745 subroutine cpu_convect_scalar_lx9(du, u, c, dx, dy, dz, Xh_GLL, &
746 coef_GLL, GLL_to_GL, nelv)
747 integer, parameter :: lx = 9
748 integer, intent(in) :: nelv
749 type(space_t), intent(in) :: Xh_GLL
750 type(coef_t), intent(in) :: coef_GLL
751 type(interpolator_t), intent(inout) :: GLL_to_GL
752 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
753 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
754 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
755 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
756 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
757 real(kind=rp) :: ud(lx*lx*lx)
758 real(kind=rp) :: tmp
759 integer :: e, i, j, k, l, idx, n_GLL
760
761 n_gll = nelv * xh_gll%lxyz
762
763 do e = 1, nelv
764 do j = 1, lx * lx
765 do i = 1, lx
766 tmp = 0.0_rp
767 do k = 1, lx
768 tmp = tmp + dx(i,k) * u(k,j,1,e)
769 end do
770 ur(i,j,1) = tmp
771 end do
772 end do
773
774 do k = 1, lx
775 do j = 1, lx
776 do i = 1, lx
777 tmp = 0.0_rp
778 do l = 1, lx
779 tmp = tmp + dy(j,l) * u(i,l,k,e)
780 end do
781 us(i,j,k) = tmp
782 end do
783 end do
784 end do
785
786 do k = 1, lx
787 do i = 1, lx*lx
788 tmp = 0.0_rp
789 do l = 1, lx
790 tmp = tmp + dz(k,l) * u(i,1,l,e)
791 end do
792 ut(i,1,k) = tmp
793 end do
794 end do
795
796 do i = 1, lx * lx * lx
797 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
798 + c(i,e,3) * ut(i,1,1)
799 end do
800 idx = (e-1) * xh_gll%lxyz+1
801 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
802 end do
803 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
804 call col2(du, coef_gll%Binv, n_gll)
805
806 end subroutine cpu_convect_scalar_lx9
807
808 subroutine cpu_convect_scalar_lx8(du, u, c, dx, dy, dz, Xh_GLL, &
809 coef_GLL, GLL_to_GL, nelv)
810 integer, parameter :: lx = 8
811 integer, intent(in) :: nelv
812 type(space_t), intent(in) :: Xh_GLL
813 type(coef_t), intent(in) :: coef_GLL
814 type(interpolator_t), intent(inout) :: GLL_to_GL
815 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
816 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
817 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
818 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
819 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
820 real(kind=rp) :: ud(lx*lx*lx)
821 real(kind=rp) :: tmp
822 integer :: e, i, j, k, l, idx, n_GLL
823
824 n_gll = nelv * xh_gll%lxyz
825
826 do e = 1, nelv
827 do j = 1, lx * lx
828 do i = 1, lx
829 tmp = 0.0_rp
830 do k = 1, lx
831 tmp = tmp + dx(i,k) * u(k,j,1,e)
832 end do
833 ur(i,j,1) = tmp
834 end do
835 end do
836
837 do k = 1, lx
838 do j = 1, lx
839 do i = 1, lx
840 tmp = 0.0_rp
841 do l = 1, lx
842 tmp = tmp + dy(j,l) * u(i,l,k,e)
843 end do
844 us(i,j,k) = tmp
845 end do
846 end do
847 end do
848
849 do k = 1, lx
850 do i = 1, lx*lx
851 tmp = 0.0_rp
852 do l = 1, lx
853 tmp = tmp + dz(k,l) * u(i,1,l,e)
854 end do
855 ut(i,1,k) = tmp
856 end do
857 end do
858
859 do i = 1, lx * lx * lx
860 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
861 + c(i,e,3) * ut(i,1,1)
862 end do
863 idx = (e-1) * xh_gll%lxyz+1
864 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
865 end do
866 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
867 call col2(du, coef_gll%Binv, n_gll)
868
869 end subroutine cpu_convect_scalar_lx8
870
871 subroutine cpu_convect_scalar_lx7(du, u, c, dx, dy, dz, Xh_GLL, &
872 coef_GLL, GLL_to_GL, nelv)
873 integer, parameter :: lx = 7
874 integer, intent(in) :: nelv
875 type(space_t), intent(in) :: Xh_GLL
876 type(coef_t), intent(in) :: coef_GLL
877 type(interpolator_t), intent(inout) :: GLL_to_GL
878 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
879 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
880 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
881 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
882 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
883 real(kind=rp) :: ud(lx*lx*lx)
884 real(kind=rp) :: tmp
885 integer :: e, i, j, k, l, idx, n_GLL
886
887 n_gll = nelv * xh_gll%lxyz
888
889 do e = 1, nelv
890 do j = 1, lx * lx
891 do i = 1, lx
892 tmp = 0.0_rp
893 do k = 1, lx
894 tmp = tmp + dx(i,k) * u(k,j,1,e)
895 end do
896 ur(i,j,1) = tmp
897 end do
898 end do
899
900 do k = 1, lx
901 do j = 1, lx
902 do i = 1, lx
903 tmp = 0.0_rp
904 do l = 1, lx
905 tmp = tmp + dy(j,l) * u(i,l,k,e)
906 end do
907 us(i,j,k) = tmp
908 end do
909 end do
910 end do
911
912 do k = 1, lx
913 do i = 1, lx*lx
914 tmp = 0.0_rp
915 do l = 1, lx
916 tmp = tmp + dz(k,l) * u(i,1,l,e)
917 end do
918 ut(i,1,k) = tmp
919 end do
920 end do
921
922 do i = 1, lx * lx * lx
923 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
924 + c(i,e,3) * ut(i,1,1)
925 end do
926 idx = (e-1) * xh_gll%lxyz+1
927 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
928 end do
929 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
930 call col2(du, coef_gll%Binv, n_gll)
931
932 end subroutine cpu_convect_scalar_lx7
933
934 subroutine cpu_convect_scalar_lx6(du, u, c, dx, dy, dz, Xh_GLL, &
935 coef_GLL, GLL_to_GL, nelv)
936 integer, parameter :: lx = 6
937 integer, intent(in) :: nelv
938 type(space_t), intent(in) :: Xh_GLL
939 type(coef_t), intent(in) :: coef_GLL
940 type(interpolator_t), intent(inout) :: GLL_to_GL
941 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
942 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
943 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
944 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
945 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
946 real(kind=rp) :: ud(lx*lx*lx)
947 real(kind=rp) :: tmp
948 integer :: e, i, j, k, l, idx, n_GLL
949
950 n_gll = nelv * xh_gll%lxyz
951
952 do e = 1, nelv
953 do j = 1, lx * lx
954 do i = 1, lx
955 tmp = 0.0_rp
956 do k = 1, lx
957 tmp = tmp + dx(i,k) * u(k,j,1,e)
958 end do
959 ur(i,j,1) = tmp
960 end do
961 end do
962
963 do k = 1, lx
964 do j = 1, lx
965 do i = 1, lx
966 tmp = 0.0_rp
967 do l = 1, lx
968 tmp = tmp + dy(j,l) * u(i,l,k,e)
969 end do
970 us(i,j,k) = tmp
971 end do
972 end do
973 end do
974
975 do k = 1, lx
976 do i = 1, lx*lx
977 tmp = 0.0_rp
978 do l = 1, lx
979 tmp = tmp + dz(k,l) * u(i,1,l,e)
980 end do
981 ut(i,1,k) = tmp
982 end do
983 end do
984
985 do i = 1, lx * lx * lx
986 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
987 + c(i,e,3) * ut(i,1,1)
988 end do
989 idx = (e-1) * xh_gll%lxyz+1
990 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
991 end do
992 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
993 call col2(du, coef_gll%Binv, n_gll)
994
995 end subroutine cpu_convect_scalar_lx6
996
997 subroutine cpu_convect_scalar_lx5(du, u, c, dx, dy, dz, Xh_GLL, &
998 coef_GLL, GLL_to_GL, nelv)
999 integer, parameter :: lx = 5
1000 integer, intent(in) :: nelv
1001 type(space_t), intent(in) :: Xh_GLL
1002 type(coef_t), intent(in) :: coef_GLL
1003 type(interpolator_t), intent(inout) :: GLL_to_GL
1004 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1005 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
1006 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
1007 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1008 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
1009 real(kind=rp) :: ud(lx*lx*lx)
1010 real(kind=rp) :: tmp
1011 integer :: e, i, j, k, l, idx, n_GLL
1012
1013 n_gll = nelv * xh_gll%lxyz
1014
1015 do e = 1, nelv
1016 do j = 1, lx * lx
1017 do i = 1, lx
1018 tmp = 0.0_rp
1019 do k = 1, lx
1020 tmp = tmp + dx(i,k) * u(k,j,1,e)
1021 end do
1022 ur(i,j,1) = tmp
1023 end do
1024 end do
1025
1026 do k = 1, lx
1027 do j = 1, lx
1028 do i = 1, lx
1029 tmp = 0.0_rp
1030 do l = 1, lx
1031 tmp = tmp + dy(j,l) * u(i,l,k,e)
1032 end do
1033 us(i,j,k) = tmp
1034 end do
1035 end do
1036 end do
1037
1038 do k = 1, lx
1039 do i = 1, lx*lx
1040 tmp = 0.0_rp
1041 do l = 1, lx
1042 tmp = tmp + dz(k,l) * u(i,1,l,e)
1043 end do
1044 ut(i,1,k) = tmp
1045 end do
1046 end do
1047
1048 do i = 1, lx * lx * lx
1049 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
1050 + c(i,e,3) * ut(i,1,1)
1051 end do
1052 idx = (e-1) * xh_gll%lxyz+1
1053 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1054 end do
1055 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1056 call col2(du, coef_gll%Binv, n_gll)
1057
1058 end subroutine cpu_convect_scalar_lx5
1059
1060 subroutine cpu_convect_scalar_lx4(du, u, c, dx, dy, dz, Xh_GLL, &
1061 coef_GLL, GLL_to_GL, nelv)
1062 integer, parameter :: lx = 4
1063 integer, intent(in) :: nelv
1064 type(space_t), intent(in) :: Xh_GLL
1065 type(coef_t), intent(in) :: coef_GLL
1066 type(interpolator_t), intent(inout) :: GLL_to_GL
1067 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1068 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
1069 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
1070 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1071 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
1072 real(kind=rp) :: ud(lx*lx*lx)
1073 real(kind=rp) :: tmp
1074 integer :: e, i, j, k, l, idx, n_GLL
1075
1076 n_gll = nelv * xh_gll%lxyz
1077
1078 do e = 1, nelv
1079 do j = 1, lx * lx
1080 do i = 1, lx
1081 tmp = 0.0_rp
1082 do k = 1, lx
1083 tmp = tmp + dx(i,k) * u(k,j,1,e)
1084 end do
1085 ur(i,j,1) = tmp
1086 end do
1087 end do
1088
1089 do k = 1, lx
1090 do j = 1, lx
1091 do i = 1, lx
1092 tmp = 0.0_rp
1093 do l = 1, lx
1094 tmp = tmp + dy(j,l) * u(i,l,k,e)
1095 end do
1096 us(i,j,k) = tmp
1097 end do
1098 end do
1099 end do
1100
1101 do k = 1, lx
1102 do i = 1, lx*lx
1103 tmp = 0.0_rp
1104 do l = 1, lx
1105 tmp = tmp + dz(k,l) * u(i,1,l,e)
1106 end do
1107 ut(i,1,k) = tmp
1108 end do
1109 end do
1110
1111 do i = 1, lx * lx * lx
1112 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
1113 + c(i,e,3) * ut(i,1,1)
1114 end do
1115 idx = (e-1) * xh_gll%lxyz+1
1116 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1117 end do
1118 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1119 call col2(du, coef_gll%Binv, n_gll)
1120
1121 end subroutine cpu_convect_scalar_lx4
1122
1123 subroutine cpu_convect_scalar_lx3(du, u, c, dx, dy, dz, Xh_GLL, &
1124 coef_GLL, GLL_to_GL, nelv)
1125 integer, parameter :: lx = 3
1126 integer, intent(in) :: nelv
1127 type(space_t), intent(in) :: Xh_GLL
1128 type(coef_t), intent(in) :: coef_GLL
1129 type(interpolator_t), intent(inout) :: GLL_to_GL
1130 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1131 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
1132 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
1133 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1134 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
1135 real(kind=rp) :: ud(lx*lx*lx)
1136 real(kind=rp) :: tmp
1137 integer :: e, i, j, k, l, idx, n_GLL
1138
1139 n_gll = nelv * xh_gll%lxyz
1140
1141 do e = 1, nelv
1142 do j = 1, lx * lx
1143 do i = 1, lx
1144 tmp = 0.0_rp
1145 do k = 1, lx
1146 tmp = tmp + dx(i,k) * u(k,j,1,e)
1147 end do
1148 ur(i,j,1) = tmp
1149 end do
1150 end do
1151
1152 do k = 1, lx
1153 do j = 1, lx
1154 do i = 1, lx
1155 tmp = 0.0_rp
1156 do l = 1, lx
1157 tmp = tmp + dy(j,l) * u(i,l,k,e)
1158 end do
1159 us(i,j,k) = tmp
1160 end do
1161 end do
1162 end do
1163
1164 do k = 1, lx
1165 do i = 1, lx*lx
1166 tmp = 0.0_rp
1167 do l = 1, lx
1168 tmp = tmp + dz(k,l) * u(i,1,l,e)
1169 end do
1170 ut(i,1,k) = tmp
1171 end do
1172 end do
1173
1174 do i = 1, lx * lx * lx
1175 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
1176 + c(i,e,3) * ut(i,1,1)
1177 end do
1178 idx = (e-1) * xh_gll%lxyz+1
1179 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1180 end do
1181 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1182 call col2(du, coef_gll%Binv, n_gll)
1183
1184 end subroutine cpu_convect_scalar_lx3
1185
1186 subroutine cpu_convect_scalar_lx2(du, u, c, dx, dy, dz, Xh_GLL, &
1187 coef_GLL, GLL_to_GL, nelv)
1188 integer, parameter :: lx = 2
1189 integer, intent(in) :: nelv
1190 type(space_t), intent(in) :: Xh_GLL
1191 type(coef_t), intent(in) :: coef_GLL
1192 type(interpolator_t), intent(inout) :: GLL_to_GL
1193 real(kind=rp), intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1194 real(kind=rp), intent(in) :: u(lx, lx, lx, nelv)
1195 real(kind=rp), intent(in) :: c(lx*lx*lx, nelv, 3)
1196 real(kind=rp), dimension(lx, lx), intent(in) :: dx, dy, dz
1197 real(kind=rp), dimension(lx, lx, lx) :: ur, us, ut
1198 real(kind=rp) :: ud(lx*lx*lx)
1199 real(kind=rp) :: tmp
1200 integer :: e, i, j, k, l, idx, n_GLL
1201
1202 n_gll = nelv * xh_gll%lxyz
1203
1204 do e = 1, nelv
1205 do j = 1, lx * lx
1206 do i = 1, lx
1207 tmp = 0.0_rp
1208 do k = 1, lx
1209 tmp = tmp + dx(i,k) * u(k,j,1,e)
1210 end do
1211 ur(i,j,1) = tmp
1212 end do
1213 end do
1214
1215 do k = 1, lx
1216 do j = 1, lx
1217 do i = 1, lx
1218 tmp = 0.0_rp
1219 do l = 1, lx
1220 tmp = tmp + dy(j,l) * u(i,l,k,e)
1221 end do
1222 us(i,j,k) = tmp
1223 end do
1224 end do
1225 end do
1226
1227 do k = 1, lx
1228 do i = 1, lx*lx
1229 tmp = 0.0_rp
1230 do l = 1, lx
1231 tmp = tmp + dz(k,l) * u(i,1,l,e)
1232 end do
1233 ut(i,1,k) = tmp
1234 end do
1235 end do
1236
1237 do i = 1, lx * lx * lx
1238 ud(i) = c(i,e,1) * ur(i,1,1) + c(i,e,2) * us(i,1,1) &
1239 + c(i,e,3) * ut(i,1,1)
1240 end do
1241 idx = (e-1) * xh_gll%lxyz+1
1242 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1243 end do
1244 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1245 call col2(du, coef_gll%Binv, n_gll)
1246
1247 end subroutine cpu_convect_scalar_lx2
1248
1249end submodule cpu_convect_scalar
Definition math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:728
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators CPU backend.
Definition opr_cpu.f90:34