IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
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 </home/runner/work/IBAMR/IBAMR/include/ibamr/IBHierarchyIntegrator.h>
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... | |
Friends | |
class | IBStrategy |
class | IBEulerianForceFunction |
class | IBEulerianSourceFunction |
Additional Inherited Members | |
![]() | |
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) |
![]() | |
static const std::string | SYNCH_CURRENT_DATA_ALG = "SYNCH_CURRENT_DATA" |
static const std::string | SYNCH_NEW_DATA_ALG = "SYNCH_NEW_DATA" |
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.
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):
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.
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.
|
default |
The destructor for class IBHierarchyIntegrator unregisters the integrator object with the restart manager when the object is so registered.
|
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.
|
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.
|
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.
|
overrideprotectedvirtual |
Function to determine whether regridding should occur at the current time step.
Reimplemented from IBTK::HierarchyIntegrator.
Return a pointer to the body force variable.
Return a pointer to the source strength variable.
Pointer< IBStrategy > IBAMR::IBHierarchyIntegrator::getIBStrategy | ( | ) | const |
Return a pointer to the IBStrategy object registered with this integrator.
Return a pointer to the fluid pressure state variable.
TimeSteppingType IBAMR::IBHierarchyIntegrator::getTimeSteppingType | ( | ) | const |
Return the time stepping scheme.
IBTK::RobinPhysBdryPatchStrategy * IBAMR::IBHierarchyIntegrator::getVelocityPhysBdryOp | ( | ) | const |
Return a pointer to the velocity physical boundary conditions
Return a pointer to the fluid velocity variable.
|
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.
|
overrideprotectedvirtual |
Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm.
Reimplemented from IBTK::HierarchyIntegrator.
|
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.
|
overridevirtual |
Clean up data following call(s) to integrateHierarchy().
Reimplemented from IBTK::HierarchyIntegrator.
Reimplemented in IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::IBInterpolantHierarchyIntegrator.
|
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.
|
overrideprotectedvirtual |
Write out specialized object state to the given database.
Reimplemented from IBTK::HierarchyIntegrator.
Reimplemented in IBAMR::IBImplicitStaggeredHierarchyIntegrator, and IBAMR::IBInterpolantHierarchyIntegrator.
void IBAMR::IBHierarchyIntegrator::registerBodyForceFunction | ( | SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | F_fcn | ) |
Supply a body force (optional).
|
overridevirtual |
Register a load balancer for non-uniform load balancing.
Reimplemented from IBTK::HierarchyIntegrator.
|
overrideprotectedvirtual |
Perform necessary data movement, workload estimation, and logging prior to regridding.
Reimplemented from IBTK::HierarchyIntegrator.
|
overrideprotectedvirtual |
Perform necessary data movement and logging after regridding.
Reimplemented from IBTK::HierarchyIntegrator.
|
overrideprotectedvirtual |
Reset cached hierarchy dependent data.
Reimplemented from IBTK::HierarchyIntegrator.
|
protected |
Flags to determine whether warnings or error messages should be emitted when time step size changes are encountered.
|
protected |
Context containing all patch data indices relevant to IB operations.
|
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).
|
protected |
Estimation on the maximum fraction of fluid cells the structure has moved based on the maximum fluid velocity.
|
protected |
The regrid CFL interval indicates the number of meshwidths a particle may move in any coordinate direction between invocations of the regridding process.
|
protected |
Estimation on the maximum fraction of fluid cells the structure has moved based on the infinity norm of the structure's displacement.
|
protected |
Enum indicating the time integration employed for the IB equations.
|
protected |
Flag indicating whether to use an explicit predictor for the structure configuration in the time stepping scheme.