34submodule(
opr_cpu) cpu_convect_scalar
40 module subroutine opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, coef_gll, &
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)
57 call cpu_convect_scalar_lx18(du, u, c, dx, dy, dz, &
58 xh_gll, coef_gll, gll_to_gl, nelv)
60 call cpu_convect_scalar_lx17(du, u, c, dx, dy, dz, &
61 xh_gll, coef_gll, gll_to_gl, nelv)
63 call cpu_convect_scalar_lx16(du, u, c, dx, dy, dz, &
64 xh_gll, coef_gll, gll_to_gl, nelv)
66 call cpu_convect_scalar_lx15(du, u, c, dx, dy, dz, &
67 xh_gll, coef_gll, gll_to_gl, nelv)
69 call cpu_convect_scalar_lx14(du, u, c, dx, dy, dz, &
70 xh_gll, coef_gll, gll_to_gl, nelv)
72 call cpu_convect_scalar_lx13(du, u, c, dx, dy, dz, &
73 xh_gll, coef_gll, gll_to_gl, nelv)
75 call cpu_convect_scalar_lx12(du, u, c, dx, dy, dz, &
76 xh_gll, coef_gll, gll_to_gl, nelv)
78 call cpu_convect_scalar_lx11(du, u, c, dx, dy, dz, &
79 xh_gll, coef_gll, gll_to_gl, nelv)
81 call cpu_convect_scalar_lx10(du, u, c, dx, dy, dz, &
82 xh_gll, coef_gll, gll_to_gl, nelv)
84 call cpu_convect_scalar_lx9(du, u, c, dx, dy, dz, &
85 xh_gll, coef_gll, gll_to_gl, nelv)
87 call cpu_convect_scalar_lx8(du, u, c, dx, dy, dz, &
88 xh_gll, coef_gll, gll_to_gl, nelv)
90 call cpu_convect_scalar_lx7(du, u, c, dx, dy, dz, &
91 xh_gll, coef_gll, gll_to_gl, nelv)
93 call cpu_convect_scalar_lx6(du, u, c, dx, dy, dz, &
94 xh_gll, coef_gll, gll_to_gl, nelv)
96 call cpu_convect_scalar_lx5(du, u, c, dx, dy, dz, &
97 xh_gll, coef_gll, gll_to_gl, nelv)
99 call cpu_convect_scalar_lx4(du, u, c, dx, dy, dz, &
100 xh_gll, coef_gll, gll_to_gl, nelv)
102 call cpu_convect_scalar_lx3(du, u, c, dx, dy, dz, &
103 xh_gll, coef_gll, gll_to_gl, nelv)
105 call cpu_convect_scalar_lx2(du, u, c, dx, dy, dz, &
106 xh_gll, coef_gll, gll_to_gl, nelv)
108 call cpu_convect_scalar_lx(du, u, c, dx, dy, dz, &
109 xh_gll, coef_gll, gll_to_gl, nelv, lx)
113 end subroutine opr_cpu_convect_scalar
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)
128 integer :: e, i, j, k, l, idx, n_GLL
130 n_gll = nelv * xh_gll%lxyz
137 tmp = tmp + dx(i,k) * u(k,j,1,e)
148 tmp = tmp + dy(j,l) * u(i,l,k,e)
159 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
169 idx = (e-1) * xh_gll%lxyz+1
170 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
172 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
173 call col2(du, coef_gll%Binv, n_gll)
175 end subroutine cpu_convect_scalar_lx
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)
191 integer :: e, i, j, k, l, idx, n_GLL
193 n_gll = nelv * xh_gll%lxyz
200 tmp = tmp + dx(i,k) * u(k,j,1,e)
211 tmp = tmp + dy(j,l) * u(i,l,k,e)
222 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
232 idx = (e-1) * xh_gll%lxyz+1
233 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
235 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
236 call col2(du, coef_gll%Binv, n_gll)
238 end subroutine cpu_convect_scalar_lx18
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)
254 integer :: e, i, j, k, l, idx, n_GLL
256 n_gll = nelv * xh_gll%lxyz
263 tmp = tmp + dx(i,k) * u(k,j,1,e)
274 tmp = tmp + dy(j,l) * u(i,l,k,e)
285 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
295 idx = (e-1) * xh_gll%lxyz+1
296 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
298 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
299 call col2(du, coef_gll%Binv, n_gll)
301 end subroutine cpu_convect_scalar_lx17
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)
317 integer :: e, i, j, k, l, idx, n_GLL
319 n_gll = nelv * xh_gll%lxyz
326 tmp = tmp + dx(i,k) * u(k,j,1,e)
337 tmp = tmp + dy(j,l) * u(i,l,k,e)
348 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
358 idx = (e-1) * xh_gll%lxyz+1
359 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
361 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
362 call col2(du, coef_gll%Binv, n_gll)
364 end subroutine cpu_convect_scalar_lx16
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)
380 integer :: e, i, j, k, l, idx, n_GLL
382 n_gll = nelv * xh_gll%lxyz
389 tmp = tmp + dx(i,k) * u(k,j,1,e)
400 tmp = tmp + dy(j,l) * u(i,l,k,e)
411 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
421 idx = (e-1) * xh_gll%lxyz+1
422 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
424 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
425 call col2(du, coef_gll%Binv, n_gll)
427 end subroutine cpu_convect_scalar_lx15
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)
443 integer :: e, i, j, k, l, idx, n_GLL
445 n_gll = nelv * xh_gll%lxyz
452 tmp = tmp + dx(i,k) * u(k,j,1,e)
463 tmp = tmp + dy(j,l) * u(i,l,k,e)
474 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
484 idx = (e-1) * xh_gll%lxyz+1
485 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
487 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
488 call col2(du, coef_gll%Binv, n_gll)
490 end subroutine cpu_convect_scalar_lx14
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)
506 integer :: e, i, j, k, l, idx, n_GLL
508 n_gll = nelv * xh_gll%lxyz
515 tmp = tmp + dx(i,k) * u(k,j,1,e)
526 tmp = tmp + dy(j,l) * u(i,l,k,e)
537 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
547 idx = (e-1) * xh_gll%lxyz+1
548 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
550 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
551 call col2(du, coef_gll%Binv, n_gll)
553 end subroutine cpu_convect_scalar_lx13
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)
569 integer :: e, i, j, k, l, idx, n_GLL
571 n_gll = nelv * xh_gll%lxyz
578 tmp = tmp + dx(i,k) * u(k,j,1,e)
589 tmp = tmp + dy(j,l) * u(i,l,k,e)
600 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
610 idx = (e-1) * xh_gll%lxyz+1
611 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
613 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
614 call col2(du, coef_gll%Binv, n_gll)
616 end subroutine cpu_convect_scalar_lx12
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)
632 integer :: e, i, j, k, l, idx, n_GLL
634 n_gll = nelv * xh_gll%lxyz
641 tmp = tmp + dx(i,k) * u(k,j,1,e)
652 tmp = tmp + dy(j,l) * u(i,l,k,e)
663 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
673 idx = (e-1) * xh_gll%lxyz+1
674 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
676 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
677 call col2(du, coef_gll%Binv, n_gll)
679 end subroutine cpu_convect_scalar_lx11
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)
695 integer :: e, i, j, k, l, idx, n_GLL
697 n_gll = nelv * xh_gll%lxyz
704 tmp = tmp + dx(i,k) * u(k,j,1,e)
715 tmp = tmp + dy(j,l) * u(i,l,k,e)
726 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
736 idx = (e-1) * xh_gll%lxyz+1
737 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
739 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
740 call col2(du, coef_gll%Binv, n_gll)
742 end subroutine cpu_convect_scalar_lx10
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)
758 integer :: e, i, j, k, l, idx, n_GLL
760 n_gll = nelv * xh_gll%lxyz
767 tmp = tmp + dx(i,k) * u(k,j,1,e)
778 tmp = tmp + dy(j,l) * u(i,l,k,e)
789 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
799 idx = (e-1) * xh_gll%lxyz+1
800 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
802 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
803 call col2(du, coef_gll%Binv, n_gll)
805 end subroutine cpu_convect_scalar_lx9
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)
821 integer :: e, i, j, k, l, idx, n_GLL
823 n_gll = nelv * xh_gll%lxyz
830 tmp = tmp + dx(i,k) * u(k,j,1,e)
841 tmp = tmp + dy(j,l) * u(i,l,k,e)
852 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
862 idx = (e-1) * xh_gll%lxyz+1
863 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
865 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
866 call col2(du, coef_gll%Binv, n_gll)
868 end subroutine cpu_convect_scalar_lx8
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)
884 integer :: e, i, j, k, l, idx, n_GLL
886 n_gll = nelv * xh_gll%lxyz
893 tmp = tmp + dx(i,k) * u(k,j,1,e)
904 tmp = tmp + dy(j,l) * u(i,l,k,e)
915 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
925 idx = (e-1) * xh_gll%lxyz+1
926 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
928 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
929 call col2(du, coef_gll%Binv, n_gll)
931 end subroutine cpu_convect_scalar_lx7
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)
947 integer :: e, i, j, k, l, idx, n_GLL
949 n_gll = nelv * xh_gll%lxyz
956 tmp = tmp + dx(i,k) * u(k,j,1,e)
967 tmp = tmp + dy(j,l) * u(i,l,k,e)
978 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
988 idx = (e-1) * xh_gll%lxyz+1
989 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
991 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
992 call col2(du, coef_gll%Binv, n_gll)
994 end subroutine cpu_convect_scalar_lx6
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
1012 n_gll = nelv * xh_gll%lxyz
1019 tmp = tmp + dx(i,k) * u(k,j,1,e)
1030 tmp = tmp + dy(j,l) * u(i,l,k,e)
1041 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1051 idx = (e-1) * xh_gll%lxyz+1
1052 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1054 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1055 call col2(du, coef_gll%Binv, n_gll)
1057 end subroutine cpu_convect_scalar_lx5
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
1075 n_gll = nelv * xh_gll%lxyz
1082 tmp = tmp + dx(i,k) * u(k,j,1,e)
1093 tmp = tmp + dy(j,l) * u(i,l,k,e)
1104 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1114 idx = (e-1) * xh_gll%lxyz+1
1115 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1117 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1118 call col2(du, coef_gll%Binv, n_gll)
1120 end subroutine cpu_convect_scalar_lx4
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
1138 n_gll = nelv * xh_gll%lxyz
1145 tmp = tmp + dx(i,k) * u(k,j,1,e)
1156 tmp = tmp + dy(j,l) * u(i,l,k,e)
1167 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1177 idx = (e-1) * xh_gll%lxyz+1
1178 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1180 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1181 call col2(du, coef_gll%Binv, n_gll)
1183 end subroutine cpu_convect_scalar_lx3
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
1201 n_gll = nelv * xh_gll%lxyz
1208 tmp = tmp + dx(i,k) * u(k,j,1,e)
1219 tmp = tmp + dy(j,l) * u(i,l,k,e)
1230 tmp = tmp + dz(k,l) * u(i,1,l,e)
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)
1240 idx = (e-1) * xh_gll%lxyz+1
1241 call gll_to_gl%map(du(idx,1,1,1), ud, 1, xh_gll)
1243 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1244 call col2(du, coef_gll%Binv, n_gll)
1246 end subroutine cpu_convect_scalar_lx2
1248end submodule cpu_convect_scalar
subroutine, public col2(a, b, n)
Vector multiplication .