PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singleboundarybase.cpp
Go to the documentation of this file.
1
8// here goes headers from our PetIBM
9#include <petibm/misc.h>
11
12namespace petibm
13{
14namespace boundary
15{
17 const type::BCLoc &bcLoc,
18 const type::Field &inField,
19 const type::BCType &inType,
20 const PetscReal &inValue)
21{
22 init(inMesh, bcLoc, inField, inType, inValue);
23} // SingleBoundaryBase
24
26{
27 PetscFunctionBeginUser;
28 PetscErrorCode ierr;
29 PetscBool finalized;
30
31 ierr = PetscFinalized(&finalized); CHKERRV(ierr);
32 if (finalized) return;
33
34 comm = MPI_COMM_NULL;
35} // ~SingleBoundaryBase
36
38{
39 PetscFunctionBeginUser;
40
41 dim = -1;
42 loc = type::BCLoc(0);
43 field = type::Field(0);
44 value = normal = 0.0;
46 onThisProc = PETSC_FALSE;
47
48 comm = MPI_COMM_NULL;
49 mpiSize = mpiRank = 0;
50 mesh.reset();
51
52 PetscFunctionReturn(0);
53} // destroy
54
55PetscErrorCode SingleBoundaryBase::init(const type::Mesh &inMesh,
56 const type::BCLoc &bcLoc,
57 const type::Field &inField,
58 const type::BCType &inType,
59 const PetscReal &inValue)
60{
61 PetscFunctionBeginUser;
62
63 PetscErrorCode ierr;
64
65 // create a shared pointer to mesh; bad practice...
66 mesh = inMesh;
67
68 // obtain MPI information from CartesianMesh object
69 comm = mesh->comm;
70 mpiSize = mesh->mpiSize;
71 mpiRank = mesh->mpiRank;
72
73 // set dim
74 dim = mesh->dim;
75
76 // set the location
77 loc = bcLoc;
78
79 // set field
80 field = inField;
81
82 // set type
83 type = inType;
84
85 // set value
86 value = inValue;
87
88 // set normal
89 normal = ((int(loc) % 2) == 0) ? -1.0 : 1.0;
90
91 // set onThisProc
92 ierr = misc::checkBoundaryProc(mesh->da[int(field)], mesh->n[int(field)],
93 loc, onThisProc); CHKERRQ(ierr);
94
95 // for processes on this boundary, set up IDs and values
96 if (onThisProc)
97 {
99 CHKERRQ(ierr);
100 }
101
102 ierr = MPI_Barrier(comm); CHKERRQ(ierr);
103
104 PetscFunctionReturn(0);
105} // init
106
107PetscErrorCode SingleBoundaryBase::setGhostICs(const Vec &vec)
108{
109 PetscFunctionBeginUser;
110
111 PetscErrorCode ierr;
112
113 PetscReal targetValue;
114
115 for (auto &it : points)
116 {
117 ierr = VecGetValues(vec, 1, &(it.second.targetPackedId), &targetValue);
118 CHKERRQ(ierr);
119
120 ierr = setGhostICsKernel(targetValue, it.second); CHKERRQ(ierr);
121 }
122
123 PetscFunctionReturn(0);
124} // setGhostICs
125
126PetscErrorCode SingleBoundaryBase::updateEqs(const Vec &vec,
127 const PetscReal &dt)
128{
129 PetscFunctionBeginUser;
130
131 PetscErrorCode ierr;
132
133 PetscReal targetValue;
134
135 for (auto &it : points)
136 {
137 ierr = VecGetValues(vec, 1, &(it.second.targetPackedId), &targetValue);
138 CHKERRQ(ierr);
139
140 ierr = updateEqsKernel(targetValue, dt, it.second); CHKERRQ(ierr);
141 }
142
143 PetscFunctionReturn(0);
144} // updateEqs
145
146PetscErrorCode SingleBoundaryBase::updateGhostValues(const Vec &vec)
147{
148 PetscFunctionBeginUser;
149
150 PetscErrorCode ierr;
151
152 PetscReal targetValue;
153
154 for (auto &it : points)
155 {
156 ierr = VecGetValues(vec, 1, &(it.second.targetPackedId), &targetValue);
157 CHKERRQ(ierr);
158
159 it.second.value = it.second.a0 * targetValue + it.second.a1;
160 }
161
162 PetscFunctionReturn(0);
163} // updateGhostValues
164
166{
167 PetscFunctionBeginUser;
168
169 PetscErrorCode ierr;
170
171 for (auto &it : points)
172 {
173 ierr = VecSetValue(lclVec, it.second.lclId, it.second.value,
174 INSERT_VALUES); CHKERRQ(ierr);
175 }
176
177 PetscFunctionReturn(0);
178} // copyValues2LocalVec
179
180} // end of namespace boundary
181} // end of namespace petibm
type::BCType type
The type of boundary conditions.
PetscErrorCode updateEqs(const Vec &vec, const PetscReal &dt)
Modify the coefficients in the equation of ghost points.
virtual PetscErrorCode destroy()
Manually destroy data.
type::GhostPointsList points
The list of ghost points on this boundary and at this field.
virtual PetscErrorCode updateEqsKernel(const PetscReal &targetValue, const PetscReal &dt, type::GhostPointInfo &p)=0
Underlying kernel for updating the coefficients of the equation.
PetscBool onThisProc
Indicate if this process holds part of this boundary.
type::Field field
The field of which the ghost points represent.
SingleBoundaryBase()=default
Default constructor.
PetscErrorCode init(const type::Mesh &mesh, const type::BCLoc &loc, const type::Field &field, const type::BCType &type, const PetscReal &bcValue)
Underlying initialization function.
PetscErrorCode setGhostICs(const Vec &vec)
Set up the initial values of the ghost points.
PetscErrorCode updateGhostValues(const Vec &vec)
Update the values of ghost points using the equation.
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.
PetscMPIInt mpiSize
The size of the MPI communicator.
virtual PetscErrorCode setGhostICsKernel(const PetscReal &targetValue, type::GhostPointInfo &p)=0
The underlying kernel for setting initial values and equations.
MPI_Comm comm
MPI communicator.
virtual ~SingleBoundaryBase()
Default destructor.
PetscErrorCode copyValues2LocalVec(Vec &lclVec)
Copy the values of ghost points to a local Vec.
PetscMPIInt mpiRank
The rank of this process.
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
PetscErrorCode checkBoundaryProc(const DM &da, const type::IntVec1D &n, const type::BCLoc &loc, PetscBool &onThisProc)
Check if a boundary is on this process.
Definition: misc.cpp:85
PetscErrorCode getGhostPointList(const type::Mesh &mesh, const type::Field &field, const type::BCLoc &loc, type::GhostPointsList &points)
Get a list of ghost points on a desired boundary.
Definition: misc.cpp:129
Field
Legal fields.
Definition: type.h:80
BCType
Type of boundary conditions.
Definition: type.h:94
std::map< MatStencil, GhostPointInfo > GhostPointsList
A map between MatStencil and GhostPointInfo.
Definition: type.h:177
BCLoc
Location of a boundary.
Definition: type.h:108
Prototypes of some miscellaneous functions.
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of the class SingleBoundaryBase.