16#include <yaml-cpp/yaml.h>
86 std::vector<type::SingleBody>
bodies;
105 const YAML::Node &node);
111 virtual PetscErrorCode
destroy();
130 PetscErrorCode
findProc(
const PetscInt &bIdx,
const PetscInt &ptIdx,
131 PetscMPIInt &proc)
const;
143 PetscErrorCode
getGlobalIndex(
const PetscInt &bIdx,
const PetscInt &ptIdx,
144 const PetscInt &dof, PetscInt &idx)
const;
155 PetscErrorCode
getGlobalIndex(
const PetscInt &bIdx,
const MatStencil &s,
156 PetscInt &idx)
const;
169 const PetscInt &ptIdx,
171 PetscInt &idx)
const;
184 PetscInt &idx)
const;
218 PetscErrorCode
init(
const MPI_Comm &
comm,
const PetscInt &
dim,
219 const YAML::Node &node);
262typedef std::shared_ptr<body::BodyPackBase>
BodyPack;
281PetscErrorCode
createBodyPack(
const MPI_Comm &comm,
const PetscInt &dim,
Base (abstract) class for a pack of bodies.
PetscInt nLclPts
Total number of local Lagrangian points.
PetscMPIInt mpiSize
Total number of processes.
std::string info
String used to print information.
PetscErrorCode printInfo() const
Print information about the bodies.
BodyPackBase()=default
Default constructor.
type::IntVec1D nLclAllProcs
Number of local packed variables of all processes.
virtual PetscErrorCode destroy()
Manually destroy data.
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.
PetscErrorCode createInfoString()
Create a string used to print information.
MPI_Comm comm
Reference to the MPI communicator.
PetscErrorCode calculateAvgForces(const Vec &f, type::RealVec2D &fAvg)
Calculate the averaged force of each body.
virtual ~BodyPackBase()
Default destructor.
PetscErrorCode createDmPack()
Create DMComposite object for all bodies.
PetscMPIInt mpiRank
Rank of the local process.
std::vector< type::SingleBody > bodies
Vector of SingleBody objects.
type::IntVec1D offsetsAllProcs
Offsets of packed variables of all processes.
PetscErrorCode updateMeshIdx(const type::Mesh &mesh)
Get the index of closest Eulerian mesh cell for each local Lagrangian point.
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.
DM dmPack
DMComposite of DMDA objects of all SingleBody objects.
PetscErrorCode findProc(const PetscInt &bIdx, const PetscInt &ptIdx, PetscMPIInt &proc) const
Find which process owns the target Lagrangian point of target body.
PetscInt nPts
Total number of Lagrangian points.
PetscErrorCode init(const MPI_Comm &comm, const PetscInt &dim, const YAML::Node &node)
Initialize the pack of bodies.
PetscInt nBodies
Number of bodies in this pack.
PetscErrorCode createBodyPack(const MPI_Comm &comm, const PetscInt &dim, const YAML::Node &node, type::BodyPack &bodies)
Factory function to create a pack of bodies.
std::shared_ptr< body::BodyPackBase > BodyPack
Definition of type::BodyPack.
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
std::vector< RealVec1D > RealVec2D
2D std::vector holding PetscReal.
std::vector< PetscInt > IntVec1D
1D std::vector holding PetscInt.
A toolbox for building flow solvers.
body::SingleBodyBase, type::SingleBody factory function.