IBAMR  IBAMR version 0.19.
Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
IBAMR::IBHierarchyIntegrator Class Referenceabstract

Class IBHierarchyIntegrator provides an abstract interface for a time integrator for various versions of the immersed boundary method on an AMR grid hierarchy, along with basic data management for variables defined on that hierarchy. More...

#include <ibamr/IBHierarchyIntegrator.h>

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

Classes

class  IBEulerianForceFunction
 A class to communicate the Eulerian body force computed by class IBHierarchyIntegrator to the incompressible Navier-Stokes solver. More...
 
class  IBEulerianSourceFunction
 A class to communicate the Eulerian fluid source-sink distribution computed by class IBHierarchyIntegrator to the incompressible Navier-Stokes solver. More...
 

Public Types

using PreprocessIntegrateHierarchyCallbackFcnPtr = void(*)(double current_time, double new_time, int num_cycles, void *ctx)
 
using IntegrateHierarchyCallbackFcnPtr = void(*)(double current_time, double new_time, int cycle_num, void *ctx)
 
using PostprocessIntegrateHierarchyCallbackFcnPtr = void(*)(double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles, void *ctx)
 
using ApplyGradientDetectorCallbackFcnPtr = void(*)(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, void *ctx)
 
using RegridHierarchyCallbackFcnPtr = void(*)(SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, double data_time, bool initial_time, void *ctx)
 

Public Member Functions

 ~IBHierarchyIntegrator ()=default
 
TimeSteppingType getTimeSteppingType () const
 
SAMRAI::tbox::Pointer< IBStrategygetIBStrategy () const
 
void registerBodyForceFunction (SAMRAI::tbox::Pointer< IBTK::CartGridFunction > F_fcn)
 
virtual void registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer) override
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > getVelocityVariable () const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > getPressureVariable () const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > getBodyForceVariable () const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > getFluidSourceVariable () const
 
IBTK::RobinPhysBdryPatchStrategygetVelocityPhysBdryOp () const
 
void preprocessIntegrateHierarchy (double current_time, double new_time, int num_cycles=1) override
 
void postprocessIntegrateHierarchy (double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles=1) override
 
void initializeHierarchyIntegrator (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
void initializePatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
const std::string & getName () const
 
virtual void advanceHierarchy (double dt)
 
double getMinimumTimeStepSize ()
 
double getMaximumTimeStepSize ()
 
void synchronizeHierarchyData (VariableContextType ctx_type)
 
void resetTimeDependentHierarchyData (double new_time)
 
void resetIntegratorToPreadvanceState ()
 
virtual void regridHierarchy ()
 
bool atRegridPoint () const
 
double getIntegratorTime () const
 
double getStartTime () const
 
double getEndTime () const
 
int getIntegratorStep () const
 
int getMaxIntegratorSteps () const
 
bool stepsRemaining () const
 
void updateWorkloadEstimates ()
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > getPatchHierarchy () const
 
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > getGriddingAlgorithm () const
 
int getWorkloadDataIndex () const
 
void registerVisItDataWriter (SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > visit_writer)
 
SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > getVisItDataWriter () const
 
void setupPlotData ()
 
virtual int getNumberOfCycles () const
 
virtual int getCurrentCycleNumber () const
 
virtual double getCurrentTimeStepSize () const
 
void integrateHierarchy (double current_time, double new_time, int cycle_num=0)
 
void skipCycle (double current_time, double new_time, int cycle_num=0)
 
void registerPreprocessIntegrateHierarchyCallback (PreprocessIntegrateHierarchyCallbackFcnPtr callback, void *ctx=nullptr)
 
void registerIntegrateHierarchyCallback (IntegrateHierarchyCallbackFcnPtr callback, void *ctx=nullptr)
 
void registerPostprocessIntegrateHierarchyCallback (PostprocessIntegrateHierarchyCallbackFcnPtr callback, void *ctx=nullptr)
 
void registerApplyGradientDetectorCallback (ApplyGradientDetectorCallbackFcnPtr callback, void *ctx=nullptr)
 
void registerRegridHierarchyCallback (RegridHierarchyCallbackFcnPtr, void *ctx=nullptr)
 
void initializeCompositeHierarchyData (double init_data_time, bool initial_time)
 
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=SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > >(nullptr), bool allocate_data=true) override
 
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
 
void resetHierarchyConfiguration (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_level, int finest_level) override
 
virtual void resetHierarchyConfiguration (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level)=0
 
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
 
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)
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextgetContext (VariableContextType ctx_type) const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextgetCurrentContext () const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextgetNewContext () const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextgetScratchContext () const
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextgetPlotContext () const
 
bool isAllocatedPatchData (int data_idx, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) const
 
void allocatePatchData (int data_idx, double data_time, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) const
 
void deallocatePatchData (int data_idx, int coarsest_ln=invalid_level_number, int finest_ln=invalid_level_number) const
 
SAMRAI::tbox::Pointer< HierarchyMathOps > getHierarchyMathOps () 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< 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 putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 
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

virtual void computeFluidSources (const int data_idx, const double data_time)
 
void regridHierarchyBeginSpecialized () override
 
void regridHierarchyEndSpecialized () override
 
 IBHierarchyIntegrator (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, SAMRAI::tbox::Pointer< IBStrategy > ib_method_ops, SAMRAI::tbox::Pointer< INSHierarchyIntegrator > ins_hier_integrator, bool register_for_restart=true)
 
bool atRegridPointSpecialized () const override
 
void initializeLevelDataSpecialized (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 resetHierarchyConfigurationSpecialized (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_level, int finest_level) override
 
void applyGradientDetectorSpecialized (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 putToDatabaseSpecialized (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override
 
virtual void addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx) override
 
virtual void integrateHierarchySpecialized (double current_time, double new_time, int cycle_num=0)=0
 
virtual double getMinimumTimeStepSizeSpecialized ()
 
virtual double getMaximumTimeStepSizeSpecialized ()
 
virtual void synchronizeHierarchyDataSpecialized (VariableContextType ctx_type)
 
virtual void resetTimeDependentHierarchyDataSpecialized (double new_time)
 
virtual void resetIntegratorToPreadvanceStateSpecialized ()
 
virtual void setupPlotDataSpecialized ()
 
virtual void initializeCompositeHierarchyDataSpecialized (double init_data_time, bool initial_time)
 
virtual void executePreprocessIntegrateHierarchyCallbackFcns (double current_time, double new_time, int num_cycles)
 
virtual void executeIntegrateHierarchyCallbackFcns (double current_time, double new_time, int cycle_num)
 
virtual void executePostprocessIntegrateHierarchyCallbackFcns (double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles)
 
virtual void executeApplyGradientDetectorCallbackFcns (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)
 
virtual void executeRegridHierarchyCallbackFcns (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, double data_time, bool initial_time)
 
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
 
void registerChildHierarchyIntegrator (HierarchyIntegrator *child_integrator)
 
void registerParentHierarchyIntegrator (HierarchyIntegrator *parent_integrator)
 
SAMRAI::tbox::Pointer< HierarchyMathOps > buildHierarchyMathOps (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy)
 
void setupTagBuffer (SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg)
 
bool regriddingHierarchy () const
 
bool atRegridTimeStep () const
 

Protected Attributes

bool d_integrator_is_initialized = false
 
TimeSteppingType d_time_stepping_type = MIDPOINT_RULE
 
bool d_use_structure_predictor
 
bool d_error_on_dt_change = true
 
bool d_warn_on_dt_change = false
 
SAMRAI::tbox::Pointer< INSHierarchyIntegratord_ins_hier_integrator
 
double d_regrid_fluid_cfl_interval = -1.0
 
double d_regrid_structure_cfl_interval = -1.0
 
double d_regrid_fluid_cfl_estimate = 0.0
 
double d_regrid_structure_cfl_estimate = 0.0
 
SAMRAI::tbox::Pointer< IBStrategyd_ib_method_ops
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > d_hier_velocity_data_ops
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > d_hier_pressure_data_ops
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyCellDataOpsReal< NDIM, double > > d_hier_cc_data_ops
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_u_var
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_p_var
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_f_var
 
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > d_q_var
 
int d_u_idx = IBTK::invalid_index
 
int d_p_idx = IBTK::invalid_index
 
int d_f_idx = IBTK::invalid_index
 
int d_f_current_idx = IBTK::invalid_index
 
int d_q_idx = IBTK::invalid_index
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextd_ib_context
 
SAMRAI::hier::ComponentSelector d_ib_data
 
IBTK::RobinPhysBdryPatchStrategyd_u_phys_bdry_op
 
IBTK::RobinPhysBdryPatchStrategyd_p_phys_bdry_op
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > d_u_ghostfill_alg
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > d_f_prolong_alg
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > d_p_ghostfill_alg
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > d_q_prolong_alg
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineOperator< NDIM > > d_u_ghostfill_op
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineOperator< NDIM > > d_f_prolong_op
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineOperator< NDIM > > d_p_ghostfill_op
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineOperator< NDIM > > d_q_prolong_op
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > > d_u_coarsen_alg
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > > d_p_coarsen_alg
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenOperator< NDIM > > d_u_coarsen_op
 
SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenOperator< NDIM > > d_p_coarsen_op
 
SAMRAI::tbox::Pointer< IBTK::CartGridFunctiond_body_force_fcn
 
std::string d_object_name
 
bool d_registered_for_restart = false
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > d_hierarchy
 
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > d_gridding_alg
 
SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > d_load_balancer
 
SAMRAI::tbox::Pointer< SAMRAI::pdat::CellVariable< NDIM, double > > d_workload_var
 
int d_workload_idx = IBTK::invalid_index
 
bool d_hierarchy_is_initialized = false
 
bool d_may_need_to_reset_hierarchy_configuration = false
 
HierarchyIntegratord_parent_integrator = nullptr
 
std::set< HierarchyIntegrator * > d_child_integrators
 
SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > d_visit_writer
 
double d_integrator_time = std::numeric_limits<double>::quiet_NaN()
 
double d_start_time = 0.0
 
double d_end_time = std::numeric_limits<double>::max()
 
double d_dt_init = std::numeric_limits<double>::max()
 
double d_dt_min = 0.0
 
double d_dt_max = std::numeric_limits<double>::max()
 
double d_dt_growth_factor = 2.0
 
int d_integrator_step = 0
 
int d_max_integrator_steps = std::numeric_limits<int>::max()
 
std::deque< doubled_dt_previous
 
int d_num_cycles = 1
 
int d_current_num_cycles = -1
 
int d_current_cycle_num = -1
 
double d_current_dt = std::numeric_limits<double>::quiet_NaN()
 
int d_regrid_interval = 1
 
RegridMode d_regrid_mode = STANDARD
 
bool d_enable_logging = false
 
bool d_enable_logging_solver_iterations = false
 
std::string d_bdry_extrap_type = "LINEAR"
 
SAMRAI::tbox::Array< intd_tag_buffer = { 0 }
 
SAMRAI::tbox::Pointer< HierarchyMathOps > d_hier_math_ops
 
bool d_manage_hier_math_ops = true
 
std::list< SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > > d_state_variables
 
std::list< SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > > d_scratch_variables
 
std::list< SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > > d_copy_scratch_to_current_fast
 
std::list< SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > > d_copy_scratch_to_current_slow
 
SAMRAI::hier::ComponentSelector d_current_data
 
SAMRAI::hier::ComponentSelector d_new_data
 
SAMRAI::hier::ComponentSelector d_scratch_data
 
SAMRAI::hier::ComponentSelector d_plot_data
 
std::map< SAMRAI::hier::Variable< NDIM > *, SAMRAI::tbox::Pointer< CartGridFunction > > d_state_var_init_fcns
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextd_current_context
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextd_new_context
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextd_scratch_context
 
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContextd_plot_context
 
SAMRAI::hier::ComponentSelector d_fill_after_regrid_bc_idxs
 
SAMRAI::xfer::RefineAlgorithm< NDIM > d_fill_after_regrid_prolong_alg
 
std::unique_ptr< SAMRAI::xfer::RefinePatchStrategy< NDIM > > d_fill_after_regrid_phys_bdry_bc_op
 
std::vector< PreprocessIntegrateHierarchyCallbackFcnPtrd_preprocess_integrate_hierarchy_callbacks
 
std::vector< void * > d_preprocess_integrate_hierarchy_callback_ctxs
 
std::vector< IntegrateHierarchyCallbackFcnPtrd_integrate_hierarchy_callbacks
 
std::vector< void * > d_integrate_hierarchy_callback_ctxs
 
std::vector< PostprocessIntegrateHierarchyCallbackFcnPtrd_postprocess_integrate_hierarchy_callbacks
 
std::vector< void * > d_postprocess_integrate_hierarchy_callback_ctxs
 
std::vector< ApplyGradientDetectorCallbackFcnPtrd_apply_gradient_detector_callbacks
 
std::vector< void * > d_apply_gradient_detector_callback_ctxs
 
std::vector< RegridHierarchyCallbackFcnPtrd_regrid_hierarchy_callbacks
 
std::vector< void * > d_regrid_hierarchy_callback_ctxs
 

Static Protected Attributes

static const std::string SYNCH_CURRENT_DATA_ALG
 
static const std::string SYNCH_NEW_DATA_ALG
 

Private Types

using RefineAlgorithmMap = std::map< std::string, SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineAlgorithm< NDIM > > >
 
using RefinePatchStrategyMap = std::map< std::string, std::unique_ptr< SAMRAI::xfer::RefinePatchStrategy< NDIM > > >
 
using RefineScheduleMap = std::map< std::string, std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > >
 
using CoarsenAlgorithmMap = std::map< std::string, SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenAlgorithm< NDIM > > >
 
using CoarsenPatchStrategyMap = std::map< std::string, std::unique_ptr< SAMRAI::xfer::CoarsenPatchStrategy< NDIM > > >
 
using CoarsenScheduleMap = std::map< std::string, std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > >
 

Private Member Functions

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

Private Attributes

bool d_regridding_hierarchy = false
 
bool d_at_regrid_time_step = false
 
RefineAlgorithmMap d_ghostfill_algs
 
RefinePatchStrategyMap d_ghostfill_strategies
 
RefineScheduleMap d_ghostfill_scheds
 
RefineAlgorithmMap d_prolong_algs
 
RefinePatchStrategyMap d_prolong_strategies
 
RefineScheduleMap d_prolong_scheds
 
CoarsenAlgorithmMap d_coarsen_algs
 
CoarsenPatchStrategyMap d_coarsen_strategies
 
CoarsenScheduleMap d_coarsen_scheds
 

Friends

class IBStrategy
 
class IBEulerianForceFunction
 
class IBEulerianSourceFunction
 

Detailed Description

Options Controlling Regridding

Most IBAMR applications involve structure meshes defined on the finest level of the patch hierarchy. These structures move: in particular, they can potentially move off the finest grid level, causing the interaction routines to no longer work.

This class offers two different strategies for calculating how much the structure has moved (which it then uses to specify that a regrid is required):

  1. Estimation based on the fluid: The structure is assumed to move at the maximum velocity of the fluid at each time step. This class will regrid once a single fluid point has moved at least regrid_fluid_cfl_interval cell widths since the last regrid.

  2. Estimation based on the structure: The structure's displacement is calculated directly from its position vector and that value is used to determine the number of cell widths the structure has moved. In this case we will regrid once a point on the structure has moved at least regrid_structure_cfl_interval cell widths since the last regrid.

Both regrid_fluid_cfl_interval and regrid_structure_cfl_interval can be specified in the input database. For backwards compatibility the value regrid_cfl_interval is equivalent to regrid_fluid_cfl_interval. At the present time regrid_structure_cfl_interval is not implemented for all IBStrategy classes.

Alternatively, one can request that the solver (regardless of any computed displacement or velocity) regrid every time a fixed number of timesteps have occurred by specifying regrid_interval in the input database. If either regrid_structure_cfl_interval or regrid_fluid_cfl_interval are provided in the input database then regrid_interval is ignored.

Member Typedef Documentation

◆ PreprocessIntegrateHierarchyCallbackFcnPtr

using IBTK::HierarchyIntegrator::PreprocessIntegrateHierarchyCallbackFcnPtr = void (*)(double current_time, double new_time, int num_cycles, void* ctx)
inherited

Callback function specification to enable further specialization of preprocessIntegrateHierarchy().

◆ IntegrateHierarchyCallbackFcnPtr

using IBTK::HierarchyIntegrator::IntegrateHierarchyCallbackFcnPtr = void (*)(double current_time, double new_time, int cycle_num, void* ctx)
inherited

Callback function specification to enable further specialization of integrateHierarchy().

◆ PostprocessIntegrateHierarchyCallbackFcnPtr

using IBTK::HierarchyIntegrator::PostprocessIntegrateHierarchyCallbackFcnPtr = void (*)(double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles, void* ctx)
inherited

Callback function specification to enable further specialization of postprocessIntegrateHierarchy().

◆ ApplyGradientDetectorCallbackFcnPtr

using IBTK::HierarchyIntegrator::ApplyGradientDetectorCallbackFcnPtr = void (*)(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, void* ctx)
inherited

Callback function specification to enable further specialization of applyGradientDetector().

◆ RegridHierarchyCallbackFcnPtr

using IBTK::HierarchyIntegrator::RegridHierarchyCallbackFcnPtr = void (*)(SAMRAI::tbox::Pointer<SAMRAI::hier::BasePatchHierarchy<NDIM> > hierarchy, double data_time, bool initial_time, void* ctx)
inherited

Callback function specification to enable further specialization of regridHierarchy().

◆ RefineAlgorithmMap

◆ RefinePatchStrategyMap

using IBTK::HierarchyIntegrator::RefinePatchStrategyMap = std::map<std::string, std::unique_ptr<SAMRAI::xfer::RefinePatchStrategy<NDIM> > >
privateinherited

◆ RefineScheduleMap

using IBTK::HierarchyIntegrator::RefineScheduleMap = std::map<std::string, std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > > >
privateinherited

◆ CoarsenAlgorithmMap

◆ CoarsenPatchStrategyMap

using IBTK::HierarchyIntegrator::CoarsenPatchStrategyMap = std::map<std::string, std::unique_ptr<SAMRAI::xfer::CoarsenPatchStrategy<NDIM> > >
privateinherited

◆ CoarsenScheduleMap

using IBTK::HierarchyIntegrator::CoarsenScheduleMap = std::map<std::string, std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenSchedule<NDIM> > > >
privateinherited

Constructor & Destructor Documentation

◆ ~IBHierarchyIntegrator()

IBAMR::IBHierarchyIntegrator::~IBHierarchyIntegrator ( )
default

The destructor for class IBHierarchyIntegrator unregisters the integrator object with the restart manager when the object is so registered.

◆ IBHierarchyIntegrator() [1/3]

IBAMR::IBHierarchyIntegrator::IBHierarchyIntegrator ( const std::string &  object_name,
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database input_db,
SAMRAI::tbox::Pointer< IBStrategy ib_method_ops,
SAMRAI::tbox::Pointer< INSHierarchyIntegrator ins_hier_integrator,
bool  register_for_restart = true 
)
protected

The constructor for class IBHierarchyIntegrator sets some default values, reads in configuration information from input and restart databases, and registers the integrator object with the restart manager when requested.

◆ IBHierarchyIntegrator() [2/3]

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

◆ IBHierarchyIntegrator() [3/3]

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

Member Function Documentation

◆ getTimeSteppingType()

TimeSteppingType IBAMR::IBHierarchyIntegrator::getTimeSteppingType ( ) const

Return the time stepping scheme.

◆ getIBStrategy()

SAMRAI::tbox::Pointer<IBStrategy> IBAMR::IBHierarchyIntegrator::getIBStrategy ( ) const

Return a pointer to the IBStrategy object registered with this integrator.

◆ registerBodyForceFunction()

void IBAMR::IBHierarchyIntegrator::registerBodyForceFunction ( SAMRAI::tbox::Pointer< IBTK::CartGridFunction F_fcn)

Supply a body force (optional).

◆ registerLoadBalancer()

virtual void IBAMR::IBHierarchyIntegrator::registerLoadBalancer ( SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > >  load_balancer)
overridevirtual

Register a load balancer for non-uniform load balancing.

Reimplemented from IBTK::HierarchyIntegrator.

◆ getVelocityVariable()

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::getVelocityVariable ( ) const

Return a pointer to the fluid velocity variable.

◆ getPressureVariable()

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::getPressureVariable ( ) const

Return a pointer to the fluid pressure state variable.

◆ getBodyForceVariable()

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::getBodyForceVariable ( ) const

Return a pointer to the body force variable.

◆ getFluidSourceVariable()

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::getFluidSourceVariable ( ) const

Return a pointer to the source strength variable.

◆ getVelocityPhysBdryOp()

IBTK::RobinPhysBdryPatchStrategy* IBAMR::IBHierarchyIntegrator::getVelocityPhysBdryOp ( ) const

Return a pointer to the velocity physical boundary conditions

◆ preprocessIntegrateHierarchy()

void IBAMR::IBHierarchyIntegrator::preprocessIntegrateHierarchy ( double  current_time,
double  new_time,
int  num_cycles = 1 
)
overridevirtual

Basic functions to prepare to advance data from current_time to new_time.

A default implementation is provided that sets the current values of num_cycles and the time step size and checks to see if the time step size has changed.

Reimplemented from IBTK::HierarchyIntegrator.

Reimplemented in IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::IBInterpolantHierarchyIntegrator.

◆ postprocessIntegrateHierarchy()

void IBAMR::IBHierarchyIntegrator::postprocessIntegrateHierarchy ( double  current_time,
double  new_time,
bool  skip_synchronize_new_state_data,
int  num_cycles = 1 
)
overridevirtual

◆ initializeHierarchyIntegrator()

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

Initialize the variables, basic communications algorithms, solvers, and other data structures used by this time integrator object.

This method is called automatically by initializePatchHierarchy() prior to the construction of the patch hierarchy. It is also possible for users to make an explicit call to initializeHierarchyIntegrator() prior to calling initializePatchHierarchy().

Implements IBTK::HierarchyIntegrator.

Reimplemented in IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::IBInterpolantHierarchyIntegrator.

◆ initializePatchHierarchy()

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

Initialize the AMR patch hierarchy and data defined on the hierarchy at the start of a computation. If the computation is begun from a restart file, the patch hierarchy and patch data are read from the hierarchy database. Otherwise, the patch hierarchy and patch data are initialized by the gridding algorithm associated with the integrator object.

The implementation of this function assumes that the hierarchy exists upon entry to the function, but that it contains no patch levels. On return from this function, the state of the integrator object will be such that it is possible to step through time via the advanceHierarchy() function.

Reimplemented from IBTK::HierarchyIntegrator.

◆ computeFluidSources()

virtual void IBAMR::IBHierarchyIntegrator::computeFluidSources ( const int  data_idx,
const double  data_time 
)
protectedvirtual

Compute the values of the fluid sources which may be provided by d_ib_method_ops.

◆ regridHierarchyBeginSpecialized()

void IBAMR::IBHierarchyIntegrator::regridHierarchyBeginSpecialized ( )
overrideprotectedvirtual

Perform necessary data movement, workload estimation, and logging prior to regridding.

Reimplemented from IBTK::HierarchyIntegrator.

◆ regridHierarchyEndSpecialized()

void IBAMR::IBHierarchyIntegrator::regridHierarchyEndSpecialized ( )
overrideprotectedvirtual

Perform necessary data movement and logging after regridding.

Reimplemented from IBTK::HierarchyIntegrator.

◆ atRegridPointSpecialized()

bool IBAMR::IBHierarchyIntegrator::atRegridPointSpecialized ( ) const
overrideprotectedvirtual

Function to determine whether regridding should occur at the current time step.

Reimplemented from IBTK::HierarchyIntegrator.

◆ initializeLevelDataSpecialized()

void IBAMR::IBHierarchyIntegrator::initializeLevelDataSpecialized ( 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 
)
overrideprotectedvirtual

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

Reimplemented from IBTK::HierarchyIntegrator.

◆ resetHierarchyConfigurationSpecialized()

void IBAMR::IBHierarchyIntegrator::resetHierarchyConfigurationSpecialized ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  coarsest_level,
int  finest_level 
)
overrideprotectedvirtual

Reset cached hierarchy dependent data.

Reimplemented from IBTK::HierarchyIntegrator.

◆ applyGradientDetectorSpecialized()

void IBAMR::IBHierarchyIntegrator::applyGradientDetectorSpecialized ( 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 
)
overrideprotectedvirtual

Set integer tags to "one" in cells where refinement of the given level should occur according to the magnitude of the fluid vorticity.

Reimplemented from IBTK::HierarchyIntegrator.

◆ putToDatabaseSpecialized()

void IBAMR::IBHierarchyIntegrator::putToDatabaseSpecialized ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overrideprotectedvirtual

Write out specialized object state to the given database.

Reimplemented from IBTK::HierarchyIntegrator.

Reimplemented in IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::IBInterpolantHierarchyIntegrator.

◆ addWorkloadEstimate()

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

Add the work contributions (excluding the background grid) for the current hierarchy into the variable with index workload_data_idx. The only direct workload contribution of this hierarchy manager is usually the work done by the IBStrategy object.

Reimplemented from IBTK::HierarchyIntegrator.

◆ operator=()

IBHierarchyIntegrator& IBAMR::IBHierarchyIntegrator::operator= ( const IBHierarchyIntegrator 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::IBHierarchyIntegrator::getFromInput ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db,
bool  is_from_restart 
)
private

Read input values from a given database.

◆ getFromRestart()

void IBAMR::IBHierarchyIntegrator::getFromRestart ( )
private

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

◆ getName()

const std::string& IBTK::HierarchyIntegrator::getName ( ) const
inherited

Return the name of the hierarchy integrator object.

◆ advanceHierarchy()

virtual void IBTK::HierarchyIntegrator::advanceHierarchy ( double  dt)
virtualinherited

Integrate data on all patches on all levels of the patch hierarchy over the specified time increment.

◆ getMinimumTimeStepSize()

double IBTK::HierarchyIntegrator::getMinimumTimeStepSize ( )
inherited

Return the current value of the minimum time step size for the integrator object.

Subclasses can control the method used to determined the time step size by overriding the protected virtual member function getMinimumTimeStepSizeSpecialized().

◆ getMaximumTimeStepSize()

double IBTK::HierarchyIntegrator::getMaximumTimeStepSize ( )
inherited

Return the current value of the maximum time step size for the integrator object.

Subclasses can control the method used to determined the time step size by overriding the protected virtual member function getMaximumTimeStepSizeSpecialized().

◆ synchronizeHierarchyData()

void IBTK::HierarchyIntegrator::synchronizeHierarchyData ( VariableContextType  ctx_type)
inherited

Synchronize data defined on the grid hierarchy.

Subclasses can control the method used to synchronize data on the grid hierarchy by overriding the protected virtual member function synchronizeHierarchyDataSpecialized().

◆ resetTimeDependentHierarchyData()

void IBTK::HierarchyIntegrator::resetTimeDependentHierarchyData ( double  new_time)
inherited

Reset the current data to equal the new data, update the time level of the current data, and deallocate the scratch and new data.

Subclasses can control the method used to reset data on the grid hierarchy by overriding the protected virtual member function resetTimeDependentHierarchyDataSpecialized().

◆ resetIntegratorToPreadvanceState()

void IBTK::HierarchyIntegrator::resetIntegratorToPreadvanceState ( )
inherited

Reset the hierarchy integrator to the state at the beginning of the current time step.

Subclasses can control the method used to reset data on the grid hierarchy by overriding the protected virtual member function resetTimeDependentHierarchyDataSpecialized().

◆ regridHierarchy()

virtual void IBTK::HierarchyIntegrator::regridHierarchy ( )
virtualinherited

Virtual method to regrid the patch hierarchy.

A default implementation is provided that calls GriddingAlgorithm::regidAllFinerLevels() to regrid the patch hierarchy. Subclasses can control the method used to regrid the patch hierarchy by overriding this public virtual member function.

Before regridding, this method calls regridHierarchyBeginSpecialized() on the current integrator and all child integrators.

After regridding and (optionally) checking the new grid volume, this method calls regridHierarchyEndSpecialized() on the current integrator and all child integrators. It then calls the following methods (and, therefore, the specialized methods on the current and all child integrators) in the following order:

  1. initializeCompositeHierarchyData
  2. synchronizeHierarchyData
Warning
This class assumes, but does not enforce, that this method is only called on the parent integrator. A future release of IBAMR will enforce this assumption.

◆ atRegridPoint()

bool IBTK::HierarchyIntegrator::atRegridPoint ( ) const
inherited

Return a boolean value that indicates whether regridding should occur at the current time step.

Subclasses can control the method used to trigger adaptive regridding by overriding the protected virtual member function atRegridPointSpecialized().

◆ getIntegratorTime()

double IBTK::HierarchyIntegrator::getIntegratorTime ( ) const
inherited

Return the current integration time.

◆ getStartTime()

double IBTK::HierarchyIntegrator::getStartTime ( ) const
inherited

Return the initial integration time.

◆ getEndTime()

double IBTK::HierarchyIntegrator::getEndTime ( ) const
inherited

Return the final integration time.

◆ getIntegratorStep()

int IBTK::HierarchyIntegrator::getIntegratorStep ( ) const
inherited

Return the number of time steps taken by the integrator.

◆ getMaxIntegratorSteps()

int IBTK::HierarchyIntegrator::getMaxIntegratorSteps ( ) const
inherited

Return the maximum number of time steps allowed by the integrator.

◆ stepsRemaining()

bool IBTK::HierarchyIntegrator::stepsRemaining ( ) const
inherited

Return true if any time steps remain, false otherwise.

◆ updateWorkloadEstimates()

void IBTK::HierarchyIntegrator::updateWorkloadEstimates ( )
inherited

Update workload estimates on each level of the patch hierarchy. Does nothing unless a load balancer has previously been attached via HierarchyIntegrator::register_load_balancer. This function is usually called immediately before regridding so that, should a load balancer be registered with the current class, that load balancer can distribute patches in a way that gives each processor a roughly equal amount of work (instead of simply an equal number of cells).

If you want to visualize the workload then you will have to ensure that workload estimates are recomputed after regridding (by default they are not). This is because workload estimates are only used when moving data between processors: they are computed immediately before the domain is repartitioned and, therefore, their patch data is always invalid at points where output is written. The easiest way to get around this is to enable logging (which will result in workload estimates being updated after regridding) or to manually call updateWorkloadEstimates at points where output is written. The former is more efficient since most applications regrid less frequently than they write visualization files.

Once the data is available, it can be attached to the visit data writer in the usual way. Here is one way to do set this up at the same time the visit data writer is registered:

Pointer<VisItDataWriter<NDIM> > visit_data_writer = app_initializer->getVisItDataWriter();
if (uses_visit)
{
time_integrator->registerVisItDataWriter(visit_data_writer);
// This is only valid if a load balancer has previously been attached
// to the time integrator or, should it exist, its parent integrator;
// otherwise no workload estimates are computed
visit_data_writer->registerPlotQuantity("workload", "SCALAR",
time_integrator->getWorkloadDataIndex());
}
Note
This function computes workloads by setting the estimated work value on all cells to 1 and then calling addWorkloadEstimate on the current object and on each child hierarchy integrator.
The workload estimate itself is stored in the variable with index HierarchyIntegrator::d_workload_idx.

@seealso HierarchyIntegrator::getWorkloadDataIndex()

◆ getPatchHierarchy()

SAMRAI::tbox::Pointer<SAMRAI::hier::PatchHierarchy<NDIM> > IBTK::HierarchyIntegrator::getPatchHierarchy ( ) const
inherited

Return a pointer to the patch hierarchy managed by the integrator.

◆ getGriddingAlgorithm()

SAMRAI::tbox::Pointer<SAMRAI::mesh::GriddingAlgorithm<NDIM> > IBTK::HierarchyIntegrator::getGriddingAlgorithm ( ) const
inherited

Return a pointer to the gridding algorithm object managed by the integrator.

◆ getWorkloadDataIndex()

int IBTK::HierarchyIntegrator::getWorkloadDataIndex ( ) const
inherited

Return the workload variable's patch data index. If the workload variable has not yet been set up by this object (or, should it exist, this object's parent hierarchy integrator) then IBTK::invalid_index is returned instead.

Note
this is only implemented since this information is not readily available from the variable database and there is no other way to access this variable. The main use of this variable is for optionally adding the workload estimates to the VisIt data writer.

◆ registerVisItDataWriter()

void IBTK::HierarchyIntegrator::registerVisItDataWriter ( SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > >  visit_writer)
inherited

Register a VisIt data writer so the integrator can output data files that may be postprocessed with the VisIt visualization tool.

◆ getVisItDataWriter()

SAMRAI::tbox::Pointer<SAMRAI::appu::VisItDataWriter<NDIM> > IBTK::HierarchyIntegrator::getVisItDataWriter ( ) const
inherited

Get a pointer to the VisIt data writer registered with the solver.

◆ setupPlotData()

void IBTK::HierarchyIntegrator::setupPlotData ( )
inherited

Prepare variables for plotting.

Subclasses can control the method used to setup plot data by overriding the protected virtual member function setupPlotData().

Note
Subclasses are allowed to require that this function be called immediately before writing visualization data.

◆ getNumberOfCycles()

virtual int IBTK::HierarchyIntegrator::getNumberOfCycles ( ) const
virtualinherited

Routines to implement the time integration scheme.

Virtual method to return the number of cycles to perform for the present time step.

Reimplemented in IBAMR::INSHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, and IBAMR::IBImplicitStaggeredHierarchyIntegrator.

◆ getCurrentCycleNumber()

virtual int IBTK::HierarchyIntegrator::getCurrentCycleNumber ( ) const
virtualinherited

Virtual method to return the current cycle number within the present time step.

The default implementation returns a value of -1 when it is not advancing the hierarchy.

◆ getCurrentTimeStepSize()

virtual double IBTK::HierarchyIntegrator::getCurrentTimeStepSize ( ) const
virtualinherited

Virtual method to return the current time step size.

The default implementation returns the value numeric_limits<>::quiet_NaN() when it is not advancing the hierarchy.

◆ integrateHierarchy()

void IBTK::HierarchyIntegrator::integrateHierarchy ( double  current_time,
double  new_time,
int  cycle_num = 0 
)
inherited

Advance data from current_time to new_time. The current implementation calls doIntegrateHierarchy() followed by executeIntegrateHierarchyCallbackFcns().

◆ skipCycle()

void IBTK::HierarchyIntegrator::skipCycle ( double  current_time,
double  new_time,
int  cycle_num = 0 
)
inherited

Method to skip a cycle of the time integration scheme (e.g. for cases in which time integration is handled by another class).

◆ registerPreprocessIntegrateHierarchyCallback()

void IBTK::HierarchyIntegrator::registerPreprocessIntegrateHierarchyCallback ( PreprocessIntegrateHierarchyCallbackFcnPtr  callback,
void *  ctx = nullptr 
)
inherited

Register a callback function to enable further specialization of preprocessIntegrateHierarchy().

◆ registerIntegrateHierarchyCallback()

void IBTK::HierarchyIntegrator::registerIntegrateHierarchyCallback ( IntegrateHierarchyCallbackFcnPtr  callback,
void *  ctx = nullptr 
)
inherited

Register a callback function to enable further specialization of integrateHierarchy().

◆ registerPostprocessIntegrateHierarchyCallback()

void IBTK::HierarchyIntegrator::registerPostprocessIntegrateHierarchyCallback ( PostprocessIntegrateHierarchyCallbackFcnPtr  callback,
void *  ctx = nullptr 
)
inherited

Register a callback function to enable further specialization of postprocessIntegrateHierarchy().

◆ registerApplyGradientDetectorCallback()

void IBTK::HierarchyIntegrator::registerApplyGradientDetectorCallback ( ApplyGradientDetectorCallbackFcnPtr  callback,
void *  ctx = nullptr 
)
inherited

Register a callback function to enable further specialization of applyGradientDetector().

◆ registerRegridHierarchyCallback()

void IBTK::HierarchyIntegrator::registerRegridHierarchyCallback ( RegridHierarchyCallbackFcnPtr  ,
void *  ctx = nullptr 
)
inherited

Register a callback function to enable further specialization of regridHierarchy().

Note
This function will be called after all data owned by the hierarchy is initialized and set, but before the data is synchronized. If data owned by the integrator is used in this function, it should be synchronized prior to use by an appropriate coarsening operation.

◆ initializeCompositeHierarchyData()

void IBTK::HierarchyIntegrator::initializeCompositeHierarchyData ( double  init_data_time,
bool  initial_time 
)
inherited

Perform data initialization after the entire hierarchy has been constructed.

◆ initializeLevelData() [1/2]

void IBTK::HierarchyIntegrator::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 = SAMRAI::tbox::PointerSAMRAI::hier::BasePatchLevel< NDIM > >(nullptr),
bool  allocate_data = true 
)
overrideinherited

Implementations of functions declared in the SAMRAI::mesh::StandardTagAndInitStrategy abstract base class.

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

Note
Subclasses should not override the implementation of this function provided by class HierarchyIntegrator. Instead, they should override the protected virtual member function initializeLevelDataSpecialized().
See also
SAMRAI::mesh::StandardTagAndInitStrategy::initializeLevelData

◆ 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() [1/2]

void IBTK::HierarchyIntegrator::resetHierarchyConfiguration ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  coarsest_level,
int  finest_level 
)
overrideinherited

Reset cached hierarchy dependent data.

Note
Subclasses should not override the implementation of this function provided by class HierarchyIntegrator. Instead, they should override the protected virtual member function resetHierarchyConfigurationSpecialized().
See also
SAMRAI::mesh::StandardTagAndInitStrategy::resetHierarchyConfiguration

◆ 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() [1/2]

void IBTK::HierarchyIntegrator::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 
)
overrideinherited

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

Note
Subclasses should not override the implementation of this function provided by class HierarchyIntegrator. Instead, they should override the protected virtual member function applyGradientDetectorSpecialized().
See also
SAMRAI::mesh::StandardTagAndInitStrategy::applyGradientDetector

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

◆ getContext()

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::getContext ( VariableContextType  ctx_type) const
inherited

Routines to access to the variable contexts maintained by the integrator.

Return a pointer the variable context corresponding to the specified variable context type.

◆ getCurrentContext()

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::getCurrentContext ( ) const
inherited

Return a pointer to the "current" variable context used by integrator. Current data corresponds to state data at the beginning of a time step, or when a new level is initialized.

◆ getNewContext()

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::getNewContext ( ) const
inherited

Return a pointer to the "new" variable context used by integrator. New data corresponds to advanced state data at the end of a time step. The data is one time step later than the "current" data.

◆ getScratchContext()

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::getScratchContext ( ) const
inherited

Return a pointer to the "scratch" variable context used by integrator. Scratch data typically corresponds to storage that user-routines in the concrete GodunovAdvector object manipulate; in particular, scratch data contains ghost cells.

◆ getPlotContext()

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::getPlotContext ( ) const
inherited

Return a pointer to the "plot" variable context used by integrator. Plot data is only read from by the VisItDataWriter and is allocated by setupPlotData() and deallocated by regridHierarchyBegin().

◆ isAllocatedPatchData()

bool IBTK::HierarchyIntegrator::isAllocatedPatchData ( int  data_idx,
int  coarsest_ln = invalid_level_number,
int  finest_ln = invalid_level_number 
) const
inherited

Check whether a patch data index corresponds to allocated data over the specified range of patch level numbers.

NOTE: This method will return "false" without error for invalid (i.e., negative) patch data indices.

◆ allocatePatchData()

void IBTK::HierarchyIntegrator::allocatePatchData ( int  data_idx,
double  data_time,
int  coarsest_ln = invalid_level_number,
int  finest_ln = invalid_level_number 
) const
inherited

Allocate a patch data index over the specified range of patch level numbers.

◆ deallocatePatchData()

void IBTK::HierarchyIntegrator::deallocatePatchData ( int  data_idx,
int  coarsest_ln = invalid_level_number,
int  finest_ln = invalid_level_number 
) const
inherited

Deallocate a patch data index over the specified range of patch level numbers.

◆ getHierarchyMathOps()

SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::HierarchyIntegrator::getHierarchyMathOps ( ) const
inherited

Routines to access utility classeses managed by the integrator.

◆ registerVariable() [1/2]

void IBTK::HierarchyIntegrator::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< CartGridFunction init_fcn = nullptr,
const bool  register_for_restart = true 
)
inherited

Routines to register new variables with the integrator.

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 IBTK::HierarchyIntegrator::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 
)
inherited

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.

◆ putToDatabase()

void IBTK::HierarchyIntegrator::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overridevirtualinherited

Implementations of functions declared in the SAMRAI::tbox::Serializable abstract base class.

Write out object state to the given database.

Note
Subclasses should not override the implementation of this function provided by class HierarchyIntegrator. Instead, they should override the protected virtual member function putToDatabaseSpecialized().

Implements SAMRAI::tbox::Serializable.

◆ integrateHierarchySpecialized()

virtual void IBTK::HierarchyIntegrator::integrateHierarchySpecialized ( double  current_time,
double  new_time,
int  cycle_num = 0 
)
protectedpure virtualinherited

Pure virtual function that integrates the hierarchy from current_time to new_time.

Implementations of this virtual function are not required to synchronize data on the patch hierarchy. Data synchronization may be done (optionally) in a specialization of the public virtual member function postprocessIntegrateHierarchy().

This function is called by integrateHierarchy() prior to the execution of callbacks.

Note
Inheriting classes should call their base class versions of this method.

Implemented in IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::BrinkmanAdvDiffSemiImplicitHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::IBExplicitHierarchyIntegrator, IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::IBInterpolantHierarchyIntegrator.

◆ getMinimumTimeStepSizeSpecialized()

