An implementation of SingleBoundaryBase for convective BC.
More...
#include <singleboundaryconvective.h>
◆ SingleBoundaryConvective()
petibm::boundary::SingleBoundaryConvective::SingleBoundaryConvective |
( |
const type::Mesh & |
mesh, |
|
|
const type::BCLoc & |
loc, |
|
|
const type::Field & |
field, |
|
|
const PetscReal & |
value |
|
) |
| |
Constructor.
- Parameters
-
mesh | [in] a Mesh instance. |
loc | [in] the location of the target boundary. |
field | [in] the target field. |
value | [in] BC value. |
Definition at line 46 of file singleboundaryconvective.cpp.
◆ ~SingleBoundaryConvective()
virtual petibm::boundary::SingleBoundaryConvective::~SingleBoundaryConvective |
( |
| ) |
|
|
virtualdefault |
◆ copyValues2LocalVec()
PetscErrorCode petibm::boundary::SingleBoundaryBase::copyValues2LocalVec |
( |
Vec & |
lclVec | ) |
|
|
inherited |
Copy the values of ghost points to a local Vec.
- Parameters
-
lclVec | [in, out] a local Vec with memory allocations for ghost points. |
- Returns
- PetscErrorCode.
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.
◆ destroy()
PetscErrorCode petibm::boundary::SingleBoundaryConvective::destroy |
( |
| ) |
|
|
virtual |
◆ init()
Underlying initialization function.
- Parameters
-
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. |
- Returns
- PetscErrorCode.
Definition at line 55 of file singleboundarybase.cpp.
◆ setGhostICs()
PetscErrorCode petibm::boundary::SingleBoundaryBase::setGhostICs |
( |
const Vec & |
vec | ) |
|
|
inherited |
Set up the initial values of the ghost points.
- Parameters
-
vec | [in] a packed solution Vec containing initial values. |
- Returns
- PetscErrorCode.
Definition at line 107 of file singleboundarybase.cpp.
◆ setGhostICsKernel()
PetscErrorCode petibm::boundary::SingleBoundaryConvective::setGhostICsKernel |
( |
const PetscReal & |
targetValue, |
|
|
type::GhostPointInfo & |
p |
|
) |
| |
|
protectedvirtual |
◆ updateEqs()
PetscErrorCode petibm::boundary::SingleBoundaryBase::updateEqs |
( |
const Vec & |
vec, |
|
|
const PetscReal & |
dt |
|
) |
| |
|
inherited |
Modify the coefficients in the equation of ghost points.
- Parameters
-
vec | [in] a packed solution Vec at current time step. |
dt | [in] a PetscReal representing time-step size. |
- Returns
- PetscErrorCode.
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.
◆ updateEqsKernel()
PetscErrorCode petibm::boundary::SingleBoundaryConvective::updateEqsKernel |
( |
const PetscReal & |
targetValue, |
|
|
const PetscReal & |
dt, |
|
|
type::GhostPointInfo & |
p |
|
) |
| |
|
protectedvirtual |
Underlying kernel for updating the coefficients of the equation.
- Parameters
-
targetValue | [in] the value of corresponding boundary point. |
dt | [in] the size of a time step. |
p | [in, out] the target ghost point. |
- Returns
- PetscErrorCode.
Implements petibm::boundary::SingleBoundaryBase.
Definition at line 87 of file singleboundaryconvective.cpp.
◆ updateGhostValues()
PetscErrorCode petibm::boundary::SingleBoundaryBase::updateGhostValues |
( |
const Vec & |
vec | ) |
|
|
inherited |
Update the values of ghost points using the equation.
- Parameters
-
vec | [in] a packed solution Vec at current time step. |
- Returns
- PetscErrorCode.
Definition at line 146 of file singleboundarybase.cpp.
◆ comm
MPI_Comm petibm::boundary::SingleBoundaryBase::comm |
|
protectedinherited |
◆ dim
PetscInt petibm::boundary::SingleBoundaryBase::dim |
|
inherited |
◆ field
type::Field petibm::boundary::SingleBoundaryBase::field |
|
inherited |
The field of which the ghost points represent.
Definition at line 40 of file singleboundary.h.
◆ kernel
PetscErrorCode(* petibm::boundary::SingleBoundaryConvective::kernel) (const PetscReal &targetValue, const PetscReal &dt, const PetscReal &normal, const PetscReal &value, type::GhostPointInfo &p) |
|
protected |
Underlying kernel that will determined during runtime according to the location of the boundary and the target field in a staggered grid.
- Parameters
-
targetValue | [in] The value at time of the corresponding boundary point. |
dt | [in] Time-step size. |
normal | [in] Normal of the boundary: represent *PLUS faces, while means *MINUS faces. |
value | [in] The value of the boundary condition. |
p | [in, out] A petibm::type::GhostPointInfo. |
- Returns
- PetscErrorCode.
Definition at line 61 of file singleboundaryconvective.h.
◆ loc
◆ mesh
type::Mesh petibm::boundary::SingleBoundaryBase::mesh |
|
protectedinherited |
◆ mpiRank
PetscMPIInt petibm::boundary::SingleBoundaryBase::mpiRank |
|
protectedinherited |
◆ mpiSize
PetscMPIInt petibm::boundary::SingleBoundaryBase::mpiSize |
|
protectedinherited |
◆ normal
PetscReal petibm::boundary::SingleBoundaryBase::normal |
|
inherited |
◆ onThisProc
PetscBool petibm::boundary::SingleBoundaryBase::onThisProc |
|
inherited |
Indicate if this process holds part of this boundary.
Definition at line 55 of file singleboundary.h.
◆ points
The list of ghost points on this boundary and at this field.
Definition at line 52 of file singleboundary.h.
◆ type
◆ value
PetscReal petibm::boundary::SingleBoundaryBase::value |
|
inherited |
The documentation for this class was generated from the following files: