PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
|
An implementation of body::SingleBodyBase that uses point data as input. More...
#include <singlebodypoints.h>
Public Member Functions | |
SingleBodyPoints (const MPI_Comm &comm, const PetscInt &dim, const std::string &name, const std::string &filePath) | |
Constructor. Initialize a single body. More... | |
virtual | ~SingleBodyPoints ()=default |
Default destructor. More... | |
virtual PetscErrorCode | findProc (const PetscInt &i, PetscMPIInt &p) const |
Get the process id owning a given Lagrangian point. More... | |
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. More... | |
virtual PetscErrorCode | getGlobalIndex (const MatStencil &s, PetscInt &idx) const |
Get the global index of a Lagrangian point in a DMDA object given a MatStencil. More... | |
virtual PetscErrorCode | calculateAvgForces (const Vec &f, type::RealVec1D &fAvg) const |
Calculate the averaged force of this body. More... | |
virtual PetscErrorCode | updateMeshIdx (const type::Mesh &mesh) |
Get the index of closest Eulerian mesh cell for each local Lagrangian point. More... | |
virtual PetscErrorCode | readBody (const std::string &filepath) |
virtual PetscErrorCode | writeBody (const std::string &filepath) |
virtual PetscErrorCode | destroy () |
Manually destroy data. More... | |
PetscErrorCode | printInfo () const |
Print information about the body. More... | |
Public Attributes | |
PetscInt | dim |
Number of dimensions. More... | |
std::string | name |
Name of the body. More... | |
std::string | filePath |
Path of the file with the body coordinates. More... | |
PetscInt | nPts |
Total number of Lagrangian points. More... | |
type::RealVec2D | coords |
Coordinates of ALL Lagrangian points. More... | |
type::RealVec2D | coords0 |
Initial coordinates of ALL Lagrangian points. More... | |
PetscInt | nLclPts |
Local number of Lagrangian points. More... | |
type::IntVec2D | meshIdx |
Index of the closest Eulerian mesh cell for each local Lagrangian point. More... | |
PetscInt | bgPt |
Global index of the first local Lagrangian point. More... | |
PetscInt | edPt |
Global index of the last local Lagrangian point. More... | |
DM | da |
Parallel layout of the body (as a 1D DMDA object). More... | |
std::string | info |
String with information about the body. More... | |
Protected Member Functions | |
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. More... | |
PetscErrorCode | createDMDA () |
Create a parallel layout (1D DMDA object) of the body. More... | |
PetscErrorCode | createInfoString () |
Create a string used to print information about the body. More... | |
Protected Attributes | |
MPI_Comm | comm |
MPI communicator. More... | |
PetscMPIInt | mpiSize |
Total number of processes. More... | |
PetscMPIInt | mpiRank |
Rank of the local process. More... | |
type::IntVec1D | nLclAllProcs |
Vector with the number of local unknowns on each process. More... | |
type::IntVec1D | offsetsAllProcs |
Offset on each process. More... | |
An implementation of body::SingleBodyBase that uses point data as input.
This implementation uses an ASCII file of coordinates of Lagrangian points as its input.
Users should not initialize an instance of this class directly. They should use petibm::body::createSingleBody. Actually, the design of PetIBM is to handle multiple bodies, users should use petibm::type::BodyPack to handle their bodies even if there is only one body.
Definition at line 31 of file singlebodypoints.h.
petibm::body::SingleBodyPoints::SingleBodyPoints | ( | const MPI_Comm & | comm, |
const PetscInt & | dim, | ||
const std::string & | name, | ||
const std::string & | filePath | ||
) |
Constructor. Initialize a single body.
comm | [in] MPI communicator. |
dim | [in] Number of dimensions. |
name | [in] Name of the body. |
filePath | [in] Path of the file with the coordinates. |
Definition at line 19 of file singlebodypoints.cpp.
|
virtualdefault |
Default destructor.
|
virtual |
Calculate the averaged force of this body.
f | [in] Vec of forces on each Lagrangian point of this body. |
fAvg | [out] return averaged force with length equal to dimension. |
Note: fAvg should have the correct size. This function does not check if fAvg has been allocated correctly.
Implements petibm::body::SingleBodyBase.
Definition at line 207 of file singlebodypoints.cpp.
|
protected |
Create a parallel layout (1D DMDA object) of the body.
Definition at line 61 of file singlebodypoints.cpp.
|
protected |
Create a string used to print information about the body.
Definition at line 121 of file singlebodypoints.cpp.
|
virtualinherited |
Manually destroy data.
Definition at line 53 of file singlebody.cpp.
|
virtual |
Get the process id owning a given Lagrangian point.
i | [in] Index of the Lagrangian point. |
p | [out] Process id. |
Note: there is no need to provide which degree of freedom is requested as all degrees for a Lagrangian point are on the same process.
Implements petibm::body::SingleBodyBase.
Definition at line 156 of file singlebodypoints.cpp.
|
virtual |
Get the global index of a Lagrangian point in a DMDA object given a MatStencil.
s | [in] MatStencil of the Lagrangian point. |
idx | [out] Global index of the Lagrangian point. |
Implements petibm::body::SingleBodyBase.
Definition at line 195 of file singlebodypoints.cpp.
|
virtual |
Get the global index of a Lagrangian point in a DMDA object given its degree of freedom.
i | [in] Index of the Lagrangian point. |
dof | [in] Degree of freedom. |
idx | [out] Global index of the Lagrangian point. |
Implements petibm::body::SingleBodyBase.
Definition at line 174 of file singlebodypoints.cpp.
|
protected |
Initialize the body, reading coordinates from given file.
comm | [in] MPI communicator. |
dim | [in] Number of dimensions. |
name | [in] Name of the body. |
filePath | [in] Path of the file with the coordinates. |
Definition at line 27 of file singlebodypoints.cpp.
|
inherited |
Print information about the body.
Definition at line 72 of file singlebody.cpp.
|
virtual |
Implements petibm::body::SingleBodyBase.
Definition at line 238 of file singlebodypoints.cpp.
|
virtual |
Get the index of closest Eulerian mesh cell for each local Lagrangian point.
mesh | [in] Structured Cartesian mesh. |
Implements petibm::body::SingleBodyBase.
Definition at line 95 of file singlebodypoints.cpp.
|
virtual |
Implements petibm::body::SingleBodyBase.
Definition at line 252 of file singlebodypoints.cpp.
|
inherited |
Global index of the first local Lagrangian point.
Definition at line 65 of file singlebody.h.
|
protectedinherited |
MPI communicator.
Definition at line 180 of file singlebody.h.
|
inherited |
Coordinates of ALL Lagrangian points.
Definition at line 50 of file singlebody.h.
|
inherited |
Initial coordinates of ALL Lagrangian points.
Definition at line 53 of file singlebody.h.
|
inherited |
Parallel layout of the body (as a 1D DMDA object).
Definition at line 71 of file singlebody.h.
|
inherited |
Number of dimensions.
Definition at line 38 of file singlebody.h.
|
inherited |
Global index of the last local Lagrangian point.
Definition at line 68 of file singlebody.h.
|
inherited |
Path of the file with the body coordinates.
Definition at line 44 of file singlebody.h.
|
inherited |
String with information about the body.
Definition at line 74 of file singlebody.h.
|
inherited |
Index of the closest Eulerian mesh cell for each local Lagrangian point.
Definition at line 62 of file singlebody.h.
|
protectedinherited |
Rank of the local process.
Definition at line 186 of file singlebody.h.
|
protectedinherited |
Total number of processes.
Definition at line 183 of file singlebody.h.
|
inherited |
Name of the body.
Definition at line 41 of file singlebody.h.
|
protectedinherited |
Vector with the number of local unknowns on each process.
Definition at line 189 of file singlebody.h.
|
inherited |
Local number of Lagrangian points.
Definition at line 56 of file singlebody.h.
|
inherited |
Total number of Lagrangian points.
Definition at line 47 of file singlebody.h.
|
protectedinherited |
Offset on each process.
Definition at line 192 of file singlebody.h.