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

Class IBInterpolantMethod is an implementation of the abstract base class IBStrategy that provides a simple functionality of interpolating and spreading scalar or vector-valued quantities to and from the moving IB structure on the background Cartesian grid. More...

#include <ibamr/IBInterpolantMethod.h>

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

Public Member Functions

 IBInterpolantMethod (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, int no_structures=1, bool register_for_restart=true)
 Constructor. More...
 
 ~IBInterpolantMethod ()
 Destructor. More...
 
void registerEulerianVariables () override
 
void registerEulerianCommunicationAlgorithms () override
 
virtual void registerVariableAndHierarchyIntegrator (const std::string &var_name, const int var_depth, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > variable, SAMRAI::tbox::Pointer< IBTK::HierarchyIntegrator > hier_integrator)
 Register a variable and the HierarchyIntegrator managing the variable with this class. More...
 
void registerLInitStrategy (SAMRAI::tbox::Pointer< IBTK::LInitStrategy > l_initializer)
 
void freeLInitStrategy ()
 
IBTK::LDataManagergetLDataManager () const
 
int getStructuresLevelNumber () const
 Get the level on which the structures reside. More...
 
int getStructureHandle (const int lag_idx) const
 Get the structure handle to which this Lagrangian index belongs. More...
 
unsigned int getNumberOfNodes (const unsigned int struct_no) const
 Get the number of nodes for this structure. More...
 
void registerLSiloDataWriter (SAMRAI::tbox::Pointer< IBTK::LSiloDataWriter > silo_writer)
 
const SAMRAI::hier::IntVector< NDIM > & getMinimumGhostCellWidth () const override
 
void setupTagBuffer (SAMRAI::tbox::Array< int > &tag_buffer, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) const override
 
void preprocessIntegrateData (double current_time, double new_time, int num_cycles) override
 
void postprocessIntegrateData (double current_time, double new_time, int num_cycles) override
 
void interpolateVelocity (int u_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &u_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &u_ghost_fill_scheds, double data_time) override
 
void interpolateQ (double data_time)
 
void interpolateQ ()
 
void forwardEulerStep (double current_time, double new_time) override
 
void backwardEulerStep (double current_time, double new_time) override
 
void midpointStep (double current_time, double new_time) override
 
void trapezoidalStep (double current_time, double new_time) override
 
void updateMeshPosition (double current_time, double new_time, const std::vector< Eigen::Vector3d > &U, const std::vector< Eigen::Vector3d > &W)
 
void computeLagrangianForce (double data_time) override
 
void spreadForce (int f_data_idx, IBTK::RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds, double data_time) override
 
void spreadQ (double data_time)
 
int getEulerianQPatchDataIndex (const std::string &var_name)
 
void copyEulerianDataToIntegrator (double data_time)
 
void initializePatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg, int u_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &u_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &u_ghost_fill_scheds, int integrator_step, double init_data_time, bool initial_time) override
 
void addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx) override
 
void beginDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
void endDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
void initializeLevelData (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double init_data_time, bool can_be_refined, bool initial_time, SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > > old_level, bool allocate_data) override
 
void resetHierarchyConfiguration (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_level, int finest_level) override
 
void applyGradientDetector (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool uses_richardson_extrapolation_too) override
 
void putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 
virtual void registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver)
 
virtual void inactivateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max())
 
virtual void activateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max())
 
virtual bool getLagrangianStructureIsActivated (int structure_number=0, int level_number=std::numeric_limits< int >::max()) const
 
virtual double getMaxPointDisplacement () const
 
void setUseFixedLEOperators (bool use_fixed_coupling_ops=true)
 
virtual void updateFixedLEOperators ()
 
virtual void setUseMultistepTimeStepping (unsigned int n_previous_steps=1)
 
virtual void AB2Step (double current_time, double new_time)
 
virtual bool hasFluidSources () const
 
virtual void computeLagrangianFluidSource (double data_time)
 
virtual void spreadFluidSource (int q_data_idx, IBTK::RobinPhysBdryPatchStrategy *q_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &q_prolongation_scheds, double data_time)
 
virtual void interpolatePressure (int p_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &p_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &p_ghost_fill_scheds, double data_time)
 
virtual void preprocessSolveFluidEquations (double current_time, double new_time, int cycle_num)
 
virtual void postprocessSolveFluidEquations (double current_time, double new_time, int cycle_num)
 
virtual void postprocessData ()
 
virtual void registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx)
 
virtual void initializeLevelData (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double init_data_time, const bool can_be_refined, const bool initial_time, const tbox::Pointer< hier::BasePatchLevel< DIM > > old_level=tbox::Pointer< hier::BasePatchLevel< DIM > >(NULL), const bool allocate_data=true)=0
 
virtual void resetHierarchyConfiguration (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level)=0
 
virtual void applyGradientDetector (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double error_data_time, const int tag_index, const bool initial_time, const bool uses_richardson_extrapolation_too)
 
virtual double getLevelDt (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double dt_time, const bool initial_time)
 
virtual double advanceLevel (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const double current_time, const double new_time, const bool first_step, const bool last_step, const bool regrid_advance=false)
 
virtual void resetTimeDependentData (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double new_time, const bool can_be_refined)
 
virtual void resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< DIM > > level)
 
virtual void applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< DIM > > level, const double error_data_time, const int tag_index, const double deltat, const int error_coarsen_ratio, const bool initial_time, const bool uses_gradient_detector_too)
 
virtual void coarsenDataForRichardsonExtrapolation (const tbox::Pointer< hier::PatchHierarchy< DIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< DIM > > coarser_level, const double coarsen_data_time, const bool before_advance)
 

Protected Member Functions

void getPositionData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **X_data, double data_time)
 
void copyEulerianDataFromIntegrator (const std::string &var_name, int q_data_idx, double data_time)
 
