PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singleboundaryconvective.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 ~SingleBoundaryConvective() = default;
37
38 // re-implementation of SingleBoundaryBase::destroy
39 virtual PetscErrorCode destroy();
40
41protected:
42 // implementation of SingleBoundaryBase::setGhostICsKernel
43 virtual PetscErrorCode setGhostICsKernel(const PetscReal &targetValue,
45
46 // implementation of SingleBoundaryBase::updateEqsKernel
47 virtual PetscErrorCode updateEqsKernel(const PetscReal &targetValue,
48 const PetscReal &dt,
50
61 PetscErrorCode (*kernel)(const PetscReal &targetValue, const PetscReal &dt,
62 const PetscReal &normal, const PetscReal &value,
64}; // SingleBoundaryConvective
65
66} // end of namespace boundary
67} // 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.
PetscReal normal
The direction of normal vector.
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 convective BC.
virtual PetscErrorCode setGhostICsKernel(const PetscReal &targetValue, type::GhostPointInfo &p)
The underlying kernel for setting initial values and equations.
PetscErrorCode(* kernel)(const PetscReal &targetValue, const PetscReal &dt, const PetscReal &normal, const PetscReal &value, type::GhostPointInfo &p)
Underlying kernel that will determined during runtime according to the location of the boundary and t...
virtual PetscErrorCode updateEqsKernel(const PetscReal &targetValue, const PetscReal &dt, type::GhostPointInfo &p)
Underlying kernel for updating the coefficients of the equation.
SingleBoundaryConvective(const type::Mesh &mesh, const type::BCLoc &loc, const type::Field &field, const PetscReal &value)
Constructor.
virtual ~SingleBoundaryConvective()=default
Default destructor.
virtual PetscErrorCode destroy()
Manually destroy data.
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