33submodule(
opr_sx) sx_convect_scalar
40 module subroutine opr_sx_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 sx_convect_scalar_lx18(du, u, c, dx, dy, dz, &
58 xh_gll, coef_gll, gll_to_gl, nelv)
60 call sx_convect_scalar_lx17(du, u, c, dx, dy, dz, &
61 xh_gll, coef_gll, gll_to_gl, nelv)
63 call sx_convect_scalar_lx16(du, u, c, dx, dy, dz, &
64 xh_gll, coef_gll, gll_to_gl, nelv)
66 call sx_convect_scalar_lx15(du, u, c, dx, dy, dz, &
67 xh_gll, coef_gll, gll_to_gl, nelv)
69 call sx_convect_scalar_lx14(du, u, c, dx, dy, dz, &
70 xh_gll, coef_gll, gll_to_gl, nelv)
72 call sx_convect_scalar_lx13(du, u, c, dx, dy, dz, &
73 xh_gll, coef_gll, gll_to_gl, nelv)
75 call sx_convect_scalar_lx12(du, u, c, dx, dy, dz, &
76 xh_gll, coef_gll, gll_to_gl, nelv)
78 call sx_convect_scalar_lx11(du, u, c, dx, dy, dz, &
79 xh_gll, coef_gll, gll_to_gl, nelv)
81 call sx_convect_scalar_lx10(du, u, c, dx, dy, dz, &
82 xh_gll, coef_gll, gll_to_gl, nelv)
84 call sx_convect_scalar_lx9(du, u, c, dx, dy, dz, &
85 xh_gll, coef_gll, gll_to_gl, nelv)
87 call sx_convect_scalar_lx8(du, u, c, dx, dy, dz, &
88 xh_gll, coef_gll, gll_to_gl, nelv)
90 call sx_convect_scalar_lx7(du, u, c, dx, dy, dz, &
91 xh_gll, coef_gll, gll_to_gl, nelv)
93 call sx_convect_scalar_lx6(du, u, c, dx, dy, dz, &
94 xh_gll, coef_gll, gll_to_gl, nelv)
96 call sx_convect_scalar_lx5(du, u, c, dx, dy, dz, &
97 xh_gll, coef_gll, gll_to_gl, nelv)
99 call sx_convect_scalar_lx4(du, u, c, dx, dy, dz, &
100 xh_gll, coef_gll, gll_to_gl, nelv)
102 call sx_convect_scalar_lx3(du, u, c, dx, dy, dz, &
103 xh_gll, coef_gll, gll_to_gl, nelv)
105 call sx_convect_scalar_lx2(du, u, c, dx, dy, dz, &
106 xh_gll, coef_gll, gll_to_gl, nelv)
108 call sx_convect_scalar_lx(du, u, c, dx, dy, dz, &
109 xh_gll, coef_gll, gll_to_gl, nelv, lx)
113 end subroutine opr_sx_convect_scalar
115 subroutine sx_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) :: ur(lx, lx, lx, nelv)
126 real(kind=
rp) :: us(lx, lx, lx, nelv)
127 real(kind=
rp) :: ut(lx, lx, lx, nelv)
128 real(kind=
rp) :: ud(lx, lx, lx, nelv)
129 real(kind=
rp) :: wr, ws, wt
130 integer :: e, i, j, k, jj, kk, n_GLL
131 n_gll = nelv * xh_gll%lxyz
133 do jj = 1, lx * lx * nelv
136 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
149 ws = ws + dy(j, kk) * u(i, kk, k, e)
164 wt = wt + dz(k, kk) * u(i, j, kk, e)
172 do i = 1, lx * lx * lx
174 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
175 + c(i,e,2) * us(i,1,1,e) &
176 + c(i,e,3) * ut(i,1,1,e) )
179 call gll_to_gl%map(du, ud, nelv, xh_gll)
180 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
181 call col2(du, coef_gll%Binv, n_gll)
183 end subroutine sx_convect_scalar_lx
185 subroutine sx_convect_scalar_lx18(du, u, c, dx, dy, dz, Xh_GLL, &
186 coef_GLL, GLL_to_GL, nelv)
187 integer,
parameter :: lx = 18
188 integer,
intent(in) :: nelv
189 type(space_t),
intent(in) :: Xh_GLL
190 type(coef_t),
intent(in) :: coef_GLL
191 type(interpolator_t),
intent(inout) :: GLL_to_GL
192 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
193 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
194 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
195 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
196 real(kind=
rp) :: ur(lx, lx, lx, nelv)
197 real(kind=
rp) :: us(lx, lx, lx, nelv)
198 real(kind=
rp) :: ut(lx, lx, lx, nelv)
199 real(kind=
rp) :: ud(lx, lx, lx, nelv)
200 real(kind=
rp) :: wr, ws, wt
201 integer :: e, i, j, k, jj, kk, n_GLL
202 n_gll = nelv * xh_gll%lxyz
204 do jj = 1, lx * lx * nelv
207 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
220 ws = ws + dy(j, kk) * u(i, kk, k, e)
235 wt = wt + dz(k, kk) * u(i, j, kk, e)
243 do i = 1, lx * lx * lx
245 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
246 + c(i,e,2) * us(i,1,1,e) &
247 + c(i,e,3) * ut(i,1,1,e) )
250 call gll_to_gl%map(du, ud, nelv, xh_gll)
251 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
252 call col2(du, coef_gll%Binv, n_gll)
254 end subroutine sx_convect_scalar_lx18
256 subroutine sx_convect_scalar_lx17(du, u, c, dx, dy, dz, Xh_GLL, &
257 coef_GLL, GLL_to_GL, nelv)
258 integer,
parameter :: lx = 17
259 integer,
intent(in) :: nelv
260 type(space_t),
intent(in) :: Xh_GLL
261 type(coef_t),
intent(in) :: coef_GLL
262 type(interpolator_t),
intent(inout) :: GLL_to_GL
263 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
264 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
265 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
266 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
267 real(kind=
rp) :: ur(lx, lx, lx, nelv)
268 real(kind=
rp) :: us(lx, lx, lx, nelv)
269 real(kind=
rp) :: ut(lx, lx, lx, nelv)
270 real(kind=
rp) :: ud(lx, lx, lx, nelv)
271 real(kind=
rp) :: wr, ws, wt
272 integer :: e, i, j, k, jj, kk, n_GLL
273 n_gll = nelv * xh_gll%lxyz
275 do jj = 1, lx * lx * nelv
278 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
291 ws = ws + dy(j, kk) * u(i, kk, k, e)
306 wt = wt + dz(k, kk) * u(i, j, kk, e)
314 do i = 1, lx * lx * lx
316 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
317 + c(i,e,2) * us(i,1,1,e) &
318 + c(i,e,3) * ut(i,1,1,e) )
321 call gll_to_gl%map(du, ud, nelv, xh_gll)
322 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
323 call col2(du, coef_gll%Binv, n_gll)
325 end subroutine sx_convect_scalar_lx17
327 subroutine sx_convect_scalar_lx16(du, u, c, dx, dy, dz, Xh_GLL, &
328 coef_GLL, GLL_to_GL, nelv)
329 integer,
parameter :: lx = 16
330 integer,
intent(in) :: nelv
331 type(space_t),
intent(in) :: Xh_GLL
332 type(coef_t),
intent(in) :: coef_GLL
333 type(interpolator_t),
intent(inout) :: GLL_to_GL
334 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
335 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
336 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
337 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
338 real(kind=
rp) :: ur(lx, lx, lx, nelv)
339 real(kind=
rp) :: us(lx, lx, lx, nelv)
340 real(kind=
rp) :: ut(lx, lx, lx, nelv)
341 real(kind=
rp) :: ud(lx, lx, lx, nelv)
342 real(kind=
rp) :: wr, ws, wt
343 integer :: e, i, j, k, jj, kk, n_GLL
344 n_gll = nelv * xh_gll%lxyz
346 do jj = 1, lx * lx * nelv
349 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
362 ws = ws + dy(j, kk) * u(i, kk, k, e)
377 wt = wt + dz(k, kk) * u(i, j, kk, e)
385 do i = 1, lx * lx * lx
387 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
388 + c(i,e,2) * us(i,1,1,e) &
389 + c(i,e,3) * ut(i,1,1,e) )
392 call gll_to_gl%map(du, ud, nelv, xh_gll)
393 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
394 call col2(du, coef_gll%Binv, n_gll)
396 end subroutine sx_convect_scalar_lx16
398 subroutine sx_convect_scalar_lx15(du, u, c, dx, dy, dz, Xh_GLL, &
399 coef_GLL, GLL_to_GL, nelv)
400 integer,
parameter :: lx = 15
401 integer,
intent(in) :: nelv
402 type(space_t),
intent(in) :: Xh_GLL
403 type(coef_t),
intent(in) :: coef_GLL
404 type(interpolator_t),
intent(inout) :: GLL_to_GL
405 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
406 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
407 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
408 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
409 real(kind=
rp) :: ur(lx, lx, lx, nelv)
410 real(kind=
rp) :: us(lx, lx, lx, nelv)
411 real(kind=
rp) :: ut(lx, lx, lx, nelv)
412 real(kind=
rp) :: ud(lx, lx, lx, nelv)
413 real(kind=
rp) :: wr, ws, wt
414 integer :: e, i, j, k, jj, kk, n_GLL
415 n_gll = nelv * xh_gll%lxyz
417 do jj = 1, lx * lx * nelv
420 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
433 ws = ws + dy(j, kk) * u(i, kk, k, e)
448 wt = wt + dz(k, kk) * u(i, j, kk, e)
456 do i = 1, lx * lx * lx
458 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
459 + c(i,e,2) * us(i,1,1,e) &
460 + c(i,e,3) * ut(i,1,1,e) )
463 call gll_to_gl%map(du, ud, nelv, xh_gll)
464 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
465 call col2(du, coef_gll%Binv, n_gll)
467 end subroutine sx_convect_scalar_lx15
469 subroutine sx_convect_scalar_lx14(du, u, c, dx, dy, dz, Xh_GLL, &
470 coef_GLL, GLL_to_GL, nelv)
471 integer,
parameter :: lx = 14
472 integer,
intent(in) :: nelv
473 type(space_t),
intent(in) :: Xh_GLL
474 type(coef_t),
intent(in) :: coef_GLL
475 type(interpolator_t),
intent(inout) :: GLL_to_GL
476 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
477 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
478 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
479 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
480 real(kind=
rp) :: ur(lx, lx, lx, nelv)
481 real(kind=
rp) :: us(lx, lx, lx, nelv)
482 real(kind=
rp) :: ut(lx, lx, lx, nelv)
483 real(kind=
rp) :: ud(lx, lx, lx, nelv)
484 real(kind=
rp) :: wr, ws, wt
485 integer :: e, i, j, k, jj, kk, n_GLL
486 n_gll = nelv * xh_gll%lxyz
488 do jj = 1, lx * lx * nelv
491 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
504 ws = ws + dy(j, kk) * u(i, kk, k, e)
519 wt = wt + dz(k, kk) * u(i, j, kk, e)
527 do i = 1, lx * lx * lx
529 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
530 + c(i,e,2) * us(i,1,1,e) &
531 + c(i,e,3) * ut(i,1,1,e) )
534 call gll_to_gl%map(du, ud, nelv, xh_gll)
535 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
536 call col2(du, coef_gll%Binv, n_gll)
538 end subroutine sx_convect_scalar_lx14
540 subroutine sx_convect_scalar_lx13(du, u, c, dx, dy, dz, Xh_GLL, &
541 coef_GLL, GLL_to_GL, nelv)
542 integer,
parameter :: lx = 13
543 integer,
intent(in) :: nelv
544 type(space_t),
intent(in) :: Xh_GLL
545 type(coef_t),
intent(in) :: coef_GLL
546 type(interpolator_t),
intent(inout) :: GLL_to_GL
547 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
548 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
549 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
550 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
551 real(kind=
rp) :: ur(lx, lx, lx, nelv)
552 real(kind=
rp) :: us(lx, lx, lx, nelv)
553 real(kind=
rp) :: ut(lx, lx, lx, nelv)
554 real(kind=
rp) :: ud(lx, lx, lx, nelv)
555 real(kind=
rp) :: wr, ws, wt
556 integer :: e, i, j, k, jj, kk, n_GLL
557 n_gll = nelv * xh_gll%lxyz
559 do jj = 1, lx * lx * nelv
562 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
575 ws = ws + dy(j, kk) * u(i, kk, k, e)
590 wt = wt + dz(k, kk) * u(i, j, kk, e)
598 do i = 1, lx * lx * lx
600 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
601 + c(i,e,2) * us(i,1,1,e) &
602 + c(i,e,3) * ut(i,1,1,e) )
605 call gll_to_gl%map(du, ud, nelv, xh_gll)
606 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
607 call col2(du, coef_gll%Binv, n_gll)
609 end subroutine sx_convect_scalar_lx13
611 subroutine sx_convect_scalar_lx12(du, u, c, dx, dy, dz, Xh_GLL, &
612 coef_GLL, GLL_to_GL, nelv)
613 integer,
parameter :: lx = 12
614 integer,
intent(in) :: nelv
615 type(space_t),
intent(in) :: Xh_GLL
616 type(coef_t),
intent(in) :: coef_GLL
617 type(interpolator_t),
intent(inout) :: GLL_to_GL
618 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
619 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
620 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
621 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
622 real(kind=
rp) :: ur(lx, lx, lx, nelv)
623 real(kind=
rp) :: us(lx, lx, lx, nelv)
624 real(kind=
rp) :: ut(lx, lx, lx, nelv)
625 real(kind=
rp) :: ud(lx, lx, lx, nelv)
626 real(kind=
rp) :: wr, ws, wt
627 integer :: e, i, j, k, jj, kk, n_GLL
628 n_gll = nelv * xh_gll%lxyz
630 do jj = 1, lx * lx * nelv
633 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
646 ws = ws + dy(j, kk) * u(i, kk, k, e)
661 wt = wt + dz(k, kk) * u(i, j, kk, e)
669 do i = 1, lx * lx * lx
671 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
672 + c(i,e,2) * us(i,1,1,e) &
673 + c(i,e,3) * ut(i,1,1,e) )
676 call gll_to_gl%map(du, ud, nelv, 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 sx_convect_scalar_lx12
682 subroutine sx_convect_scalar_lx11(du, u, c, dx, dy, dz, Xh_GLL, &
683 coef_GLL, GLL_to_GL, nelv)
684 integer,
parameter :: lx = 11
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) :: ur(lx, lx, lx, nelv)
694 real(kind=
rp) :: us(lx, lx, lx, nelv)
695 real(kind=
rp) :: ut(lx, lx, lx, nelv)
696 real(kind=
rp) :: ud(lx, lx, lx, nelv)
697 real(kind=
rp) :: wr, ws, wt
698 integer :: e, i, j, k, jj, kk, n_GLL
699 n_gll = nelv * xh_gll%lxyz
701 do jj = 1, lx * lx * nelv
704 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
717 ws = ws + dy(j, kk) * u(i, kk, k, e)
732 wt = wt + dz(k, kk) * u(i, j, kk, e)
740 do i = 1, lx * lx * lx
742 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
743 + c(i,e,2) * us(i,1,1,e) &
744 + c(i,e,3) * ut(i,1,1,e) )
747 call gll_to_gl%map(du, ud, nelv, xh_gll)
748 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
749 call col2(du, coef_gll%Binv, n_gll)
751 end subroutine sx_convect_scalar_lx11
753 subroutine sx_convect_scalar_lx10(du, u, c, dx, dy, dz, Xh_GLL, &
754 coef_GLL, GLL_to_GL, nelv)
755 integer,
parameter :: lx = 10
756 integer,
intent(in) :: nelv
757 type(space_t),
intent(in) :: Xh_GLL
758 type(coef_t),
intent(in) :: coef_GLL
759 type(interpolator_t),
intent(inout) :: GLL_to_GL
760 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
761 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
762 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
763 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
764 real(kind=
rp) :: ur(lx, lx, lx, nelv)
765 real(kind=
rp) :: us(lx, lx, lx, nelv)
766 real(kind=
rp) :: ut(lx, lx, lx, nelv)
767 real(kind=
rp) :: ud(lx, lx, lx, nelv)
768 real(kind=
rp) :: wr, ws, wt
769 integer :: e, i, j, k, jj, kk, n_GLL
770 n_gll = nelv * xh_gll%lxyz
772 do jj = 1, lx * lx * nelv
775 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
788 ws = ws + dy(j, kk) * u(i, kk, k, e)
803 wt = wt + dz(k, kk) * u(i, j, kk, e)
811 do i = 1, lx * lx * lx
813 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
814 + c(i,e,2) * us(i,1,1,e) &
815 + c(i,e,3) * ut(i,1,1,e) )
818 call gll_to_gl%map(du, ud, nelv, xh_gll)
819 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
820 call col2(du, coef_gll%Binv, n_gll)
822 end subroutine sx_convect_scalar_lx10
824 subroutine sx_convect_scalar_lx9(du, u, c, dx, dy, dz, Xh_GLL, &
825 coef_GLL, GLL_to_GL, nelv)
826 integer,
parameter :: lx = 9
827 integer,
intent(in) :: nelv
828 type(space_t),
intent(in) :: Xh_GLL
829 type(coef_t),
intent(in) :: coef_GLL
830 type(interpolator_t),
intent(inout) :: GLL_to_GL
831 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
832 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
833 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
834 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
835 real(kind=
rp) :: ur(lx, lx, lx, nelv)
836 real(kind=
rp) :: us(lx, lx, lx, nelv)
837 real(kind=
rp) :: ut(lx, lx, lx, nelv)
838 real(kind=
rp) :: ud(lx, lx, lx, nelv)
839 real(kind=
rp) :: wr, ws, wt
840 integer :: e, i, j, k, jj, kk, n_GLL
841 n_gll = nelv * xh_gll%lxyz
843 do jj = 1, lx * lx * nelv
846 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
859 ws = ws + dy(j, kk) * u(i, kk, k, e)
874 wt = wt + dz(k, kk) * u(i, j, kk, e)
882 do i = 1, lx * lx * lx
884 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
885 + c(i,e,2) * us(i,1,1,e) &
886 + c(i,e,3) * ut(i,1,1,e) )
889 call gll_to_gl%map(du, ud, nelv, xh_gll)
890 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
891 call col2(du, coef_gll%Binv, n_gll)
893 end subroutine sx_convect_scalar_lx9
895 subroutine sx_convect_scalar_lx8(du, u, c, dx, dy, dz, Xh_GLL, &
896 coef_GLL, GLL_to_GL, nelv)
897 integer,
parameter :: lx = 8
898 integer,
intent(in) :: nelv
899 type(space_t),
intent(in) :: Xh_GLL
900 type(coef_t),
intent(in) :: coef_GLL
901 type(interpolator_t),
intent(inout) :: GLL_to_GL
902 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
903 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
904 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
905 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
906 real(kind=
rp) :: ur(lx, lx, lx, nelv)
907 real(kind=
rp) :: us(lx, lx, lx, nelv)
908 real(kind=
rp) :: ut(lx, lx, lx, nelv)
909 real(kind=
rp) :: ud(lx, lx, lx, nelv)
910 real(kind=
rp) :: wr, ws, wt
911 integer :: e, i, j, k, jj, kk, n_GLL
912 n_gll = nelv * xh_gll%lxyz
914 do jj = 1, lx * lx * nelv
917 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
930 ws = ws + dy(j, kk) * u(i, kk, k, e)
945 wt = wt + dz(k, kk) * u(i, j, kk, e)
953 do i = 1, lx * lx * lx
955 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
956 + c(i,e,2) * us(i,1,1,e) &
957 + c(i,e,3) * ut(i,1,1,e) )
960 call gll_to_gl%map(du, ud, nelv, xh_gll)
961 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
962 call col2(du, coef_gll%Binv, n_gll)
964 end subroutine sx_convect_scalar_lx8
966 subroutine sx_convect_scalar_lx7(du, u, c, dx, dy, dz, Xh_GLL, &
967 coef_GLL, GLL_to_GL, nelv)
968 integer,
parameter :: lx = 7
969 integer,
intent(in) :: nelv
970 type(space_t),
intent(in) :: Xh_GLL
971 type(coef_t),
intent(in) :: coef_GLL
972 type(interpolator_t),
intent(inout) :: GLL_to_GL
973 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
974 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
975 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
976 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
977 real(kind=
rp) :: ur(lx, lx, lx, nelv)
978 real(kind=
rp) :: us(lx, lx, lx, nelv)
979 real(kind=
rp) :: ut(lx, lx, lx, nelv)
980 real(kind=
rp) :: ud(lx, lx, lx, nelv)
981 real(kind=
rp) :: wr, ws, wt
982 integer :: e, i, j, k, jj, kk, n_GLL
983 n_gll = nelv * xh_gll%lxyz
985 do jj = 1, lx * lx * nelv
988 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
1001 ws = ws + dy(j, kk) * u(i, kk, k, e)
1016 wt = wt + dz(k, kk) * u(i, j, kk, e)
1024 do i = 1, lx * lx * lx
1026 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
1027 + c(i,e,2) * us(i,1,1,e) &
1028 + c(i,e,3) * ut(i,1,1,e) )
1031 call gll_to_gl%map(du, ud, nelv, xh_gll)
1032 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1033 call col2(du, coef_gll%Binv, n_gll)
1035 end subroutine sx_convect_scalar_lx7
1037 subroutine sx_convect_scalar_lx6(du, u, c, dx, dy, dz, Xh_GLL, &
1038 coef_GLL, GLL_to_GL, nelv)
1039 integer,
parameter :: lx = 6
1040 integer,
intent(in) :: nelv
1041 type(space_t),
intent(in) :: Xh_GLL
1042 type(coef_t),
intent(in) :: coef_GLL
1043 type(interpolator_t),
intent(inout) :: GLL_to_GL
1044 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1045 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
1046 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
1047 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1048 real(kind=
rp) :: ur(lx, lx, lx, nelv)
1049 real(kind=
rp) :: us(lx, lx, lx, nelv)
1050 real(kind=
rp) :: ut(lx, lx, lx, nelv)
1051 real(kind=
rp) :: ud(lx, lx, lx, nelv)
1052 real(kind=
rp) :: wr, ws, wt
1053 integer :: e, i, j, k, jj, kk, n_GLL
1054 n_gll = nelv * xh_gll%lxyz
1056 do jj = 1, lx * lx * nelv
1059 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
1061 ur(i, jj, 1, 1) = wr
1072 ws = ws + dy(j, kk) * u(i, kk, k, e)
1087 wt = wt + dz(k, kk) * u(i, j, kk, e)
1095 do i = 1, lx * lx * lx
1097 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
1098 + c(i,e,2) * us(i,1,1,e) &
1099 + c(i,e,3) * ut(i,1,1,e) )
1102 call gll_to_gl%map(du, ud, nelv, xh_gll)
1103 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1104 call col2(du, coef_gll%Binv, n_gll)
1106 end subroutine sx_convect_scalar_lx6
1108 subroutine sx_convect_scalar_lx5(du, u, c, dx, dy, dz, Xh_GLL, &
1109 coef_GLL, GLL_to_GL, nelv)
1110 integer,
parameter :: lx = 5
1111 integer,
intent(in) :: nelv
1112 type(space_t),
intent(in) :: Xh_GLL
1113 type(coef_t),
intent(in) :: coef_GLL
1114 type(interpolator_t),
intent(inout) :: GLL_to_GL
1115 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1116 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
1117 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
1118 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1119 real(kind=
rp) :: ur(lx, lx, lx, nelv)
1120 real(kind=
rp) :: us(lx, lx, lx, nelv)
1121 real(kind=
rp) :: ut(lx, lx, lx, nelv)
1122 real(kind=
rp) :: ud(lx, lx, lx, nelv)
1123 real(kind=
rp) :: wr, ws, wt
1124 integer :: e, i, j, k, jj, kk, n_GLL
1125 n_gll = nelv * xh_gll%lxyz
1127 do jj = 1, lx * lx * nelv
1130 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
1132 ur(i, jj, 1, 1) = wr
1143 ws = ws + dy(j, kk) * u(i, kk, k, e)
1158 wt = wt + dz(k, kk) * u(i, j, kk, e)
1166 do i = 1, lx * lx * lx
1168 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
1169 + c(i,e,2) * us(i,1,1,e) &
1170 + c(i,e,3) * ut(i,1,1,e) )
1173 call gll_to_gl%map(du, ud, nelv, xh_gll)
1174 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1175 call col2(du, coef_gll%Binv, n_gll)
1177 end subroutine sx_convect_scalar_lx5
1179 subroutine sx_convect_scalar_lx4(du, u, c, dx, dy, dz, Xh_GLL, &
1180 coef_GLL, GLL_to_GL, nelv)
1181 integer,
parameter :: lx = 4
1182 integer,
intent(in) :: nelv
1183 type(space_t),
intent(in) :: Xh_GLL
1184 type(coef_t),
intent(in) :: coef_GLL
1185 type(interpolator_t),
intent(inout) :: GLL_to_GL
1186 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1187 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
1188 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
1189 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1190 real(kind=
rp) :: ur(lx, lx, lx, nelv)
1191 real(kind=
rp) :: us(lx, lx, lx, nelv)
1192 real(kind=
rp) :: ut(lx, lx, lx, nelv)
1193 real(kind=
rp) :: ud(lx, lx, lx, nelv)
1194 real(kind=
rp) :: wr, ws, wt
1195 integer :: e, i, j, k, jj, kk, n_GLL
1196 n_gll = nelv * xh_gll%lxyz
1198 do jj = 1, lx * lx * nelv
1201 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
1203 ur(i, jj, 1, 1) = wr
1214 ws = ws + dy(j, kk) * u(i, kk, k, e)
1229 wt = wt + dz(k, kk) * u(i, j, kk, e)
1237 do i = 1, lx * lx * lx
1239 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
1240 + c(i,e,2) * us(i,1,1,e) &
1241 + c(i,e,3) * ut(i,1,1,e) )
1244 call gll_to_gl%map(du, ud, nelv, xh_gll)
1245 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1246 call col2(du, coef_gll%Binv, n_gll)
1248 end subroutine sx_convect_scalar_lx4
1250 subroutine sx_convect_scalar_lx3(du, u, c, dx, dy, dz, Xh_GLL, &
1251 coef_GLL, GLL_to_GL, nelv)
1252 integer,
parameter :: lx = 3
1253 integer,
intent(in) :: nelv
1254 type(space_t),
intent(in) :: Xh_GLL
1255 type(coef_t),
intent(in) :: coef_GLL
1256 type(interpolator_t),
intent(inout) :: GLL_to_GL
1257 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1258 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
1259 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
1260 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1261 real(kind=
rp) :: ur(lx, lx, lx, nelv)
1262 real(kind=
rp) :: us(lx, lx, lx, nelv)
1263 real(kind=
rp) :: ut(lx, lx, lx, nelv)
1264 real(kind=
rp) :: ud(lx, lx, lx, nelv)
1265 real(kind=
rp) :: wr, ws, wt
1266 integer :: e, i, j, k, jj, kk, n_GLL
1267 n_gll = nelv * xh_gll%lxyz
1269 do jj = 1, lx * lx * nelv
1272 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
1274 ur(i, jj, 1, 1) = wr
1285 ws = ws + dy(j, kk) * u(i, kk, k, e)
1300 wt = wt + dz(k, kk) * u(i, j, kk, e)
1308 do i = 1, lx * lx * lx
1310 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
1311 + c(i,e,2) * us(i,1,1,e) &
1312 + c(i,e,3) * ut(i,1,1,e) )
1315 call gll_to_gl%map(du, ud, nelv, xh_gll)
1316 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1317 call col2(du, coef_gll%Binv, n_gll)
1319 end subroutine sx_convect_scalar_lx3
1321 subroutine sx_convect_scalar_lx2(du, u, c, dx, dy, dz, Xh_GLL, &
1322 coef_GLL, GLL_to_GL, nelv)
1323 integer,
parameter :: lx = 2
1324 integer,
intent(in) :: nelv
1325 type(space_t),
intent(in) :: Xh_GLL
1326 type(coef_t),
intent(in) :: coef_GLL
1327 type(interpolator_t),
intent(inout) :: GLL_to_GL
1328 real(kind=
rp),
intent(inout) :: du(xh_gll%lx, xh_gll%lx, xh_gll%lx, nelv)
1329 real(kind=
rp),
intent(in) :: u(lx, lx, lx, nelv)
1330 real(kind=
rp),
intent(in) :: c(lx*lx*lx, nelv, 3)
1331 real(kind=
rp),
dimension(lx, lx),
intent(in) :: dx, dy, dz
1332 real(kind=
rp) :: ur(lx, lx, lx, nelv)
1333 real(kind=
rp) :: us(lx, lx, lx, nelv)
1334 real(kind=
rp) :: ut(lx, lx, lx, nelv)
1335 real(kind=
rp) :: ud(lx, lx, lx, nelv)
1336 real(kind=
rp) :: wr, ws, wt
1337 integer :: e, i, j, k, jj, kk, n_GLL
1338 n_gll = nelv * xh_gll%lxyz
1340 do jj = 1, lx * lx * nelv
1343 wr = wr + dx(i, kk) * u(kk, jj, 1, 1)
1345 ur(i, jj, 1, 1) = wr
1356 ws = ws + dy(j, kk) * u(i, kk, k, e)
1371 wt = wt + dz(k, kk) * u(i, j, kk, e)
1379 do i = 1, lx * lx * lx
1381 ud(i,1,1,e) = ( c(i,e,1) * ur(i,1,1,e) &
1382 + c(i,e,2) * us(i,1,1,e) &
1383 + c(i,e,3) * ut(i,1,1,e) )
1386 call gll_to_gl%map(du, ud, nelv, xh_gll)
1387 call coef_gll%gs_h%op(du, n_gll, gs_op_add)
1388 call col2(du, coef_gll%Binv, n_gll)
1390 end subroutine sx_convect_scalar_lx2
1392end submodule sx_convect_scalar
subroutine, public col2(a, b, n)
Vector multiplication .
integer, parameter, public rp
Global precision used in computations.
Operators SX-Aurora backend.