void zeroOutEulerianData (const std::string &var_name, int q_data_idx)
 
void copyEulerianDataToIntegrator (const std::string &name, int q_data_idx, double data_time)
 
void getQData (const std::string &var_name, std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **Q_data, double data_time)
 
void computeCenterOfMass (std::vector< Eigen::Vector3d > &center_of_mass, std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &X_data)
 Compute the center of mass of the mesh. More...
 
INSHierarchyIntegratorgetINSHierarchyIntegrator () const
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > getVelocityHierarchyDataOps () const
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > getPressureHierarchyDataOps () const
 
SAMRAI::tbox::Pointer< IBTK::HierarchyMathOpsgetHierarchyMathOps () const
 
void registerVariable (int &current_idx, int &new_idx, int &scratch_idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > variable, const SAMRAI::hier::IntVector< NDIM > &scratch_ghosts=SAMRAI::hier::IntVector< NDIM >(0), const std::string &coarsen_name="NO_COARSEN", const std::string &refine_name="NO_REFINE", SAMRAI::tbox::Pointer< IBTK::CartGridFunction > init_fcn=nullptr, const bool register_for_restart=true)
 
void registerVariable (int &idx, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > variable, const SAMRAI::hier::IntVector< NDIM > &ghosts=SAMRAI::hier::IntVector< NDIM >(0), SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > ctx=SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext >(nullptr), const bool register_for_restart=true)
 
void registerGhostfillRefineAlgorithm (const std::string &name, SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > ghostfill_alg, std::unique_ptr< SAMRAI::xfer::RefinePatchStrategy< NDIM > > ghostfill_patch_strategy=nullptr)
 
void registerProlongRefineAlgorithm (const std::string &name, SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > prolong_alg, std::unique_ptr< SAMRAI::xfer::RefinePatchStrategy< NDIM > > prolong_patch_strategy=nullptr)
 
void registerCoarsenAlgorithm (const std::string &name, SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > > coarsen_alg, std::unique_ptr< SAMRAI::xfer::CoarsenPatchStrategy< NDIM > > coarsen_patch_strategy=nullptr)
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > getGhostfillRefineAlgorithm (const std::string &name) const
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > getProlongRefineAlgorithm (const std::string &name) const
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > > getCoarsenAlgorithm (const std::string &name) const
 
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & getGhostfillRefineSchedules (const std::string &name) const
 
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > & getProlongRefineSchedules (const std::string &name) const
 
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > & getCoarsenSchedules (const std::string &name) const
 

Protected Attributes

bool d_do_log = false
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > d_hierarchy
 
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > d_gridding_alg
 
double d_current_time = std::numeric_limits<double>::quiet_NaN()
 
double d_new_time = std::numeric_limits<double>::quiet_NaN()
 
double d_half_time = std::numeric_limits<double>::quiet_NaN()
 
std::vector< std::pair< int, int > > d_struct_lag_idx_range
 
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_new
 
std::vector< Eigen::Quaterniond > d_quaternion_current
 
std::vector< Eigen::Quaterniond > d_quaternion_new
 
IBTK::LDataManagerd_l_data_manager
 
std::string d_interp_kernel_fcn = "IB_4"
 
std::string d_spread_kernel_fcn = "IB_4"
 
bool d_error_if_points_leave_domain = false
 
SAMRAI::hier::IntVector< NDIM > d_ghosts
 
std::map< std::string, SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > > d_q_var
 Eulerian variables and their HierarchyIntegrators. More...
 
std::map< std::string, intd_q_depth
 
std::map< std::string, intd_q_interp_idx
 
std::map< std::string, SAMRAI::tbox::Pointer< IBTK::HierarchyIntegrator > > d_q_hier_integrator
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_current_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_new_data
 
std::map< std::string, std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > > d_Q_current_data
 
std::map< std::string, std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > > d_Q_new_data
 
SAMRAI::tbox::Pointer< IBTK::LInitStrategyd_l_initializer
 
SAMRAI::tbox::Pointer< IBTK::LSiloDataWriterd_silo_writer
 
std::string d_object_name
 
bool d_registered_for_restart
 
IBHierarchyIntegratord_ib_solver = nullptr
 
bool d_use_fixed_coupling_ops = false
 

Private Member Functions

 IBInterpolantMethod ()=delete
 Default constructor. More...
 
 IBInterpolantMethod (const IBInterpolantMethod &from)=delete
 Copy constructor. More...
 
IBInterpolantMethodoperator= (const IBInterpolantMethod &that)=delete
 Assignment operator. More...
 
void getFromInput (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db, bool is_from_restart)
 
void getFromRestart ()
 

Detailed Description

This class does not provide direct support for fluid-structure interaction, but rather facilitates Eulerian-Lagrangian data management for LMesh.

Constructor & Destructor Documentation

◆ IBInterpolantMethod() [1/3]

IBAMR::IBInterpolantMethod::IBInterpolantMethod ( std::string  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
int  no_structures = 1,
bool  register_for_restart = true 
)

◆ ~IBInterpolantMethod()

IBAMR::IBInterpolantMethod::~IBInterpolantMethod ( )

◆ IBInterpolantMethod() [2/3]

IBAMR::IBInterpolantMethod::IBInterpolantMethod ( )
privatedelete
Note
This constructor is not implemented and should not be used.

◆ IBInterpolantMethod() [3/3]

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

Member Function Documentation

◆ registerEulerianVariables()

void IBAMR::IBInterpolantMethod::registerEulerianVariables ( )
overridevirtual

Register Eulerian variables with the parent IBHierarchyIntegrator with the VariableDatabase. They are used for interpolating the Eulerian quantities on the Lagrangian mesh.

Reimplemented from IBAMR::IBStrategy.

