PetIBM 0.5.4
Toolbox and applications of the immersed-boundary method for distributed-memory architectures
petibm::body::SingleBodyBase Class Referenceabstract

Base (abstract) class for a single body. More...

#include <singlebody.h>

Inheritance diagram for petibm::body::SingleBodyBase:
[legend]

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
 

Detailed Description

Base (abstract) class for a single body.

See also
Immersed-boundary bodies, petibm::type::SingleBody, petibm::body::createSingleBody

There is currently only one implementation for this abstract class: body::SingleBodyPoints.

Definition at line 32 of file singlebody.h.

Constructor & Destructor Documentation

◆ SingleBodyBase() [1/2]

petibm::body::SingleBodyBase::SingleBodyBase ( const MPI_Comm &  comm,
const PetscInt &  dim,
const std::string &  name,
const std::string &  filePath 
)

Constructor. Initialize the body.

Parameters
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.

◆ SingleBodyBase() [2/2]

petibm::body::SingleBodyBase::SingleBodyBase ( )
default

Default constructor.

◆ ~SingleBodyBase()

petibm::body::SingleBodyBase::~SingleBodyBase ( )
virtual

Default destructor.

Definition at line 39 of file singlebody.cpp.

Member Function Documentation

◆ calculateAvgForces()

virtual PetscErrorCode petibm::body::SingleBodyBase::calculateAvgForces ( const Vec &  f,
type::RealVec1D fAvg 
) const
pure virtual

Calculate the averaged force of this body.

Parameters
f[in] Vec of forces on each Lagrangian point of this body.
fAvg[out] return averaged force with length equal to dimension.
Returns
PetscErrorCode.

Note: fAvg should have the correct size. This function does not check if fAvg has been allocated correctly.

Implemented in petibm::body::SingleBodyPoints.

◆ destroy()

PetscErrorCode petibm::body::SingleBodyBase::destroy ( )
virtual

Manually destroy data.

Definition at line 53 of file singlebody.cpp.

◆ findProc()

virtual PetscErrorCode petibm::body::SingleBodyBase::findProc ( const PetscInt &  i,
PetscMPIInt &  p 
) const
pure virtual

Get the process id owning a given Lagrangian point.

Parameters
i[in] Index of the Lagrangian point.
p[out] Process id.
Returns
PetscErrorCode.

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.

◆ getGlobalIndex() [1/2]

virtual PetscErrorCode petibm::body::SingleBodyBase::getGlobalIndex ( const MatStencil &  s,
PetscInt &  idx 
) const
pure virtual

Get the global index of a Lagrangian point in a DMDA object given a MatStencil.

Parameters
s[in] MatStencil of the Lagrangian point.
idx[out] Global index of the Lagrangian point.
Returns
PetscErrorCode.

Implemented in petibm::body::SingleBodyPoints.

◆ getGlobalIndex() [2/2]

virtual PetscErrorCode petibm::body::SingleBodyBase::getGlobalIndex ( const PetscInt &  i,
const PetscInt &  dof,
PetscInt &  idx 
) const
pure virtual

Get the global index of a Lagrangian point in a DMDA object given its degree of freedom.

Parameters
i[in] Index of the Lagrangian point.
dof[in] Degree of freedom.
idx[out] Global index of the Lagrangian point.
Returns
PetscErrorCode.

Implemented in petibm::body::SingleBodyPoints.

◆ printInfo()

PetscErrorCode petibm::body::SingleBodyBase::printInfo ( ) const

Print information about the body.

Returns
PetscErrorCode.

Definition at line 72 of file singlebody.cpp.

◆ readBody()

virtual PetscErrorCode petibm::body::SingleBodyBase::readBody ( const std::string &  filepath)
pure virtual

◆ updateMeshIdx()

virtual PetscErrorCode petibm::body::SingleBodyBase::updateMeshIdx ( const type::Mesh mesh)
pure virtual

Get the index of closest Eulerian mesh cell for each local Lagrangian point.

Parameters
mesh[in] Structured Cartesian mesh.
Returns
PetscErrorCode.

Implemented in petibm::body::SingleBodyPoints.

◆ writeBody()

virtual PetscErrorCode petibm::body::SingleBodyBase::writeBody ( const std::string &  filepath)
pure virtual

Friends And Related Function Documentation

◆ BodyPackBase

friend class BodyPackBase
friend

Definition at line 34 of file singlebody.h.

Member Data Documentation

◆ bgPt

PetscInt petibm::body::SingleBodyBase::bgPt

Global index of the first local Lagrangian point.

Definition at line 65 of file singlebody.h.

◆ comm

MPI_Comm petibm::body::SingleBodyBase::comm
protected

MPI communicator.

Definition at line 180 of file singlebody.h.

◆ coords

type::RealVec2D petibm::body::SingleBodyBase::coords

Coordinates of ALL Lagrangian points.

Definition at line 50 of file singlebody.h.

◆ coords0

type::RealVec2D petibm::body::SingleBodyBase::coords0

Initial coordinates of ALL Lagrangian points.

Definition at line 53 of file singlebody.h.

◆ da

DM petibm::body::SingleBodyBase::da

Parallel layout of the body (as a 1D DMDA object).

Definition at line 71 of file singlebody.h.

◆ dim

PetscInt petibm::body::SingleBodyBase::dim

Number of dimensions.

Definition at line 38 of file singlebody.h.

◆ edPt

PetscInt petibm::body::SingleBodyBase::edPt

Global index of the last local Lagrangian point.

Definition at line 68 of file singlebody.h.

◆ filePath

std::string petibm::body::SingleBodyBase::filePath

Path of the file with the body coordinates.

Definition at line 44 of file singlebody.h.

◆ info

std::string petibm::body::SingleBodyBase::info

String with information about the body.

Definition at line 74 of file singlebody.h.

◆ meshIdx

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.

◆ mpiRank

PetscMPIInt petibm::body::SingleBodyBase::mpiRank
protected

Rank of the local process.

Definition at line 186 of file singlebody.h.

◆ mpiSize

PetscMPIInt petibm::body::SingleBodyBase::mpiSize
protected

Total number of processes.

Definition at line 183 of file singlebody.h.

◆ name

std::string petibm::body::SingleBodyBase::name

Name of the body.

Definition at line 41 of file singlebody.h.

◆ nLclAllProcs

type::IntVec1D petibm::body::SingleBodyBase::nLclAllProcs
protected

Vector with the number of local unknowns on each process.

Definition at line 189 of file singlebody.h.

◆ nLclPts

PetscInt petibm::body::SingleBodyBase::nLclPts

Local number of Lagrangian points.

Definition at line 56 of file singlebody.h.

◆ nPts

PetscInt petibm::body::SingleBodyBase::nPts

Total number of Lagrangian points.

Definition at line 47 of file singlebody.h.

◆ offsetsAllProcs

type::IntVec1D petibm::body::SingleBodyBase::offsetsAllProcs
protected

Offset on each process.

Definition at line 192 of file singlebody.h.


The documentation for this class was generated from the following files: