IBAMR  IBAMR version 0.19.
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
IBAMR::CIBStrategy Class Referenceabstract

Class CIBStrategy is a lightweight abstract strategy class which provides support for constraint based IB methods for rigid bodies.

#include <ibamr/CIBStrategy.h>

Inheritance diagram for IBAMR::CIBStrategy:
Inheritance graph
[legend]

Public Member Functions

 CIBStrategy (const unsigned int parts)
 Constructor of the class. More...
 
virtual ~CIBStrategy ()
 Destructor of the class. More...
 
virtual void setConstraintForce (Vec L, const double data_time, const double scale=1.0)=0
 Prepare the implementation class for sprading constraint force. In particular, set the constraint Lagrangian force in the internal data structure of the class. More...
 
virtual void getConstraintForce (Vec *L, const double data_time)=0
 Get the constraint rigid body force at the specified time within the current time interval. Generally, the implementation class maintains and stores the constraint force. This routine is called by constraint solver to update the contraint force after the (converged) solution is obtained. More...
 
virtual void getFreeRigidVelocities (Vec *U, const double data_time)
 Get the free rigid velocities (DOFs) at the specified time within the current time interval. This routine is called by constraint solver to update the free rigid DOFs after the (converged) solution is obtained. More...
 
virtual void getNetExternalForceTorque (Vec *F, const double data_time)
 Get net external force and torque at the specified time within the current time interval. This routine is called by constraint solver to form the appropriate RHS. More...
 
virtual void subtractMeanConstraintForce (Vec L, int f_data_idx, const double scale=1.0)=0
 Subtract the mean of constraint force from the background Eulerian grid. This is required for certain cases like periodic steady Stokes. More...
 
virtual void setInterpolatedVelocityVector (Vec V, const double data_time)
 Prepare the implementation class for getting the interpolated fluid velocity on the Lagrangian vector V. More...
 
virtual void getInterpolatedVelocity (Vec V, const double data_time, const double scale=1.0)=0
 Get the interpolated velocity from the Eulerian grid at the specified time. More...
 
virtual void computeMobilityRegularization (Vec D, Vec L, const double scale=1.0)=0
 Compute regularization vector for the mobility problem. More...
 
unsigned int getNumberOfRigidStructures () const
 Get number of rigid structures registered with this class. More...
 
virtual unsigned int getNumberOfNodes (const unsigned int part) const =0
 Get number of nodes associated with the particular structure. More...
 
void setInitialCenterOfMass (const unsigned int part, const Eigen::Vector3d &XCOM_0)
 Set the initial center of mass location for the structures. More...
 
void setSolveRigidBodyVelocity (const unsigned int part, const IBTK::FreeRigidDOFVector &solve_rigid_dofs)
 Set what rigid DOFs need to be solved for this particular structure. More...
 
const IBTK::FreeRigidDOFVectorgetSolveRigidBodyVelocity (const unsigned int part, int &num_free_dofs) const
 Query what rigid DOFs need to be solved for. More...
 
virtual void setRigidBodyVelocity (const unsigned int part, const IBTK::RigidDOFVector &U, Vec V)=0
 Set the rigid body velocity at the nodal/marker points contained in the Vec V. More...
 
virtual void setRigidBodyVelocity (const unsigned int part, Vec U, Vec V)
 Set the rigid body velocity at the nodal/marker points contained in the Vec V. More...
 
virtual void setRigidBodyVelocity (Vec U, Vec V, const bool only_free_dofs, const bool only_imposed_dofs, const bool all_dofs=false)
 Set the rigid body velocity at the nodal/marker points contained in the Vec V. More...
 
virtual void computeNetRigidGeneralizedForce (const unsigned int part, Vec L, IBTK::RigidDOFVector &F)=0
 Compute total force and torque on the structure. More...
 
virtual void computeNetRigidGeneralizedForce (const unsigned int part, Vec L, Vec F)
 Compute total force and torque on the structure. More...
 
virtual void computeNetRigidGeneralizedForce (Vec L, Vec F, const bool only_free_dofs, const bool only_imposed_dofs, const bool all_dofs=false)
 Compute total force and torque on the structure. More...
 
const IBTK::RigidDOFVectorgetNetRigidGeneralizedForce (const unsigned int part)
 Get total torque and force on the structure at new_time within the current time interval. More...
 
void updateFreeDOFsMapping ()
 Update the mapping of free DOFs for all structures if they are collected in a global vector. More...
 
void updateNewRigidBodyVelocity (const unsigned int part, const IBTK::RigidDOFVector &U)
 Update the rigid body velocity obtained from the constraint Stokes solver for free-moving case. More...
 
void updateNewRigidBodyVelocity (const unsigned int part, Vec U)
 Update the rigid body velocity obtained from the constraint Stokes solver for free-moving case. More...
 
void updateNewRigidBodyVelocity (Vec U, const bool only_free_dofs, const bool only_imposed_dofs, const bool all_dofs=false)
 Update the rigid body velocity obtained from the constraint Stokes solver for free-moving case. More...
 
virtual void copyVecToArray (Vec b, double *array, const std::vector< unsigned > &struct_ids, const int data_depth, const int array_rank)
 Copy data from distributed PETSc Vec for specified stucture indices to an array defined on a single processor. A default empty implementation is provided. More...
 
virtual void copyFreeDOFsVecToArray (Vec b, double *array, const std::vector< unsigned > &struct_ids, const int array_rank)
 Copy data from distributed PETSc Vec for specified stucture indices to an array defined on a single processor. A default implementation is provided. More...
 
virtual void copyArrayToVec (Vec b, double *array, const std::vector< unsigned > &struct_ids, const int data_depth, const int array_rank)
 Copy data from array defined on a single processor for specified stucture indices to distributed PETScVec. A default empty implementation is provided. More...
 
virtual void copyFreeDOFsArrayToVec (Vec b, double *array, const std::vector< unsigned > &struct_ids, const int array_rank)
 Copy data from array defined on a single processor for specified stucture indices to distributed PETScVec. A default implementation is provided. More...
 
void getCurrentRigidBodyVelocity (const unsigned int part, IBTK::RigidDOFVector &U)
 Get the rigid body translational velocity at the beginning of the timestep. More...
 
void getNewRigidBodyVelocity (const unsigned int part, IBTK::RigidDOFVector &U)
 Get the rigid body translational velocity at the end of the timestep. More...
 
const Eigen::Vector3d & getCurrentBodyCenterOfMass (const unsigned int part)
 Get body center of mass at the current time step. More...
 
const Eigen::Vector3d & getNewBodyCenterOfMass (const unsigned int part)
 Get body center of mass at half time step. More...
 
virtual void constructMobilityMatrix (const std::string &mat_name, MobilityMatrixType mat_type, Mat &mobility_mat, const std::vector< unsigned > &prototype_struct_ids, const double *grid_dx, const double *domain_extents, const bool initial_time, double rho, double mu, const std::pair< double, double > &scale, double f_periodic_corr, const int managing_rank)
 Construct dense mobility matrix for the prototypical structures identified by their indices. More...
 
virtual void constructGeometricMatrix (const std::string &mat_name, Mat &geometric_mat, const std::vector< unsigned > &prototype_struct_ids, const bool initial_time, const int managing_rank)
 Construct a geometric matrix for the prototypical structures identified by their indices. A geometric matrix maps center of mass rigid body velocity to nodal velocities. Geometric matrix is generally used with a dense mobility matrices to construct an associated body-mobility matrix algebrically. More...
 
virtual void rotateArray (double *array, const std::vector< unsigned > &struct_ids, const bool use_transpose, const int managing_rank, const int depth)
 Rotate vector using rotation matrix to/from the reference frame of the structures (which is at the initial time of the simulation). More...
 

Static Public Member Functions

static void vecToRDV (Vec U, IBTK::RigidDOFVector &Ur)
 Set the DOFs from PETSc Vec U to RigidDOFVector Ur. More...
 
static void rdvToVec (const IBTK::RigidDOFVector &Ur, Vec &U)
 Set the DOFs from RigidDOFVector Ur to PETSc Vec U. More...
 
static void eigenToRDV (const Eigen::Vector3d &U, const Eigen::Vector3d &W, IBTK::RigidDOFVector &UW)
 Set the DOFs from Eigen::Vector3d U and W to RigidDOFVector UW. More...
 
static void rdvToEigen (const IBTK::RigidDOFVector &UW, Eigen::Vector3d &U, Eigen::Vector3d &W)
 Set the DOFs from RigidDOFVector UW to Eigen::Vector3d U and W. More...
 

Protected Member Functions

