| 
    PetIBM 0.5.4
    
   Toolbox and applications of the immersed-boundary method for distributed-memory architectures 
   | 
 
Base (abstract) class for a single body. More...
#include <singlebody.h>
Public Member Functions | |
| SingleBodyBase (const MPI_Comm &comm, const PetscInt &dim, const std::string &name, const std::string &filePath) | |
| Constructor. Initialize the body.  More... | |
| SingleBodyBase ()=default | |
| Default constructor.  More... | |
| virtual | ~SingleBodyBase () | 
| Default destructor.  More... | |
| virtual PetscErrorCode | destroy () | 
| Manually destroy data.  More... | |
| PetscErrorCode | printInfo () const | 
| Print information about the body.  More... | |
| virtual PetscErrorCode | findProc (const PetscInt &i, PetscMPIInt &p) const =0 | 
| Get the process id owning a given Lagrangian point.  More... | |
| 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.  More... | |
| 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.  More... | |
| virtual PetscErrorCode | calculateAvgForces (const Vec &f, type::RealVec1D &fAvg) const =0 | 
| Calculate the averaged force of this body.  More... | |
| virtual PetscErrorCode | updateMeshIdx (const type::Mesh &mesh)=0 | 
| Get the index of closest Eulerian mesh cell for each local Lagrangian point.  More... | |
| virtual PetscErrorCode | readBody (const std::string &filepath)=0 | 
| virtual PetscErrorCode | writeBody (const std::string &filepath)=0 | 
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 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... | |
Friends | |
| class | BodyPackBase | 
Base (abstract) class for a single body.
There is currently only one implementation for this abstract class: body::SingleBodyPoints.
Definition at line 32 of file singlebody.h.
| petibm::body::SingleBodyBase::SingleBodyBase | ( | const MPI_Comm & | comm, | 
| const PetscInt & | dim, | ||
| const std::string & | name, | ||
| const std::string & | filePath | ||
| ) | 
Constructor. Initialize the body.
| comm | [in] MPI communicator. | 
| dim | [in] Number of dimensions. | 
| name | [in] Name of body. | 
| filePath | [in] Path of the file with data of the body. | 
Data in filePath depend on the class used. For example, with the class body::SingleBodyPoints, the filePath should contain the number of Lagrangian points followed by the list of coordinates (placed in columns, the first column being the x coordinate). Note: body::SingleBodyPoints is currently the only sub-class implemented. 
Definition at line 16 of file singlebody.cpp.
      
  | 
  default | 
Default constructor.
      
  | 
  virtual | 
Default destructor.
Definition at line 39 of file singlebody.cpp.
      
  | 
  pure 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.
Implemented in petibm::body::SingleBodyPoints.
      
  | 
  virtual | 
Manually destroy data.
Definition at line 53 of file singlebody.cpp.
      
  | 
  pure 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.
Implemented in petibm::body::SingleBodyPoints.
      
  | 
  pure 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. | 
Implemented in petibm::body::SingleBodyPoints.
      
  | 
  pure 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. | 
Implemented in petibm::body::SingleBodyPoints.
| PetscErrorCode petibm::body::SingleBodyBase::printInfo | ( | ) | const | 
Print information about the body.
Definition at line 72 of file singlebody.cpp.
      
  | 
  pure virtual | 
Implemented in petibm::body::SingleBodyPoints.
      
  | 
  pure virtual | 
Get the index of closest Eulerian mesh cell for each local Lagrangian point.
| mesh | [in] Structured Cartesian mesh. | 
Implemented in petibm::body::SingleBodyPoints.
      
  | 
  pure virtual | 
Implemented in petibm::body::SingleBodyPoints.
      
  | 
  friend | 
Definition at line 34 of file singlebody.h.
| PetscInt petibm::body::SingleBodyBase::bgPt | 
Global index of the first local Lagrangian point.
Definition at line 65 of file singlebody.h.
      
  | 
  protected | 
MPI communicator.
Definition at line 180 of file singlebody.h.
| type::RealVec2D petibm::body::SingleBodyBase::coords | 
Coordinates of ALL Lagrangian points.
Definition at line 50 of file singlebody.h.
| type::RealVec2D petibm::body::SingleBodyBase::coords0 | 
Initial coordinates of ALL Lagrangian points.
Definition at line 53 of file singlebody.h.
| DM petibm::body::SingleBodyBase::da | 
Parallel layout of the body (as a 1D DMDA object).
Definition at line 71 of file singlebody.h.
| PetscInt petibm::body::SingleBodyBase::dim | 
Number of dimensions.
Definition at line 38 of file singlebody.h.
| PetscInt petibm::body::SingleBodyBase::edPt | 
Global index of the last local Lagrangian point.
Definition at line 68 of file singlebody.h.
| std::string petibm::body::SingleBodyBase::filePath | 
Path of the file with the body coordinates.
Definition at line 44 of file singlebody.h.
| std::string petibm::body::SingleBodyBase::info | 
String with information about the body.
Definition at line 74 of file singlebody.h.
| type::IntVec2D petibm::body::SingleBodyBase::meshIdx | 
Index of the closest Eulerian mesh cell for each local Lagrangian point.
Definition at line 62 of file singlebody.h.
      
  | 
  protected | 
Rank of the local process.
Definition at line 186 of file singlebody.h.
      
  | 
  protected | 
Total number of processes.
Definition at line 183 of file singlebody.h.
| std::string petibm::body::SingleBodyBase::name | 
Name of the body.
Definition at line 41 of file singlebody.h.
      
  | 
  protected | 
Vector with the number of local unknowns on each process.
Definition at line 189 of file singlebody.h.
| PetscInt petibm::body::SingleBodyBase::nLclPts | 
Local number of Lagrangian points.
Definition at line 56 of file singlebody.h.
| PetscInt petibm::body::SingleBodyBase::nPts | 
Total number of Lagrangian points.
Definition at line 47 of file singlebody.h.
      
  | 
  protected | 
Offset on each process.
Definition at line 192 of file singlebody.h.