Neko  0.8.99
A portable framework for high-order spectral element flow simulations
operators Module Reference

Operators.

Functions/Subroutines

subroutine, public dudxyz (du, u, dr, ds, dt, coef)
 Compute derivative of a scalar field along a single direction. More...
 
subroutine, public div (res, ux, uy, uz, coef)
 Compute the divergence of a vector field. More...
 
subroutine, public grad (ux, uy, uz, u, coef)
 Compute the gradient of a scalar field. More...
 
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. More...
 
subroutine, public ortho (x, n, glb_n)
 Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T. More...
 
subroutine, public cdtp (dtx, x, dr, ds, dt, coef, es, ee)
 Apply D^T to a scalar field, where D is the derivative matrix. More...
 
subroutine, public conv1 (du, u, vx, vy, vz, Xh, coef, es, ee)
 Compute the advection term. More...
 
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. More...
 
subroutine, public curl (w1, w2, w3, u1, u2, u3, work1, work2, coef)
 
real(kind=rp) function, public cfl (dt, u, v, w, Xh, coef, nelv, gdim)
 
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. More...
 
subroutine, public lambda2op (lambda2, u, v, w, coef)
 Compute the Lambda2 field for a given velocity field. More...
 
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. More...
 
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. More...
 

Function/Subroutine Documentation

◆ cdtp()

subroutine, public operators::cdtp ( real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  dtx,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  x,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(in)  dr,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(in)  ds,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(in)  dt,
type(coef_t), intent(in)  coef,
integer, optional  es,
integer, optional  ee 
)
Parameters
dtxWill store the result.
xThe values of the field.
drThe derivative of r with respect to the chosen direction.
dsThe derivative of s with respect to the chosen direction.
dtThe derivative of t with respect to the chosen direction.
coefThe SEM coefficients.
esStarting element index, optional, defaults to 1.
eeEnding element index, optional, defaults to nelv.
Note
This needs to be revised... the loop over n1,n2 is probably unesccssary

Definition at line 227 of file operators.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cfl()

real(kind=rp) function, public operators::cfl ( real(kind=rp), intent(in)  dt,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, nelv), intent(in)  u,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, nelv), intent(in)  v,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, nelv), intent(in)  w,
type(space_t), intent(in)  Xh,
type(coef_t), intent(in)  coef,
integer, intent(in)  nelv,
integer, intent(in)  gdim 
)

Definition at line 391 of file operators.f90.

Here is the caller graph for this function:

◆ conv1()

subroutine, public operators::conv1 ( real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(inout)  du,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, coef%msh%nelv), intent(inout)  u,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, coef%msh%nelv), intent(inout)  vx,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, coef%msh%nelv), intent(inout)  vy,
real(kind=rp), dimension(xh%lx, xh%ly, xh%lz, coef%msh%nelv), intent(inout)  vz,
type(space_t), intent(inout)  Xh,
type(coef_t), intent(inout)  coef,
integer, optional  es,
integer, optional  ee 
)
Parameters
duHolds the result.
uThe advected field.
vxThe x component of the advecting velocity.
vyThe y component of the advecting velocity.
vzThe z component of the advecting velocity.
XhThe function space for the fields involved.
coefThe SEM coefficients.
esStarting element index, defaults to 1.
eeLast element index, defaults to mesh size.

Definition at line 271 of file operators.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ convect_scalar()

subroutine operators::convect_scalar ( real(kind=rp), dimension(xh_gll%lx, xh_gll%ly, xh_gll%lz, coef_gl%msh%nelv), intent(inout)  du,
real(kind=rp), dimension(xh_gl%lx, xh_gl%lx, xh_gl%lx, coef_gl%msh%nelv), intent(inout)  u,
real(kind=rp), dimension(xh_gl%lxyz, coef_gl%msh%nelv, 3), intent(inout)  c,
type(space_t), intent(in)  Xh_GLL,
type(space_t), intent(in)  Xh_GL,
type(coef_t), intent(in)  coef_GLL,
type(coef_t), intent(in)  coef_GL,
type(interpolator_t), intent(inout)  GLL_to_GL 
)
private
Parameters
duHolds the result
cThe convecting velocity
uThe convected scalar field
Xh_GLLThe GLL space used in simulation
Xh_GLThe GL space used for dealiasing
coefThe coefficients of the original space in simulation
coef_GLThe coefficients of the GL space used for dealiasing
GLL_to_GLthe interpolator between the GLL and GL spaces
Note
This subroutine is equal to the convop_fst_3d of the NEK5000.
This subroutine is used specifically in the OIFS scheme, calculateing eq(17) in https://publications.anl.gov/anlpubs/2017/12/140626.pdf. The convecting term is calculated in the rst format and the GL grid. Then converted back to the GLL grid, going through an ADD gatter scatter operation at the element boundaries, before being multiplied by inverse of mass matrix.