void setRotationMatrix (const std::vector< Eigen::Vector3d > &rot_vel, const std::vector< Eigen::Quaterniond > &q_old, std::vector< Eigen::Quaterniond > &q_new, std::vector< Eigen::Matrix3d > &rot_mat, const double dt)
 Fill the rotation matrix. More...
 

Protected Attributes

unsigned int d_num_rigid_parts
 
std::vector< Eigen::Vector3d > d_center_of_mass_initial
 
std::vector< Eigen::Vector3d > d_center_of_mass_current
 
std::vector< Eigen::Vector3d > d_center_of_mass_half
 
std::vector< Eigen::Vector3d > d_center_of_mass_new
 
std::vector< boold_compute_center_of_mass_initial
 
std::vector< Eigen::Quaterniond > d_quaternion_current
 
std::vector< Eigen::Quaterniond > d_quaternion_half
 
std::vector< Eigen::Quaterniond > d_quaternion_new
 
std::vector< IBTK::FRDVd_solve_rigid_vel
 
Vec d_U
 
Vec d_F
 
std::vector< std::pair< int, int > > d_free_dofs_map
 
bool d_free_dofs_map_updated
 
std::vector< Eigen::Vector3d > d_trans_vel_current
 
std::vector< Eigen::Vector3d > d_trans_vel_half
 
std::vector< Eigen::Vector3d > d_trans_vel_new
 
std::vector< Eigen::Vector3d > d_rot_vel_current
 
std::vector< Eigen::Vector3d > d_rot_vel_half
 
std::vector< Eigen::Vector3d > d_rot_vel_new
 
std::vector< IBTK::RigidDOFVectord_net_rigid_generalized_force
 

Private Member Functions

 CIBStrategy (const CIBStrategy &from)=delete
 Copy constructor. More...
 
CIBStrategyoperator= (const CIBStrategy &that)=delete
 Assignment operator. More...
 

Constructor & Destructor Documentation

◆ CIBStrategy() [1/2]

IBAMR::CIBStrategy::CIBStrategy ( const unsigned int  parts)

◆ ~CIBStrategy()

virtual IBAMR::CIBStrategy::~CIBStrategy ( )
virtual

◆ CIBStrategy() [2/2]

IBAMR::CIBStrategy::CIBStrategy ( const CIBStrategy from)
privatedelete
Note
This constructor is not implemented and should not be used.
Parameters
fromThe value to copy to this object.

Member Function Documentation

◆ setConstraintForce()

virtual void IBAMR::CIBStrategy::setConstraintForce ( Vec  L,
const double  data_time,
const double  scale = 1.0 
)
pure virtual
Parameters
LVec containing the constraint force for all structures.
data_timeTime at which constraint force is to be spread.
scaleScales the constraint force before spreading.

Implemented in IBAMR::CIBMethod.

◆ getConstraintForce()

virtual void IBAMR::CIBStrategy::getConstraintForce ( Vec *  L,
const double  data_time 
)
pure virtual
Parameters
timeTime (current_time or new_time) at which constraint force is required.

Implemented in IBAMR::CIBMethod.

◆ getFreeRigidVelocities()

virtual void IBAMR::CIBStrategy::getFreeRigidVelocities ( Vec *  U,
const double  data_time 
)
virtual
Note
A default implementation is provided that returns the vector of free DOFs.
Parameters
timeTime (current_time or new_time) at which constraint force is required.

Reimplemented in IBAMR::CIBMethod.

◆ getNetExternalForceTorque()

virtual void IBAMR::CIBStrategy::getNetExternalForceTorque ( Vec *  F,
const double  data_time 
)
virtual
Note
A default implementation is provided that returns the vector of net external force and torque for the corresponding free DOFs.
Parameters
timeTime (current_time or new_time) at which external force and torque is required.

Reimplemented in IBAMR::CIBMethod.

◆ subtractMeanConstraintForce()

virtual void IBAMR::CIBStrategy::subtractMeanConstraintForce ( Vec  L,
int  f_data_idx,
const double  scale = 1.0 
)
pure virtual
Parameters
LVec containing the constraint force.
f_data_idxPatch data index of Eulerian body force.
scaleFactor by which L is scaled.

Implemented in IBAMR::CIBMethod.

◆ setInterpolatedVelocityVector()

virtual void IBAMR::CIBStrategy::setInterpolatedVelocityVector ( Vec  V,
const double  data_time 
)
virtual
Parameters
VVector that should contain the interpolated velocity.
data_timeTime at which Eulerian velocity is to be interpolated.
Note
A default implementation is provided that does nothing.

Reimplemented in IBAMR::CIBMethod.

◆ getInterpolatedVelocity()

virtual void IBAMR::CIBStrategy::getInterpolatedVelocity ( Vec  V,
const double  data_time,
const double  scale = 1.0 
)
pure virtual
Parameters
VVector that should contain the interpolated velocity.
data_timeTime at which Eulerian velocity is to be interpolated.
scaleScales the velocity vector after interpolating from the Eulerian grid.

Implemented in IBAMR::CIBMethod.

◆ computeMobilityRegularization()

virtual void IBAMR::CIBStrategy::computeMobilityRegularization ( Vec  D,
Vec  L,
const double  scale = 1.0 
)
pure virtual
Parameters
DVector containing the regularization for the mobility problem.
LVector from which regularization is to be computed.

Implemented in IBAMR::CIBMethod.

◆ getNumberOfRigidStructures()

unsigned int IBAMR::CIBStrategy::getNumberOfRigidStructures ( ) const

◆ getNumberOfNodes()

virtual unsigned int IBAMR::CIBStrategy::getNumberOfNodes ( const unsigned int  part) const
pure virtual

Implemented in IBAMR::CIBMethod.

◆ setInitialCenterOfMass()

void IBAMR::CIBStrategy::setInitialCenterOfMass ( const unsigned int  part,
const Eigen::Vector3d &  XCOM_0 
)
Parameters
partThe rigid body for which we are setting the initial center of mass position.

◆ setSolveRigidBodyVelocity()

void IBAMR::CIBStrategy::setSolveRigidBodyVelocity ( const unsigned int  part,
const IBTK::FreeRigidDOFVector solve_rigid_dofs 
)
Parameters
partThe rigid body for which we are setting the free DOFs.

◆ getSolveRigidBodyVelocity()

const IBTK::FreeRigidDOFVector& IBAMR::CIBStrategy::getSolveRigidBodyVelocity ( const unsigned int  part,
int num_free_dofs 
) const

◆ setRigidBodyVelocity() [1/3]

virtual void IBAMR::CIBStrategy::setRigidBodyVelocity ( const unsigned int  part,
const IBTK::RigidDOFVector U,
Vec  V 
)
pure virtual
Parameters
partThe rigid part for which velocity needs to be set.
URDV contains the rigid component of velocities. For two-dimensions the vector contains the values \([u,v,\omega_z]\) and for three-dimensions the vector values are \([u,v,w,\omega_x,\omega_y,\omega_z]\).

Implemented in IBAMR::CIBMethod.

◆ setRigidBodyVelocity() [2/3]

virtual void IBAMR::CIBStrategy::setRigidBodyVelocity ( const unsigned int  part,
Vec  U,
Vec  V 
)
virtual
Parameters
partThe rigid part for which velocity needs to be set.
UVec contains the rigid component of velocities. For two-dimensions the vector contains the values \([u,v,\omega_z]\) and for three-dimensions the vector values are \([u,v,w,\omega_x,\omega_y,\omega_z]\).

◆ setRigidBodyVelocity() [3/3]

virtual void IBAMR::CIBStrategy::setRigidBodyVelocity ( Vec  U,
Vec  V,
const bool  only_free_dofs,
const bool  only_imposed_dofs,
const bool  all_dofs = false 
)
virtual
Parameters
UVec that contains the rigid component of velocities for the required components. For two-dimensions each sub Vec contains the values \([u,v,\omega_z]\) and for three-dimensions the vector values are \([u,v,w,\omega_x,\omega_y,\omega_z]\).
only_free_dofsBoolean indicating if the rigid body velocity is to be set only for free DOFS for all parts. The corresponding size of U_sub would be \( U_{sub} \leq NDIM * (NDIM + 1) / 2 \).
only_imposed_dofsBoolean indicating if the rigid body velocity is to be set only for prescribed kinematics dofs for all parts. The corresponding size of U_sub would be \( U_{sub} \leq NDIM * (NDIM + 1) / 2 \).
all_dofsBoolean indicating if the rigid body velocity is to be set for all parts. The corresponding size of U_sub would be \( U_{sub} = NDIM * (NDIM + 1) / 2 \).
Note
User is responsible for setting correct number of subvecs in U that corresponds to the particular combination of booleans.

