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::IBStrategy Class Referenceabstract

Class IBStrategy provides a generic interface for specifying the implementation details of a particular version of the IB method. More...

#include </home/runner/work/IBAMR/IBAMR/include/ibamr/IBStrategy.h>

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

Public Member Functions

 IBStrategy ()=default
 Constructor.
 
virtual void registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver)
 
virtual void registerEulerianVariables ()
 
virtual void registerEulerianCommunicationAlgorithms ()
 
virtual const SAMRAI::hier::IntVector< NDIM > & getMinimumGhostCellWidth () const =0
 
virtual void setupTagBuffer (SAMRAI::tbox::Array< int > &tag_buffer, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) const
 
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
 
virtual void preprocessIntegrateData (double current_time, double new_time, int num_cycles)
 
virtual void postprocessIntegrateData (double current_time, double new_time, int num_cycles)
 
void setUseFixedLEOperators (bool use_fixed_coupling_ops=true)
 
virtual void updateFixedLEOperators ()
 
virtual 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)=0
 
virtual void forwardEulerStep (double current_time, double new_time)=0
 
virtual void backwardEulerStep (double current_time, double new_time)
 
virtual void midpointStep (double current_time, double new_time)=0
 
virtual void trapezoidalStep (double current_time, double new_time)=0
 
virtual void computeLagrangianForce (double data_time)=0
 
virtual 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)=0
 
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 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)
 
virtual void registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx)
 
virtual void addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx)
 
virtual void beginDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg)
 
virtual void endDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg)
 
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

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

IBHierarchyIntegratord_ib_solver = nullptr
 
bool d_use_fixed_coupling_ops = false
 

Detailed Description

Class IBStrategy provides a generic interface for specifying the implementation details of a particular version of the IB method.

Options Controlling Interpolation and Spreading

All classes inheriting from IBStrategy read the following values from the input database to determine which interpolation and spreading delta functions should be used:

Member Function Documentation

◆ activateLagrangianStructure()

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

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

◆ addWorkloadEstimate()

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

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

An empty default implementation is provided.

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

◆ applyGradientDetector()

void IBAMR::IBStrategy::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

An empty default implementation is provided.

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

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

◆ backwardEulerStep()

void IBAMR::IBStrategy::backwardEulerStep ( double  current_time,
double  new_time 
)
virtual

Advance the positions of the Lagrangian structure using the (explicit) backward Euler method.

A default implementation is provided that emits an unrecoverable exception.

Reimplemented in IBAMR::IBMethod, IBAMR::IBInterpolantMethod, and IBAMR::CIBMethod.

◆ beginDataRedistribution()

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

Begin redistributing Lagrangian data prior to regridding the patch hierarchy.

An empty default implementation is provided.

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

◆ computeLagrangianFluidSource()

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

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

An empty default implementation is provided.

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

◆ computeLagrangianForce()

virtual void IBAMR::IBStrategy::computeLagrangianForce ( double  data_time)
pure virtual

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

Implemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, and IBAMR::GeneralizedIBMethod.

◆ endDataRedistribution()

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

Complete redistributing Lagrangian data following regridding the patch hierarchy.

An empty default implementation is provided.

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

◆ forwardEulerStep()

virtual void IBAMR::IBStrategy::forwardEulerStep ( double  current_time,
double  new_time 
)
pure virtual

◆ getCoarsenAlgorithm()

Pointer< CoarsenAlgorithm< NDIM > > IBAMR::IBStrategy::getCoarsenAlgorithm ( const std::string name) const
protected

Get coarsen algorithm.

◆ getCoarsenSchedules()

const std::vector< Pointer< CoarsenSchedule< NDIM > > > & IBAMR::IBStrategy::getCoarsenSchedules ( const std::string name) const
protected

Get coarsen schedules.

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

◆ getGhostfillRefineAlgorithm()

Pointer< RefineAlgorithm< NDIM > > IBAMR::IBStrategy::getGhostfillRefineAlgorithm ( const std::string name) const
protected

Get ghost cell-filling refine algorithm.

◆ getGhostfillRefineSchedules()

const std::vector< Pointer< RefineSchedule< NDIM > > > & IBAMR::IBStrategy::getGhostfillRefineSchedules ( const std::string name) const
protected

Get ghost cell-filling refine schedules.

◆ getHierarchyMathOps()

Pointer< HierarchyMathOps > IBAMR::IBStrategy::getHierarchyMathOps ( ) const
protected

Return a pointer to a HierarchyMathOps object.

◆ getINSHierarchyIntegrator()

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

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

◆ getLagrangianStructureIsActivated()

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

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

◆ getMaxPointDisplacement()

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

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

◆ getMinimumGhostCellWidth()

virtual const SAMRAI::hier::IntVector<NDIM>& IBAMR::IBStrategy::getMinimumGhostCellWidth ( ) const
pure virtual

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

Implemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBInterpolantMethod.

◆ getPressureHierarchyDataOps()

Pointer< HierarchyDataOpsReal< NDIM, double > > IBAMR::IBStrategy::getPressureHierarchyDataOps ( ) const
protected

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

