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