PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singleboundaryconvective.cpp
Go to the documentation of this file.
1
9
10namespace // anonymous namespace for internal linkage
11{
12PetscErrorCode kernelConvectiveDiffDir(const PetscReal &targetValue,
13 const PetscReal &dt,
14 const PetscReal &normal,
15 const PetscReal &value,
17{
18 PetscFunctionBeginUser;
19
20 p.a0 = -1.0;
21 p.a1 = p.value + targetValue -
22 2.0 * normal * dt * value * (p.value - targetValue) / p.dL;
23
24 PetscFunctionReturn(0);
25} // kernelConvectiveDiffDir
26
27PetscErrorCode kernelConvectiveSameDir(const PetscReal &targetValue,
28 const PetscReal &dt,
29 const PetscReal &normal,
30 const PetscReal &value,
32{
33 PetscFunctionBeginUser;
34
35 p.a0 = 0.0;
36 p.a1 = p.value - normal * dt * value * (p.value - targetValue) / p.dL;
37
38 PetscFunctionReturn(0);
39} // kernelConvectiveSameDir
40} // end of anonymous namespace
41
42namespace petibm
43{
44namespace boundary
45{
47 const type::BCLoc &inLoc,
48 const type::Field &inField,
49 const PetscReal &inValue)
50 : SingleBoundaryBase(inMesh, inLoc, inField, type::CONVECTIVE, inValue)
51{
52 PetscInt dir = int(loc) / 2;
53
54 if (dir == int(field))
55 kernel = &kernelConvectiveSameDir;
56 else
57 kernel = &kernelConvectiveDiffDir;
58} // SingleBoundaryConvective
59
61{
62 PetscFunctionBeginUser;
63 PetscErrorCode ierr;
64 kernel = nullptr;
65 ierr = SingleBoundaryBase::destroy(); CHKERRQ(ierr);
66 PetscFunctionReturn(0);
67} // destroy
68
70 const PetscReal &targetValue, type::GhostPointInfo &p)
71{
72 PetscFunctionBeginUser;
73
74 PetscErrorCode ierr;
75
76 // at beginning (t=0), due to lack of information of previous solution,
77 // we just assume the ghost point has the same value as target boundary
78 // point does
79 p.value = targetValue;
80
81 // due to p.value=targetValue, it doesn't matter how big is time-step size
82 ierr = kernel(targetValue, 0.0, normal, value, p); CHKERRQ(ierr);
83
84 PetscFunctionReturn(0);
85} // setGhostICsKernel
86
88 const PetscReal &targetValue, const PetscReal &dt, type::GhostPointInfo &p)
89{
90 PetscFunctionBeginUser;
91
92 PetscErrorCode ierr;
93
94 ierr = kernel(targetValue, dt, normal, value, p); CHKERRQ(ierr);
95
96 PetscFunctionReturn(0);
97} // updateEqsKernel
98
99} // end of namespace boundary
100} // end of namespace petibm
Base (abstract) class for ghost points & BC on a single boundary.
virtual PetscErrorCode destroy()
Manually destroy data.
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.
PetscReal value
A constant value representing BC value.
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 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
@ CONVECTIVE
Definition: type.h:99
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of the class SingleBoundaryConvective.
A data structure for a single ghost point.
Definition: type.h:161