IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class HierarchyIntegrator provides an abstract interface for a time integrator for a system of equations defined on an AMR grid hierarchy, along with basic data management for variables defined on that hierarchy. More...
#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/HierarchyIntegrator.h>
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 | |
HierarchyIntegrator (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, bool register_for_restart) | |
~HierarchyIntegrator () | |
const std::string & | getName () const |
virtual void | initializeHierarchyIntegrator (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg)=0 |
virtual void | initializePatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) |
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 |
virtual void | registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer) |
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 |
virtual void | preprocessIntegrateHierarchy (double current_time, double new_time, int num_cycles=1) |
void | integrateHierarchy (double current_time, double new_time, int cycle_num=0) |
void | skipCycle (double current_time, double new_time, int cycle_num=0) |
virtual void | postprocessIntegrateHierarchy (double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles=1) |
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 |
void | resetHierarchyConfiguration (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_level, int finest_level) override |
void | applyGradientDetector (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool uses_richardson_extrapolation_too) override |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | getContext (VariableContextType ctx_type) const |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | getCurrentContext () const |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | getNewContext () const |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | getScratchContext () const |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | getPlotContext () 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 ¤t_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 |
Public Member Functions inherited from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM > | |
virtual double | getLevelDt (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double dt_time, const bool initial_time) |
virtual double | advanceLevel (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< NDIM > > hierarchy, const double current_time, const double new_time, const bool first_step, const bool last_step, const bool regrid_advance=false) |
virtual void | resetTimeDependentData (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double new_time, const bool can_be_refined) |
virtual void | resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level) |
virtual void | applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< NDIM > > level, const double error_data_time, const int tag_index, const double deltat, const int error_coarsen_ratio, const bool initial_time, const bool uses_gradient_detector_too) |
virtual void | coarsenDataForRichardsonExtrapolation (const tbox::Pointer< hier::PatchHierarchy< NDIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< NDIM > > coarser_level, const double coarsen_data_time, const bool before_advance) |
Protected Attributes | |
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 |
HierarchyIntegrator * | d_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< double > | d_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< int > | d_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::VariableContext > | d_current_context |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | d_new_context |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | d_scratch_context |
SAMRAI::tbox::Pointer< SAMRAI::hier::VariableContext > | d_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< PreprocessIntegrateHierarchyCallbackFcnPtr > | d_preprocess_integrate_hierarchy_callbacks |
std::vector< void * > | d_preprocess_integrate_hierarchy_callback_ctxs |
std::vector< IntegrateHierarchyCallbackFcnPtr > | d_integrate_hierarchy_callbacks |
std::vector< void * > | d_integrate_hierarchy_callback_ctxs |
std::vector< PostprocessIntegrateHierarchyCallbackFcnPtr > | d_postprocess_integrate_hierarchy_callbacks |
std::vector< void * > | d_postprocess_integrate_hierarchy_callback_ctxs |
std::vector< ApplyGradientDetectorCallbackFcnPtr > | d_apply_gradient_detector_callbacks |
std::vector< void * > | d_apply_gradient_detector_callback_ctxs |
std::vector< RegridHierarchyCallbackFcnPtr > | d_regrid_hierarchy_callbacks |
std::vector< void * > | d_regrid_hierarchy_callback_ctxs |
Static Protected Attributes | |
static const std::string | SYNCH_CURRENT_DATA_ALG = "SYNCH_CURRENT_DATA" |
static const std::string | SYNCH_NEW_DATA_ALG = "SYNCH_NEW_DATA" |
Class HierarchyIntegrator provides an abstract interface for a time integrator for a system of equations defined on an AMR grid hierarchy, along with basic data management for variables defined on that hierarchy.
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) |
Callback function specification to enable further specialization of applyGradientDetector().
using IBTK::HierarchyIntegrator::IntegrateHierarchyCallbackFcnPtr = void (*)(double current_time, double new_time, int cycle_num, void* ctx) |
Callback function specification to enable further specialization of integrateHierarchy().
using IBTK::HierarchyIntegrator::PostprocessIntegrateHierarchyCallbackFcnPtr = void (*)(double current_time, double new_time, bool skip_synchronize_new_state_data, int num_cycles, void* ctx) |
Callback function specification to enable further specialization of postprocessIntegrateHierarchy().
using IBTK::HierarchyIntegrator::PreprocessIntegrateHierarchyCallbackFcnPtr = void (*)(double current_time, double new_time, int num_cycles, void* ctx) |
Callback function specification to enable further specialization of preprocessIntegrateHierarchy().
using IBTK::HierarchyIntegrator::RegridHierarchyCallbackFcnPtr = void (*)(SAMRAI::tbox::Pointer<SAMRAI::hier::BasePatchHierarchy<NDIM> > hierarchy, double data_time, bool initial_time, void* ctx) |
Callback function specification to enable further specialization of regridHierarchy().
IBTK::HierarchyIntegrator::HierarchyIntegrator | ( | std::string | object_name, |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | input_db, | ||
bool | register_for_restart | ||
) |
The constructor for class HierarchyIntegrator sets some default values, reads in configuration information from input and restart databases, and registers the integrator object with the restart manager when requested.
IBTK::HierarchyIntegrator::~HierarchyIntegrator | ( | ) |
The destructor for class HierarchyIntegrator unregisters the integrator object with the restart manager when the object is so registered.
|
protectedvirtual |
Virtual method to provide implementation-specific workload estimate calculations. This method will be called on each registered child integrator. The intent of this function is that any HierarchyIntegrator object that manages data whose work varies from SAMRAI cell to SAMRAI cell will represent its additional workload by summing estimates into the workload_data_idx
variable.
Reimplemented in IBAMR::IBHierarchyIntegrator.
|
virtual |
Integrate data on all patches on all levels of the patch hierarchy over the specified time increment.
void IBTK::HierarchyIntegrator::allocatePatchData | ( | int | data_idx, |
double | data_time, | ||
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) | const |
Allocate a patch data index over the specified range of patch level numbers.
|
overridevirtual |
Set integer tags to "one" in cells where refinement of the given level should occur according to user-supplied feature detection criteria.
Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
|
protectedvirtual |
Virtual method to perform implementation-specific cell tagging operations.
An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, and IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.
bool IBTK::HierarchyIntegrator::atRegridPoint | ( | ) | const |
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().
|
protectedvirtual |
Virtual method to provide implementation-specific function to determine whether regridding should occur at the current time step.
A default implementation is provided that indicates that the hierarchy should be regridded at a fixed integer interval of time steps unless a parent integrator has been registered with this integrator. If a parent integrator has been registered with this integrator, atRegridPointSpecialized() returns false, in order to allow the parent integrator to control the timing of regridding.
Reimplemented in IBAMR::IBHierarchyIntegrator.
|
inlineprotected |
Returns true when we are executing a time step in which a regridding operation was performed.
|
protected |
Build the HierarchyMathOps object.
void IBTK::HierarchyIntegrator::deallocatePatchData | ( | int | data_idx, |
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) | const |
Deallocate a patch data index over the specified range of patch level numbers.
|
protectedvirtual |
Execute any user-specified applyGradientDetector callback functions.
|
protectedvirtual |
Execute any user-specified integrateHierarchy callback functions.
|
protectedvirtual |
Execute any user-specified postprocessIntegrateHierarchy callback functions.
|
protectedvirtual |
Execute any user-specified preprocessIntegrateHierarchy callback functions.
|
protectedvirtual |
Execute any user-specified regridHierarchy callback functions.
|
protected |
Get coarsen algorithm.
|
protected |
Get coarsen schedules.
SAMRAI::tbox::Pointer<SAMRAI::hier::VariableContext> IBTK::HierarchyIntegrator::getContext | ( | VariableContextType | ctx_type | ) | const |
Routines to access to the variable contexts maintained by the integrator.
Return a pointer the variable context corresponding to the specified variable context type.
Pointer< VariableContext > IBTK::HierarchyIntegrator::getCurrentContext | ( | ) | const |
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.
|
virtual |
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.
|
virtual |
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.
double IBTK::HierarchyIntegrator::getEndTime | ( | ) | const |
Return the final integration time.
|
protected |
Get ghost cell-filling refine algorithm.
|
protected |
Get ghost cell-filling refine schedules.
Pointer< GriddingAlgorithm< NDIM > > IBTK::HierarchyIntegrator::getGriddingAlgorithm | ( | ) | const |
Return a pointer to the gridding algorithm object managed by the integrator.
Pointer< HierarchyMathOps > IBTK::HierarchyIntegrator::getHierarchyMathOps | ( | ) | const |
Routines to access utility classeses managed by the integrator.
int IBTK::HierarchyIntegrator::getIntegratorStep | ( | ) | const |
Return the number of time steps taken by the integrator.
double IBTK::HierarchyIntegrator::getIntegratorTime | ( | ) | const |
Return the current integration time.
double IBTK::HierarchyIntegrator::getMaximumTimeStepSize | ( | ) |
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().
|
protectedvirtual |
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::INSHierarchyIntegrator, IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator, and IBAMR::AdvDiffHierarchyIntegrator.
int IBTK::HierarchyIntegrator::getMaxIntegratorSteps | ( | ) | const |
Return the maximum number of time steps allowed by the integrator.
double IBTK::HierarchyIntegrator::getMinimumTimeStepSize | ( | ) |
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().
|
protectedvirtual |
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.
const std::string & IBTK::HierarchyIntegrator::getName | ( | ) | const |
Return the name of the hierarchy integrator object.
Pointer< VariableContext > IBTK::HierarchyIntegrator::getNewContext | ( | ) | const |
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.
|
virtual |
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::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::AdvDiffSemiImplicitHierarchyIntegrator.
Pointer< PatchHierarchy< NDIM > > IBTK::HierarchyIntegrator::getPatchHierarchy | ( | ) | const |
Return a pointer to the patch hierarchy managed by the integrator.
Pointer< VariableContext > IBTK::HierarchyIntegrator::getPlotContext | ( | ) | const |
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().
|
protected |
Get data-prolonging refine algorithm.
|
protected |
Get data-prolonging refine schedules.
Pointer< VariableContext > IBTK::HierarchyIntegrator::getScratchContext | ( | ) | const |
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.
double IBTK::HierarchyIntegrator::getStartTime | ( | ) | const |
Return the initial integration time.
Pointer< VisItDataWriter< NDIM > > IBTK::HierarchyIntegrator::getVisItDataWriter | ( | ) | const |
Get a pointer to the VisIt data writer registered with the solver.
int IBTK::HierarchyIntegrator::getWorkloadDataIndex | ( | ) | const |
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.
void IBTK::HierarchyIntegrator::initializeCompositeHierarchyData | ( | double | init_data_time, |
bool | initial_time | ||
) |
Perform data initialization after the entire hierarchy has been constructed.
|
protectedvirtual |
Virtual method to perform implementation-specific data initialization after the entire hierarchy has been constructed.
An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, and IBAMR::AdvDiffHierarchyIntegrator.
|
pure virtual |
Virtual method to initialize the variables, basic communications algorithms, solvers, and other data structures used by a concrete time integrator object.
This method is called automatically by initializePatchHierarchy() prior to the construction of the patch hierarchy.
Implemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBInterpolantHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, IBAMR::IBExplicitHierarchyIntegrator, IBAMR::BrinkmanAdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator, and IBAMR::AdvDiffHierarchyIntegrator.
|
overridevirtual |
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.
Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
|
protectedvirtual |
Virtual method to perform implementation-specific data initialization on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm.
An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, and IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.
|
virtual |
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 in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, and IBAMR::IBHierarchyIntegrator.
void IBTK::HierarchyIntegrator::integrateHierarchy | ( | double | current_time, |
double | new_time, | ||
int | cycle_num = 0 |
||
) |
Advance data from current_time to new_time. The current implementation calls doIntegrateHierarchy() followed by executeIntegrateHierarchyCallbackFcns().
|
protectedpure virtual |
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.
Implemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBInterpolantHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, IBAMR::IBExplicitHierarchyIntegrator, IBAMR::BrinkmanAdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, and IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.
bool IBTK::HierarchyIntegrator::isAllocatedPatchData | ( | int | data_idx, |
int | coarsest_ln = invalid_level_number , |
||
int | finest_ln = invalid_level_number |
||
) | const |
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.
|
virtual |
Virtual method to clean up data following call(s) to integrateHierarchy().
A default implementation is provided that resets the current values of num_cycles and the time step size.
Reimplemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBInterpolantHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, IBAMR::IBExplicitHierarchyIntegrator, IBAMR::BrinkmanAdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, and IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator.
|
virtual |
Virtual method 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.
Reimplemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBInterpolantHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, IBAMR::IBExplicitHierarchyIntegrator, IBAMR::BrinkmanAdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator, and IBAMR::AdvDiffHierarchyIntegrator.
|
overridevirtual |
Implementations of functions declared in the SAMRAI::tbox::Serializable abstract base class.
Write out object state to the given database.
Implements SAMRAI::tbox::Serializable.
|
protectedvirtual |
Protecethod to write implementation-specific object state to a database.
An empty default implementation is provided.
Reimplemented in IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSHierarchyIntegrator, IBAMR::IBInterpolantHierarchyIntegrator, IBAMR::IBImplicitStaggeredHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, IBAMR::IBExplicitHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, and IBAMR::AdvDiffHierarchyIntegrator.
void IBTK::HierarchyIntegrator::registerApplyGradientDetectorCallback | ( | ApplyGradientDetectorCallbackFcnPtr | callback, |
void * | ctx = nullptr |
||
) |
Register a callback function to enable further specialization of applyGradientDetector().
|
protected |
Register a "child" integrator object with this integrator object.
|
protected |
Register a coarsen algorithm.
|
protected |
Register a ghost cell-filling refine algorithm.
void IBTK::HierarchyIntegrator::registerIntegrateHierarchyCallback | ( | IntegrateHierarchyCallbackFcnPtr | callback, |
void * | ctx = nullptr |
||
) |
Register a callback function to enable further specialization of integrateHierarchy().
|
virtual |
Register a load balancer for non-uniform load balancing.
load_balancer
is registered with another object; e.g., IBHierarchyIntegrator passes load_balancer
to its IBStrategy object. It will still usually be necessary to call HierarchyIntegrator::registerLoadBalancer at the beginning top of such overloaded functions.d_workload_idx
. All inheriting classes assume that this is the cell variable associated with workload estimates and will write their own data into these arrays. Reimplemented in IBAMR::IBHierarchyIntegrator.
|
protected |
Register a "parent" integrator object with this integrator object.
void IBTK::HierarchyIntegrator::registerPostprocessIntegrateHierarchyCallback | ( | PostprocessIntegrateHierarchyCallbackFcnPtr | callback, |
void * | ctx = nullptr |
||
) |
Register a callback function to enable further specialization of postprocessIntegrateHierarchy().
void IBTK::HierarchyIntegrator::registerPreprocessIntegrateHierarchyCallback | ( | PreprocessIntegrateHierarchyCallbackFcnPtr | callback, |
void * | ctx = nullptr |
||
) |
Register a callback function to enable further specialization of preprocessIntegrateHierarchy().
|
protected |
Register a data-prolonging refine algorithm.
void IBTK::HierarchyIntegrator::registerRegridHierarchyCallback | ( | RegridHierarchyCallbackFcnPtr | callback, |
void * | ctx = nullptr |
||
) |
Register a callback function to enable further specialization of regridHierarchy().
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 |
||
) |
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.
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::Pointer<SAMRAI::hier::VariableContext>(nullptr) , |
||
const bool | register_for_restart = true |
||
) |
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.
void IBTK::HierarchyIntegrator::registerVisItDataWriter | ( | SAMRAI::tbox::Pointer< SAMRAI::appu::VisItDataWriter< NDIM > > | visit_writer | ) |
Register a VisIt data writer so the integrator can output data files that may be postprocessed with the VisIt visualization tool.
|
inlineprotected |
Returns true when we are regridding the patch hierarchy.
|
virtual |
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:
|
protectedvirtual |
Perform any necessary work relevant to data owned by the current integrator prior to regridding (e.g., calculating divergences). An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, and IBAMR::IBExplicitHierarchyIntegrator.
|
protectedvirtual |
Perform any necessary work relevant to data owned by the current integrator after regridding (e.g., calculating divergences). An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, and IBAMR::IBExplicitHierarchyIntegrator.
|
overridevirtual |
Reset cached hierarchy dependent data.
Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
|
protectedvirtual |
Virtual method to perform implementation-specific data reset operations.
An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::IBHierarchyIntegrator, IBAMR::AdvDiffSemiImplicitHierarchyIntegrator, IBAMR::AdvDiffPredictorCorrectorHierarchyIntegrator, and IBAMR::AdvDiffHierarchyIntegrator.
void IBTK::HierarchyIntegrator::resetIntegratorToPreadvanceState | ( | ) |
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().
|
protectedvirtual |
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.
void IBTK::HierarchyIntegrator::resetTimeDependentHierarchyData | ( | double | new_time | ) |
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().
|
protectedvirtual |
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.
void IBTK::HierarchyIntegrator::setupPlotData | ( | ) |
Prepare variables for plotting.
Subclasses can control the method used to setup plot data by overriding the protected virtual member function setupPlotData().
|
protectedvirtual |
Virtual method to perform implementation-specific visualization setup.
An empty default implementation is provided.
Reimplemented in IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, and IBAMR::INSCollocatedHierarchyIntegrator.
|
protected |
Setup the tag buffer.
void IBTK::HierarchyIntegrator::skipCycle | ( | double | current_time, |
double | new_time, | ||
int | cycle_num = 0 |
||
) |
Method to skip a cycle of the time integration scheme (e.g. for cases in which time integration is handled by another class).
bool IBTK::HierarchyIntegrator::stepsRemaining | ( | ) | const |
Return true if any time steps remain, false otherwise.
void IBTK::HierarchyIntegrator::synchronizeHierarchyData | ( | VariableContextType | ctx_type | ) |
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().
|
protectedvirtual |
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().
void IBTK::HierarchyIntegrator::updateWorkloadEstimates | ( | ) |
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:
|
protected |
Variable contexts.
|
protected |
Regridding-related communications algorithms and other data structures.
|
protected |
Callback functions and callback function contexts.
|
protected |
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.
|
staticprotected |
Names of special coarsen algorithms/schedules.