|
Neko 0.9.1
A portable framework for high-order spectral element flow simulations
|
Operators.
Functions/Subroutines | |
| subroutine, public | dudxyz (du, u, dr, ds, dt, coef) |
| Compute derivative of a scalar field along a single direction. | |
| subroutine, public | div (res, ux, uy, uz, coef) |
| Compute the divergence of a vector field. | |
| subroutine, public | grad (ux, uy, uz, u, coef) |
| Compute the gradient of a scalar field. | |
| 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 | 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 | conv1 (du, u, vx, vy, vz, xh, coef, es, ee) |
| Compute the advection term. | |
| 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 | 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. | |
| subroutine, public | lambda2op (lambda2, u, v, w, coef) |
| Compute the Lambda2 field for a given velocity field. | |
| 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 | 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 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 | ||
| ) |
| 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. |
| es | Starting element index, optional, defaults to 1. |
| ee | Ending element index, optional, defaults to nelv. |
Definition at line 229 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 | ||
| ) |
Definition at line 393 of file operators.f90.


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


|
private |
| du | Holds the result |
| c | The convecting velocity |
| u | The convected scalar field |
| Xh_GLL | The GLL space used in simulation |
| Xh_GL | The GL space used for dealiasing |
| coef | The coefficients of the original space in simulation |
| coef_GL | The coefficients of the GL space used for dealiasing |
| GLL_to_GL | the interpolator between the GLL and GL spaces |
Definition at line 325 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::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 | ||
| ) |
| res | Holds the resulting divergence values. |
| ux | The x component of the vector field. |
| uy | The y component of the vector field. |
| uz | The z component of the vector field. |
| coef | The SEM coefficients. |
Definition at line 100 of file operators.f90.

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


| 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 | ||
| ) |
| 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. |
Definition at line 145 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 | ||
| ) |
| 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 502 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 | ||
| ) |
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 170 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 | ||
| ) |
| 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 207 of file operators.f90.

| 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 | ||
| ) |
| phi | The iterpolated field |
| c_r1 | The covecting velocity for the first stage |
| c_r23 | The convecting velocity for the second and third stage |
| c_r4 | The convecting velocity for the fourth stage |
| Xh_GLL | The GLL space used in simulation |
| Xh_GL | The GL space used for dealiasing |
| coef | The coefficients of the original space in simulation |
| coef_GL | The coefficients of the GL space used for dealiasing |
| GLL_to_GL | the interpolator between the GLL and GL spaces |
| tau | The the starting time |
| dtau | The time step used for the Runge Kutta scheme |
| n | size of phi |
| nel | Total number of elements |
| n_GL | the size in the GL space |
Definition at line 560 of file operators.f90.


| 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 | ||
| ) |
| cr | convecting velocity in r-direction |
| cs | convecting velocity in s-direction |
| ct | convecting velocity in t-direction |
| cx | convecting velocity in x-direction |
| cy | convecting velocity in y-direction |
| cz | convecting velocity in z-direction |
| Xh | The GL space used for dealiasing |
| coef | The coeffiecients of the GL space used for dealiasing |
Definition at line 527 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 | ||
| ) |
| 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 429 of file operators.f90.

