34 submodule(
opr_cpu) cpu_convect_scalar
41 module subroutine opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, coef_gll, &
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)
58 call cpu_convect_scalar_lx18(du, u, c, dx, dy, dz, &
59 xh_gll, coef_gll, gll_to_gl, nelv)
61 call cpu_convect_scalar_lx17(du, u, c, dx, dy, dz, &
62 xh_gll, coef_gll, gll_to_gl, nelv)
64 call cpu_convect_scalar_lx16(du, u, c, dx, dy, dz, &
65 xh_gll, coef_gll, gll_to_gl, nelv)
67 call cpu_convect_scalar_lx15(du, u, c, dx, dy, dz, &
68 xh_gll, coef_gll, gll_to_gl, nelv)
70 call cpu_convect_scalar_lx14(du, u, c, dx, dy, dz, &
71 xh_gll, coef_gll, gll_to_gl, nelv)
73 call cpu_convect_scalar_lx13(du, u, c, dx, dy, dz, &
74 xh_gll, coef_gll, gll_to_gl, nelv)
76 call cpu_convect_scalar_lx12(du, u, c, dx, dy, dz, &
77 xh_gll, coef_gll, gll_to_gl, nelv)
79 call cpu_convect_scalar_lx11(du, u, c, dx, dy, dz, &
80 xh_gll, coef_gll, gll_to_gl, nelv)
82 call cpu_convect_scalar_lx10(du, u, c, dx, dy, dz, &
83 xh_gll, coef_gll, gll_to_gl, nelv)
85 call cpu_convect_scalar_lx9(du, u, c, dx, dy, dz, &
86 xh_gll, coef_gll, gll_to_gl, nelv)
88 call cpu_convect_scalar_lx8(du, u, c, dx, dy, dz, &
89 xh_gll, coef_gll, gll_to_gl, nelv)
91 call cpu_convect_scalar_lx7(du, u, c, dx, dy, dz, &
92 xh_gll, coef_gll, gll_to_gl, nelv)
94 call cpu_convect_scalar_lx6(du, u, c, dx, dy, dz, &
95 xh_gll, coef_gll, gll_to_gl, nelv)
97 call cpu_convect_scalar_lx5(du, u, c, dx, dy, dz, &
98 xh_gll, coef_gll, gll_to_gl, nelv)
100 call cpu_convect_scalar_lx4(du, u, c, dx, dy, dz, &
101 xh_gll, coef_gll, gll_to_gl, nelv)
103 call cpu_convect_scalar_lx3(du, u, c, dx, dy, dz, &
104 xh_gll, coef_gll, gll_to_gl, nelv)
106 call cpu_convect_scalar_lx2(du, u, c, dx, dy, dz, &
107 xh_gll, coef_gll, gll_to_gl, nelv)
109 call cpu_convect_scalar_lx(du, u, c, dx, dy, dz, &
110 xh_gll, coef_gll, gll_to_gl, nelv, lx)
114 end subroutine opr_cpu_convect_scalar
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)
129 integer :: e, i, j, k, l, idx, n_GLL
131 n_gll = nelv * xh_gll%lxyz
138 tmp = tmp + dx(i,k) * u(k,j,1,e)
149 tmp = tmp + dy(j,l) * u(i,l,k,e)
160 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
170 idx = (e-1) * xh_gll%lxyz+1
171 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
173 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
174 call col2(du, coef_gll%Binv, n_gll)
176 end subroutine cpu_convect_scalar_lx
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)
192 integer :: e, i, j, k, l, idx, n_GLL
194 n_gll = nelv * xh_gll%lxyz
201 tmp = tmp + dx(i,k) * u(k,j,1,e)
212 tmp = tmp + dy(j,l) * u(i,l,k,e)
223 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
233 idx = (e-1) * xh_gll%lxyz+1
234 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
236 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
237 call col2(du, coef_gll%Binv, n_gll)
239 end subroutine cpu_convect_scalar_lx18
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)
255 integer :: e, i, j, k, l, idx, n_GLL
257 n_gll = nelv * xh_gll%lxyz
264 tmp = tmp + dx(i,k) * u(k,j,1,e)
275 tmp = tmp + dy(j,l) * u(i,l,k,e)
286 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
296 idx = (e-1) * xh_gll%lxyz+1
297 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
299 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
300 call col2(du, coef_gll%Binv, n_gll)
302 end subroutine cpu_convect_scalar_lx17
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)
318 integer :: e, i, j, k, l, idx, n_GLL
320 n_gll = nelv * xh_gll%lxyz
327 tmp = tmp + dx(i,k) * u(k,j,1,e)
338 tmp = tmp + dy(j,l) * u(i,l,k,e)
349 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
359 idx = (e-1) * xh_gll%lxyz+1
360 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
362 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
363 call col2(du, coef_gll%Binv, n_gll)
365 end subroutine cpu_convect_scalar_lx16
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)
381 integer :: e, i, j, k, l, idx, n_GLL
383 n_gll = nelv * xh_gll%lxyz
390 tmp = tmp + dx(i,k) * u(k,j,1,e)
401 tmp = tmp + dy(j,l) * u(i,l,k,e)
412 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
422 idx = (e-1) * xh_gll%lxyz+1
423 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
425 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
426 call col2(du, coef_gll%Binv, n_gll)
428 end subroutine cpu_convect_scalar_lx15
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)
444 integer :: e, i, j, k, l, idx, n_GLL
446 n_gll = nelv * xh_gll%lxyz
453 tmp = tmp + dx(i,k) * u(k,j,1,e)
464 tmp = tmp + dy(j,l) * u(i,l,k,e)
475 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
485 idx = (e-1) * xh_gll%lxyz+1
486 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
488 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
489 call col2(du, coef_gll%Binv, n_gll)
491 end subroutine cpu_convect_scalar_lx14
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)
507 integer :: e, i, j, k, l, idx, n_GLL
509 n_gll = nelv * xh_gll%lxyz
516 tmp = tmp + dx(i,k) * u(k,j,1,e)
527 tmp = tmp + dy(j,l) * u(i,l,k,e)
538 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
548 idx = (e-1) * xh_gll%lxyz+1
549 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
551 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
552 call col2(du, coef_gll%Binv, n_gll)
554 end subroutine cpu_convect_scalar_lx13
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)
570 integer :: e, i, j, k, l, idx, n_GLL
572 n_gll = nelv * xh_gll%lxyz
579 tmp = tmp + dx(i,k) * u(k,j,1,e)
590 tmp = tmp + dy(j,l) * u(i,l,k,e)
601 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
611 idx = (e-1) * xh_gll%lxyz+1
612 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
614 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
615 call col2(du, coef_gll%Binv, n_gll)
617 end subroutine cpu_convect_scalar_lx12
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)
633 integer :: e, i, j, k, l, idx, n_GLL
635 n_gll = nelv * xh_gll%lxyz
642 tmp = tmp + dx(i,k) * u(k,j,1,e)
653 tmp = tmp + dy(j,l) * u(i,l,k,e)
664 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
674 idx = (e-1) * xh_gll%lxyz+1
675 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
677 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
678 call col2(du, coef_gll%Binv, n_gll)
680 end subroutine cpu_convect_scalar_lx11
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)
696 integer :: e, i, j, k, l, idx, n_GLL
698 n_gll = nelv * xh_gll%lxyz
705 tmp = tmp + dx(i,k) * u(k,j,1,e)
716 tmp = tmp + dy(j,l) * u(i,l,k,e)
727 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
737 idx = (e-1) * xh_gll%lxyz+1
738 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
740 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
741 call col2(du, coef_gll%Binv, n_gll)
743 end subroutine cpu_convect_scalar_lx10
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)
759 integer :: e, i, j, k, l, idx, n_GLL
761 n_gll = nelv * xh_gll%lxyz
768 tmp = tmp + dx(i,k) * u(k,j,1,e)
779 tmp = tmp + dy(j,l) * u(i,l,k,e)
790 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
800 idx = (e-1) * xh_gll%lxyz+1
801 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
803 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
804 call col2(du, coef_gll%Binv, n_gll)
806 end subroutine cpu_convect_scalar_lx9
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)
822 integer :: e, i, j, k, l, idx, n_GLL
824 n_gll = nelv * xh_gll%lxyz
831 tmp = tmp + dx(i,k) * u(k,j,1,e)
842 tmp = tmp + dy(j,l) * u(i,l,k,e)
853 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
863 idx = (e-1) * xh_gll%lxyz+1
864 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
866 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
867 call col2(du, coef_gll%Binv, n_gll)
869 end subroutine cpu_convect_scalar_lx8
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)
885 integer :: e, i, j, k, l, idx, n_GLL
887 n_gll = nelv * xh_gll%lxyz
894 tmp = tmp + dx(i,k) * u(k,j,1,e)
905 tmp = tmp + dy(j,l) * u(i,l,k,e)
916 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
926 idx = (e-1) * xh_gll%lxyz+1
927 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
929 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
930 call col2(du, coef_gll%Binv, n_gll)
932 end subroutine cpu_convect_scalar_lx7
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)
948 integer :: e, i, j, k, l, idx, n_GLL
950 n_gll = nelv * xh_gll%lxyz
957 tmp = tmp + dx(i,k) * u(k,j,1,e)
968 tmp = tmp + dy(j,l) * u(i,l,k,e)
979 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
989 idx = (e-1) * xh_gll%lxyz+1
990 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
992 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
993 call col2(du, coef_gll%Binv, n_gll)
995 end subroutine cpu_convect_scalar_lx6
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
1013 n_gll = nelv * xh_gll%lxyz
1020 tmp = tmp + dx(i,k) * u(k,j,1,e)
1031 tmp = tmp + dy(j,l) * u(i,l,k,e)
1042 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1052 idx = (e-1) * xh_gll%lxyz+1
1053 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1055 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1056 call col2(du, coef_gll%Binv, n_gll)
1058 end subroutine cpu_convect_scalar_lx5
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
1076 n_gll = nelv * xh_gll%lxyz
1083 tmp = tmp + dx(i,k) * u(k,j,1,e)
1094 tmp = tmp + dy(j,l) * u(i,l,k,e)
1105 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1115 idx = (e-1) * xh_gll%lxyz+1
1116 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1118 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1119 call col2(du, coef_gll%Binv, n_gll)
1121 end subroutine cpu_convect_scalar_lx4
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
1139 n_gll = nelv * xh_gll%lxyz
1146 tmp = tmp + dx(i,k) * u(k,j,1,e)
1157 tmp = tmp + dy(j,l) * u(i,l,k,e)
1168 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1178 idx = (e-1) * xh_gll%lxyz+1
1179 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1181 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1182 call col2(du, coef_gll%Binv, n_gll)
1184 end subroutine cpu_convect_scalar_lx3
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
1202 n_gll = nelv * xh_gll%lxyz
1209 tmp = tmp + dx(i,k) * u(k,j,1,e)
1220 tmp = tmp + dy(j,l) * u(i,l,k,e)
1231 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1241 idx = (e-1) * xh_gll%lxyz+1
1242 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1244 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1245 call col2(du, coef_gll%Binv, n_gll)
1247 end subroutine cpu_convect_scalar_lx2
1249 end submodule cpu_convect_scalar
subroutine, public col2(a, b, n)
Vector multiplication .
integer, parameter, public rp
Global precision used in computations.