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

Operators. More...

Functions/Subroutines

subroutine, public dudxyz (du, u, dr, ds, dt, coef)
 Compute derivative of a scalar field along a single direction. More...
 
subroutine, public opgrad (ux, uy, uz, u, coef, es, ee)
 Compute the gradient of a scalar field. 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)
 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, 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 double the strain rate tensor, i.e 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...
 

Detailed Description

Operators.

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 
)

Apply D^T to a scalar field, where D is the derivative matrix.

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.
Note
This needs to be revised... the loop over n1,n2 is probably unesccssary

Definition at line 156 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 265 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 
)

Compute the advection term.

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 186 of file operators.f90.

Here is the call 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 233 of file operators.f90.

Here is the caller 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 
)

Compute derivative of a scalar field along a single direction.

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 69 of file operators.f90.

Here is the call graph for this function:
Here is the caller 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 
)

Compute the Lambda2 field for a given velocity field.

Parameters
lambda2Holds the computed Lambda2 field.
uThe x-velocity.
vThe y-velocity.
wthe z-velocity.
coefThe SEM coefficients.

Definition at line 374 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 
)

Compute the gradient of a scalar field.

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 99 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 
)

Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.

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 136 of file operators.f90.

Here is the call 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 
)

Compute double the strain rate tensor, i.e du_i/dx_j + du_j/dx_i.

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 301 of file operators.f90.

Here is the call graph for this function: