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

An implementation of body::SingleBodyBase that uses point data as input. More...

#include <singlebodypoints.h>

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

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

Detailed Description

An implementation of body::SingleBodyBase that uses point data as input.

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

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.

Constructor & Destructor Documentation

◆ SingleBodyPoints()

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

Constructor. Initialize a single body.

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

◆ ~SingleBodyPoints()

virtual petibm::body::SingleBodyPoints::~SingleBodyPoints ( )
virtualdefault

Default destructor.

Member Function Documentation

◆ calculateAvgForces()

PetscErrorCode petibm::body::SingleBodyPoints::calculateAvgForces ( const Vec &  f,
type::RealVec1D fAvg 
) const
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.

Implements petibm::body::SingleBodyBase.

Definition at line 207 of file singlebodypoints.cpp.

◆ createDMDA()

PetscErrorCode petibm::body::SingleBodyPoints::createDMDA ( )
protected

Create a parallel layout (1D DMDA object) of the body.

Returns
PetscErrorCode.

Definition at line 61 of file singlebodypoints.cpp.

◆ createInfoString()

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

Create a string used to print information about the body.

Returns
PetscErrorCode.

Definition at line 121 of file singlebodypoints.cpp.

◆ destroy()

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

Manually destroy data.

Definition at line 53 of file singlebody.cpp.

◆ findProc()

PetscErrorCode petibm::body::SingleBodyPoints::findProc ( const PetscInt &  i,
PetscMPIInt &  p 
) const
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.

Implements petibm::body::SingleBodyBase.

Definition at line 156 of file singlebodypoints.cpp.

◆ getGlobalIndex() [1/2]

PetscErrorCode petibm::body::SingleBodyPoints::getGlobalIndex ( const MatStencil &  s,
PetscInt &  idx 
) const
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.

Implements petibm::body::SingleBodyBase.

Definition at line 195 of file singlebodypoints.cpp.

◆ getGlobalIndex() [2/2]

PetscErrorCode petibm::body::SingleBodyPoints::getGlobalIndex ( const PetscInt &  i,
const PetscInt &  dof,
PetscInt &  idx 
) const
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.

Implements petibm::body::SingleBodyBase.

Definition at line 174 of file singlebodypoints.cpp.

◆ init()

PetscErrorCode petibm::body::SingleBodyPoints::init ( const MPI_Comm &  comm,
const PetscInt &  dim,
const std::string &  name,
const std::string &  filePath 
)
protected

Initialize the body, reading coordinates from given file.

Parameters
comm[in] MPI communicator.
dim[in] Number of dimensions.
name[in] Name of the body.
filePath[in] Path of the file with the coordinates.
Returns
PetscErrorCode.

Definition at line 27 of file singlebodypoints.cpp.

◆ printInfo()

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

Print information about the body.

Returns
PetscErrorCode.

Definition at line 72 of file singlebody.cpp.

◆ readBody()

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

Implements petibm::body::SingleBodyBase.

Definition at line 238 of file singlebodypoints.cpp.

◆ updateMeshIdx()

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

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

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

Implements petibm::body::SingleBodyBase.

Definition at line 95 of file singlebodypoints.cpp.

◆ writeBody()

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

Implements petibm::body::SingleBodyBase.

Definition at line 252 of file singlebodypoints.cpp.

Member Data Documentation

◆ bgPt

PetscInt petibm::body::SingleBodyBase::bgPt
inherited

Global index of the first local Lagrangian point.

Definition at line 65 of file singlebody.h.

◆ comm

MPI_Comm petibm::body::SingleBodyBase::comm
protectedinherited

MPI communicator.

Definition at line 180 of file singlebody.h.

◆ coords

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

Coordinates of ALL Lagrangian points.

Definition at line 50 of file singlebody.h.

◆ coords0

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

Initial coordinates of ALL Lagrangian points.

Definition at line 53 of file singlebody.h.

◆ da

DM petibm::body::SingleBodyBase::da
inherited

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

Definition at line 71 of file singlebody.h.

◆ dim

PetscInt petibm::body::SingleBodyBase::dim
inherited

Number of dimensions.

Definition at line 38 of file singlebody.h.

◆ edPt

PetscInt petibm::body::SingleBodyBase::edPt
inherited

Global index of the last local Lagrangian point.

Definition at line 68 of file singlebody.h.

◆ filePath

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

Path of the file with the body coordinates.

Definition at line 44 of file singlebody.h.

◆ info

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

String with information about the body.

Definition at line 74 of file singlebody.h.

◆ meshIdx

type::IntVec2D petibm::body::SingleBodyBase::meshIdx
inherited

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
protectedinherited

Rank of the local process.

Definition at line 186 of file singlebody.h.

◆ mpiSize

PetscMPIInt petibm::body::SingleBodyBase::mpiSize
protectedinherited

Total number of processes.

Definition at line 183 of file singlebody.h.

◆ name

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

Name of the body.

Definition at line 41 of file singlebody.h.

◆ nLclAllProcs

type::IntVec1D petibm::body::SingleBodyBase::nLclAllProcs
protectedinherited

Vector with the number of local unknowns on each process.

Definition at line 189 of file singlebody.h.

◆ nLclPts

PetscInt petibm::body::SingleBodyBase::nLclPts
inherited

Local number of Lagrangian points.

Definition at line 56 of file singlebody.h.

◆ nPts

PetscInt petibm::body::SingleBodyBase::nPts
inherited

Total number of Lagrangian points.

Definition at line 47 of file singlebody.h.

◆ offsetsAllProcs

type::IntVec1D petibm::body::SingleBodyBase::offsetsAllProcs
protectedinherited

Offset on each process.

Definition at line 192 of file singlebody.h.


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