39 opr_cpu_conv1, opr_cpu_convect_scalar, opr_cpu_cdtp, &
42 opr_sx_conv1, opr_sx_convect_scalar, opr_sx_cdtp, &
43 opr_sx_dudxyz, opr_sx_lambda2, opr_sx_set_convect_rst
75 subroutine dudxyz (du, u, dr, ds, dt, coef)
76 type(
coef_t),
intent(in),
target :: coef
77 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(inout) :: du
78 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(in) :: u, dr, ds, dt
81 call opr_sx_dudxyz(du, u, dr, ds, dt, coef)
87 call opr_cpu_dudxyz(du, u, dr, ds, dt, coef)
98 subroutine div(res, ux, uy, uz, coef)
99 type(
coef_t),
intent(in),
target :: coef
100 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(inout) :: res
101 real(kind=
rp),
dimension(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
coef%msh%nelv),
intent(in) :: ux, uy, uz
113 call dudxyz(res, ux, coef%drdx, coef%dsdx, coef%dtdx, coef)
116 call dudxyz(work%x, uy, coef%drdy, coef%dsdy, coef%dtdy, coef)
120 call add2(res, work%x, work%size())
124 call dudxyz(work%x, uz, coef%drdz, coef%dsdz, coef%dtdz, coef)
128 call add2(res, work%x, work%size())
141 subroutine grad(ux, uy, uz, u, coef)
142 type(
coef_t),
intent(in) :: coef
143 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: ux
144 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uy
145 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uz
146 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: u
148 call dudxyz(ux, u, coef%drdx, coef%dsdx, coef%dtdx, coef)
149 call dudxyz(uy, u, coef%drdy, coef%dsdy, coef%dtdy, coef)
150 call dudxyz(uz, u, coef%drdz, coef%dsdz, coef%dtdz, coef)
166 subroutine opgrad(ux, uy, uz, u, coef, es, ee)
167 type(
coef_t),
intent(in) :: coef
168 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: ux
169 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uy
170 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: uz
171 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: u
172 integer,
optional :: es, ee
173 integer :: eblk_start, eblk_end
175 if (
present(es))
then
181 if (
present(ee))
then
184 eblk_end = coef%msh%nelv
188 call opr_sx_opgrad(ux, uy, uz, u, coef)
194 call opr_cpu_opgrad(ux, uy, uz, u, coef, eblk_start, eblk_end)
203 subroutine ortho(x, n, glb_n)
204 integer,
intent(in) :: n
205 integer,
intent(in) :: glb_n
206 real(kind=
rp),
dimension(n),
intent(inout) :: x
207 real(kind=
rp) :: rlam
209 rlam =
glsum(x, n)/glb_n
210 call cadd(x, -rlam, n)
225 subroutine cdtp (dtx, x, dr, ds, dt, coef, es, ee)
226 type(
coef_t),
intent(in) :: coef
227 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: dtx
228 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(inout) :: x
229 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: dr
230 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: ds
231 real(kind=
rp),
dimension(coef%Xh%lxyz, coef%msh%nelv),
intent(in) :: dt
232 integer,
optional :: es, ee
233 integer :: eblk_start, eblk_end
235 if (
present(es))
then
241 if (
present(ee))
then
244 eblk_end = coef%msh%nelv
248 call opr_sx_cdtp(dtx, x, dr, ds, dt, coef)
254 call opr_cpu_cdtp(dtx, x, dr, ds, dt, coef, eblk_start, eblk_end)
269 subroutine conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
270 type(
space_t),
intent(inout) :: xh
271 type(
coef_t),
intent(inout) :: coef
272 real(kind=
rp),
intent(inout) :: du(xh%lxyz, coef%msh%nelv)
273 real(kind=
rp),
intent(inout) :: u(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
274 real(kind=
rp),
intent(inout) :: vx(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
275 real(kind=
rp),
intent(inout) :: vy(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
276 real(kind=
rp),
intent(inout) :: vz(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
277 integer,
optional :: es, ee
278 integer :: eblk_end, eblk_start
280 associate(nelv => coef%msh%nelv, gdim => coef%msh%gdim)
281 if (
present(es))
then
287 if (
present(ee))
then
290 eblk_end = coef%msh%nelv
294 call opr_sx_conv1(du, u, vx, vy, vz, xh, coef, nelv)
300 call opr_cpu_conv1(du, u, vx, vy, vz, xh, coef, eblk_start, eblk_end)
323 type(space_t),
intent(in) :: xh_gl
324 type(space_t),
intent(in) :: xh_gll
325 type(coef_t),
intent(in) :: coef_GLL
326 type(coef_t),
intent(in) :: coef_GL
327 type(interpolator_t),
intent(inout) :: GLL_to_GL
328 real(kind=rp),
intent(inout) :: &
329 du(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv)
330 real(kind=rp),
intent(inout) :: &
331 u(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv)
332 real(kind=rp),
intent(inout) :: c(xh_gl%lxyz, coef_gl%msh%nelv, 3)
334 if (neko_bcknd_sx .eq. 1)
then
335 call opr_sx_convect_scalar(du, u, c, xh_gll, xh_gl, &
336 coef_gll, coef_gl, gll_to_gl)
337 else if (neko_bcknd_xsmm .eq. 1)
then
338 call opr_xsmm_convect_scalar(du, u, c, xh_gll, xh_gl, &
339 coef_gll, coef_gl, gll_to_gl)
341 call opr_cpu_convect_scalar(du, u, c, xh_gll, xh_gl, &
342 coef_gll, coef_gl, gll_to_gl)
357 subroutine curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
358 type(field_t),
intent(inout) :: w1
359 type(field_t),
intent(inout) :: w2
360 type(field_t),
intent(inout) :: w3
361 type(field_t),
intent(inout) :: u1
362 type(field_t),
intent(inout) :: u2
363 type(field_t),
intent(inout) :: u3
364 type(field_t),
intent(inout) :: work1
365 type(field_t),
intent(inout) :: work2
366 type(coef_t),
intent(in) :: coef
368 if (neko_bcknd_sx .eq. 1)
then
369 call opr_sx_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
370 else if (neko_bcknd_xsmm .eq. 1)
then
371 call opr_xsmm_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
372 else if (neko_bcknd_device .eq. 1)
then
373 call opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
375 call opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, coef)
389 function cfl(dt, u, v, w, Xh, coef, nelv, gdim)
390 type(space_t),
intent(in) :: xh
391 type(coef_t),
intent(in) :: coef
392 integer,
intent(in) :: nelv, gdim
393 real(kind=rp),
intent(in) :: dt
394 real(kind=rp),
dimension(Xh%lx, Xh%ly, Xh%lz, nelv),
intent(in) :: u, v, w
398 if (neko_bcknd_sx .eq. 1)
then
399 cfl = opr_sx_cfl(dt, u, v, w, xh, coef, nelv)
400 else if (neko_bcknd_device .eq. 1)
then
401 cfl = opr_device_cfl(dt, u, v, w, xh, coef, nelv, gdim)
403 cfl = opr_cpu_cfl(dt, u, v, w, xh, coef, nelv, gdim)
406 if (.not. neko_device_mpi)
then
407 call mpi_allreduce(mpi_in_place,
cfl, 1, &
408 mpi_real_precision, mpi_max, neko_comm, ierr)
425 subroutine strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
426 type(field_t),
intent(in) :: u, v, w
427 type(coef_t),
intent(in) :: coef
428 real(kind=rp),
intent(inout) :: s11(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
429 real(kind=rp),
intent(inout) :: s22(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
430 real(kind=rp),
intent(inout) :: s33(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
431 real(kind=rp),
intent(inout) :: s12(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
432 real(kind=rp),
intent(inout) :: s13(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
433 real(kind=rp),
intent(inout) :: s23(u%Xh%lx, u%Xh%ly, u%Xh%lz, u%msh%nelv)
435 type(c_ptr) :: s11_d, s22_d, s33_d, s12_d, s23_d, s13_d
437 integer :: nelv, lxyz
439 if (neko_bcknd_device .eq. 1)
then
440 s11_d = device_get_ptr(s11)
441 s22_d = device_get_ptr(s22)
442 s33_d = device_get_ptr(s33)
443 s12_d = device_get_ptr(s12)
444 s23_d = device_get_ptr(s23)
445 s13_d = device_get_ptr(s13)
452 call dudxyz (s12, u%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
453 call dudxyz (s11, v%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
454 if (neko_bcknd_device .eq. 1)
then
455 call device_add2(s12_d, s11_d, nelv*lxyz)
457 call add2(s12, s11, nelv*lxyz)
460 call dudxyz (s13, u%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
461 call dudxyz (s11, w%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
462 if (neko_bcknd_device .eq. 1)
then
463 call device_add2(s13_d, s11_d, nelv*lxyz)
465 call add2(s13, s11, nelv*lxyz)
468 call dudxyz (s23, v%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
469 call dudxyz (s11, w%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
470 if (neko_bcknd_device .eq. 1)
then
471 call device_add2(s23_d, s11_d, nelv*lxyz)
473 call add2(s23, s11, nelv*lxyz)
476 call dudxyz (s11, u%x, coef%drdx, coef%dsdx, coef%dtdx, coef)
477 call dudxyz (s22, v%x, coef%drdy, coef%dsdy, coef%dtdy, coef)
478 call dudxyz (s33, w%x, coef%drdz, coef%dsdz, coef%dtdz, coef)
480 if (neko_bcknd_device .eq. 1)
then
481 call device_cmult(s12_d, 0.5_rp, nelv*lxyz)
482 call device_cmult(s13_d, 0.5_rp, nelv*lxyz)
483 call device_cmult(s23_d, 0.5_rp, nelv*lxyz)
485 call cmult(s12, 0.5_rp, nelv*lxyz)
486 call cmult(s13, 0.5_rp, nelv*lxyz)
487 call cmult(s23, 0.5_rp, nelv*lxyz)
498 subroutine lambda2op(lambda2, u, v, w, coef)
499 type(coef_t),
intent(in) :: coef
500 type(field_t),
intent(inout) ::
lambda2
501 type(field_t),
intent(in) :: u, v, w
503 if (neko_bcknd_sx .eq. 1)
then
504 call opr_sx_lambda2(
lambda2, u, v, w, coef)
505 else if (neko_bcknd_device .eq. 1)
then
506 call opr_device_lambda2(
lambda2, u, v, w, coef)
508 call opr_cpu_lambda2(
lambda2, u, v, w, coef)
524 type(space_t),
intent(inout) :: xh
525 type(coef_t),
intent(inout) :: coef
526 real(kind=rp),
dimension(Xh%lxyz, coef%msh%nelv), &
527 intent(inout) :: cr, cs, ct
528 real(kind=rp),
dimension(Xh%lxyz, coef%msh%nelv), &
529 intent(in) :: cx, cy, cz
531 if (neko_bcknd_sx .eq. 1)
then
532 call opr_sx_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
533 else if (neko_bcknd_xsmm .eq. 1)
then
534 call opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
536 call opr_cpu_set_convect_rst(cr, cs, ct, cx, cy, cz, xh, coef)
556 subroutine runge_kutta(phi, c_r1, c_r23, c_r4, Xh_GLL, Xh_GL, coef, &
557 coef_GL, GLL_to_GL, tau, dtau, n, nel, n_GL)
558 type(space_t),
intent(inout) :: xh_gll
559 type(space_t),
intent(inout) :: xh_gl
560 type(coef_t),
intent(inout) :: coef
561 type(coef_t),
intent(inout) :: coef_gl
562 type(interpolator_t) :: gll_to_gl
563 real(kind=rp),
intent(inout) :: tau, dtau
564 integer,
intent(in) :: n, nel, n_gl
565 real(kind=rp),
dimension(n),
intent(inout) :: phi
566 real(kind=rp),
dimension(3 * n_GL),
intent(inout) :: c_r1, c_r23, c_r4
567 real(kind=rp) :: c1, c2, c3
569 real(kind=rp),
dimension(n) :: u1, r1, r2, r3, r4
570 real(kind=rp),
dimension(n_GL) :: u1_gl
578 call invcol3 (u1, phi, coef%B, n)
579 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
582 call col2(r1, coef%B, n)
585 call add3s2 (u1, phi, r1, c1, c2, n)
586 call invcol2 (u1, coef%B, n)
587 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
590 call col2(r2, coef%B, n)
593 call add3s2 (u1, phi, r2, c1, c2, n)
594 call invcol2 (u1, coef%B, n)
595 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
598 call col2(r3, coef%B, n)
601 call add3s2 (u1, phi, r3, c1, c3, n)
602 call invcol2 (u1, coef%B, n)
603 call gll_to_gl%map(u1_gl, u1, nel, xh_gl)
606 call col2(r4, coef%B, n)
611 phi(i) = phi(i) + c1 * (r1(i) + r4(i)) + c2 * (r2(i) + r3(i))
617 Return the device pointer for an associated Fortran array.
subroutine, public device_add2(a_d, b_d, n)
Vector addition .
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
Device abstraction, common interface for various accelerators.
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 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 rp
Global precision used in computations.
subroutine, public runge_kutta(phi, c_r1, c_r23, c_r4, 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 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, 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 div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
subroutine convect_scalar(du, u, c, 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, 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)
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 ortho(x, n, glb_n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
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.
real(kind=rp) function, public cfl(dt, u, v, w, Xh, coef, nelv, gdim)
subroutine, public conv1(du, u, vx, vy, vz, Xh, coef, es, ee)
Compute the advection term.
subroutine, public opr_cpu_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
subroutine, public opr_cpu_lambda2(lambda2, u, v, w, coef)
real(kind=rp) function, public opr_cpu_cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Operators accelerator backends.
real(kind=rp) function, public opr_device_cfl(dt, u, v, w, Xh, coef, nelv, gdim)
subroutine, public opr_device_cdtp(dtx, x, dr, ds, dt, coef)
subroutine, public opr_device_conv1(du, u, vx, vy, vz, Xh, coef, nelv, gdim)
subroutine, public opr_device_dudxyz(du, u, dr, ds, dt, coef)
subroutine, public opr_device_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
subroutine, public opr_device_opgrad(ux, uy, uz, u, coef)
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, c, Xh_GLL, Xh_GL, coef_GLL, coef_GL, GLL_to_GL)
subroutine, public opr_xsmm_cdtp(dtx, x, dr, ds, dt, coef)
subroutine, public opr_xsmm_set_convect_rst(cr, cs, ct, cx, cy, cz, Xh, coef)
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_curl(w1, w2, w3, u1, u2, u3, work1, work2, c_Xh)
subroutine, public opr_xsmm_conv1(du, u, vx, vy, vz, Xh, coef, nelv, gdim)
Defines a registry for storing and requesting temporary fields This can be used when you have a funct...
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,...
Interpolation between two space::space_t.
The function space for the SEM solution fields.