◆ computeNetRigidGeneralizedForce() [1/3]

virtual void IBAMR::CIBStrategy::computeNetRigidGeneralizedForce ( const unsigned int  part,
Vec  L,
IBTK::RigidDOFVector F 
)
pure virtual
Parameters
LThe Lagrange multiplier vector.
FVector RDV storing the net generalized force.
FRDV storing the net generalized force.

Implemented in IBAMR::CIBMethod.

◆ computeNetRigidGeneralizedForce() [2/3]

virtual void IBAMR::CIBStrategy::computeNetRigidGeneralizedForce ( const unsigned int  part,
Vec  L,
Vec  F 
)
virtual
Parameters
partThe structure index.
LThe Lagrange multiplier vector.
FVec storing the net generalized force.

◆ computeNetRigidGeneralizedForce() [3/3]

virtual void IBAMR::CIBStrategy::computeNetRigidGeneralizedForce ( Vec  L,
Vec  F,
const bool  only_free_dofs,
const bool  only_imposed_dofs,
const bool  all_dofs = false 
)
virtual
Parameters
LThe Lagrange multiplier vector.
FVec storing the net generalized force.
only_free_dofsBoolean indicating if the net generalized force and torque is to be computed only for free dofs of all bodies.
only_imposed_dofsBoolean indicating if the net generalized force and torque is to be computed for imposed dofs of all bodies.
all_dofsBoolean indicating if the net generalized force and torque is to be computed for all dofs of all bodies.
Note
User is responsible for setting correct number of subvecs in F that corresponds to the particular combination of booleans.

◆ getNetRigidGeneralizedForce()

const IBTK::RigidDOFVector& IBAMR::CIBStrategy::getNetRigidGeneralizedForce ( const unsigned int  part)
Parameters
partThe rigid part.

◆ updateFreeDOFsMapping()

void IBAMR::CIBStrategy::updateFreeDOFsMapping ( )

◆ updateNewRigidBodyVelocity() [1/3]

void IBAMR::CIBStrategy::updateNewRigidBodyVelocity ( const unsigned int  part,
const IBTK::RigidDOFVector U 
)

◆ updateNewRigidBodyVelocity() [2/3]

void IBAMR::CIBStrategy::updateNewRigidBodyVelocity ( const unsigned int  part,
Vec  U 
)

◆ updateNewRigidBodyVelocity() [3/3]

void IBAMR::CIBStrategy::updateNewRigidBodyVelocity ( Vec  U,
const bool  only_free_dofs,
const bool  only_imposed_dofs,
const bool  all_dofs = false 
)

◆ copyVecToArray()

virtual void IBAMR::CIBStrategy::copyVecToArray ( Vec  b,
double array,
const std::vector< unsigned > &  struct_ids,
const int  data_depth,
const int  array_rank 
)
virtual
Parameters
bPETSc Vec to copy from. The Vec stores data for nodal/marker points.
arrayData pointer to copy to.
struct_idsVector of structure indices.
data_depthDepth of the data stored at each Lagrangian node.
array_rankRank of the processor on which the array is located.
Note
The size of array is assummed to be sum of nodes of all the structures given in struct_ids times the data_depth.

Reimplemented in IBAMR::CIBMethod.

◆ copyFreeDOFsVecToArray()

virtual void IBAMR::CIBStrategy::copyFreeDOFsVecToArray ( Vec  b,
double array,
const std::vector< unsigned > &  struct_ids,
const int  array_rank 
)
virtual
Parameters
bPETSc Vec to copy from. The Vec stores only free DOFs of all the structures.
arrayData pointer to copy to. It is a linear array of maximum free DOFs of the passed structure IDs.
struct_idsVector of structure indices.
array_rankRank of the processor on which the array is located.
Note
The size of array is assummed to be sum of maximum number of free degrees of freedom of all the structures given in struct_ids. The caller is responsible for allocating and destroying array memory outside of this routine.

◆ copyArrayToVec()

virtual void IBAMR::CIBStrategy::copyArrayToVec ( Vec  b,
double array,
const std::vector< unsigned > &  struct_ids,
const int  data_depth,
const int  array_rank 
)
virtual
Parameters
bCopy to PETSc Vec.
arrayCopy from data pointer.
struct_idsVector of structure indices.
data_depthDepth of the data stored at each Lagrangian node.
array_rankRank of the processor on which the array is located.
Note
The size of array is assummed to be sum of nodes of all the structures given in struct_ids times the data_depth. The caller is responsible for allocating and destroying array memory outside of this routine.

