PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
bodypack.h
Go to the documentation of this file.
1
8#pragma once
9
10#include <memory>
11#include <string>
12#include <vector>
13
14#include <petscdm.h>
15#include <petscsys.h>
16#include <yaml-cpp/yaml.h>
17
18#include <petibm/singlebody.h>
19
51namespace petibm
52{
58namespace body
59{
71{
72public:
74 PetscInt dim;
75
77 PetscInt nBodies;
78
80 PetscInt nPts;
81
83 PetscInt nLclPts;
84
86 std::vector<type::SingleBody> bodies;
87
90
92 std::string info;
93
95 BodyPackBase() = default;
96
104 BodyPackBase(const MPI_Comm &comm, const PetscInt &dim,
105 const YAML::Node &node);
106
108 virtual ~BodyPackBase();
109
111 virtual PetscErrorCode destroy();
112
118 PetscErrorCode printInfo() const;
119
130 PetscErrorCode findProc(const PetscInt &bIdx, const PetscInt &ptIdx,
131 PetscMPIInt &proc) const;
132
143 PetscErrorCode getGlobalIndex(const PetscInt &bIdx, const PetscInt &ptIdx,
144 const PetscInt &dof, PetscInt &idx) const;
145
155 PetscErrorCode getGlobalIndex(const PetscInt &bIdx, const MatStencil &s,
156 PetscInt &idx) const;
157
168 PetscErrorCode getPackedGlobalIndex(const PetscInt &bIdx,
169 const PetscInt &ptIdx,
170 const PetscInt &dof,
171 PetscInt &idx) const;
172
182 PetscErrorCode getPackedGlobalIndex(const PetscInt &bIdx,
183 const MatStencil &s,
184 PetscInt &idx) const;
185
197 PetscErrorCode calculateAvgForces(const Vec &f, type::RealVec2D &fAvg);
198
206 PetscErrorCode updateMeshIdx(const type::Mesh &mesh);
207
208protected:
218 PetscErrorCode init(const MPI_Comm &comm, const PetscInt &dim,
219 const YAML::Node &node);
220
222 MPI_Comm comm;
223
225 PetscMPIInt mpiSize;
226
228 PetscMPIInt mpiRank;
229
232
235
241 PetscErrorCode createDmPack();
242
248 PetscErrorCode createInfoString();
249
250}; // BodyPackBase
251
252} // end of namespace body
253
254namespace type
255{
262typedef std::shared_ptr<body::BodyPackBase> BodyPack;
263
264} // end of namespace type
265
266namespace body
267{
281PetscErrorCode createBodyPack(const MPI_Comm &comm, const PetscInt &dim,
282 const YAML::Node &node, type::BodyPack &bodies);
283
284} // end of namespace body
285
286} // end of namespace petibm
Base (abstract) class for a pack of bodies.
Definition: bodypack.h:71
PetscInt nLclPts
Total number of local Lagrangian points.
Definition: bodypack.h:83
PetscMPIInt mpiSize
Total number of processes.
Definition: bodypack.h:225
std::string info
String used to print information.
Definition: bodypack.h:92
PetscErrorCode printInfo() const
Print information about the bodies.
Definition: bodypack.cpp:195
BodyPackBase()=default
Default constructor.
type::IntVec1D nLclAllProcs
Number of local packed variables of all processes.
Definition: bodypack.h:231
virtual PetscErrorCode destroy()
Manually destroy data.
Definition: bodypack.cpp:35
PetscErrorCode getPackedGlobalIndex(const PetscInt &bIdx, const PetscInt &ptIdx, const PetscInt &dof, PetscInt &idx) const
Find packed global index of a DoF of Lagrangian point of a body.
Definition: bodypack.cpp:261
PetscErrorCode createInfoString()
Create a string used to print information.
Definition: bodypack.cpp:174
MPI_Comm comm
Reference to the MPI communicator.
Definition: bodypack.h:222
PetscInt dim
Dimension.
Definition: bodypack.h:74
PetscErrorCode calculateAvgForces(const Vec &f, type::RealVec2D &fAvg)
Calculate the averaged force of each body.
Definition: bodypack.cpp:298
virtual ~BodyPackBase()
Default destructor.
Definition: bodypack.cpp:21
PetscErrorCode createDmPack()
Create DMComposite object for all bodies.
Definition: bodypack.cpp:144
PetscMPIInt mpiRank
Rank of the local process.
Definition: bodypack.h:228
std::vector< type::SingleBody > bodies
Vector of SingleBody objects.
Definition: bodypack.h:86
type::IntVec1D offsetsAllProcs
Offsets of packed variables of all processes.
Definition: bodypack.h:234
PetscErrorCode updateMeshIdx(const type::Mesh &mesh)
Get the index of closest Eulerian mesh cell for each local Lagrangian point.
Definition: bodypack.cpp:160
PetscErrorCode getGlobalIndex(const PetscInt &bIdx, const PetscInt &ptIdx, const PetscInt &dof, PetscInt &idx) const
Find unpacked global index of a DoF of Lagrangian point of a body.
Definition: bodypack.cpp:229
DM dmPack
DMComposite of DMDA objects of all SingleBody objects.
Definition: bodypack.h:89
PetscErrorCode findProc(const PetscInt &bIdx, const PetscInt &ptIdx, PetscMPIInt &proc) const
Find which process owns the target Lagrangian point of target body.
Definition: bodypack.cpp:211
PetscInt nPts
Total number of Lagrangian points.
Definition: bodypack.h:80
PetscErrorCode init(const MPI_Comm &comm, const PetscInt &dim, const YAML::Node &node)
Initialize the pack of bodies.
Definition: bodypack.cpp:55
PetscInt nBodies
Number of bodies in this pack.
Definition: bodypack.h:77
PetscErrorCode createBodyPack(const MPI_Comm &comm, const PetscInt &dim, const YAML::Node &node, type::BodyPack &bodies)
Factory function to create a pack of bodies.
Definition: bodypack.cpp:326
std::shared_ptr< body::BodyPackBase > BodyPack
Definition of type::BodyPack.
Definition: bodypack.h:262
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
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
A toolbox for building flow solvers.
Definition: bodypack.h:52
body::SingleBodyBase, type::SingleBody factory function.