Neko
0.8.1
A portable framework for high-order spectral element flow simulations
|
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... | |
Operators.
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.
dtx | Will store the result. |
x | The values of the field. |
dr | The derivative of r with respect to the chosen direction. |
ds | The derivative of s with respect to the chosen direction. |
dt | The derivative of t with respect to the chosen direction. |
coef | The SEM coefficients. |
Definition at line 156 of file operators.f90.
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 | ||
) |
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.
du | Holds the result. |
u | The advected field. |
vx | The x component of the advecting velocity. |
vy | The y component of the advecting velocity. |
vz | The z component of the advecting velocity. |
Xh | The function space for the fields involved. |
coef | The SEM coefficients. |
es | Starting element index, defaults to 1. |
ee | Last element index, defaults to mesh size. |
Definition at line 186 of file operators.f90.
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 | ||
) |
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.
du | Holds the resulting derivative values. |
u | The values of the field. |
dr | The derivative of r with respect to the chosen direction. |
ds | The derivative of s with respect to the chosen direction. |
dt | The derivative of t with respect to the chosen direction. |
coef | The SEM coefficients. |
Definition at line 69 of file operators.f90.
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.
lambda2 | Holds the computed Lambda2 field. |
u | The x-velocity. |
v | The y-velocity. |
w | the z-velocity. |
coef | The SEM coefficients. |
Definition at line 374 of file operators.f90.
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.
ux | Will store the x component of the gradient. |
uy | Will store the y component of the gradient. |
uz | Will store the z component of the gradient. |
u | The values of the field. |
coef | The SEM coefficients. |
es | Starting element index, optional, defaults to 1. |
ee | Ending element index, optional, defaults to nelv . |
Definition at line 99 of file operators.f90.
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.
x | The vector to orthogonolize. |
n | The size of x . |
glb_n | The global number of elements of x across all MPI ranks. Be careful with overflow! |
Definition at line 136 of file operators.f90.
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.
s11 | Will hold the 1,1 component of the strain rate tensor. |
s22 | Will hold the 2,2 component of the strain rate tensor. |
s33 | Will hold the 3,3 component of the strain rate tensor. |
s12 | Will hold the 1,2 component of the strain rate tensor. |
s13 | Will hold the 1,3 component of the strain rate tensor. |
s23 | Will hold the 2,3 component of the strain rate tensor. |
u | The x component of velocity. |
v | The y component of velocity. |
w | The z component of velocity. |
coef | The SEM coefficients. |
[in] | w | velocity components |
Definition at line 301 of file operators.f90.