Reimplemented in IBAMR::CIBMethod.

◆ copyFreeDOFsArrayToVec()

virtual void IBAMR::CIBStrategy::copyFreeDOFsArrayToVec ( Vec  b,
double array,
const std::vector< unsigned > &  struct_ids,
const int  array_rank 
)
virtual
Parameters
bCopy to PETSc Vec. The Vec stores only free DOFs of all the structures.
arrayCopy from data pointer. It is a linear array of maximum free DOFs of the passed structure IDs.
struct_idsVector of structure indices.
array_rankRank of the processor on which the array is located.
Note
The size of array is assummed to be sum of maximum number of free degrees of freedom of all the structures given in struct_ids. The caller is responsible for allocating and destroying array memory outside of this routine.

◆ vecToRDV()

static void IBAMR::CIBStrategy::vecToRDV ( Vec  U,
IBTK::RigidDOFVector Ur 
)
static

◆ rdvToVec()

static void IBAMR::CIBStrategy::rdvToVec ( const IBTK::RigidDOFVector Ur,
Vec &  U 
)
static

◆ eigenToRDV()

static void IBAMR::CIBStrategy::eigenToRDV ( const Eigen::Vector3d &  U,
const Eigen::Vector3d &  W,
IBTK::RigidDOFVector UW 
)
static

◆ rdvToEigen()

static void IBAMR::CIBStrategy::rdvToEigen ( const IBTK::RigidDOFVector UW,
Eigen::Vector3d &  U,
Eigen::Vector3d &  W 
)
static

◆ getCurrentRigidBodyVelocity()

void IBAMR::CIBStrategy::getCurrentRigidBodyVelocity ( const unsigned int  part,
IBTK::RigidDOFVector U 
)

◆ getNewRigidBodyVelocity()

void IBAMR::CIBStrategy::getNewRigidBodyVelocity ( const unsigned int  part,
IBTK::RigidDOFVector U 
)

◆ getCurrentBodyCenterOfMass()

const Eigen::Vector3d& IBAMR::CIBStrategy::getCurrentBodyCenterOfMass ( const unsigned int  part)

◆ getNewBodyCenterOfMass()

const Eigen::Vector3d& IBAMR::CIBStrategy::getNewBodyCenterOfMass ( const unsigned int  part)

◆ constructMobilityMatrix()

virtual void IBAMR::CIBStrategy::constructMobilityMatrix ( const std::string &  mat_name,
MobilityMatrixType  mat_type,
Mat &  mobility_mat,
const std::vector< unsigned > &  prototype_struct_ids,
const double grid_dx,
const double domain_extents,
const bool  initial_time,
double  rho,
double  mu,
const std::pair< double, double > &  scale,
double  f_periodic_corr,
const int  managing_rank 
)
virtual
Note
A default empty implementation is provided in this class. The derived class provides the actual implementation.
Parameters
mat_nameMatrix handle.
mat_typeMobility matrix type, e.g., RPY, EMPIRICAL, etc.
mobility_matDense sequential mobility matrix. The matrix is stored in column-major(FORTRAN) order.
Note
Must be allocated prior to entering this routine.
Parameters
prototype_struct_idsIndices of the structures as registered with
See also
IBAMR::IBStrategy class. A combined dense mobility matrix will formed for multiple structures.
Parameters
grid_dxNDIM vector of grid spacing of structure level.
domain_extentsNDIM vector of domain length.
initial_timeBoolean to indicate if the mobility matrix is to be generated for the initial position of the structures.
rhoFluid density
muFluid viscosity.
scaleScale for improving the conditioning number of dense mobility matrix. The matrix is scaled as: \( = \alpha * mobility_mat + \beta * identity_mat. \)
managing_rankRank of the processor managing this dense matrix.

Reimplemented in IBAMR::CIBMethod.

◆ constructGeometricMatrix()

