194 const PetscBool &normalize = PETSC_TRUE);
224 const PetscBool &normalize = PETSC_TRUE);
246 Mat &L, Mat &LCorrection);
294PetscErrorCode
createBnHead(
const Mat &Op,
const PetscReal &dt,
295 const PetscReal &coeff,
const PetscInt &n,
317PetscErrorCode
createBn(
const Mat &Op,
const Mat &R,
const Mat &MHead,
318 const PetscReal &dt,
const PetscReal &coeff,
319 const PetscInt &n, Mat &Bn);
336PetscErrorCode
createBn(
const Mat &Op,
const Mat &M,
const PetscReal &dt,
337 const PetscReal &coeff,
const PetscInt &n, Mat &Bn);
364 const PetscInt &kernelSize,
body::BodyPackBase, type::BodyPack, and factory function.
boundary::BoundaryBase, type::Boundary, and factory function.
Prototype of Delta functions.
std::shared_ptr< body::BodyPackBase > BodyPack
Definition of type::BodyPack.
std::shared_ptr< boundary::BoundaryBase > Boundary
Type definition of petibm::type::Boundary.
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
PetscErrorCode createRInv(const type::Mesh &mesh, Mat &RInv)
Create a matrix of inversed flux areas, .
PetscErrorCode createMHead(const type::Mesh &mesh, Mat &MHead)
Create a matrix of cell widths of velocity points, .
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 , .
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.
PetscErrorCode createIdentity(const type::Mesh &mesh, Mat &I)
Create an identity matrix for velocity fields.
PetscErrorCode createGradient(const type::Mesh &mesh, Mat &G, const PetscBool &normalize=PETSC_TRUE)
Create a gradient operator, , for pressure field.
PetscErrorCode createLaplacian(const type::Mesh &mesh, const type::Boundary &bc, Mat &L, Mat &LCorrection)
Create a Laplacian operator, , and boundary correction, for velocity fields.
PetscErrorCode createM(const type::Mesh &mesh, Mat &M)
Create a matrix .
PetscErrorCode createBnHead(const Mat &Op, const PetscReal &dt, const PetscReal &coeff, const PetscInt &n, Mat &BnHead)
Create non-normalized matrix of approximated , .
PetscErrorCode createConvection(const type::Mesh &mesh, const type::Boundary &bc, Mat &H)
Create a matrix-free Mat for convection operator, .
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, .
PetscErrorCode createR(const type::Mesh &mesh, Mat &R)
Create a matrix of flux areas, .
Prototype of mesh::MeshBase, type::Mesh, and factory function.
std::function< PetscReal(const PetscReal &r, const PetscReal &dr)> DeltaKernel
Typedef to choose the regularized delta kernel to use.
A toolbox for building flow solvers.
Definition of user-defined types for convenience.