PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singleboundarydirichlet.h
Go to the documentation of this file.
1
8#pragma once
9
10// here goes headers from our PetIBM
12
13namespace petibm
14{
15namespace boundary
16{
23{
24public:
33 const type::Field &field, const PetscReal &value);
34
36 virtual ~SingleBoundaryDirichlet() = default;
37
38protected:
39 // implementation of SingleBoundaryBase::setGhostICsKernel
40 virtual PetscErrorCode setGhostICsKernel(const PetscReal &targetValue,
42
43 // implementation of SingleBoundaryBase::updateEqsKernel
44 virtual PetscErrorCode updateEqsKernel(const PetscReal &targetValue,
45 const PetscReal &dt,
47
48}; // SingleBoundaryDirichlet
49
50} // end of namespace boundary
51} // end of namespace petibm
Base (abstract) class for ghost points & BC on a single boundary.
type::Field field
The field of which the ghost points represent.
type::BCLoc loc
The location of this boundary.
type::Mesh mesh
The corresponding Mesh object.
PetscReal value
A constant value representing BC value.
An implementation of SingleBoundaryBase for Dirichlet BC.
virtual PetscErrorCode setGhostICsKernel(const PetscReal &targetValue, type::GhostPointInfo &p)
The underlying kernel for setting initial values and equations.
virtual ~SingleBoundaryDirichlet()=default
Default destructor.
SingleBoundaryDirichlet(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.
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
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of the class SingleBoundaryBase.
A data structure for a single ghost point.
Definition: type.h:161