virtual double IBTK::HierarchyIntegrator::getMinimumTimeStepSizeSpecialized ( )
protectedvirtualinherited

Virtual method to compute an implementation-specific minimum stable time step size. Implementations should ensure that the returned time step is consistent across all processors.

A default implementation is provided that returns min(dt_max,dt_growth_factor*dt_current). The growth condition prevents excessive changes in the time step size as the computation progresses.

◆ getMaximumTimeStepSizeSpecialized()

virtual double IBTK::HierarchyIntegrator::getMaximumTimeStepSizeSpecialized ( )
protectedvirtualinherited

Virtual method to compute an implementation-specific maximum stable time step size. Implementations should ensure that the returned time step is consistent across all processors.

A default implementation is provided that returns min(dt_max,dt_growth_factor*dt_current). The growth condition prevents excessive changes in the time step size as the computation progresses.

Reimplemented in IBAMR::AdvDiffHierarchyIntegrator, IBAMR::INSHierarchyIntegrator, and IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.

◆ synchronizeHierarchyDataSpecialized()

virtual void IBTK::HierarchyIntegrator::synchronizeHierarchyDataSpecialized ( VariableContextType  ctx_type)
protectedvirtualinherited

Virtual method to perform implementation-specific data synchronization.

A default implementation is provided that synchronizes state data registered with the HierarchyIntegrator object using the coarsen operations specified by calls to registerVariable().

◆ resetTimeDependentHierarchyDataSpecialized()

virtual void IBTK::HierarchyIntegrator::resetTimeDependentHierarchyDataSpecialized ( double  new_time)
protectedvirtualinherited

Virtual method to perform implementation-specific data reset operations.

A default implementation is provided that first swaps the current and new PatchData pointers, and then deallocates the new and scratch data contexts.

Reimplemented in IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.

◆ resetIntegratorToPreadvanceStateSpecialized()

virtual void IBTK::HierarchyIntegrator::resetIntegratorToPreadvanceStateSpecialized ( )
protectedvirtualinherited

Virtual method to perform implementation-specific data reset operations.

A default implementation is provided that deallocates the new and scratch data contexts when data associated with those contexts have been allocated.

Reimplemented in IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.

◆ setupPlotDataSpecialized()

virtual void IBTK::HierarchyIntegrator::setupPlotDataSpecialized ( )
protectedvirtualinherited

◆ initializeCompositeHierarchyDataSpecialized()

virtual void IBTK::HierarchyIntegrator::initializeCompositeHierarchyDataSpecialized ( double  init_data_time,
bool  initial_time 
)
protectedvirtualinherited

Virtual method to perform implementation-specific data initialization after the entire hierarchy has been constructed.

An empty default implementation is provided.

Reimplemented in IBAMR::AdvDiffHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, and IBAMR::INSCollocatedHierarchyIntegrator.

◆ executePreprocessIntegrateHierarchyCallbackFcns()

virtual void IBTK::HierarchyIntegrator::executePreprocessIntegrateHierarchyCallbackFcns ( double  current_time,
double  new_time,
int  num_cycles 
)
protectedvirtualinherited

Execute any user-specified preprocessIntegrateHierarchy callback functions.

◆ executeIntegrateHierarchyCallbackFcns()

virtual void IBTK::HierarchyIntegrator::executeIntegrateHierarchyCallbackFcns ( double  current_time,
double  new_time,
int  cycle_num 
)
protectedvirtualinherited

Execute any user-specified integrateHierarchy callback functions.

◆ executePostprocessIntegrateHierarchyCallbackFcns()

virtual void IBTK::HierarchyIntegrator::executePostprocessIntegrateHierarchyCallbackFcns ( double  current_time,
double  new_time,
bool  skip_synchronize_new_state_data,
int  num_cycles 
)
protectedvirtualinherited

Execute any user-specified postprocessIntegrateHierarchy callback functions.

◆ executeApplyGradientDetectorCallbackFcns()

virtual void IBTK::HierarchyIntegrator::executeApplyGradientDetectorCallbackFcns ( 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 
)
protectedvirtualinherited

Execute any user-specified applyGradientDetector callback functions.

◆ executeRegridHierarchyCallbackFcns()

virtual void IBTK::HierarchyIntegrator::executeRegridHierarchyCallbackFcns ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
double  data_time,
bool  initial_time 
)
protectedvirtualinherited

Execute any user-specified regridHierarchy callback functions.

◆ registerGhostfillRefineAlgorithm()

void IBTK::HierarchyIntegrator::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 IBTK::HierarchyIntegrator::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 IBTK::HierarchyIntegrator::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> > IBTK::HierarchyIntegrator::getGhostfillRefineAlgorithm ( const std::string &  name) const
protectedinherited

Get ghost cell-filling refine algorithm.

◆ getProlongRefineAlgorithm()

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

Get data-prolonging refine algorithm.

◆ getCoarsenAlgorithm()

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

Get coarsen algorithm.

◆ getGhostfillRefineSchedules()

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

Get ghost cell-filling refine schedules.

◆ getProlongRefineSchedules()

const std::vector<SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineSchedule<NDIM> > >& IBTK::HierarchyIntegrator::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> > >& IBTK::HierarchyIntegrator::getCoarsenSchedules ( const std::string &  name) const
protectedinherited

Get coarsen schedules.

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

◆ registerChildHierarchyIntegrator()

void IBTK::HierarchyIntegrator::registerChildHierarchyIntegrator ( HierarchyIntegrator child_integrator)
protectedinherited

Register a "child" integrator object with this integrator object.

Note
Multiple child integrator objects may be registered with a single parent integrator object.

◆ registerParentHierarchyIntegrator()

void IBTK::HierarchyIntegrator::registerParentHierarchyIntegrator ( HierarchyIntegrator parent_integrator)
protectedinherited

Register a "parent" integrator object with this integrator object.

Note
Only a single parent integrator object may be registered with a particular child integrator object.

◆ buildHierarchyMathOps()

SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::HierarchyIntegrator::buildHierarchyMathOps ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy)
protectedinherited

Build the HierarchyMathOps object.

◆ setupTagBuffer()

void IBTK::HierarchyIntegrator::setupTagBuffer ( SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg)
protectedinherited

Setup the tag buffer.

◆ regriddingHierarchy()

bool IBTK::HierarchyIntegrator::regriddingHierarchy ( ) const
inlineprotectedinherited

Returns true when we are regridding the patch hierarchy.

◆ atRegridTimeStep()

bool IBTK::HierarchyIntegrator::atRegridTimeStep ( ) const
inlineprotectedinherited

Returns true when we are executing a time step in which a regridding operation was performed.

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

Friends And Related Function Documentation

◆ IBStrategy

friend class IBStrategy
friend

◆ IBEulerianForceFunction

friend class IBEulerianForceFunction
friend

◆ IBEulerianSourceFunction

friend class IBEulerianSourceFunction
friend

Member Data Documentation

◆ d_integrator_is_initialized

bool IBAMR::IBHierarchyIntegrator::d_integrator_is_initialized = false
protected

◆ d_time_stepping_type

TimeSteppingType IBAMR::IBHierarchyIntegrator::d_time_stepping_type = MIDPOINT_RULE
protected

Enum indicating the time integration employed for the IB equations.

◆ d_use_structure_predictor

bool IBAMR::IBHierarchyIntegrator::d_use_structure_predictor
protected

Flag indicating whether to use an explicit predictor for the structure configuration in the time stepping scheme.

◆ d_error_on_dt_change

bool IBAMR::IBHierarchyIntegrator::d_error_on_dt_change = true
protected

Flags to determine whether warnings or error messages should be emitted when time step size changes are encountered.

◆ d_warn_on_dt_change

bool IBAMR::IBHierarchyIntegrator::d_warn_on_dt_change = false
protected

◆ d_ins_hier_integrator

SAMRAI::tbox::Pointer<INSHierarchyIntegrator> IBAMR::IBHierarchyIntegrator::d_ins_hier_integrator
protected

◆ d_regrid_fluid_cfl_interval

double IBAMR::IBHierarchyIntegrator::d_regrid_fluid_cfl_interval = -1.0
protected

The regrid CFL interval indicates the number of meshwidths a particle may move in any coordinate direction between invocations of the regridding process.

Note
Currently, when the CFL-based regrid interval is specified, it is always used instead of the fixed-step regrid interval.

◆ d_regrid_structure_cfl_interval

double IBAMR::IBHierarchyIntegrator::d_regrid_structure_cfl_interval = -1.0
protected

◆ d_regrid_fluid_cfl_estimate

