39 opr_cpu_conv1, opr_cpu_convect_scalar, opr_cpu_cdtp, &
43 opr_sx_conv1, opr_sx_convect_scalar, opr_sx_cdtp, &
44 opr_sx_dudxyz, opr_sx_lambda2, opr_sx_set_convect_rst
67 use mpi_f08,
only : mpi_allreduce, mpi_in_place, mpi_max, mpi_sum
68 use,
intrinsic :: iso_c_binding, only : c_ptr
89 subroutine dudxyz (du, u, dr, ds, dt, coef)
90 type(
coef_t),
intent(in),
target :: coef
91 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(inout) :: du
92 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(in) :: u, dr, ds, dt
95 call opr_sx_dudxyz(du, u, dr, ds, dt, coef)
101 call opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
112 subroutine div(res, ux, uy, uz, coef)
113 type(
coef_t),
intent(in),
target :: coef
114 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(inout) :: res
115 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(in) :: ux, uy, uz
127 call dudxyz(res, ux, coef%drdx, coef%dsdx, coef%dtdx, coef)
130 call dudxyz(work%x, uy, coef%drdy, coef%dsdy, coef%dtdy, coef)
134 call add2(res, work%x, work%size())
138 call dudxyz(work%x, uz, coef%drdz, coef%dsdz, coef%dtdz, coef)
142 call add2(res, work%x, work%size())
155 subroutine grad(ux, uy, uz, u, coef)
156 type(
coef_t),
intent(in) :: coef
157 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: ux
158 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uy
159 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uz
160 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: u
162 call dudxyz(ux, u, coef%drdx, coef%dsdx, coef%dtdx, coef)
163 call dudxyz(uy, u, coef%drdy, coef%dsdy, coef%dtdy, coef)
164 call dudxyz(uz, u, coef%drdz, coef%dsdz, coef%dtdz, coef)
180 subroutine opgrad(ux, uy, uz, u, coef, es, ee)
181 type(
coef_t),
intent(in) :: coef
182 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: ux
183 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uy
184 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uz
185 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: u
186 integer,
optional :: es, ee
187 integer :: eblk_start, eblk_end
189 if (
present(es))
then
195 if (
present(ee))
then
198 eblk_end = coef%msh%nelv
202 call opr_sx_opgrad(ux, uy, uz, u, coef)
208 call opr_cpu_opgrad(ux, uy, uz, u, coef, eblk_start, eblk_end)
217 subroutine ortho(x, glb_n_points, n)
218 integer,
intent(in) :: n
219 integer(kind=i8),
intent(in) :: glb_n_points
220 real(kind=
rp),
dimension(n),
intent(inout) :: x
228 c =
glsum(x, n)/glb_n_points
245 subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
246 type(
coef_t),
intent(in) :: coef
247 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: dtx
248 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: x
249 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: dr
250 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: ds
251 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: dt
252 integer,
optional :: es, ee
253 integer :: eblk_start, eblk_end
255 if (
present(es))
then
261 if (
present(ee))
then
264 eblk_end = coef%msh%nelv
268 call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
274 call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
289 subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
290 type(
space_t),
intent(in) :: xh
291 type(
coef_t),
intent(in) :: coef
292 real(kind=
rp),
intent(inout) :: du(xh%lxyz, coef%msh%nelv)
293 real(kind=
rp),
intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
294 real(kind=
rp),
intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
295 real(kind=
rp),
intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
296 real(kind=
rp),
intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
297 integer,
optional :: es, ee
298 integer :: eblk_end, eblk_start
300 associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
301 if (
present(es))
then
307 if (
present(ee))
then
310 eblk_end = coef%msh%nelv
314 call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
320 call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
343 subroutine convect_scalar(du, u, cr, cs, ct, Xh_GLL, Xh_GL, coef_GLL, &
345 type(space_t),
intent(in) :: xh_gl
346 type(space_t),
intent(in) :: xh_gll
347 type(coef_t),
intent(in) :: coef_GLL
348 type(coef_t),
intent(in) :: coef_GL
349 type(interpolator_t),
intent(inout) :: GLL_to_GL
350 real(kind=rp),
intent(inout) :: &
351 du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
352 real(kind=rp),
intent(inout) :: &
353 u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
354 type(field_t),
intent(inout) :: cr, cs, ct
357 if (neko_bcknd_sx .eq. 1)
then
358 call opr_sx_convect_scalar(du, u, cr%x, cs%x, ct%x, &
359 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
360 else if (neko_bcknd_xsmm .eq. 1)
then
361 call opr_xsmm_convect_scalar(du, u, cr%x, cs%x, ct%x, &
362 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
363 else if (neko_bcknd_device .eq. 1)
then
364 u_d = device_get_ptr(u)
365 call opr_device_convect_scalar(du, u_d, cr%x_d, cs%x_d, ct%x_d, &
366 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
368 call opr_cpu_convect_scalar(du, u, cr%x, cs%x, ct%x, &
369 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
384 subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef, event)
385 type(field_t),
intent(inout) :: w1
386 type(field_t),
intent(inout) :: w2
387 type(field_t),
intent(inout) :: w3
388 type(field_t),
intent(in) :: u1
389 type(field_t),
intent(in) :: u2
390 type(field_t),
intent(in) :: u3
391 type(field_t),
intent(inout) :: work1
392 type(field_t),
intent(inout) :: work2
393 type(coef_t),
intent(in) :: coef
394 type(c_ptr),
optional,
intent(inout) :: event
396 if (neko_bcknd_sx .eq. 1)
then
397 call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
398 else if (neko_bcknd_xsmm .eq. 1)
then
399 call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
400 else if (neko_bcknd_device .eq. 1)
then
401 if (
present(event))
then
402 call opr_device_curl(w1, w2, w3, u1, u2, u3, &
403 work1, work2, coef, event)
405 call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
408 call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
422 function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
423 type(space_t),
intent(in) :: xh
424 type(coef_t),
intent(in) :: coef
425 integer,
intent(in) :: nelv, gdim
426 real(kind=rp),
intent(in) :: dt
427 real(kind=rp),
dimension(Xh%lx, Xh%ly, Xh%lz, nelv),
intent(in) :: u, v, w
431 if (neko_bcknd_sx .eq. 1)
then
432 cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
433 else if (neko_bcknd_device .eq. 1)
then
434 cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
436 cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
439 if (.not. neko_device_mpi)
then
440 call mpi_allreduce(mpi_in_place,
cfl, 1, &
441 mpi_real_precision, mpi_max, neko_comm, ierr)
454 type(space_t),
intent(in) :: xh
455 type(coef_t),
intent(in) :: coef
456 integer,
intent(in) :: nelv, gdim
457 real(kind=rp),
intent(in) :: dt
458 real(kind=rp),
dimension(Xh%lx, Xh%ly, Xh%lz, nelv),
intent(in) :: max_wave_speed
462 xh, coef, nelv, gdim)
478 subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
479 type(field_t),
intent(in) :: u, v, w
480 type(coef_t),
intent(in) :: coef
481 real(kind=rp),
intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
482 real(kind=rp),
intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
483 real(kind=rp),
intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
484 real(kind=rp),
intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
485 real(kind=rp),
intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
486 real(kind=rp),
intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
488 type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
490 integer :: nelv, lxyz
492 if (neko_bcknd_device .eq. 1)
then
493 s11_d = device_get_ptr(s11)
494 s22_d = device_get_ptr(s22)
495 s33_d = device_get_ptr(s33)
496 s12_d = device_get_ptr(s12)
497 s23_d = device_get_ptr(s23)
498 s13_d = device_get_ptr(s13)
505 call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
506 call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
507 if (neko_bcknd_device .eq. 1)
then
508 call device_add2(s12_d, s11_d, nelv*lxyz)
510 call add2(s12, s11, nelv*lxyz)
513 call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
514 call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
515 if (neko_bcknd_device .eq. 1)
then
516 call device_add2(s13_d, s11_d, nelv*lxyz)
518 call add2(s13, s11, nelv*lxyz)
521 call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
522 call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
523 if (neko_bcknd_device .eq. 1)
then
524 call device_add2(s23_d, s11_d, nelv*lxyz)
526 call add2(s23, s11, nelv*lxyz)
529 call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
530 call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
531 call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
533 if (neko_bcknd_device .eq. 1)
then
534 call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
535 call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
536 call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
538 call cmult(s12, 0.5_rp, nelv*lxyz)
539 call cmult(s13, 0.5_rp, nelv*lxyz)
540 call cmult(s23, 0.5_rp, nelv*lxyz)
551 subroutine lambda2op(lambda2, u, v, w, coef)
552 type(coef_t),
intent(in) :: coef
553 type(field_t),
intent(inout) ::
lambda2
554 type(field_t),
intent(in) :: u, v, w
556 if (neko_bcknd_sx .eq. 1)
then
557 call opr_sx_lambda2(
lambda2, u, v, w, coef)
558 else if (neko_bcknd_device .eq. 1)
then
559 call opr_device_lambda2(
lambda2, u, v, w, coef)
561 call opr_cpu_lambda2(
lambda2, u, v, w, coef)
577 type(space_t),
intent(inout) :: xh
578 type(coef_t),
intent(inout) :: coef
579 type(field_t),
intent(inout) :: cr, cs, ct
580 real(kind=rp),
dimension(Xh%lxyz, coef%msh%nelv), &
581 intent(in) :: cx, cy, cz
582 type(c_ptr) :: cx_d, cy_d, cz_d
584 if (neko_bcknd_sx .eq. 1)
then
585 call opr_sx_set_convect_rst(cr%x, cs%x, ct%x, cx, cy, cz, xh, coef)
586 else if (neko_bcknd_xsmm .eq. 1)
then
587 call opr_xsmm_set_convect_rst(cr%x, cs%x, ct%x, cx, cy, cz, xh, coef)
588 else if (neko_bcknd_device .eq. 1)
then
589 cx_d = device_get_ptr(cx)
590 cy_d = device_get_ptr(cy)
591 cz_d = device_get_ptr(cz)
592 call opr_device_set_convect_rst(cr%x_d, cs%x_d, ct%x_d, &
593 cx_d, cy_d, cz_d, xh, coef)
595 call opr_cpu_set_convect_rst(cr%x, cs%x, ct%x, cx, cy, cz, xh, coef)
615 subroutine runge_kutta(phi, conv_k1, conv_k23, conv_k4, Xh_GLL, Xh_GL, &
616 coef, coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
617 type(space_t),
intent(in) :: xh_gll
618 type(space_t),
intent(inout) :: xh_gl
619 type(coef_t),
intent(in) :: coef
620 type(coef_t),
intent(inout) :: coef_gl
621 type(interpolator_t) :: gll_to_gl
622 real(kind=rp),
intent(inout) :: tau, dtau
623 integer,
intent(in) :: n, nel, n_gl
624 type(field_t),
intent(inout) :: phi
625 type(field_list_t) :: conv_k1, conv_k23, conv_k4
626 real(kind=rp) :: c1, c2, c3
627 type(field_t),
pointer :: u1, k1, k2, k3, k4
628 type(vector_t),
pointer :: u1_gl
629 integer :: ind(6), i, e
631 call neko_scratch_registry%request_field(u1, ind(1), .false.)
632 call neko_scratch_registry%request_field(k1, ind(2), .false.)
633 call neko_scratch_registry%request_field(k2, ind(3), .false.)
634 call neko_scratch_registry%request_field(k3, ind(4), .false.)
635 call neko_scratch_registry%request_field(k4, ind(5), .false.)
636 call neko_scratch_registry%request_vector(u1_gl, ind(6), n_gl, .false.)
642 if (neko_bcknd_device .eq. 1)
then
645 call device_invcol3(u1%x_d, phi%x_d, coef%B_d, n)
646 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
648 conv_k1%items(2)%ptr, conv_k1%items(3)%ptr, &
649 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
650 call device_col2(k1%x_d, coef%B_d, n)
653 call device_add3s2(u1%x_d, phi%x_d, k1%x_d, c1, c2, n)
654 call device_invcol2(u1%x_d, coef%B_d, n)
655 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
657 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
658 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
659 call device_col2(k2%x_d, coef%B_d, n)
662 call device_add3s2(u1%x_d, phi%x_d, k2%x_d, c1, c2, n)
663 call device_invcol2(u1%x_d, coef%B_d, n)
664 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
666 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
667 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
668 call device_col2(k3%x_d, coef%B_d, n)
671 call device_add3s2(u1%x_d, phi%x_d, k3%x_d, c1, c3, n)
672 call device_invcol2(u1%x_d, coef%B_d, n)
673 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
675 conv_k4%items(2)%ptr, conv_k4%items(3)%ptr, &
676 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
677 call device_col2(k4%x_d, coef%B_d, n)
682 call device_add5s4(phi%x_d, k1%x_d, k2%x_d, k3%x_d, k4%x_d, &
688 call invcol3(u1%x, phi%x, coef%B, n)
689 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
691 conv_k1%items(2)%ptr, conv_k1%items(3)%ptr, &
692 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
693 call col2(k1%x, coef%B, n)
696 call add3s2(u1%x, phi%x, k1%x, c1, c2, n)
697 call invcol2(u1%x, coef%B, n)
698 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
700 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
701 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
702 call col2(k2%x, coef%B, n)
705 call add3s2(u1%x, phi%x, k2%x, c1, c2, n)
706 call invcol2(u1%x, coef%B, n)
707 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
709 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
710 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
711 call col2(k3%x, coef%B, n)
714 call add3s2(u1%x, phi%x, k3%x, c1, c3, n)
715 call invcol2(u1%x, coef%B, n)
716 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
718 conv_k4%items(2)%ptr, conv_k4%items(3)%ptr, &
719 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
720 call col2(k4%x, coef%B, n)
724 call add5s4(phi%x, k1%x, k2%x, k3%x, k4%x, c1, c2, c2, c1, n)
727 call neko_scratch_registry%relinquish_field(ind)
732 real(kind=rp),
dimension(:),
intent(inout) :: vx, vy, vz
733 integer,
intent(in) :: idir
734 type(coef_t),
intent(in) :: coef
735 if (coef%cyclic)
then
736 if (neko_bcknd_device .eq. 1)
then
737 call opr_device_rotate_cyc_r1(vx, vy, vz, idir, coef)
739 call opr_cpu_rotate_cyc_r1(vx, vy, vz, idir, coef)
745 real(kind=rp),
dimension(:,:,:,:),
intent(inout) :: vx, vy, vz
746 integer,
intent(in) :: idir
747 type(coef_t),
intent(in) :: coef
748 if (coef%cyclic)
then
749 if (neko_bcknd_device .eq. 1)
then
750 call opr_device_rotate_cyc_r4(vx, vy, vz, idir, coef)
752 call opr_cpu_rotate_cyc_r4(vx, vy, vz, idir, coef)
Return the device pointer for an associated Fortran array.
Map a Fortran array to a device (allocate and associate)
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
type(mpi_comm), public neko_comm
MPI communicator.
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n, strm)
Returns .
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_invcol3(a_d, b_d, c_d, n, strm)
Vector division .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_add5s4(a_d, b_d, c_d, d_d, e_d, c1, c2, c3, c4, n, strm)
Returns .
subroutine, public device_invcol2(a_d, b_d, n, strm)
Vector division .
Device abstraction, common interface for various accelerators.
subroutine, public device_free(x_d)
Deallocate memory on the device.
subroutine, public field_rzero(a, n)
Zero a real vector.
Routines to interpolate between different spaces.
A simulation component that computes lambda2 The values are stored in the field registry under the na...
subroutine, public cmult(a, c, n)
Multiplication by constant c .
subroutine, public invcol2(a, b, n)
Vector division .
subroutine, public cadd(a, s, n)
Add a scalar to vector .
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
subroutine, public add2(a, b, n)
Vector addition .
subroutine, public invcol3(a, b, c, n)
Invert a vector .
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public add5s4(a, b, c, d, e, c1, c2, c3, c4, n)
Returns .
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter neko_bcknd_sx
integer, parameter neko_bcknd_device
logical, parameter neko_device_mpi
integer, parameter neko_bcknd_xsmm
integer, parameter, public i8
integer, parameter, public rp
Global precision used in computations.
subroutine, public set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
Transforms the convecting velocity field to the rst form of the GL space.
subroutine, public ortho(x, glb_n_points, n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
real(kind=rp) function, public cfl(dt, u, v, w, xh, coef, nelv, gdim)
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
subroutine rotate_cyc_r1(vx, vy, vz, idir, coef)
subroutine convect_scalar(du, u, cr, cs, ct, xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
Apply the convecting velocity c to the to the scalar field u, used in the OIFS scheme.
subroutine rotate_cyc_r4(vx, vy, vz, idir, coef)
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
real(kind=rp) function, public cfl_compressible(dt, max_wave_speed, xh, coef, nelv, gdim)
subroutine, public conv1(du, u, vx, vy, vz, xh, coef, es, ee)
Compute the advection term.
subroutine, public grad(ux, uy, uz, u, coef)
Compute the gradient of a scalar field.
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef, event)
subroutine, public strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
subroutine, public lambda2op(lambda2, u, v, w, coef)
Compute the Lambda2 field for a given velocity field.
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
subroutine, public dudxyz(du, u, dr, ds, dt, coef)
Compute derivative of a scalar field along a single direction.
subroutine, public runge_kutta(phi, conv_k1, conv_k23, conv_k4, xh_gll, xh_gl, coef, coef_gl, gll_to_gl, tau, dtau, n, nel, n_gl)
Compute one step of Runge Kutta time interpolation for OIFS scheme.
subroutine, public opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
subroutine, public opr_cpu_rotate_cyc_r1(vx, vy, vz, idir, coef)
real(kind=rp) function, public opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
subroutine, public opr_cpu_lambda2(lambda2, u, v, w, coef)
subroutine, public opr_cpu_rotate_cyc_r4(vx, vy, vz, idir, coef)
Operators accelerator backends.
subroutine, public opr_device_convect_scalar(du, u_d, cr_d, cs_d, ct_d, xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
subroutine, public opr_device_rotate_cyc_r1(vx, vy, vz, idir, coef)
subroutine, public opr_device_cdtp(dtx, x, dr, ds, dt, coef)
subroutine, public opr_device_rotate_cyc_r4(vx, vy, vz, idir, coef)
real(kind=rp) function, public opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
subroutine, public opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh, event)
subroutine, public opr_device_set_convect_rst(cr_d, cs_d, ct_d, cx_d, cy_d, cz_d, xh, coef)
subroutine, public opr_device_dudxyz(du, u, dr, ds, dt, coef)
subroutine, public opr_device_opgrad(ux, uy, uz, u, coef)
subroutine, public opr_device_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
subroutine, public opr_device_lambda2(lambda2, u, v, w, coef)
Operators SX-Aurora backend.
subroutine, public opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
Operators libxsmm backend.
subroutine, public opr_xsmm_convect_scalar(du, u, cr, cs, ct, xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
subroutine, public opr_xsmm_conv1(du, u, vx, vy, vz, xh, coef, nelv, gdim)
subroutine, public opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
subroutine, public opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_xh)
subroutine, public opr_xsmm_dudxyz(du, u, dr, ds, dt, coef)
subroutine, public opr_xsmm_opgrad(ux, uy, uz, u, coef)
subroutine, public opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a function space.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
field_list_t, To be able to group fields together
Interpolation between two space::space_t.
The function space for the SEM solution fields.