◆ registerEulerianCommunicationAlgorithms()

void IBAMR::IBInterpolantMethod::registerEulerianCommunicationAlgorithms ( )
overridevirtual

Register ghost cell filling algrorithms for the Eulerian variables.

Reimplemented from IBAMR::IBStrategy.

◆ registerVariableAndHierarchyIntegrator()

virtual void IBAMR::IBInterpolantMethod::registerVariableAndHierarchyIntegrator ( const std::string &  var_name,
const int  var_depth,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  variable,
SAMRAI::tbox::Pointer< IBTK::HierarchyIntegrator hier_integrator 
)
virtual

◆ registerLInitStrategy()

void IBAMR::IBInterpolantMethod::registerLInitStrategy ( SAMRAI::tbox::Pointer< IBTK::LInitStrategy l_initializer)

Supply a Lagrangian initialization object.

◆ freeLInitStrategy()

void IBAMR::IBInterpolantMethod::freeLInitStrategy ( )

Free references to Lagrangian initialization objects.

◆ getLDataManager()

IBTK::LDataManager* IBAMR::IBInterpolantMethod::getLDataManager ( ) const

Return a pointer to the Lagrangian data manager object.

◆ getStructuresLevelNumber()

int IBAMR::IBInterpolantMethod::getStructuresLevelNumber ( ) const

◆ getStructureHandle()

int IBAMR::IBInterpolantMethod::getStructureHandle ( const int  lag_idx) const

◆ getNumberOfNodes()

unsigned int IBAMR::IBInterpolantMethod::getNumberOfNodes ( const unsigned int  struct_no) const
inline

◆ registerLSiloDataWriter()

void IBAMR::IBInterpolantMethod::registerLSiloDataWriter ( SAMRAI::tbox::Pointer< IBTK::LSiloDataWriter silo_writer)

Register a Lagrangian Silo data writer so this class will write plot files that may be postprocessed with the VisIt visualization tool.

◆ getMinimumGhostCellWidth()

const SAMRAI::hier::IntVector<NDIM>& IBAMR::IBInterpolantMethod::getMinimumGhostCellWidth ( ) const
overridevirtual

Return the number of ghost cells required by the Lagrangian-Eulerian interaction routines.

Implements IBAMR::IBStrategy.

◆ setupTagBuffer()

void IBAMR::IBInterpolantMethod::setupTagBuffer ( SAMRAI::tbox::Array< int > &  tag_buffer,
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg 
) const
overridevirtual

Setup the tag buffer.

Reimplemented from IBAMR::IBStrategy.

◆ preprocessIntegrateData()

void IBAMR::IBInterpolantMethod::preprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
overridevirtual

Method to prepare to advance data from current_time to new_time.

Reimplemented from IBAMR::IBStrategy.

◆ postprocessIntegrateData()

void IBAMR::IBInterpolantMethod::postprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
overridevirtual

Method to clean up data following call(s) to integrateHierarchy().

Reimplemented from IBAMR::IBStrategy.

◆ interpolateVelocity()

void IBAMR::IBInterpolantMethod::interpolateVelocity ( int  u_data_idx,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &  u_synch_scheds,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &  u_ghost_fill_scheds,
double  data_time 
)
overridevirtual

Interpolate the Eulerian velocity to the curvilinear mesh at the specified time within the current time interval.

Implements IBAMR::IBStrategy.

◆ interpolateQ() [1/2]

void IBAMR::IBInterpolantMethod::interpolateQ ( double  data_time)

Interpolate the Eulerian quantities to the curvilinear mesh at the specified time within the current time interval.

◆ interpolateQ() [2/2]

void IBAMR::IBInterpolantMethod::interpolateQ ( )

Interpolate the Eulerian quantities to the curvilinear mesh at both current and new time.

◆ forwardEulerStep()

void IBAMR::IBInterpolantMethod::forwardEulerStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the forward Euler method.

Implements IBAMR::IBStrategy.

◆ backwardEulerStep()

void IBAMR::IBInterpolantMethod::backwardEulerStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the backward Euler method.

Reimplemented from IBAMR::IBStrategy.

◆ midpointStep()

void IBAMR::IBInterpolantMethod::midpointStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the midpoint rule.

Implements IBAMR::IBStrategy.

◆ trapezoidalStep()

void IBAMR::IBInterpolantMethod::trapezoidalStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the trapezoidal rule.

Implements IBAMR::IBStrategy.

◆ updateMeshPosition()

void IBAMR::IBInterpolantMethod::updateMeshPosition ( double  current_time,
double  new_time,
const std::vector< Eigen::Vector3d > &  U,
const std::vector< Eigen::Vector3d > &  W 
)

Update the position of the mesh at new time.

X^n+1 = R[theta]*X^n + U*dt

◆ computeLagrangianForce()

void IBAMR::IBInterpolantMethod::computeLagrangianForce ( double  data_time)
overridevirtual

Compute the Lagrangian force at the specified time within the current time interval.

Implements IBAMR::IBStrategy.

◆ spreadForce()

void IBAMR::IBInterpolantMethod::spreadForce ( int  f_data_idx,
IBTK::RobinPhysBdryPatchStrategy f_phys_bdry_op,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &  f_prolongation_scheds,
double  data_time 
)
overridevirtual

Spread the Lagrangian force to the Cartesian grid at the specified time within the current time interval.

Implements IBAMR::IBStrategy.

◆ spreadQ()

void IBAMR::IBInterpolantMethod::spreadQ ( double  data_time)

Spread Q on the Eulerian grid.

◆ getEulerianQPatchDataIndex()

int IBAMR::IBInterpolantMethod::getEulerianQPatchDataIndex ( const std::string &  var_name)
inline

Get the patch data index of Q.

◆ copyEulerianDataToIntegrator() [1/2]

void IBAMR::IBInterpolantMethod::copyEulerianDataToIntegrator ( double  data_time)

Copy the spreaded data back to the integrator.

◆ initializePatchHierarchy()

void IBAMR::IBInterpolantMethod::initializePatchHierarchy ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg,
int  u_data_idx,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &  u_synch_scheds,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &  u_ghost_fill_scheds,
int  integrator_step,
double  init_data_time,
bool  initial_time 
)
overridevirtual

Initialize Lagrangian data corresponding to the given AMR patch hierarchy at the start of a computation. If the computation is begun from a restart file, data may be read from the restart databases.

A patch data descriptor is provided for the Eulerian velocity in case initialization requires interpolating Eulerian data. Ghost cells for Eulerian data will be filled upon entry to this function.

Reimplemented from IBAMR::IBStrategy.

◆ addWorkloadEstimate()

void IBAMR::IBInterpolantMethod::addWorkloadEstimate ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
const int  workload_data_idx 
)
overridevirtual

Add the estimated computational work from the current object per cell into the specified workload_data_idx.

Reimplemented from IBAMR::IBStrategy.

◆ beginDataRedistribution()

void IBAMR::IBInterpolantMethod::beginDataRedistribution ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg 
)
overridevirtual

Begin redistributing Lagrangian data prior to regridding the patch hierarchy.

Reimplemented from IBAMR::IBStrategy.

◆ endDataRedistribution()

void IBAMR::IBInterpolantMethod::endDataRedistribution ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg 
)
overridevirtual

Complete redistributing Lagrangian data following regridding the patch hierarchy.

Reimplemented from IBAMR::IBStrategy.

◆ initializeLevelData() [1/2]

void IBAMR::IBInterpolantMethod::initializeLevelData ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  level_number,
double  init_data_time,
bool  can_be_refined,
bool  initial_time,
SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > >  old_level,
bool  allocate_data 
)
override

Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm.

See also
SAMRAI::mesh::StandardTagAndInitStrategy::initializeLevelData

◆ resetHierarchyConfiguration() [1/2]

void IBAMR::IBInterpolantMethod::resetHierarchyConfiguration ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  coarsest_level,
int  finest_level 
)
override

Reset cached hierarchy dependent data.

See also
SAMRAI::mesh::StandardTagAndInitStrategy::resetHierarchyConfiguration

◆ applyGradientDetector() [1/2]

void IBAMR::IBInterpolantMethod::applyGradientDetector ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  level_number,
double  error_data_time,
int  tag_index,
bool  initial_time,
bool  uses_richardson_extrapolation_too 
)
override

Set integer tags to "one" in cells where refinement of the given level should occur according to user-supplied feature detection criteria.

See also
SAMRAI::mesh::StandardTagAndInitStrategy::applyGradientDetector

◆ putToDatabase()

void IBAMR::IBInterpolantMethod::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overridevirtual

Write out object state to the given database.

Implements SAMRAI::tbox::Serializable.

◆ getPositionData()

void IBAMR::IBInterpolantMethod::getPositionData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  X_data,
double  data_time 
)
protected

Get the structure position data.

◆ copyEulerianDataFromIntegrator()

void IBAMR::IBInterpolantMethod::copyEulerianDataFromIntegrator ( const std::string &  var_name,
int  q_data_idx,
double  data_time 
)
protected

Get the q data defined on Eulerian grid for interpolation.

◆ zeroOutEulerianData()

void IBAMR::IBInterpolantMethod::zeroOutEulerianData ( const std::string &  var_name,
int  q_data_idx 
)
protected

Get the q data defined on Eulerian grid for spreading.

◆ copyEulerianDataToIntegrator() [2/2]

void IBAMR::IBInterpolantMethod::copyEulerianDataToIntegrator ( const std::string &  name,
int  q_data_idx,
double  data_time 
)
protected

Copy the q_data_idx to suitable intergrator.

◆ getQData()

void IBAMR::IBInterpolantMethod::getQData ( const std::string &  var_name,
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  Q_data,
double  data_time 
)
protected

Get the Q data defined on Lagrangian mesh.

◆ computeCenterOfMass()

void IBAMR::IBInterpolantMethod::computeCenterOfMass ( std::vector< Eigen::Vector3d > &  center_of_mass,
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &  X_data 
)
protected

◆ operator=()

IBInterpolantMethod& IBAMR::IBInterpolantMethod::operator= ( const IBInterpolantMethod 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.

◆ getFromInput()

void IBAMR::IBInterpolantMethod::getFromInput ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db,
bool  is_from_restart 
)
private

Read input values from a given database.

◆ getFromRestart()

void IBAMR::IBInterpolantMethod::getFromRestart ( )
private

Read object state from the restart file and initialize class data members.

◆ registerIBHierarchyIntegrator()

virtual void IBAMR::IBStrategy::registerIBHierarchyIntegrator ( IBHierarchyIntegrator ib_solver)
virtualinherited

Register the IBHierarchyIntegrator object that is using this strategy class.

Reimplemented in IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ inactivateLagrangianStructure()

virtual void IBAMR::IBStrategy::inactivateLagrangianStructure ( int  structure_number = 0,
int  level_number = std::numeric_limits< int >::max() 
)
virtualinherited

Inactivate a structure or part. Such a structure will be ignored by FSI calculations: i.e., it will have its velocity set to zero and no forces will be spread from the structure to the Eulerian grid.

Parameters
[in]structure_numberNumber of the structure/part.
[in]level_numberLevel on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy.
Note
This method should be implemented by inheriting classes for which the notion of activating and inactivating a structure 'owned' by this class makes sense (for example, IBAMR::IBStrategySet does not directly own any structures, but IBAMR::IBMethod and IBAMR::IBFEMethod both do). The default implementation unconditionally aborts the program.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBMethod.

◆ activateLagrangianStructure()

virtual void IBAMR::IBStrategy::activateLagrangianStructure ( int  structure_number = 0,
int  level_number = std::numeric_limits< int >::max() 
)
virtualinherited

Activate a previously inactivated structure or part to be used again in FSI calculations.

Parameters
[in]structure_numberNumber of the structure/part.
[in]level_numberLevel on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy.
Note
This method should be implemented by inheriting classes for which the notion of activating and inactivating a structure 'owned' by this class makes sense (for example, IBAMR::IBStrategySet does not directly own any structures, but IBAMR::IBMethod and IBAMR::IBFEMethod both do). The default implementation unconditionally aborts the program.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBMethod.

◆ getLagrangianStructureIsActivated()

virtual bool IBAMR::IBStrategy::getLagrangianStructureIsActivated ( int  structure_number = 0,
int  level_number = std::numeric_limits< int >::max() 
) const
virtualinherited

Determine whether or not the given structure or part is currently activated.

Parameters
[in]structure_numberNumber of the structure/part.
[in]level_numberLevel on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy.
Note
This method should be implemented by inheriting classes for which the notion of activating and inactivating a structure 'owned' by this class makes sense (for example, IBAMR::IBStrategySet does not directly own any structures, but IBAMR::IBMethod and IBAMR::IBFEMethod both do). The default implementation unconditionally aborts the program.

Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBMethod.

◆ getMaxPointDisplacement()

virtual double IBAMR::IBStrategy::getMaxPointDisplacement ( ) const
virtualinherited

Get the ratio of the maximum point displacement of all the structures owned by the current class to the cell width of the grid level on which the structure is assigned. This value is useful for determining if the Eulerian patch hierarchy needs to be regridded.

Note
The process of regridding is distinct, for some IBStrategy objects (like IBFEMethod), from forming (or reforming) the association between Lagrangian structures and patches. In particular, this function computes the distance between the current position of the structure and the structure at the point of the last regrid, which may not be the same point at which we last rebuilt the structure-to-patch mappings. The reassociation check should be implemented in postprocessIntegrateData().

Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBStrategySet.

◆ setUseFixedLEOperators()

void IBAMR::IBStrategy::setUseFixedLEOperators ( bool  use_fixed_coupling_ops = true)
inherited

Indicate whether "fixed" interpolation and spreading operators should be used during Lagrangian-Eulerian interaction.

◆ updateFixedLEOperators()

virtual void IBAMR::IBStrategy::updateFixedLEOperators ( )
virtualinherited

Update the positions used for the "fixed" interpolation and spreading operators.

A default implementation is provided that emits an unrecoverable exception.

Reimplemented in IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ setUseMultistepTimeStepping()

virtual void IBAMR::IBStrategy::setUseMultistepTimeStepping ( unsigned int  n_previous_steps = 1)
virtualinherited

Indicate that multistep time stepping will be used.

A default implementation is provided that emits an unrecoverable exception.

Parameters
[in]n_previous_stepsNumber of previous solution values that can be used by the multistep scheme.

Reimplemented in IBAMR::IBFEMethod.

◆ AB2Step()

virtual void IBAMR::IBStrategy::AB2Step ( double  current_time,
double  new_time 
)
virtualinherited

Advance the positions of the Lagrangian structure using the standard 2nd-order Adams-Bashforth rule.

A default implementation is provided that emits an unrecoverable exception.

Reimplemented in IBAMR::IBFEMethod.

◆ hasFluidSources()

virtual bool IBAMR::IBStrategy::hasFluidSources ( ) const
virtualinherited

Indicate whether there are any internal fluid sources/sinks.

A default implementation is provided that returns false.

Reimplemented in IBAMR::IBFEMethod, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ computeLagrangianFluidSource()

virtual void IBAMR::IBStrategy::computeLagrangianFluidSource ( double  data_time)
virtualinherited

Compute the Lagrangian source/sink density at the specified time within the current time interval.

An empty default implementation is provided.

Reimplemented in IBAMR::IBFEMethod, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ spreadFluidSource()

virtual void IBAMR::IBStrategy::spreadFluidSource ( int  q_data_idx,
IBTK::RobinPhysBdryPatchStrategy q_phys_bdry_op,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &  q_prolongation_scheds,
double  data_time 
)
virtualinherited

Spread the Lagrangian source/sink density to the Cartesian grid at the specified time within the current time interval.

An empty default implementation is provided.

Reimplemented in IBAMR::IBFEMethod, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ interpolatePressure()

virtual void IBAMR::IBStrategy::interpolatePressure ( int  p_data_idx,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &  p_synch_scheds,
const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &  p_ghost_fill_scheds,
double  data_time 
)
virtualinherited

Compute the pressures at the positions of any distributed internal fluid sources or sinks.

An empty default implementation is provided.

Reimplemented in IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ preprocessSolveFluidEquations()

virtual void IBAMR::IBStrategy::preprocessSolveFluidEquations ( double  current_time,
double  new_time,
int  cycle_num 
)
virtualinherited

Execute user-defined routines just before solving the fluid equations.

An empty default implementation is provided.

Reimplemented in IBAMR::IBLevelSetMethod, IBAMR::IBStrategySet, IBAMR::CIBMethod, and IBAMR::ConstraintIBMethod.

◆ postprocessSolveFluidEquations()

virtual void IBAMR::IBStrategy::postprocessSolveFluidEquations ( double  current_time,
double  new_time,
int  cycle_num 
)
virtualinherited

Execute user-defined routines just after solving the fluid equations.

An empty default implementation is provided.

Reimplemented in IBAMR::IBLevelSetMethod, IBAMR::IBStrategySet, and IBAMR::ConstraintIBMethod.

◆ postprocessData()

virtual void IBAMR::IBStrategy::postprocessData ( )
virtualinherited

Execute user-defined post-processing operations.

An empty default implementation is provided.

Reimplemented in IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.

◆ registerLoadBalancer()

virtual void IBAMR::IBStrategy::registerLoadBalancer ( SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > >  load_balancer,
int  workload_data_idx 
)
virtualinherited

Register a load balancer and work load patch data index with the IB strategy object.

An empty default implementation is provided.

Deprecated:
This method is no longer necessary with the current workload estimation scheme.

Reimplemented in IBAMR::IBMethod, IBAMR::IBStrategySet, and IBAMR::IMPMethod.

◆ initializeLevelData() [2/2]

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::initializeLevelData ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const double  init_data_time,
const bool  can_be_refined,
const bool  initial_time,
const tbox::Pointer< hier::BasePatchLevel< DIM > >  old_level = tbox::Pointerhier::BasePatchLevel<DIM> >(NULL),
const bool  allocate_data = true 
)
pure virtualinherited

Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm. The level number indicates that of the new level.

Generally, when data is set, it is interpolated from coarser levels in the hierarchy. If the old level pointer in the argument list is non-null, then data is copied from the old level to the new level on regions of intersection between those levels before interpolation occurs. In this case, the level number must match that of the old level. The specific operations that occur when initializing level data are determined by the particular solution methods in use; i.e., in the subclass of this abstract base class.

The boolean argument initial_time indicates whether the level is being introduced for the first time (i.e., at initialization time), or after some regrid process during the calculation beyond the initial hierarchy construction. This information is provided since the initialization of the data may be different in each of those circumstances. The can_be_refined boolean argument indicates whether the level is the finest allowable level in the hierarchy.

◆ resetHierarchyConfiguration() [2/2]

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::resetHierarchyConfiguration ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  coarsest_level,
const int  finest_level 
)
pure virtualinherited

After hierarchy levels have changed and data has been initialized on the new levels, this routine can be used to reset any information needed by the solution method that is particular to the hierarchy configuration. For example, the solution procedure may cache communication schedules to amortize the cost of data movement on the AMR patch hierarchy. This function will be called by the gridding algorithm after the initialization occurs so that the algorithm-specific subclass can reset such things. Also, if the solution method must make the solution consistent across multiple levels after the hierarchy is changed, this process may be invoked by this routine. Of course the details of these processes are determined by the particular solution methods in use.

The level number arguments indicate the coarsest and finest levels in the current hierarchy configuration that have changed. It should be assumed that all intermediate levels have changed as well.

◆ applyGradientDetector() [2/2]

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::applyGradientDetector ( const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const double  error_data_time,
const int  tag_index,
const bool  initial_time,
const bool  uses_richardson_extrapolation_too 
)
virtualinherited

Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied gradient criteria. The double time argument is the regrid time. The integer "tag_index" argument is the patch descriptor index of the cell-centered integer tag array on each patch in the hierarchy. The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. Typically, this information is passed to the user's patch tagging routines since the error estimator or gradient detector may be different in each case.

The boolean uses_richardson_extrapolation_too is true when Richardson extrapolation error estimation is used in addition to the gradient detector, and false otherwise. This argument helps the user to manage multiple regridding criteria.

This routine is only when gradient detector is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

◆ getINSHierarchyIntegrator()

INSHierarchyIntegrator* IBAMR::IBStrategy::getINSHierarchyIntegrator ( ) const
protectedinherited

Return a pointer to the INSHierarchyIntegrator object being used with the IBHierarchyIntegrator class registered with this IBStrategy object.

◆ getVelocityHierarchyDataOps()

SAMRAI::tbox::Pointer<SAMRAI::math::HierarchyDataOpsReal<NDIM, double> > IBAMR::IBStrategy::getVelocityHierarchyDataOps ( ) const
protectedinherited

Return a pointer to the HierarchyDataOpsReal object associated with velocity-like variables.

◆ getPressureHierarchyDataOps()

SAMRAI::tbox::Pointer<SAMRAI::math::HierarchyDataOpsReal<NDIM, double> > IBAMR::IBStrategy::getPressureHierarchyDataOps ( ) const
protectedinherited

Return a pointer to the HierarchyDataOpsReal object associated with pressure-like variables.

◆ getHierarchyMathOps()

SAMRAI::tbox::Pointer<IBTK::HierarchyMathOps> IBAMR::IBStrategy::getHierarchyMathOps ( ) const
protectedinherited

Return a pointer to a HierarchyMathOps object.

◆ registerVariable() [1/2]

void IBAMR::IBStrategy::registerVariable ( int current_idx,
int new_idx,
int scratch_idx,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  variable,
const SAMRAI::hier::IntVector< NDIM > &  scratch_ghosts = SAMRAI::hier::IntVector< NDIM >(0),
const std::string &  coarsen_name = "NO_COARSEN",
const std::string &  refine_name = "NO_REFINE",
SAMRAI::tbox::Pointer< IBTK::CartGridFunction init_fcn = nullptr,
const bool  register_for_restart = true 
)
protectedinherited

Register a state variable with the integrator. When a refine operator is specified, the data for the variable are automatically maintained as the patch hierarchy evolves.

All state variables are registered with three contexts: current, new, and scratch. The current context of a state variable is maintained from time step to time step and, if the necessary coarsen and refine operators are specified, as the patch hierarchy evolves.

◆ registerVariable() [2/2]

void IBAMR::IBStrategy::registerVariable ( int idx,
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > >  variable,
const SAMRAI::hier::IntVector< NDIM > &  ghosts = SAMRAI::hier::IntVector< NDIM >(0),
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext ctx = SAMRAI::tbox::PointerSAMRAI::hier::VariableContext >(nullptr),
const bool  register_for_restart = true 
)
protectedinherited

Register a variable with the integrator that may not be maintained from time step to time step.

By default, variables are registered with the scratch context, which is deallocated after each time step.

◆ registerGhostfillRefineAlgorithm()

void IBAMR::IBStrategy::registerGhostfillRefineAlgorithm ( const std::string &  name,
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > >  ghostfill_alg,
std::unique_ptr< SAMRAI::xfer::RefinePatchStrategy< NDIM > >  ghostfill_patch_strategy = nullptr 
)
protectedinherited

Register a ghost cell-filling refine algorithm.

◆ registerProlongRefineAlgorithm()

void IBAMR::IBStrategy::registerProlongRefineAlgorithm ( const std::string &  name,
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > >  prolong_alg,
std::unique_ptr< SAMRAI::xfer::RefinePatchStrategy< NDIM > >  prolong_patch_strategy = nullptr 
)
protectedinherited

Register a data-prolonging refine algorithm.

◆ registerCoarsenAlgorithm()

void IBAMR::IBStrategy::registerCoarsenAlgorithm ( const std::string &  name,
SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > >  coarsen_alg,
std::unique_ptr< SAMRAI::xfer::CoarsenPatchStrategy< NDIM > >  coarsen_patch_strategy = nullptr 
)
protectedinherited

Register a coarsen algorithm.

◆ getGhostfillRefineAlgorithm()

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineAlgorithm<NDIM> > IBAMR::IBStrategy::getGhostfillRefineAlgorithm ( const std::string &  name) const
protectedinherited

Get ghost cell-filling refine algorithm.

◆ getProlongRefineAlgorithm()

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineAlgorithm<NDIM> > IBAMR::IBStrategy::getProlongRefineAlgorithm ( const std::string &  name) const
protectedinherited

Get data-prolonging refine algorithm.

◆ getCoarsenAlgorithm()

SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenAlgorithm<NDIM> > IBAMR::IBStrategy::getCoarsenAlgorithm ( const std::string &  name) const
protectedinherited

Get coarsen algorithm.

◆ getGhostfillRefineSchedules()

const std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >& IBAMR::IBStrategy::getGhostfillRefineSchedules ( const std::string &  name) const
protectedinherited

Get ghost cell-filling refine schedules.

◆ getProlongRefineSchedules()

const std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >& IBAMR::IBStrategy::getProlongRefineSchedules ( const std::string &  name) const
protectedinherited

Get data-prolonging refine schedules.

Note
These schedules are allocated only for level numbers >= 1.

◆ getCoarsenSchedules()

const std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenSchedule<NDIM> > >& IBAMR::IBStrategy::getCoarsenSchedules ( const std::string &  name) const
protectedinherited

Get coarsen schedules.

Note
These schedules are allocated only for level numbers >= 1.

◆ getLevelDt()

virtual double SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::getLevelDt ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level,
const double  dt_time,
const bool  initial_time 
)
virtualinherited

Determine time increment to advance data on level. The recompute_dt option specifies whether to compute the timestep using the current level data or to return the value stored by the time integrator. The default true setting means the timestep will be computed if no value is supplied.

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

◆ advanceLevel()

virtual double SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::advanceLevel ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level,
const tbox::Pointer< hier::BasePatchHierarchy< DIM > >  hierarchy,
const double  current_time,
const double  new_time,
const bool  first_step,
const bool  last_step,
const bool  regrid_advance = false 
)
virtualinherited

Advance data on all patches on specified patch level from current time (current_time) to new time (new_time). This routine is called only during time-dependent regridding procedures, such as Richardson extrapolation. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed. The boolean arguments are used to determine the state of the algorithm and the data when the advance routine is called. Note that this advance function is also used during normal time integration steps.

When this function is called, the level data required to begin the advance must be allocated and be defined appropriately. Typically, this is equivalent to what is needed to initialize a new level after regridding. Upon exiting this routine, both current and new data may exist on the level. This data is needed until level synchronization occurs, in general. Current and new data may be reset by calling the member function resetTimeDependentData().

This routine is called from two different points within the Richardson exptrapolation process: to advance a temporary level that is coarser than the hierarchy level on which error estimation is performed, and to advance the hierarchy level itself. In the first case, the values of the boolean flags are:

  • first_step = true.
  • last_step = true.
  • regrid_advance = true.

In the second case, the values of the boolean flags are:

  • first_step (when regridding during time integration sequence) = true when the level is not coarsest level to synchronize immediately before the regridding process; else, false. (when generating initial hierarchy construction) = true, even though there may be multiple advance steps.
  • last_step = true when the advance is the last in the Richardson extrapolation step sequence; else false.
  • regrid_advance = true.

◆ resetTimeDependentData()

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::resetTimeDependentData ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level,
const double  new_time,
const bool  can_be_refined 
)
virtualinherited

Reset time-dependent data storage for the specified patch level.

This routine only applies when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

◆ resetDataToPreadvanceState()

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::resetDataToPreadvanceState ( const tbox::Pointer< hier::BasePatchLevel< DIM > >  level)
virtualinherited

Reset data on the patch level by destroying all patch data other than that which is needed to initialize the solution on that level. In other words, this is the data needed to begin a time integration step on the level.

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

◆ applyRichardsonExtrapolation()

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::applyRichardsonExtrapolation ( const tbox::Pointer< hier::PatchLevel< DIM > >  level,
const double  error_data_time,
const int  tag_index,
const double  deltat,
const int  error_coarsen_ratio,
const bool  initial_time,
const bool  uses_gradient_detector_too 
)
virtualinherited

Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied Richardson extrapolation criteria. The "error_data_time" argument is the regrid time. The "deltat" argument is the time increment to advance the solution on the level to be refined. Note that that level is finer than the level in the argument list, in general. The ratio between the argument level and the actual hierarchy level is given by the integer "coarsen ratio".