double IBAMR::IBHierarchyIntegrator::d_regrid_fluid_cfl_estimate = 0.0
protected

Estimation on the maximum fraction of fluid cells the structure has moved based on the maximum fluid velocity.

◆ d_regrid_structure_cfl_estimate

double IBAMR::IBHierarchyIntegrator::d_regrid_structure_cfl_estimate = 0.0
protected

Estimation on the maximum fraction of fluid cells the structure has moved based on the infinity norm of the structure's displacement.

◆ d_ib_method_ops

SAMRAI::tbox::Pointer<IBStrategy> IBAMR::IBHierarchyIntegrator::d_ib_method_ops
protected

◆ d_hier_velocity_data_ops

SAMRAI::tbox::Pointer<SAMRAI::math::HierarchyDataOpsReal<NDIM, double> > IBAMR::IBHierarchyIntegrator::d_hier_velocity_data_ops
protected

◆ d_hier_pressure_data_ops

SAMRAI::tbox::Pointer<SAMRAI::math::HierarchyDataOpsReal<NDIM, double> > IBAMR::IBHierarchyIntegrator::d_hier_pressure_data_ops
protected

◆ d_hier_cc_data_ops

SAMRAI::tbox::Pointer<SAMRAI::math::HierarchyCellDataOpsReal<NDIM, double> > IBAMR::IBHierarchyIntegrator::d_hier_cc_data_ops
protected

◆ d_u_var

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::d_u_var
protected

◆ d_p_var

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::d_p_var
protected

◆ d_f_var

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::d_f_var
protected

◆ d_q_var

SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > IBAMR::IBHierarchyIntegrator::d_q_var
protected

◆ d_u_idx

int IBAMR::IBHierarchyIntegrator::d_u_idx = IBTK::invalid_index
protected

◆ d_p_idx

int IBAMR::IBHierarchyIntegrator::d_p_idx = IBTK::invalid_index
protected

◆ d_f_idx

int IBAMR::IBHierarchyIntegrator::d_f_idx = IBTK::invalid_index
protected

◆ d_f_current_idx

int IBAMR::IBHierarchyIntegrator::d_f_current_idx = IBTK::invalid_index
protected

◆ d_q_idx

int IBAMR::IBHierarchyIntegrator::d_q_idx = IBTK::invalid_index
protected

◆ d_ib_context

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBAMR::IBHierarchyIntegrator::d_ib_context
protected

Context containing all patch data indices relevant to IB operations.

◆ d_ib_data

SAMRAI::hier::ComponentSelector IBAMR::IBHierarchyIntegrator::d_ib_data
protected

ComponentSelector corresponding to d_ib_context. Also contains patch data indices for relevant cloned indices (which, as they are clones, cannot be placed in the Context).

◆ d_u_phys_bdry_op

IBTK::RobinPhysBdryPatchStrategy* IBAMR::IBHierarchyIntegrator::d_u_phys_bdry_op
protected

◆ d_p_phys_bdry_op

IBTK::RobinPhysBdryPatchStrategy * IBAMR::IBHierarchyIntegrator::d_p_phys_bdry_op
protected

◆ d_u_ghostfill_alg

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineAlgorithm<NDIM> > IBAMR::IBHierarchyIntegrator::d_u_ghostfill_alg
protected

◆ d_f_prolong_alg

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineAlgorithm<NDIM> > IBAMR::IBHierarchyIntegrator::d_f_prolong_alg
protected

◆ d_p_ghostfill_alg

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineAlgorithm<NDIM> > IBAMR::IBHierarchyIntegrator::d_p_ghostfill_alg
protected

◆ d_q_prolong_alg

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineAlgorithm<NDIM> > IBAMR::IBHierarchyIntegrator::d_q_prolong_alg
protected

◆ d_u_ghostfill_op

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineOperator<NDIM> > IBAMR::IBHierarchyIntegrator::d_u_ghostfill_op
protected

◆ d_f_prolong_op

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineOperator<NDIM> > IBAMR::IBHierarchyIntegrator::d_f_prolong_op
protected

◆ d_p_ghostfill_op

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineOperator<NDIM> > IBAMR::IBHierarchyIntegrator::d_p_ghostfill_op
protected

◆ d_q_prolong_op

SAMRAI::tbox::Pointer<SAMRAI::xfer::RefineOperator<NDIM> > IBAMR::IBHierarchyIntegrator::d_q_prolong_op
protected

◆ d_u_coarsen_alg

SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenAlgorithm<NDIM> > IBAMR::IBHierarchyIntegrator::d_u_coarsen_alg
protected

◆ d_p_coarsen_alg

SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenAlgorithm<NDIM> > IBAMR::IBHierarchyIntegrator::d_p_coarsen_alg
protected

◆ d_u_coarsen_op

SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenOperator<NDIM> > IBAMR::IBHierarchyIntegrator::d_u_coarsen_op
protected

◆ d_p_coarsen_op

SAMRAI::tbox::Pointer<SAMRAI::xfer::CoarsenOperator<NDIM> > IBAMR::IBHierarchyIntegrator::d_p_coarsen_op
protected

◆ d_body_force_fcn

SAMRAI::tbox::Pointer<IBTK::CartGridFunction> IBAMR::IBHierarchyIntegrator::d_body_force_fcn
protected

◆ d_object_name

std::string IBTK::HierarchyIntegrator::d_object_name
protectedinherited

◆ d_registered_for_restart

bool IBTK::HierarchyIntegrator::d_registered_for_restart = false
protectedinherited

◆ d_hierarchy

SAMRAI::tbox::Pointer<SAMRAI::hier::PatchHierarchy<NDIM> > IBTK::HierarchyIntegrator::d_hierarchy
protectedinherited

◆ d_gridding_alg

SAMRAI::tbox::Pointer<SAMRAI::mesh::GriddingAlgorithm<NDIM> > IBTK::HierarchyIntegrator::d_gridding_alg
protectedinherited

◆ d_load_balancer

SAMRAI::tbox::Pointer<SAMRAI::mesh::LoadBalancer<NDIM> > IBTK::HierarchyIntegrator::d_load_balancer
protectedinherited

◆ d_workload_var

SAMRAI::tbox::Pointer<SAMRAI::pdat::CellVariable<NDIM, double> > IBTK::HierarchyIntegrator::d_workload_var
protectedinherited

◆ d_workload_idx

int IBTK::HierarchyIntegrator::d_workload_idx = IBTK::invalid_index
protectedinherited

The index of the workload estimate variable. If the current integrator is a child integrator then this variable index may not be set to the correct variable index since the parent is assumed to manage the variable, and the correct index is always provided when calling HierarchyIntegrator::addWorkloadEstimate() from parent integrators.

If necessary, this variable can be retrieved from the variable database under the name object_name + "::workload", where object_name is the name of the parent hierarchy integrator.

◆ d_hierarchy_is_initialized

bool IBTK::HierarchyIntegrator::d_hierarchy_is_initialized = false
protectedinherited

◆ d_may_need_to_reset_hierarchy_configuration

bool IBTK::HierarchyIntegrator::d_may_need_to_reset_hierarchy_configuration = false
protectedinherited

◆ d_parent_integrator

HierarchyIntegrator* IBTK::HierarchyIntegrator::d_parent_integrator = nullptr
protectedinherited

◆ d_child_integrators

std::set<HierarchyIntegrator*> IBTK::HierarchyIntegrator::d_child_integrators
protectedinherited

◆ d_visit_writer

SAMRAI::tbox::Pointer<SAMRAI::appu::VisItDataWriter<NDIM> > IBTK::HierarchyIntegrator::d_visit_writer
protectedinherited

◆ d_integrator_time

double IBTK::HierarchyIntegrator::d_integrator_time = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_start_time

double IBTK::HierarchyIntegrator::d_start_time = 0.0
protectedinherited

◆ d_end_time

double IBTK::HierarchyIntegrator::d_end_time = std::numeric_limits<double>::max()
protectedinherited

◆ d_dt_init

double IBTK::HierarchyIntegrator::d_dt_init = std::numeric_limits<double>::max()
protectedinherited

◆ d_dt_min

double IBTK::HierarchyIntegrator::d_dt_min = 0.0
protectedinherited

◆ d_dt_max

double IBTK::HierarchyIntegrator::d_dt_max = std::numeric_limits<double>::max()
protectedinherited

◆ d_dt_growth_factor

double IBTK::HierarchyIntegrator::d_dt_growth_factor = 2.0
protectedinherited

◆ d_integrator_step

int IBTK::HierarchyIntegrator::d_integrator_step = 0
protectedinherited

