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
91 subroutine dudxyz (du, u, dr, ds, dt, coef)
92 type(
coef_t),
intent(in),
target :: coef
93 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(inout) :: du
94 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(in) :: u, dr, ds, dt
97 call opr_sx_dudxyz(du, u, dr, ds, dt, coef)
103 call opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
114 subroutine div(res, ux, uy, uz, coef)
115 type(
coef_t),
intent(in),
target :: coef
116 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(inout) :: res
117 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(in) :: ux, uy, uz
129 call dudxyz(res, ux, coef%drdx, coef%dsdx, coef%dtdx, coef)
132 call dudxyz(work%x, uy, coef%drdy, coef%dsdy, coef%dtdy, coef)
136 call add2(res, work%x, work%size())
140 call dudxyz(work%x, uz, coef%drdz, coef%dsdz, coef%dtdz, coef)
144 call add2(res, work%x, work%size())
157 subroutine grad(ux, uy, uz, u, coef)
158 type(
coef_t),
intent(in) :: coef
159 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: ux
160 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uy
161 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uz
162 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: u
164 call dudxyz(ux, u, coef%drdx, coef%dsdx, coef%dtdx, coef)
165 call dudxyz(uy, u, coef%drdy, coef%dsdy, coef%dtdy, coef)
166 call dudxyz(uz, u, coef%drdz, coef%dsdz, coef%dtdz, coef)
182 subroutine opgrad(ux, uy, uz, u, coef, es, ee)
183 type(
coef_t),
intent(in) :: coef
184 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: ux
185 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uy
186 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uz
187 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: u
188 integer,
optional :: es, ee
189 integer :: eblk_start, eblk_end
191 if (
present(es))
then
197 if (
present(ee))
then
200 eblk_end = coef%msh%nelv
204 call opr_sx_opgrad(ux, uy, uz, u, coef)
210 call opr_cpu_opgrad(ux, uy, uz, u, coef, eblk_start, eblk_end)
219 subroutine ortho(x, glb_n_points, n)
220 integer,
intent(in) :: n
221 integer(kind=i8),
intent(in) :: glb_n_points
222 real(kind=
rp),
dimension(n),
intent(inout) :: x
227 call neko_log%deprecated(
'Operator: ortho, implicit device', &
228 'v2.0.0',
'Please call device_ortho instead.')
234 c =
glsum(x, n) / glb_n_points
251 subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
252 type(
coef_t),
intent(in) :: coef
253 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: dtx
254 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: x
255 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: dr
256 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: ds
257 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: dt
258 integer,
optional :: es, ee
259 integer :: eblk_start, eblk_end
261 if (
present(es))
then
267 if (
present(ee))
then
270 eblk_end = coef%msh%nelv
274 call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
280 call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
295 subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
296 type(
space_t),
intent(in) :: xh
297 type(
coef_t),
intent(in) :: coef
298 real(kind=
rp),
intent(inout) :: du(xh%lxyz, coef%msh%nelv)
299 real(kind=
rp),
intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
300 real(kind=
rp),
intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
301 real(kind=
rp),
intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
302 real(kind=
rp),
intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
303 integer,
optional :: es, ee
304 integer :: eblk_end, eblk_start
306 associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
307 if (
present(es))
then
313 if (
present(ee))
then
316 eblk_end = coef%msh%nelv
320 call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
326 call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
349 subroutine convect_scalar(du, u, cr, cs, ct, Xh_GLL, Xh_GL, coef_GLL, &
351 type(space_t),
intent(in) :: xh_gl
352 type(space_t),
intent(in) :: xh_gll
353 type(coef_t),
intent(in) :: coef_GLL
354 type(coef_t),
intent(in) :: coef_GL
355 type(interpolator_t),
intent(inout) :: GLL_to_GL
356 real(kind=rp),
intent(inout) :: &
357 du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
358 real(kind=rp),
intent(inout) :: &
359 u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
360 type(field_t),
intent(inout) :: cr, cs, ct
363 if (neko_bcknd_sx .eq. 1)
then
364 call opr_sx_convect_scalar(du, u, cr%x, cs%x, ct%x, &
365 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
366 else if (neko_bcknd_xsmm .eq. 1)
then
367 call opr_xsmm_convect_scalar(du, u, cr%x, cs%x, ct%x, &
368 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
369 else if (neko_bcknd_device .eq. 1)
then
370 u_d = device_get_ptr(u)
371 call opr_device_convect_scalar(du, u_d, cr%x_d, cs%x_d, ct%x_d, &
372 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
374 call opr_cpu_convect_scalar(du, u, cr%x, cs%x, ct%x, &
375 xh_gll, xh_gl, coef_gll, coef_gl, gll_to_gl)
390 subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef, event)
391 type(field_t),
intent(inout) :: w1
392 type(field_t),
intent(inout) :: w2
393 type(field_t),
intent(inout) :: w3
394 type(field_t),
intent(in) :: u1
395 type(field_t),
intent(in) :: u2
396 type(field_t),
intent(in) :: u3
397 type(field_t),
intent(inout) :: work1
398 type(field_t),
intent(inout) :: work2
399 type(coef_t),
intent(in) :: coef
400 type(c_ptr),
optional,
intent(inout) :: event
402 if (neko_bcknd_sx .eq. 1)
then
403 call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
404 else if (neko_bcknd_xsmm .eq. 1)
then
405 call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
406 else if (neko_bcknd_device .eq. 1)
then
407 if (
present(event))
then
408 call opr_device_curl(w1, w2, w3, u1, u2, u3, &
409 work1, work2, coef, event)
411 call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
414 call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
428 function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
429 type(space_t),
intent(in) :: xh
430 type(coef_t),
intent(in) :: coef
431 integer,
intent(in) :: nelv, gdim
432 real(kind=rp),
intent(in) :: dt
433 real(kind=rp),
dimension(Xh%lx, Xh%ly, Xh%lz, nelv),
intent(in) :: u, v, w
437 if (neko_bcknd_sx .eq. 1)
then
438 cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
439 else if (neko_bcknd_device .eq. 1)
then
440 cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
442 cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
445 if (.not. neko_device_mpi)
then
446 call mpi_allreduce(mpi_in_place,
cfl, 1, &
447 mpi_real_precision, mpi_max, neko_comm, ierr)
460 type(space_t),
intent(in) :: xh
461 type(coef_t),
intent(in) :: coef
462 integer,
intent(in) :: nelv, gdim
463 real(kind=rp),
intent(in) :: dt
464 real(kind=rp),
dimension(Xh%lx, Xh%ly, Xh%lz, nelv),
intent(in) :: &
469 xh, coef, nelv, gdim)
485 subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
486 type(field_t),
intent(in) :: u, v, w
487 type(coef_t),
intent(in) :: coef
488 real(kind=rp),
intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
489 real(kind=rp),
intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
490 real(kind=rp),
intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
491 real(kind=rp),
intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
492 real(kind=rp),
intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
493 real(kind=rp),
intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
495 type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
497 integer :: nelv, lxyz
499 if (neko_bcknd_device .eq. 1)
then
500 s11_d = device_get_ptr(s11)
501 s22_d = device_get_ptr(s22)
502 s33_d = device_get_ptr(s33)
503 s12_d = device_get_ptr(s12)
504 s23_d = device_get_ptr(s23)
505 s13_d = device_get_ptr(s13)
512 call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
513 call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
514 if (neko_bcknd_device .eq. 1)
then
515 call device_add2(s12_d, s11_d, nelv*lxyz)
517 call add2(s12, s11, nelv*lxyz)
520 call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
521 call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
522 if (neko_bcknd_device .eq. 1)
then
523 call device_add2(s13_d, s11_d, nelv*lxyz)
525 call add2(s13, s11, nelv*lxyz)
528 call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
529 call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
530 if (neko_bcknd_device .eq. 1)
then
531 call device_add2(s23_d, s11_d, nelv*lxyz)
533 call add2(s23, s11, nelv*lxyz)
536 call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
537 call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
538 call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
540 if (neko_bcknd_device .eq. 1)
then
541 call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
542 call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
543 call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
545 call cmult(s12, 0.5_rp, nelv*lxyz)
546 call cmult(s13, 0.5_rp, nelv*lxyz)
547 call cmult(s23, 0.5_rp, nelv*lxyz)
558 subroutine lambda2op(lambda2, u, v, w, coef)
559 type(coef_t),
intent(in) :: coef
560 type(field_t),
intent(inout) ::
lambda2
561 type(field_t),
intent(in) :: u, v, w
563 if (neko_bcknd_sx .eq. 1)
then
564 call opr_sx_lambda2(
lambda2, u, v, w, coef)
565 else if (neko_bcknd_device .eq. 1)
then
566 call opr_device_lambda2(
lambda2, u, v, w, coef)
568 call opr_cpu_lambda2(
lambda2, u, v, w, coef)
584 type(space_t),
intent(inout) :: xh
585 type(coef_t),
intent(inout) :: coef
586 type(field_t),
intent(inout) :: cr, cs, ct
587 real(kind=rp),
dimension(Xh%lxyz, coef%msh%nelv), &
588 intent(in) :: cx, cy, cz
589 type(c_ptr) :: cx_d, cy_d, cz_d
591 if (neko_bcknd_sx .eq. 1)
then
592 call opr_sx_set_convect_rst(cr%x, cs%x, ct%x, cx, cy, cz, xh, coef)
593 else if (neko_bcknd_xsmm .eq. 1)
then
594 call opr_xsmm_set_convect_rst(cr%x, cs%x, ct%x, cx, cy, cz, xh, coef)
595 else if (neko_bcknd_device .eq. 1)
then
596 cx_d = device_get_ptr(cx)
597 cy_d = device_get_ptr(cy)
598 cz_d = device_get_ptr(cz)
599 call opr_device_set_convect_rst(cr%x_d, cs%x_d, ct%x_d, &
600 cx_d, cy_d, cz_d, xh, coef)
602 call opr_cpu_set_convect_rst(cr%x, cs%x, ct%x, cx, cy, cz, xh, coef)
622 subroutine runge_kutta(phi, conv_k1, conv_k23, conv_k4, Xh_GLL, Xh_GL, &
623 coef, coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
624 type(space_t),
intent(in) :: xh_gll
625 type(space_t),
intent(inout) :: xh_gl
626 type(coef_t),
intent(in) :: coef
627 type(coef_t),
intent(inout) :: coef_gl
628 type(interpolator_t) :: gll_to_gl
629 real(kind=rp),
intent(inout) :: tau, dtau
630 integer,
intent(in) :: n, nel, n_gl
631 type(field_t),
intent(inout) :: phi
632 type(field_list_t) :: conv_k1, conv_k23, conv_k4
633 real(kind=rp) :: c1, c2, c3
634 type(field_t),
pointer :: u1, k1, k2, k3, k4
635 type(vector_t),
pointer :: u1_gl
636 integer :: ind(6), i, e
638 call neko_scratch_registry%request_field(u1, ind(1), .false.)
639 call neko_scratch_registry%request_field(k1, ind(2), .false.)
640 call neko_scratch_registry%request_field(k2, ind(3), .false.)
641 call neko_scratch_registry%request_field(k3, ind(4), .false.)
642 call neko_scratch_registry%request_field(k4, ind(5), .false.)
643 call neko_scratch_registry%request_vector(u1_gl, ind(6), n_gl, .false.)
649 if (neko_bcknd_device .eq. 1)
then
652 call device_invcol3(u1%x_d, phi%x_d, coef%B_d, n)
653 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
655 conv_k1%items(2)%ptr, conv_k1%items(3)%ptr, &
656 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
657 call device_col2(k1%x_d, coef%B_d, n)
660 call device_add3s2(u1%x_d, phi%x_d, k1%x_d, c1, c2, n)
661 call device_invcol2(u1%x_d, coef%B_d, n)
662 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
664 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
665 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
666 call device_col2(k2%x_d, coef%B_d, n)
669 call device_add3s2(u1%x_d, phi%x_d, k2%x_d, c1, c2, n)
670 call device_invcol2(u1%x_d, coef%B_d, n)
671 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
673 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
674 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
675 call device_col2(k3%x_d, coef%B_d, n)
678 call device_add3s2(u1%x_d, phi%x_d, k3%x_d, c1, c3, n)
679 call device_invcol2(u1%x_d, coef%B_d, n)
680 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
682 conv_k4%items(2)%ptr, conv_k4%items(3)%ptr, &
683 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
684 call device_col2(k4%x_d, coef%B_d, n)
689 call device_add5s4(phi%x_d, k1%x_d, k2%x_d, k3%x_d, k4%x_d, &
695 call invcol3(u1%x, phi%x, coef%B, n)
696 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
698 conv_k1%items(2)%ptr, conv_k1%items(3)%ptr, &
699 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
700 call col2(k1%x, coef%B, n)
703 call add3s2(u1%x, phi%x, k1%x, c1, c2, n)
704 call invcol2(u1%x, coef%B, n)
705 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
707 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
708 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
709 call col2(k2%x, coef%B, n)
712 call add3s2(u1%x, phi%x, k2%x, c1, c2, n)
713 call invcol2(u1%x, coef%B, n)
714 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
716 conv_k23%items(2)%ptr, conv_k23%items(3)%ptr, &
717 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
718 call col2(k3%x, coef%B, n)
721 call add3s2(u1%x, phi%x, k3%x, c1, c3, n)
722 call invcol2(u1%x, coef%B, n)
723 call gll_to_gl%map(u1_gl%x, u1%x, nel, xh_gl)
725 conv_k4%items(2)%ptr, conv_k4%items(3)%ptr, &
726 xh_gll, xh_gl, coef, coef_gl, gll_to_gl)
727 call col2(k4%x, coef%B, n)
731 call add5s4(phi%x, k1%x, k2%x, k3%x, k4%x, c1, c2, c2, c1, n)
734 call neko_scratch_registry%relinquish_field(ind)
739 real(kind=rp),
dimension(:),
intent(inout) :: vx, vy, vz
740 integer,
intent(in) :: idir
741 type(coef_t),
intent(in) :: coef
742 if (coef%cyclic)
then
743 if (neko_bcknd_device .eq. 1)
then
744 call opr_device_rotate_cyc_r1(vx, vy, vz, idir, coef)
746 call opr_cpu_rotate_cyc_r1(vx, vy, vz, idir, coef)
752 real(kind=rp),
dimension(:,:,:,:),
intent(inout) :: vx, vy, vz
753 integer,
intent(in) :: idir
754 type(coef_t),
intent(in) :: coef
755 if (coef%cyclic)
then
756 if (neko_bcknd_device .eq. 1)
then
757 call opr_device_rotate_cyc_r4(vx, vy, vz, idir, coef)
759 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...
type(log_t), public neko_log
Global log stream.
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.