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

Base (abstract) class for a pack of bodies. More...

#include <bodypack.h>

Public Member Functions

 BodyPackBase ()=default
 Default constructor. More...
 
 BodyPackBase (const MPI_Comm &comm, const PetscInt &dim, const YAML::Node &node)
 Constructor. Initialize the pack of bodies. More...
 
virtual ~BodyPackBase ()
 Default destructor. More...
 
virtual PetscErrorCode destroy ()
 Manually destroy data. More...
 
PetscErrorCode printInfo () const
 Print information about the bodies. More...
 
PetscErrorCode findProc (const PetscInt &bIdx, const PetscInt &ptIdx, PetscMPIInt &proc) const
 Find which process owns the target Lagrangian point of target body. More...
 
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. More...
 
PetscErrorCode getGlobalIndex (const PetscInt &bIdx, const MatStencil &s, PetscInt &idx) const
 Find unpacked global index of a DoF of Lagrangian point of a body. More...
 
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. More...
 
PetscErrorCode getPackedGlobalIndex (const PetscInt &bIdx, const MatStencil &s, PetscInt &idx) const
 Find packed global index of a DoF of Lagrangian point of a body. More...
 
PetscErrorCode calculateAvgForces (const Vec &f, type::RealVec2D &fAvg)
 Calculate the averaged force of each body. More...
 
PetscErrorCode updateMeshIdx (const type::Mesh &mesh)
 Get the index of closest Eulerian mesh cell for each local Lagrangian point. More...
 

Public Attributes

PetscInt dim
 Dimension. More...
 
PetscInt nBodies
 Number of bodies in this pack. More...
 
PetscInt nPts
 Total number of Lagrangian points. More...
 
PetscInt nLclPts
 Total number of local Lagrangian points. More...
 
std::vector< type::SingleBodybodies
 Vector of SingleBody objects. More...
 
DM dmPack
 DMComposite of DMDA objects of all SingleBody objects. More...
 
std::string info
 String used to print information. More...
 

Protected Member Functions

PetscErrorCode init (const MPI_Comm &comm, const PetscInt &dim, const YAML::Node &node)
 Initialize the pack of bodies. More...
 
PetscErrorCode createDmPack ()
 Create DMComposite object for all bodies. More...
 
PetscErrorCode createInfoString ()
 Create a string used to print information. More...
 

Protected Attributes

MPI_Comm comm
 Reference to the MPI communicator. More...
 
PetscMPIInt mpiSize
 Total number of processes. More...
 
PetscMPIInt mpiRank
 Rank of the local process. More...
 
type::IntVec1D nLclAllProcs
 Number of local packed variables of all processes. More...
 
type::IntVec1D offsetsAllProcs
 Offsets of packed variables of all processes. More...
 

Detailed Description

Base (abstract) class for a pack of bodies.

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

This class is designed to be an abstract class though, it actually has full implementations of all members because there is currently no necessary to do otherwise. It's just simply a collection of bodies, though it also deals with indexing of all Lagrangian points in a globally packed Vec.

Definition at line 70 of file bodypack.h.

Constructor & Destructor Documentation

◆ BodyPackBase() [1/2]

petibm::body::BodyPackBase::BodyPackBase ( )
default

Default constructor.

◆ BodyPackBase() [2/2]

petibm::body::BodyPackBase::BodyPackBase ( const MPI_Comm &  comm,
const PetscInt &  dim,
const YAML::Node &  node 
)

Constructor. Initialize the pack of bodies.

Parameters
comm[in] MPI communicator.
dim[in] Number of dimensions.
node[in] YAML configuration node.

Definition at line 15 of file bodypack.cpp.

◆ ~BodyPackBase()

petibm::body::BodyPackBase::~BodyPackBase ( )
virtual

Default destructor.

Definition at line 21 of file bodypack.cpp.

Member Function Documentation

◆ calculateAvgForces()

PetscErrorCode petibm::body::BodyPackBase::calculateAvgForces ( const Vec &  f,
type::RealVec2D fAvg 
)

Calculate the averaged force of each body.

Parameters
f[in] Packed force Vec of Lagrangian points.
fAvg[out] Averaged force for each body.
Returns
PetscErrorCode.

Note: if fAvg doesn't have the correct size, the function will resize it.

Definition at line 298 of file bodypack.cpp.

◆ createDmPack()

PetscErrorCode petibm::body::BodyPackBase::createDmPack ( )
protected

Create DMComposite object for all bodies.

Returns
PetscErrorCode.

Definition at line 144 of file bodypack.cpp.

◆ createInfoString()

PetscErrorCode petibm::body::BodyPackBase::createInfoString ( )
protected

Create a string used to print information.

Returns
PetscErrorCode.

Definition at line 174 of file bodypack.cpp.

◆ destroy()

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

Manually destroy data.

Definition at line 35 of file bodypack.cpp.

◆ findProc()