Definition at line 323 of file operators.f90.

Here is the caller graph for this function:

◆ curl()

subroutine, public operators::curl ( type(field_t), intent(inout)  w1,
type(field_t), intent(inout)  w2,
type(field_t), intent(inout)  w3,
type(field_t), intent(inout)  u1,
type(field_t), intent(inout)  u2,
type(field_t), intent(inout)  u3,
type(field_t), intent(inout)  work1,
type(field_t), intent(inout)  work2,
type(coef_t), intent(in)  coef 
)

Definition at line 359 of file operators.f90.

Here is the caller graph for this function:

◆ div()

subroutine, public operators::div ( real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(inout)  res,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  ux,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  uy,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  uz,
type(coef_t), intent(in), target  coef 
)
Parameters
resHolds the resulting divergence values.
uxThe x component of the vector field.
uyThe y component of the vector field.
uzThe z component of the vector field.
coefThe SEM coefficients.

Definition at line 100 of file operators.f90.

Here is the call graph for this function:

◆ dudxyz()

subroutine, public operators::dudxyz ( real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(inout)  du,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  u,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  dr,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  ds,
real(kind=rp), dimension(coef%xh%lx, coef%xh%ly, coef%xh%lz, coef%msh%nelv), intent(in)  dt,
type(coef_t), intent(in), target  coef 
)
Parameters
duHolds the resulting derivative values.
uThe values of the field.
drThe derivative of r with respect to the chosen direction.
dsThe derivative of s with respect to the chosen direction.
dtThe derivative of t with respect to the chosen direction.
coefThe SEM coefficients.

Definition at line 75 of file operators.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ grad()

subroutine, public operators::grad ( real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  ux,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  uy,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  uz,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(in)  u,
type(coef_t), intent(in)  coef 
)
Parameters
uxWill store the x component of the gradient.
uyWill store the y component of the gradient.
uzWill store the z component of the gradient.
uThe values of the field.
coefThe SEM coefficients.

Definition at line 143 of file operators.f90.

Here is the call graph for this function:

◆ lambda2op()

subroutine, public operators::lambda2op ( type(field_t), intent(inout)  lambda2,
type(field_t), intent(in)  u,
type(field_t), intent(in)  v,
type(field_t), intent(in)  w,
type(coef_t), intent(in)  coef 
)
Parameters
lambda2Holds the computed Lambda2 field.
uThe x-velocity.
vThe y-velocity.
wthe z-velocity.
coefThe SEM coefficients.

Definition at line 500 of file operators.f90.

Here is the caller graph for this function:

◆ opgrad()

subroutine, public operators::opgrad ( real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  ux,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  uy,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(inout)  uz,
real(kind=rp), dimension(coef%xh%lxyz, coef%msh%nelv), intent(in)  u,
type(coef_t), intent(in)  coef,
integer, optional  es,
integer, optional  ee 
)

By providing es and ee, it is possible to compute only for a range of element indices.

Parameters
uxWill store the x component of the gradient.
uyWill store the y component of the gradient.
uzWill store the z component of the gradient.
uThe values of the field.
coefThe SEM coefficients.
esStarting element index, optional, defaults to 1.
eeEnding element index, optional, defaults to nelv.
Note
Equals wgradm1 in Nek5000, the weak form of the gradient.

Definition at line 168 of file operators.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ortho()

subroutine, public operators::ortho ( real(kind=rp), dimension(n), intent(inout)  x,
integer, intent(in)  n,
integer, intent(in)  glb_n 
)
Parameters
xThe vector to orthogonolize.
nThe size of x.
glb_nThe global number of elements of x across all MPI ranks. Be careful with overflow!

