PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
singlebodypoints.h
Go to the documentation of this file.
1
8#pragma once
9
10#include <petibm/singlebody.h>
11
12namespace petibm
13{
14namespace body
15{
32{
33public:
42 SingleBodyPoints(const MPI_Comm &comm, const PetscInt &dim,
43 const std::string &name, const std::string &filePath);
44
46 virtual ~SingleBodyPoints() = default;
47
48 // implementation of SingleBodyBase::findProc
49 virtual PetscErrorCode findProc(const PetscInt &i, PetscMPIInt &p) const;
50
51 // implementation of SingleBodyBase::getGlobalIndex
52 virtual PetscErrorCode getGlobalIndex(const PetscInt &i,
53 const PetscInt &dof,
54 PetscInt &idx) const;
55
56 // implementation of SingleBodyBase::getGlobalIndex
57 virtual PetscErrorCode getGlobalIndex(const MatStencil &s,
58 PetscInt &idx) const;
59
60 // implementation of SingleBodyBase::calculateAvgForces
61 virtual PetscErrorCode calculateAvgForces(const Vec &f,
62 type::RealVec1D &fAvg) const;
63
64 // implementation of SingleBodyBase::updateMeshIdx
65 virtual PetscErrorCode updateMeshIdx(const type::Mesh &mesh);
66
67 virtual PetscErrorCode readBody(const std::string &filepath);
68
69 virtual PetscErrorCode writeBody(const std::string &filepath);
70
71protected:
82 PetscErrorCode init(const MPI_Comm &comm, const PetscInt &dim,
83 const std::string &name, const std::string &filePath);
84
90 PetscErrorCode createDMDA();
91
97 PetscErrorCode createInfoString();
98
99}; // SingleBodyPoints
100
101} // end of namespace body
102
103} // end of namespace petibm
Base (abstract) class for a single body.
Definition: singlebody.h:33
std::string name
Name of the body.
Definition: singlebody.h:41
MPI_Comm comm
MPI communicator.
Definition: singlebody.h:180
PetscInt dim
Number of dimensions.
Definition: singlebody.h:38
std::string filePath
Path of the file with the body coordinates.
Definition: singlebody.h:44
An implementation of body::SingleBodyBase that uses point data as input.
virtual PetscErrorCode writeBody(const std::string &filepath)
virtual PetscErrorCode updateMeshIdx(const type::Mesh &mesh)
Get the index of closest Eulerian mesh cell for each local Lagrangian point.
PetscErrorCode createInfoString()
Create a string used to print information about the body.
virtual ~SingleBodyPoints()=default
Default destructor.
PetscErrorCode init(const MPI_Comm &comm, const PetscInt &dim, const std::string &name, const std::string &filePath)
Initialize the body, reading coordinates from given file.
virtual PetscErrorCode getGlobalIndex(const PetscInt &i, const PetscInt &dof, PetscInt &idx) const
Get the global index of a Lagrangian point in a DMDA object given its degree of freedom.
SingleBodyPoints(const MPI_Comm &comm, const PetscInt &dim, const std::string &name, const std::string &filePath)
Constructor. Initialize a single body.
PetscErrorCode createDMDA()
Create a parallel layout (1D DMDA object) of the body.
virtual PetscErrorCode findProc(const PetscInt &i, PetscMPIInt &p) const
Get the process id owning a given Lagrangian point.
virtual PetscErrorCode readBody(const std::string &filepath)
virtual PetscErrorCode calculateAvgForces(const Vec &f, type::RealVec1D &fAvg) const
Calculate the averaged force of this body.
std::shared_ptr< mesh::MeshBase > Mesh
Type definition of Mesh.
Definition: mesh.h:348
std::vector< PetscReal > RealVec1D
1D std::vector holding PetscReal.
Definition: type.h:140
A toolbox for building flow solvers.
Definition: bodypack.h:52
body::SingleBodyBase, type::SingleBody factory function.