◆ getProlongRefineAlgorithm()

Pointer< RefineAlgorithm< NDIM > > IBAMR::IBStrategy::getProlongRefineAlgorithm ( const std::string name) const
protected

Get data-prolonging refine algorithm.

◆ getProlongRefineSchedules()

const std::vector< Pointer< RefineSchedule< NDIM > > > & IBAMR::IBStrategy::getProlongRefineSchedules ( const std::string name) const
protected

Get data-prolonging refine schedules.

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

◆ getVelocityHierarchyDataOps()

Pointer< HierarchyDataOpsReal< NDIM, double > > IBAMR::IBStrategy::getVelocityHierarchyDataOps ( ) const
protected

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

◆ hasFluidSources()

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

Indicate whether there are any internal fluid sources/sinks.

A default implementation is provided that returns false.

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

◆ inactivateLagrangianStructure()

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

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

◆ initializeLevelData()

void IBAMR::IBStrategy::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

An empty default implementation is provided.

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

Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, and IBAMR::IBStrategySet.

◆ initializePatchHierarchy()

void IBAMR::IBStrategy::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 
)
virtual

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.

An empty default implementation is provided.

Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.

◆ interpolatePressure()

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

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

An empty default implementation is provided.

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

◆ interpolateVelocity()

virtual void IBAMR::IBStrategy::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 
)
pure virtual

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

Implemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.

◆ midpointStep()

virtual void IBAMR::IBStrategy::midpointStep ( double  current_time,
double  new_time 
)
pure virtual

◆ postprocessData()

void IBAMR::IBStrategy::postprocessData ( )
virtual

Execute user-defined post-processing operations.

An empty default implementation is provided.

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

◆ postprocessIntegrateData()

void IBAMR::IBStrategy::postprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
virtual

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

An empty default implementation is provided.

Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.

◆ postprocessSolveFluidEquations()

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

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

An empty default implementation is provided.

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

◆ preprocessIntegrateData()

void IBAMR::IBStrategy::preprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
virtual

Method to prepare to advance data from current_time to new_time.

An empty default implementation is provided.

Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.

◆ preprocessSolveFluidEquations()

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

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

An empty default implementation is provided.

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

◆ putToDatabase()

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

Write out object state to the given database.

An empty default implementation is provided.

Implements SAMRAI::tbox::Serializable.

Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, and IBAMR::IBStrategySet.

◆ 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 
)
protected

Register a coarsen algorithm.

◆ registerEulerianCommunicationAlgorithms()

void IBAMR::IBStrategy::registerEulerianCommunicationAlgorithms ( )
virtual

Register Eulerian refinement or coarsening algorithms with the parent IBHierarchyIntegrator using the two versions of the protected methods IBStrategy::registerGhostfillRefineAlgorithm(), IBStrategy::registerProlongRefineAlgorithm(), and IBStrategy::registerCoarsenAlgorithm().

An empty default implementation is provided.

Reimplemented in IBAMR::IBStrategySet, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.

◆ registerEulerianVariables()

void IBAMR::IBStrategy::registerEulerianVariables ( )
virtual

Register Eulerian variables with the parent IBHierarchyIntegrator with the VariableDatabase, or via the various versions of the protected method IBStrategy::registerVariable().

An empty default implementation is provided.

Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.

◆ 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 
)
protected

Register a ghost cell-filling refine algorithm.

◆ registerIBHierarchyIntegrator()

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

Register the IBHierarchyIntegrator object that is using this strategy class.

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

◆ registerLoadBalancer()

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

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::IIMethod, IBAMR::IBStrategySet, and IBAMR::IBMethod.

◆ 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 
)
protected

Register a data-prolonging refine algorithm.

◆ 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 = SAMRAI::tbox::Pointer<IBTK::CartGridFunction>(NULL) 
)
protected

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::Pointer<SAMRAI::hier::VariableContext>(NULL) 
)
protected

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.

◆ resetHierarchyConfiguration()

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

Reset cached hierarchy dependent data.

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

An empty default implementation is provided.

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

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

◆ setupTagBuffer()

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

Setup the tag buffer.

A default implementation is provided that sets the tag buffer to be at least the minimum ghost cell width.

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

◆ setUseFixedLEOperators()

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

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

◆ spreadFluidSource()

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

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::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.

◆ spreadForce()

virtual void IBAMR::IBStrategy::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 
)
pure virtual

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

Implemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.

◆ trapezoidalStep()

virtual void IBAMR::IBStrategy::trapezoidalStep ( double  current_time,
double  new_time 
)
pure virtual

Advance the positions of the Lagrangian structure using the (explicit) trapezoidal rule.

Implemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.

◆ updateFixedLEOperators()

void IBAMR::IBStrategy::updateFixedLEOperators ( )
virtual

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

A default implementation is provided that emits an unrecoverable exception.

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

Member Data Documentation

◆ d_ib_solver

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

The IBHierarchyIntegrator object that is using this strategy class.

◆ d_use_fixed_coupling_ops

bool IBAMR::IBStrategy::d_use_fixed_coupling_ops = false
protected

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


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