Definition at line 205 of file operators.f90.

Here is the call graph for this function:

◆ runge_kutta()

subroutine, public operators::runge_kutta ( real(kind=rp), dimension(n), intent(inout)  phi,
real(kind=rp), dimension(3 * n_gl), intent(inout)  c_r1,
real(kind=rp), dimension(3 * n_gl), intent(inout)  c_r23,
real(kind=rp), dimension(3 * n_gl), intent(inout)  c_r4,
type(space_t), intent(inout)  Xh_GLL,
type(space_t), intent(inout)  Xh_GL,
type(coef_t), intent(inout)  coef,
type(coef_t), intent(inout)  coef_GL,
type(interpolator_t GLL_to_GL,
real(kind=rp), intent(inout)  tau,
real(kind=rp), intent(inout)  dtau,
integer, intent(in)  n,
integer, intent(in)  nel,
integer, intent(in)  n_GL 
)
Parameters
phiThe iterpolated field
c_r1The covecting velocity for the first stage
c_r23The convecting velocity for the second and third stage
c_r4The convecting velocity for the fourth stage
Xh_GLLThe GLL space used in simulation
Xh_GLThe GL space used for dealiasing
coefThe coefficients of the original space in simulation
coef_GLThe coefficients of the GL space used for dealiasing
GLL_to_GLthe interpolator between the GLL and GL spaces
tauThe the starting time
dtauThe time step used for the Runge Kutta scheme
nsize of phi
nelTotal number of elements
n_GLthe size in the GL space

Definition at line 558 of file operators.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_convect_rst()

subroutine, public operators::set_convect_rst ( real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(inout)  cr,
real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(inout)  cs,
real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(inout)  ct,
real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(in)  cx,
real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(in)  cy,
real(kind=rp), dimension(xh%lxyz, coef%msh%nelv), intent(in)  cz,
type(space_t), intent(inout)  Xh,
type(coef_t), intent(inout)  coef 
)
Parameters
crconvecting velocity in r-direction
csconvecting velocity in s-direction
ctconvecting velocity in t-direction
cxconvecting velocity in x-direction
cyconvecting velocity in y-direction
czconvecting velocity in z-direction
XhThe GL space used for dealiasing
coefThe coeffiecients of the GL space used for dealiasing
Note
This subroutine is equal to the set_convect_new subroutine of NEK5000

Definition at line 525 of file operators.f90.

Here is the caller graph for this function:

◆ strain_rate()

subroutine, public operators::strain_rate ( real(kind=rp), dimension(u%xh%lx, u%xh%ly, u%xh%lz, u%msh%nelv), intent(inout)  s11,
real(kind=rp), dimension(u%xh%lx, u%xh%ly, u%xh%lz, u%msh%nelv), intent(inout)  s22,
real(kind=rp), dimension(u%xh%lx, u%xh%ly, u%xh%lz, u%msh%nelv), intent(inout)  s33,
real(kind=rp), dimension(u%xh%lx, u%xh%ly, u%xh%lz, u%msh%nelv), intent(inout)  s12,
real(kind=rp), dimension(u%xh%lx, u%xh%ly, u%xh%lz, u%msh%nelv), intent(inout)  s13,
real(kind=rp), dimension(u%xh%lx, u%xh%ly, u%xh%lz, u%msh%nelv), intent(inout)  s23,
type(field_t), intent(in)  u,
type(field_t), intent(in)  v,
type(field_t), intent(in)  w,
type(coef_t), intent(in)  coef 
)
Parameters
s11Will hold the 1,1 component of the strain rate tensor.
s22Will hold the 2,2 component of the strain rate tensor.
s33Will hold the 3,3 component of the strain rate tensor.
s12Will hold the 1,2 component of the strain rate tensor.
s13Will hold the 1,3 component of the strain rate tensor.
s23Will hold the 2,3 component of the strain rate tensor.
uThe x component of velocity.
vThe y component of velocity.
wThe z component of velocity.
coefThe SEM coefficients.
Note
Similar to comp_sij in Nek5000.
Parameters
[in]wvelocity components

Definition at line 427 of file operators.f90.

Here is the call graph for this function:
Here is the caller graph for this function: