IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
IBAMR::IBInterpolantMethod Class Reference

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 </home/runner/work/IBAMR/IBAMR/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.
 
 ~IBInterpolantMethod ()
 Destructor.
 
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.
 
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.
 
int getStructureHandle (const int lag_idx) const
 Get the structure handle to which this Lagrangian index belongs.
 
unsigned int getNumberOfNodes (const unsigned int struct_no) const
 Get the number of nodes for this structure.
 
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 IBTK::EigenAlignedVector< Eigen::Vector3d > &U, const IBTK::EigenAlignedVector< 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
 
- Public Member Functions inherited from IBAMR::IBStrategy
 IBStrategy ()=default
 Constructor.
 
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 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)
 
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
 
- Public Member Functions inherited from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >
virtual double getLevelDt (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double dt_time, const bool initial_time)
 
virtual double advanceLevel (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< NDIM > > 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< NDIM > > level, const double new_time, const bool can_be_refined)
 
virtual void resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level)
 
virtual void applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< NDIM > > 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< NDIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< NDIM > > 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 (IBTK::EigenAlignedVector< Eigen::Vector3d > &center_of_mass, std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &X_data)
 Compute the center of mass of the mesh.
 
- Protected Member Functions inherited from IBAMR::IBStrategy
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=SAMRAI::tbox::Pointer< IBTK::CartGridFunction >(NULL))
 
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 >(NULL))
 
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
 
IBTK::EigenAlignedVector< Eigen::Vector3d > d_center_of_mass_initial
 
IBTK::EigenAlignedVector< Eigen::Vector3d > d_center_of_mass_current
 
IBTK::EigenAlignedVector< Eigen::Vector3d > d_center_of_mass_new
 
IBTK::EigenAlignedVector< Eigen::Quaterniond > d_quaternion_current
 
IBTK::EigenAlignedVector< 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.
 
std::map< std::string, int > d_q_depth
 
std::map< std::string, int > d_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
 
- Protected Attributes inherited from IBAMR::IBStrategy
IBHierarchyIntegratord_ib_solver = nullptr
 
bool d_use_fixed_coupling_ops = false
 

Detailed Description

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.

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

Member Function Documentation

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

◆ applyGradientDetector()

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 
)
overridevirtual

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

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.

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

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

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

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

◆ copyEulerianDataToIntegrator() [1/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.

◆ copyEulerianDataToIntegrator() [2/2]

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

Copy the spreaded data back to the integrator.

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

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

◆ freeLInitStrategy()

void IBAMR::IBInterpolantMethod::freeLInitStrategy ( )

Free references to Lagrangian initialization objects.

◆ getEulerianQPatchDataIndex()

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

Get the patch data index of Q.

◆ getLDataManager()

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

Return a pointer to the Lagrangian data manager object.

◆ getMinimumGhostCellWidth()

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

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

Implements IBAMR::IBStrategy.

◆ getPositionData()

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

Get the structure position data.

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

◆ initializeLevelData()

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 
)
overridevirtual

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

Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.

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

◆ interpolateQ() [1/2]

void IBAMR::IBInterpolantMethod::interpolateQ ( )

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

◆ interpolateQ() [2/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.

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

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

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

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

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

◆ registerEulerianCommunicationAlgorithms()

void IBAMR::IBInterpolantMethod::registerEulerianCommunicationAlgorithms ( )
overridevirtual

Register ghost cell filling algrorithms for the Eulerian variables.

Reimplemented from IBAMR::IBStrategy.

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

◆ registerLInitStrategy()

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

Supply a Lagrangian initialization object.

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

◆ resetHierarchyConfiguration()

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

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

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

◆ 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 IBTK::EigenAlignedVector< Eigen::Vector3d > &  U,
const IBTK::EigenAlignedVector< Eigen::Vector3d > &  W 
)

Update the position of the mesh at new time.

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

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

Member Data Documentation

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


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