PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
|
Collections of factory functions of operators. More...
Typedefs | |
typedef std::function< PetscReal(const PetscInt &, const PetscInt &, const PetscInt &)> | KernelType |
a short name for the signatures of kernels. More... | |
typedef std::vector< MatStencil > | StencilVec |
STL vector holding MatStencil. More... | |
typedef std::function< StencilVec(const PetscInt &, const PetscInt &, const PetscInt &)> | GetNeighborFunc |
a function type for functions returning neighbors' stencils. More... | |
typedef std::function< type::RealVec1D(const PetscInt &, const PetscInt &, const PetscInt &)> | Kernel |
a function type for functions computing matrix entries' values. More... | |
Functions | |
PetscErrorCode | createR (const type::Mesh &mesh, Mat &R) |
Create a matrix of flux areas, . More... | |
PetscErrorCode | createRInv (const type::Mesh &mesh, Mat &RInv) |
Create a matrix of inversed flux areas, . More... | |
PetscErrorCode | createMHead (const type::Mesh &mesh, Mat &MHead) |
Create a matrix of cell widths of velocity points, . More... | |
PetscErrorCode | createM (const type::Mesh &mesh, Mat &M) |
Create a matrix . More... | |
PetscErrorCode | createIdentity (const type::Mesh &mesh, Mat &I) |
Create an identity matrix for velocity fields. More... | |
PetscErrorCode | createGradient (const type::Mesh &mesh, Mat &G, const PetscBool &normalize=PETSC_TRUE) |
Create a gradient operator, , for pressure field. More... | |
PetscErrorCode | createDivergence (const type::Mesh &mesh, const type::Boundary &bc, Mat &D, Mat &DCorrection, const PetscBool &normalize=PETSC_TRUE) |
Create a divergence operator, , and corresponding boundary correction, , for velocity fields. More... | |
PetscErrorCode | createLaplacian (const type::Mesh &mesh, const type::Boundary &bc, Mat &L, Mat &LCorrection) |
Create a Laplacian operator, , and boundary correction, for velocity fields. More... | |
PetscErrorCode | createConvection (const type::Mesh &mesh, const type::Boundary &bc, Mat &H) |
Create a matrix-free Mat for convection operator, . More... | |
PetscErrorCode | createBnHead (const Mat &Op, const PetscReal &dt, const PetscReal &coeff, const PetscInt &n, Mat &BnHead) |
Create non-normalized matrix of approximated , . More... | |
PetscErrorCode | createBn (const Mat &Op, const Mat &R, const Mat &MHead, const PetscReal &dt, const PetscReal &coeff, const PetscInt &n, Mat &Bn) |
Create normalized matrix of approximated , . More... | |
PetscErrorCode | createBn (const Mat &Op, const Mat &M, const PetscReal &dt, const PetscReal &coeff, const PetscInt &n, Mat &Bn) |
Create normalized matrix of approximated , . More... | |
PetscErrorCode | createDelta (const type::Mesh &mesh, const type::Boundary &bc, const type::BodyPack &bodies, const delta::DeltaKernel &kernel, const PetscInt &kernelSize, Mat &Op) |
Create a Delta operator, . More... | |
PetscErrorCode | getEulerianNeighbors (const type::Mesh &mesh, const PetscInt &dof, const type::IntVec1D &IJK, const std::vector< bool > &periodic, const PetscInt &window, type::IntVec2D &ijk, type::RealVec2D &xyz) |
PetscErrorCode | createDiagMatrix (const type::Mesh &mesh, const std::vector< KernelType > &kernels, Mat &M) |
a function to create diagonal matrix according to input kernel. More... | |
Collections of factory functions of operators.
typedef std::function<StencilVec(const PetscInt &, const PetscInt &, const PetscInt &)> petibm::operators::GetNeighborFunc |
a function type for functions returning neighbors' stencils.
Definition at line 28 of file creategradient.cpp.
typedef std::function<type::RealVec1D(const PetscInt &, const PetscInt &, const PetscInt &)> petibm::operators::Kernel |
a function type for functions computing matrix entries' values.
Definition at line 33 of file creategradient.cpp.
typedef std::function<PetscReal(const PetscInt &, const PetscInt &, const PetscInt &)> petibm::operators::KernelType |
a short name for the signatures of kernels.
A kernel is a function that can be used to calculate the diagonal values according to a given set of i, j, k index
Definition at line 29 of file creatediagmatrix.cpp.
typedef std::vector<MatStencil> petibm::operators::StencilVec |
STL vector holding MatStencil.
Definition at line 23 of file creategradient.cpp.
PetscErrorCode petibm::operators::createDiagMatrix | ( | const type::Mesh & | mesh, |
const std::vector< KernelType > & | kernels, | ||
Mat & | M | ||
) |
a function to create diagonal matrix according to input kernel.
mesh | an instance of CartesianMesh. |
kernels | a length 3 STL vector holding KernelType. |
M | the returned matrix. |
This is not designed for public use. It's only valid for the functions defined in this source file.
Definition at line 43 of file creatediagmatrix.cpp.
PetscErrorCode petibm::operators::getEulerianNeighbors | ( | const type::Mesh & | mesh, |
const PetscInt & | dof, | ||
const type::IntVec1D & | IJK, | ||
const std::vector< bool > & | periodic, | ||
const PetscInt & | window, | ||
type::IntVec2D & | ijk, | ||
type::RealVec2D & | xyz | ||
) |
Definition at line 171 of file createdelta.cpp.