The integer "tag_index" argument is the patch descriptor index of the cell-centered integer tag array on each patch in the hierarchy.

The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. Typically, this information is passed to the user's patch tagging routines since the application of the Richardson extrapolation process may be different in each case.

The boolean uses_gradient_detector_too is true when a gradient detector procedure is used in addition to Richardson extrapolation, and false otherwise. This argument helps the user to manage multiple regridding criteria.

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

◆ coarsenDataForRichardsonExtrapolation()

virtual void SAMRAI::mesh::StandardTagAndInitStrategy< DIM >::coarsenDataForRichardsonExtrapolation ( const tbox::Pointer< hier::PatchHierarchy< DIM > >  hierarchy,
const int  level_number,
const tbox::Pointer< hier::PatchLevel< DIM > >  coarser_level,
const double  coarsen_data_time,
const bool  before_advance 
)
virtualinherited

Coarsen solution data from level to coarse_level for Richardson extrapolation. Note that this routine will be called twice during the Richardson extrapolation error estimation process, once to set data on the coarser level and once to coarsen data from after advancing the fine level. The init_coarse_level boolean argument indicates whether data is set on the coarse level by coarsening the "old" time level solution or by coarsening the "new" solution on the fine level (i.e., after it has been advanced).

This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.

Member Data Documentation

◆ d_do_log

bool IBAMR::IBInterpolantMethod::d_do_log = false
protected

◆ d_hierarchy

SAMRAI::tbox::Pointer<SAMRAI::hier::PatchHierarchy<NDIM> > IBAMR::IBInterpolantMethod::d_hierarchy
protected

◆ d_gridding_alg

SAMRAI::tbox::Pointer<SAMRAI::mesh::GriddingAlgorithm<NDIM> > IBAMR::IBInterpolantMethod::d_gridding_alg
protected

◆ d_current_time

double IBAMR::IBInterpolantMethod::d_current_time = std::numeric_limits<double>::quiet_NaN()
protected

◆ d_new_time

double IBAMR::IBInterpolantMethod::d_new_time = std::numeric_limits<double>::quiet_NaN()
protected

◆ d_half_time

double IBAMR::IBInterpolantMethod::d_half_time = std::numeric_limits<double>::quiet_NaN()
protected

◆ d_struct_lag_idx_range

std::vector<std::pair<int, int> > IBAMR::IBInterpolantMethod::d_struct_lag_idx_range
protected

Vector of Lagrangian indices of all structures.

◆ d_num_rigid_parts

unsigned int IBAMR::IBInterpolantMethod::d_num_rigid_parts
protected

◆ d_center_of_mass_initial

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

◆ d_center_of_mass_current

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

◆ d_center_of_mass_new

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

◆ d_quaternion_current

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

◆ d_quaternion_new

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

◆ d_l_data_manager

IBTK::LDataManager* IBAMR::IBInterpolantMethod::d_l_data_manager
protected

◆ d_interp_kernel_fcn

std::string IBAMR::IBInterpolantMethod::d_interp_kernel_fcn = "IB_4"
protected

◆ d_spread_kernel_fcn

std::string IBAMR::IBInterpolantMethod::d_spread_kernel_fcn = "IB_4"
protected

◆ d_error_if_points_leave_domain

bool IBAMR::IBInterpolantMethod::d_error_if_points_leave_domain = false
protected

◆ d_ghosts

SAMRAI::hier::IntVector<NDIM> IBAMR::IBInterpolantMethod::d_ghosts
protected

◆ d_q_var

std::map<std::string, SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > > IBAMR::IBInterpolantMethod::d_q_var
protected

◆ d_q_depth

std::map<std::string, int> IBAMR::IBInterpolantMethod::d_q_depth
protected

◆ d_q_interp_idx

std::map<std::string, int> IBAMR::IBInterpolantMethod::d_q_interp_idx
protected

◆ d_q_hier_integrator

std::map<std::string, SAMRAI::tbox::Pointer<IBTK::HierarchyIntegrator> > IBAMR::IBInterpolantMethod::d_q_hier_integrator
protected

◆ d_X_current_data

std::vector<SAMRAI::tbox::Pointer<IBTK::LData> > IBAMR::IBInterpolantMethod::d_X_current_data
protected

◆ d_X_new_data

std::vector<SAMRAI::tbox::Pointer<IBTK::LData> > IBAMR::IBInterpolantMethod::d_X_new_data
protected

◆ d_Q_current_data

std::map<std::string, std::vector<SAMRAI::tbox::Pointer<IBTK::LData> > > IBAMR::IBInterpolantMethod::d_Q_current_data
protected

◆ d_Q_new_data

std::map<std::string, std::vector<SAMRAI::tbox::Pointer<IBTK::LData> > > IBAMR::IBInterpolantMethod::d_Q_new_data
protected

◆ d_l_initializer

SAMRAI::tbox::Pointer<IBTK::LInitStrategy> IBAMR::IBInterpolantMethod::d_l_initializer
protected

◆ d_silo_writer

SAMRAI::tbox::Pointer<IBTK::LSiloDataWriter> IBAMR::IBInterpolantMethod::d_silo_writer
protected

◆ d_object_name

std::string IBAMR::IBInterpolantMethod::d_object_name
protected

◆ d_registered_for_restart

bool IBAMR::IBInterpolantMethod::d_registered_for_restart
protected

◆ d_ib_solver

IBHierarchyIntegrator* IBAMR::IBStrategy::d_ib_solver = nullptr
protectedinherited

The IBHierarchyIntegrator object that is using this strategy class.

◆ d_use_fixed_coupling_ops

bool IBAMR::IBStrategy::d_use_fixed_coupling_ops = false
protectedinherited

Whether to use "fixed" Lagrangian-Eulerian coupling operators.


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