PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singleboundaryneumann.cpp
Go to the documentation of this file.
1
9
10namespace petibm
11{
12namespace boundary
13{
15 const type::BCLoc &inLoc,
16 const type::Field &inField,
17 const PetscReal &inValue)
18 : SingleBoundaryBase(inMesh, inLoc, inField, type::NEUMANN, inValue)
19{
20} // SingleBoundaryNeumann
21
23 const PetscReal &targetValue, type::GhostPointInfo &p)
24{
25 PetscFunctionBeginUser;
26
27 p.a0 = 1.0;
28 p.a1 = normal * p.dL * value;
29
30 p.value = p.a0 * targetValue + p.a1;
31
32 PetscFunctionReturn(0);
33} // setGhostICsKernel
34
36 const PetscReal &targetValue, const PetscReal &dt, type::GhostPointInfo &p)
37{
38 PetscFunctionBeginUser;
39 // for time-independent Neumann BC, the coefficient a0 & a1 won't change
40 PetscFunctionReturn(0);
41} // updateEqsKernel
42
43} // end of namespace boundary
44} // end of namespace petibm
Base (abstract) class for ghost points & BC on a single boundary.
PetscReal normal
The direction of normal vector.
PetscReal value
A constant value representing BC value.
SingleBoundaryNeumann(const type::Mesh &mesh, const type::BCLoc &loc, const type::Field &field, const PetscReal &value)
Constructor.
virtual PetscErrorCode updateEqsKernel(const PetscReal &targetValue, const PetscReal &dt, type::GhostPointInfo &p)
Underlying kernel for updating the coefficients of the equation.
virtual PetscErrorCode setGhostICsKernel(const PetscReal &targetValue, type::GhostPointInfo &p)
The underlying kernel for setting initial values and equations.
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
Field
Legal fields.
Definition: type.h:80
BCLoc
Location of a boundary.
Definition: type.h:108
@ NEUMANN
Definition: type.h:98
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of the class SingleBoundaryNeumann.
A data structure for a single ghost point.
Definition: type.h:161