virtual void IBAMR::CIBStrategy::constructGeometricMatrix ( const std::string &  mat_name,
Mat &  geometric_mat,
const std::vector< unsigned > &  prototype_struct_ids,
const bool  initial_time,
const int  managing_rank 
)
virtual
Note
A default empty implementation is provided in this class. The derived class provides the actual implementation.
Parameters
mat_nameMatrix handle.
geometric_matDense sequential geometric matrix. The matrix is stored in column-major(FORTRAN) order.
Note
Must be allocated prior to entering this routine.
Parameters
prototype_struct_idsIndices of the structures as registered with
See also
IBAMR::IBStrategy class. A combined block-diagonal geometric matrix will be formed for multiple structures.
Parameters
initial_timeBoolean to indicate if the corresponding geometric matrix is to be generated for the initial position of the structures.
managing_rankRank of the processor managing this dense matrix.

Reimplemented in IBAMR::CIBMethod.

◆ rotateArray()

virtual void IBAMR::CIBStrategy::rotateArray ( double array,
const std::vector< unsigned > &  struct_ids,
const bool  use_transpose,
const int  managing_rank,
const int  depth 
)
virtual
Parameters
arrayRaw data pointer containing the vector enteries.
struct_idsStructure ID indices.
use_transposeUse transpose of rotation matrix to rotate the vector.
Note
Transpose of rotation matrix is its inverse and it takes the vector back to its reference frame.
Parameters
managing_rankRank of the processor managing the matrix.
depthDepth of the data array components.

Reimplemented in IBAMR::CIBMethod.

◆ setRotationMatrix()

void IBAMR::CIBStrategy::setRotationMatrix ( const std::vector< Eigen::Vector3d > &  rot_vel,
const std::vector< Eigen::Quaterniond > &  q_old,
std::vector< Eigen::Quaterniond > &  q_new,
std::vector< Eigen::Matrix3d > &  rot_mat,
const double  dt 
)
protected
Parameters
q_oldPrevious applied quaternion.
q_newNew quaternion to set.
rot_matMatrix to set.
dtTime interval of rotation.

◆ operator=()

CIBStrategy& IBAMR::CIBStrategy::operator= ( const CIBStrategy that)
privatedelete
Note
This operator is not implemented and should not be used.
Parameters
thatThe value to assign to this object.
Returns
A reference to this object.

Member Data Documentation

◆ d_num_rigid_parts

unsigned int IBAMR::CIBStrategy::d_num_rigid_parts
protected

◆ d_center_of_mass_initial

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_center_of_mass_initial
protected

Center of mass.

◆ d_center_of_mass_current

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_center_of_mass_current
protected

◆ d_center_of_mass_half

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_center_of_mass_half
protected

◆ d_center_of_mass_new

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_center_of_mass_new
protected

◆ d_compute_center_of_mass_initial

std::vector<bool> IBAMR::CIBStrategy::d_compute_center_of_mass_initial
protected

◆ d_quaternion_current

std::vector<Eigen::Quaterniond> IBAMR::CIBStrategy::d_quaternion_current
protected

Quaternion of the body.

◆ d_quaternion_half

std::vector<Eigen::Quaterniond> IBAMR::CIBStrategy::d_quaternion_half
protected

◆ d_quaternion_new

std::vector<Eigen::Quaterniond> IBAMR::CIBStrategy::d_quaternion_new
protected

◆ d_solve_rigid_vel

std::vector<IBTK::FRDV> IBAMR::CIBStrategy::d_solve_rigid_vel
protected

◆ d_U

Vec IBAMR::CIBStrategy::d_U
protected

◆ d_F

Vec IBAMR::CIBStrategy::d_F
protected

◆ d_free_dofs_map

std::vector<std::pair<int, int> > IBAMR::CIBStrategy::d_free_dofs_map
protected

◆ d_free_dofs_map_updated

bool IBAMR::CIBStrategy::d_free_dofs_map_updated
protected

◆ d_trans_vel_current

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_trans_vel_current
protected

Rigid body velocity of the structures.

◆ d_trans_vel_half

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_trans_vel_half
protected

◆ d_trans_vel_new

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_trans_vel_new
protected

◆ d_rot_vel_current

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_rot_vel_current
protected

◆ d_rot_vel_half

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_rot_vel_half
protected

◆ d_rot_vel_new

std::vector<Eigen::Vector3d> IBAMR::CIBStrategy::d_rot_vel_new
protected

◆ d_net_rigid_generalized_force

std::vector<IBTK::RigidDOFVector> IBAMR::CIBStrategy::d_net_rigid_generalized_force
protected

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