PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
|
Base (abstract) class for ghost points & BC on a single boundary. More...
#include <singleboundary.h>
Public Member Functions | |
SingleBoundaryBase (const type::Mesh &mesh, const type::BCLoc &loc, const type::Field &field, const type::BCType &type, const PetscReal &value) | |
Constructor. More... | |
SingleBoundaryBase ()=default | |
Default constructor. More... | |
virtual | ~SingleBoundaryBase () |
Default destructor. More... | |
virtual PetscErrorCode | destroy () |
Manually destroy data. More... | |
PetscErrorCode | setGhostICs (const Vec &vec) |
Set up the initial values of the ghost points. More... | |
PetscErrorCode | updateEqs (const Vec &vec, const PetscReal &dt) |
Modify the coefficients in the equation of ghost points. More... | |
PetscErrorCode | updateGhostValues (const Vec &vec) |
Update the values of ghost points using the equation. More... | |
PetscErrorCode | copyValues2LocalVec (Vec &lclVec) |
Copy the values of ghost points to a local Vec. More... | |
Public Attributes | |
PetscInt | dim |
Dimension. More... | |
type::BCLoc | loc |
The location of this boundary. More... | |
type::Field | field |
The field of which the ghost points represent. More... | |
type::BCType | type |
The type of boundary conditions. More... | |
PetscReal | value |
A constant value representing BC value. More... | |
PetscReal | normal |
The direction of normal vector. More... | |
type::GhostPointsList | points |
The list of ghost points on this boundary and at this field. More... | |
PetscBool | onThisProc |
Indicate if this process holds part of this boundary. More... | |
Protected Member Functions | |
PetscErrorCode | init (const type::Mesh &mesh, const type::BCLoc &loc, const type::Field &field, const type::BCType &type, const PetscReal &bcValue) |
Underlying initialization function. More... | |
virtual PetscErrorCode | setGhostICsKernel (const PetscReal &targetValue, type::GhostPointInfo &p)=0 |
The underlying kernel for setting initial values and equations. More... | |
virtual PetscErrorCode | updateEqsKernel (const PetscReal &targetValue, const PetscReal &dt, type::GhostPointInfo &p)=0 |
Underlying kernel for updating the coefficients of the equation. More... | |
Protected Attributes | |
MPI_Comm | comm |
MPI communicator. More... | |
PetscMPIInt | mpiSize |
The size of the MPI communicator. More... | |
PetscMPIInt | mpiRank |
The rank of this process. More... | |
type::Mesh | mesh |
The corresponding Mesh object. More... | |
Base (abstract) class for ghost points & BC on a single boundary.
Definition at line 30 of file singleboundary.h.
petibm::boundary::SingleBoundaryBase::SingleBoundaryBase | ( | const type::Mesh & | mesh, |
const type::BCLoc & | loc, | ||
const type::Field & | field, | ||
const type::BCType & | type, | ||
const PetscReal & | value | ||
) |
Constructor.
mesh | [in] a Mesh instance. |
loc | [in] the location of the target boundary. |
field | [in] the target field. |
type | [in] the type of BC. |
value | [in] BC value. |
Definition at line 16 of file singleboundarybase.cpp.
|
default |
Default constructor.
|
virtual |
Default destructor.
Definition at line 25 of file singleboundarybase.cpp.
PetscErrorCode petibm::boundary::SingleBoundaryBase::copyValues2LocalVec | ( | Vec & | lclVec | ) |
Copy the values of ghost points to a local Vec.
lclVec | [in, out] a local Vec with memory allocations for ghost points. |
In PetIBM, we use a global packed Vec for velocity fields. But in some occasions, we may need local Vecs that have ghost points in them. So we have to copy the ghost values from this instance to the local Vecs.
Definition at line 165 of file singleboundarybase.cpp.
|
virtual |
Manually destroy data.
Reimplemented in petibm::boundary::SingleBoundaryConvective.
Definition at line 37 of file singleboundarybase.cpp.
|
protected |
Underlying initialization function.
mesh | [in] a Mesh instance. |
loc | [in] the location of the target boundary. |
field | [in] the target field. |
type | [in] the type of BC. |
bcValue | [in] BC value. |
Definition at line 55 of file singleboundarybase.cpp.
PetscErrorCode petibm::boundary::SingleBoundaryBase::setGhostICs | ( | const Vec & | vec | ) |
Set up the initial values of the ghost points.
vec | [in] a packed solution Vec containing initial values. |
Definition at line 107 of file singleboundarybase.cpp.
|
protectedpure virtual |
The underlying kernel for setting initial values and equations.
targetValue | [in] the value of the corresponding boundary point. |
p | [in, out] the target ghost point. |
Implemented in petibm::boundary::SingleBoundaryConvective, petibm::boundary::SingleBoundaryDirichlet, petibm::boundary::SingleBoundaryNeumann, and petibm::boundary::SingleBoundaryPeriodic.
PetscErrorCode petibm::boundary::SingleBoundaryBase::updateEqs | ( | const Vec & | vec, |
const PetscReal & | dt | ||
) |
Modify the coefficients in the equation of ghost points.
vec | [in] a packed solution Vec at current time step. |
dt | [in] a PetscReal representing time-step size. |
The equation of ghost points means the relationship between boundary points and ghost points. The equation has a form u_ghost = a0 * u_boundary + a1. This function changes a1 and a0 according to the type of BC.
Definition at line 126 of file singleboundarybase.cpp.
|
protectedpure virtual |
Underlying kernel for updating the coefficients of the equation.
targetValue | [in] the value of corresponding boundary point. |
dt | [in] the size of a time step. |
p | [in, out] the target ghost point. |
Implemented in petibm::boundary::SingleBoundaryConvective, petibm::boundary::SingleBoundaryDirichlet, petibm::boundary::SingleBoundaryNeumann, and petibm::boundary::SingleBoundaryPeriodic.
PetscErrorCode petibm::boundary::SingleBoundaryBase::updateGhostValues | ( | const Vec & | vec | ) |
Update the values of ghost points using the equation.
vec | [in] a packed solution Vec at current time step. |
Definition at line 146 of file singleboundarybase.cpp.
|
protected |
MPI communicator.
Definition at line 154 of file singleboundary.h.
PetscInt petibm::boundary::SingleBoundaryBase::dim |
Dimension.
Definition at line 34 of file singleboundary.h.
type::Field petibm::boundary::SingleBoundaryBase::field |
The field of which the ghost points represent.
Definition at line 40 of file singleboundary.h.
type::BCLoc petibm::boundary::SingleBoundaryBase::loc |
The location of this boundary.
Definition at line 37 of file singleboundary.h.
|
protected |
The corresponding Mesh object.
Definition at line 163 of file singleboundary.h.
|
protected |
The rank of this process.
Definition at line 160 of file singleboundary.h.
|
protected |
The size of the MPI communicator.
Definition at line 157 of file singleboundary.h.
PetscReal petibm::boundary::SingleBoundaryBase::normal |
The direction of normal vector.
Definition at line 49 of file singleboundary.h.
PetscBool petibm::boundary::SingleBoundaryBase::onThisProc |
Indicate if this process holds part of this boundary.
Definition at line 55 of file singleboundary.h.
type::GhostPointsList petibm::boundary::SingleBoundaryBase::points |
The list of ghost points on this boundary and at this field.
Definition at line 52 of file singleboundary.h.
type::BCType petibm::boundary::SingleBoundaryBase::type |
The type of boundary conditions.
Definition at line 43 of file singleboundary.h.
PetscReal petibm::boundary::SingleBoundaryBase::value |
A constant value representing BC value.
Definition at line 46 of file singleboundary.h.