PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singlebody.h
Go to the documentation of this file.
1
8#pragma once
9
10#include <memory>
11
12#include <petscdm.h>
13#include <petscsys.h>
14
15#include <petibm/mesh.h>
16#include <petibm/type.h>
17
18namespace petibm
19{
20namespace body
21{
22// TODO: the way we store coords and meshIdx is different from CartesianMesh.
23// Will this cause confusion?
33{
34 friend class BodyPackBase;
35
36public:
38 PetscInt dim;
39
41 std::string name;
42
44 std::string filePath;
45
47 PetscInt nPts;
48
51
54
56 PetscInt nLclPts;
57
63
65 PetscInt bgPt;
66
68 PetscInt edPt;
69
71 DM da;
72
74 std::string info;
75
91 SingleBodyBase(const MPI_Comm &comm, const PetscInt &dim,
92 const std::string &name, const std::string &filePath);
93
95 SingleBodyBase() = default;
96
98 virtual ~SingleBodyBase();
99
101 virtual PetscErrorCode destroy();
102
108 PetscErrorCode printInfo() const;
109
121 virtual PetscErrorCode findProc(const PetscInt &i,
122 PetscMPIInt &p) const = 0;
123
134 virtual PetscErrorCode getGlobalIndex(const PetscInt &i,
135 const PetscInt &dof,
136 PetscInt &idx) const = 0;
137
147 virtual PetscErrorCode getGlobalIndex(const MatStencil &s,
148 PetscInt &idx) const = 0;
149
161 virtual PetscErrorCode calculateAvgForces(const Vec &f,
162 type::RealVec1D &fAvg) const = 0;
163
171 virtual PetscErrorCode updateMeshIdx(const type::Mesh &mesh) = 0;
172
173 virtual PetscErrorCode readBody(const std::string &filepath) = 0;
174
175 virtual PetscErrorCode writeBody(const std::string &filepath) = 0;
176
177
178protected:
180 MPI_Comm comm;
181
183 PetscMPIInt mpiSize;
184
186 PetscMPIInt mpiRank;
187
190
193
194}; // SingleBodyBase
195
196} // end of namespace body
197
198namespace type
199{
206typedef std::shared_ptr<body::SingleBodyBase> SingleBody;
207} // end of namespace type
208
209namespace body
210{
226PetscErrorCode createSingleBody(const MPI_Comm &comm, const PetscInt &dim,
227 const std::string &type,
228 const std::string &name,
229 const std::string &filePath,
230 type::SingleBody &body);
231
232} // end of namespace body
233
234} // end of namespace petibm
Base (abstract) class for a pack of bodies.
Definition: bodypack.h:71
Base (abstract) class for a single body.
Definition: singlebody.h:33
PetscInt nPts
Total number of Lagrangian points.
Definition: singlebody.h:47
std::string name
Name of the body.
Definition: singlebody.h:41
virtual PetscErrorCode calculateAvgForces(const Vec &f, type::RealVec1D &fAvg) const =0
Calculate the averaged force of this body.
virtual PetscErrorCode writeBody(const std::string &filepath)=0
PetscInt edPt
Global index of the last local Lagrangian point.
Definition: singlebody.h:68
MPI_Comm comm
MPI communicator.
Definition: singlebody.h:180
PetscMPIInt mpiRank
Rank of the local process.
Definition: singlebody.h:186
virtual PetscErrorCode readBody(const std::string &filepath)=0
type::RealVec2D coords0
Initial coordinates of ALL Lagrangian points.
Definition: singlebody.h:53
virtual PetscErrorCode findProc(const PetscInt &i, PetscMPIInt &p) const =0
Get the process id owning a given Lagrangian point.
PetscMPIInt mpiSize
Total number of processes.
Definition: singlebody.h:183
PetscInt bgPt
Global index of the first local Lagrangian point.
Definition: singlebody.h:65
PetscInt dim
Number of dimensions.
Definition: singlebody.h:38
virtual PetscErrorCode destroy()
Manually destroy data.
Definition: singlebody.cpp:53
DM da
Parallel layout of the body (as a 1D DMDA object).
Definition: singlebody.h:71
PetscErrorCode printInfo() const
Print information about the body.
Definition: singlebody.cpp:72
type::IntVec1D offsetsAllProcs
Offset on each process.
Definition: singlebody.h:192
virtual PetscErrorCode getGlobalIndex(const PetscInt &i, const PetscInt &dof, PetscInt &idx) const =0
Get the global index of a Lagrangian point in a DMDA object given its degree of freedom.
type::IntVec2D meshIdx
Index of the closest Eulerian mesh cell for each local Lagrangian point.
Definition: singlebody.h:62
virtual PetscErrorCode getGlobalIndex(const MatStencil &s, PetscInt &idx) const =0
Get the global index of a Lagrangian point in a DMDA object given a MatStencil.
type::IntVec1D nLclAllProcs
Vector with the number of local unknowns on each process.
Definition: singlebody.h:189
std::string filePath
Path of the file with the body coordinates.
Definition: singlebody.h:44
type::RealVec2D coords
Coordinates of ALL Lagrangian points.
Definition: singlebody.h:50
SingleBodyBase()=default
Default constructor.
std::string info
String with information about the body.
Definition: singlebody.h:74
virtual PetscErrorCode updateMeshIdx(const type::Mesh &mesh)=0
Get the index of closest Eulerian mesh cell for each local Lagrangian point.
PetscInt nLclPts
Local number of Lagrangian points.
Definition: singlebody.h:56
virtual ~SingleBodyBase()
Default destructor.
Definition: singlebody.cpp:39
PetscErrorCode createSingleBody(const MPI_Comm &comm, const PetscInt &dim, const std::string &type, const std::string &name, const std::string &filePath, type::SingleBody &body)
Factory function to create a single body.
Definition: singlebody.cpp:83
std::shared_ptr< body::SingleBodyBase > SingleBody
Definition of type::SingleBody.
Definition: singlebody.h:206
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
std::vector< IntVec1D > IntVec2D
2D std::vector holding PetscInt.
Definition: type.h:135
std::vector< RealVec1D > RealVec2D
2D std::vector holding PetscReal.
Definition: type.h:142
std::vector< PetscInt > IntVec1D
1D std::vector holding PetscInt.
Definition: type.h:133
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
Prototype of mesh::MeshBase, type::Mesh, and factory function.
A toolbox for building flow solvers.
Definition: bodypack.h:52
Definition of user-defined types for convenience.