◆ d_max_integrator_steps

int IBTK::HierarchyIntegrator::d_max_integrator_steps = std::numeric_limits<int>::max()
protectedinherited

◆ d_dt_previous

std::deque<double> IBTK::HierarchyIntegrator::d_dt_previous
protectedinherited

◆ d_num_cycles

int IBTK::HierarchyIntegrator::d_num_cycles = 1
protectedinherited

◆ d_current_num_cycles

int IBTK::HierarchyIntegrator::d_current_num_cycles = -1
protectedinherited

◆ d_current_cycle_num

int IBTK::HierarchyIntegrator::d_current_cycle_num = -1
protectedinherited

◆ d_current_dt

double IBTK::HierarchyIntegrator::d_current_dt = std::numeric_limits<double>::quiet_NaN()
protectedinherited

◆ d_regrid_interval

int IBTK::HierarchyIntegrator::d_regrid_interval = 1
protectedinherited

◆ d_regrid_mode

RegridMode IBTK::HierarchyIntegrator::d_regrid_mode = STANDARD
protectedinherited

◆ d_enable_logging

bool IBTK::HierarchyIntegrator::d_enable_logging = false
protectedinherited

◆ d_enable_logging_solver_iterations

bool IBTK::HierarchyIntegrator::d_enable_logging_solver_iterations = false
protectedinherited

◆ d_bdry_extrap_type

std::string IBTK::HierarchyIntegrator::d_bdry_extrap_type = "LINEAR"
protectedinherited

◆ d_tag_buffer

SAMRAI::tbox::Array<int> IBTK::HierarchyIntegrator::d_tag_buffer = { 0 }
protectedinherited

◆ d_hier_math_ops

SAMRAI::tbox::Pointer<HierarchyMathOps> IBTK::HierarchyIntegrator::d_hier_math_ops
protectedinherited

◆ d_manage_hier_math_ops

bool IBTK::HierarchyIntegrator::d_manage_hier_math_ops = true
protectedinherited

◆ d_state_variables

std::list<SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > > IBTK::HierarchyIntegrator::d_state_variables
protectedinherited

◆ d_scratch_variables

std::list<SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > > IBTK::HierarchyIntegrator::d_scratch_variables
protectedinherited

◆ d_copy_scratch_to_current_fast

std::list<SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > > IBTK::HierarchyIntegrator::d_copy_scratch_to_current_fast
protectedinherited

◆ d_copy_scratch_to_current_slow

std::list<SAMRAI::tbox::Pointer<SAMRAI::hier::Variable<NDIM> > > IBTK::HierarchyIntegrator::d_copy_scratch_to_current_slow
protectedinherited

◆ d_current_data

SAMRAI::hier::ComponentSelector IBTK::HierarchyIntegrator::d_current_data
protectedinherited

◆ d_new_data

SAMRAI::hier::ComponentSelector IBTK::HierarchyIntegrator::d_new_data
protectedinherited

◆ d_scratch_data

SAMRAI::hier::ComponentSelector IBTK::HierarchyIntegrator::d_scratch_data
protectedinherited

◆ d_plot_data

SAMRAI::hier::ComponentSelector IBTK::HierarchyIntegrator::d_plot_data
protectedinherited

◆ d_state_var_init_fcns

std::map<SAMRAI::hier::Variable<NDIM>*, SAMRAI::tbox::Pointer<CartGridFunction> > IBTK::HierarchyIntegrator::d_state_var_init_fcns
protectedinherited

◆ d_current_context

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::d_current_context
protectedinherited

Variable contexts.

◆ d_new_context

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::d_new_context
protectedinherited

◆ d_scratch_context

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::d_scratch_context
protectedinherited

◆ d_plot_context

SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::d_plot_context
protectedinherited

◆ SYNCH_CURRENT_DATA_ALG

const std::string IBTK::HierarchyIntegrator::SYNCH_CURRENT_DATA_ALG
staticprotectedinherited

Names of special coarsen algorithms/schedules.

◆ SYNCH_NEW_DATA_ALG

const std::string IBTK::HierarchyIntegrator::SYNCH_NEW_DATA_ALG
staticprotectedinherited

◆ d_fill_after_regrid_bc_idxs

SAMRAI::hier::ComponentSelector IBTK::HierarchyIntegrator::d_fill_after_regrid_bc_idxs
protectedinherited

Regridding-related communications algorithms and other data structures.

◆ d_fill_after_regrid_prolong_alg

SAMRAI::xfer::RefineAlgorithm<NDIM> IBTK::HierarchyIntegrator::d_fill_after_regrid_prolong_alg
protectedinherited

◆ d_fill_after_regrid_phys_bdry_bc_op

std::unique_ptr<SAMRAI::xfer::RefinePatchStrategy<NDIM> > IBTK::HierarchyIntegrator::d_fill_after_regrid_phys_bdry_bc_op
protectedinherited

◆ d_preprocess_integrate_hierarchy_callbacks

std::vector<PreprocessIntegrateHierarchyCallbackFcnPtr> IBTK::HierarchyIntegrator::d_preprocess_integrate_hierarchy_callbacks
protectedinherited

Callback functions and callback function contexts.

◆ d_preprocess_integrate_hierarchy_callback_ctxs

std::vector<void*> IBTK::HierarchyIntegrator::d_preprocess_integrate_hierarchy_callback_ctxs
protectedinherited

◆ d_integrate_hierarchy_callbacks

std::vector<IntegrateHierarchyCallbackFcnPtr> IBTK::HierarchyIntegrator::d_integrate_hierarchy_callbacks
protectedinherited

◆ d_integrate_hierarchy_callback_ctxs

std::vector<void*> IBTK::HierarchyIntegrator::d_integrate_hierarchy_callback_ctxs
protectedinherited

◆ d_postprocess_integrate_hierarchy_callbacks

std::vector<PostprocessIntegrateHierarchyCallbackFcnPtr> IBTK::HierarchyIntegrator::d_postprocess_integrate_hierarchy_callbacks
protectedinherited

◆ d_postprocess_integrate_hierarchy_callback_ctxs

std::vector<void*> IBTK::HierarchyIntegrator::d_postprocess_integrate_hierarchy_callback_ctxs
protectedinherited

◆ d_apply_gradient_detector_callbacks

std::vector<ApplyGradientDetectorCallbackFcnPtr> IBTK::HierarchyIntegrator::d_apply_gradient_detector_callbacks
protectedinherited

◆ d_apply_gradient_detector_callback_ctxs

std::vector<void*> IBTK::HierarchyIntegrator::d_apply_gradient_detector_callback_ctxs
protectedinherited

◆ d_regrid_hierarchy_callbacks

std::vector<RegridHierarchyCallbackFcnPtr> IBTK::HierarchyIntegrator::d_regrid_hierarchy_callbacks
protectedinherited

◆ d_regrid_hierarchy_callback_ctxs

std::vector<void*> IBTK::HierarchyIntegrator::d_regrid_hierarchy_callback_ctxs
protectedinherited

◆ d_regridding_hierarchy

bool IBTK::HierarchyIntegrator::d_regridding_hierarchy = false
privateinherited

◆ d_at_regrid_time_step

bool IBTK::HierarchyIntegrator::d_at_regrid_time_step = false
privateinherited

◆ d_ghostfill_algs

RefineAlgorithmMap IBTK::HierarchyIntegrator::d_ghostfill_algs
privateinherited

◆ d_ghostfill_strategies

RefinePatchStrategyMap IBTK::HierarchyIntegrator::d_ghostfill_strategies
privateinherited

◆ d_ghostfill_scheds

RefineScheduleMap IBTK::HierarchyIntegrator::d_ghostfill_scheds
privateinherited

◆ d_prolong_algs

RefineAlgorithmMap IBTK::HierarchyIntegrator::d_prolong_algs
privateinherited

◆ d_prolong_strategies

RefinePatchStrategyMap IBTK::HierarchyIntegrator::d_prolong_strategies
privateinherited

◆ d_prolong_scheds

RefineScheduleMap IBTK::HierarchyIntegrator::d_prolong_scheds
privateinherited

◆ d_coarsen_algs

CoarsenAlgorithmMap IBTK::HierarchyIntegrator::d_coarsen_algs
privateinherited

◆ d_coarsen_strategies

CoarsenPatchStrategyMap IBTK::HierarchyIntegrator::d_coarsen_strategies
privateinherited

◆ d_coarsen_scheds

CoarsenScheduleMap IBTK::HierarchyIntegrator::d_coarsen_scheds
privateinherited

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