PetscErrorCode petibm::body::BodyPackBase::findProc ( const PetscInt &  bIdx,
const PetscInt &  ptIdx,
PetscMPIInt &  proc 
) const

Find which process owns the target Lagrangian point of target body.

Parameters
bIdx[in] Index of target body.
ptIdx[in] Index of target point.
proc[out] Process index.
Returns
PetscErrorCode.

Definition at line 211 of file bodypack.cpp.

◆ getGlobalIndex() [1/2]

PetscErrorCode petibm::body::BodyPackBase::getGlobalIndex ( const PetscInt &  bIdx,
const MatStencil &  s,
PetscInt &  idx 
) const

Find unpacked global index of a DoF of Lagrangian point of a body.

Parameters
bIdx[in] Index of target body.
s[in] PETSc MatStencil object of target point.
idx[out] Unpacked global index.
Returns
PetscErrorCode.

Definition at line 248 of file bodypack.cpp.

◆ getGlobalIndex() [2/2]

PetscErrorCode petibm::body::BodyPackBase::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.

Parameters
bIdx[in] Index of target body.
ptIdx[in] Index of target point.
dof[in] Index of target DoF.
idx[out] Unpacked global index.
Returns
PetscErrorCode.

Definition at line 229 of file bodypack.cpp.

◆ getPackedGlobalIndex() [1/2]

PetscErrorCode petibm::body::BodyPackBase::getPackedGlobalIndex ( const PetscInt &  bIdx,
const MatStencil &  s,
PetscInt &  idx 
) const

Find packed global index of a DoF of Lagrangian point of a body.

Parameters
bIdx[in] Index of target body.
s[in] PETSc MatStencil object of target point.
idx[in] Packed global index.
Returns
PetscErrorCode.

Definition at line 285 of file bodypack.cpp.

◆ getPackedGlobalIndex() [2/2]

PetscErrorCode petibm::body::BodyPackBase::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.

Parameters
bIdx[in] Index of target body.
ptIdx[in] Index of target point.
dof[in] Index of target DoF.
idx[out] Packed global index.
Returns
PetscErrorCode.

Definition at line 261 of file bodypack.cpp.

◆ init()

PetscErrorCode petibm::body::BodyPackBase::init ( const MPI_Comm &  comm,
const PetscInt &  dim,
const YAML::Node &  node 
)
protected

Initialize the pack of bodies.

Parameters
comm[in] MPI communicator.
dim[in] Number of dimensions.
node[in] YAML configuration node.
Returns
PetscErrorCode.

Definition at line 55 of file bodypack.cpp.

◆ printInfo()

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

Print information about the bodies.

Returns
PetscErrorCode.

Definition at line 195 of file bodypack.cpp.

◆ updateMeshIdx()

PetscErrorCode petibm::body::BodyPackBase::updateMeshIdx ( const type::Mesh mesh)

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

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

Definition at line 160 of file bodypack.cpp.

Member Data Documentation

◆ bodies

std::vector<type::SingleBody> petibm::body::BodyPackBase::bodies

Vector of SingleBody objects.

Definition at line 86 of file bodypack.h.

◆ comm

MPI_Comm petibm::body::BodyPackBase::comm
protected

Reference to the MPI communicator.

Definition at line 222 of file bodypack.h.

◆ dim

PetscInt petibm::body::BodyPackBase::dim

Dimension.

Definition at line 74 of file bodypack.h.

◆ dmPack

DM petibm::body::BodyPackBase::dmPack

DMComposite of DMDA objects of all SingleBody objects.

Definition at line 89 of file bodypack.h.

◆ info

std::string petibm::body::BodyPackBase::info

String used to print information.

Definition at line 92 of file bodypack.h.

◆ mpiRank

PetscMPIInt petibm::body::BodyPackBase::mpiRank
protected

Rank of the local process.

Definition at line 228 of file bodypack.h.

◆ mpiSize

PetscMPIInt petibm::body::BodyPackBase::mpiSize
protected

Total number of processes.

Definition at line 225 of file bodypack.h.

◆ nBodies

PetscInt petibm::body::BodyPackBase::nBodies

Number of bodies in this pack.

Definition at line 77 of file bodypack.h.

◆ nLclAllProcs

type::IntVec1D petibm::body::BodyPackBase::nLclAllProcs
protected

Number of local packed variables of all processes.

Definition at line 231 of file bodypack.h.

◆ nLclPts

PetscInt petibm::body::BodyPackBase::nLclPts

Total number of local Lagrangian points.

Definition at line 83 of file bodypack.h.

◆ nPts

PetscInt petibm::body::BodyPackBase::nPts

Total number of Lagrangian points.

Definition at line 80 of file bodypack.h.

◆ offsetsAllProcs

type::IntVec1D petibm::body::BodyPackBase::offsetsAllProcs
protected

Offsets of packed variables of all processes.

Definition at line 234 of file bodypack.h.


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