IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class IBStrategy provides a generic interface for specifying the implementation details of a particular version of the IB method. More...
#include </home/runner/work/IBAMR/IBAMR/include/ibamr/IBStrategy.h>
Public Member Functions | |
IBStrategy ()=default | |
Constructor. | |
virtual void | registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver) |
virtual void | registerEulerianVariables () |
virtual void | registerEulerianCommunicationAlgorithms () |
virtual const SAMRAI::hier::IntVector< NDIM > & | getMinimumGhostCellWidth () const =0 |
virtual void | setupTagBuffer (SAMRAI::tbox::Array< int > &tag_buffer, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) const |
virtual void | inactivateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) |
virtual void | activateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) |
virtual bool | getLagrangianStructureIsActivated (int structure_number=0, int level_number=std::numeric_limits< int >::max()) const |
virtual double | getMaxPointDisplacement () const |
virtual void | preprocessIntegrateData (double current_time, double new_time, int num_cycles) |
virtual void | postprocessIntegrateData (double current_time, double new_time, int num_cycles) |
void | setUseFixedLEOperators (bool use_fixed_coupling_ops=true) |
virtual void | updateFixedLEOperators () |
virtual void | interpolateVelocity (int u_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &u_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &u_ghost_fill_scheds, double data_time)=0 |
virtual void | setUseMultistepTimeStepping (unsigned int n_previous_steps=1) |
virtual void | forwardEulerStep (double current_time, double new_time)=0 |
virtual void | backwardEulerStep (double current_time, double new_time) |
virtual void | midpointStep (double current_time, double new_time)=0 |
virtual void | trapezoidalStep (double current_time, double new_time)=0 |
virtual void | AB2Step (double current_time, double new_time) |
virtual void | computeLagrangianForce (double data_time)=0 |
virtual void | spreadForce (int f_data_idx, IBTK::RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds, double data_time)=0 |
virtual bool | hasFluidSources () const |
virtual void | computeLagrangianFluidSource (double data_time) |
virtual void | spreadFluidSource (int q_data_idx, IBTK::RobinPhysBdryPatchStrategy *q_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &q_prolongation_scheds, double data_time) |
virtual void | interpolatePressure (int p_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &p_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &p_ghost_fill_scheds, double data_time) |
virtual void | preprocessSolveFluidEquations (double current_time, double new_time, int cycle_num) |
virtual void | postprocessSolveFluidEquations (double current_time, double new_time, int cycle_num) |
virtual void | postprocessData () |
virtual void | initializePatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg, int u_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &u_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &u_ghost_fill_scheds, int integrator_step, double init_data_time, bool initial_time) |
virtual void | registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx) |
virtual void | addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx) |
virtual void | beginDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) |
virtual void | endDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) |
void | initializeLevelData (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double init_data_time, bool can_be_refined, bool initial_time, SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > > old_level, bool allocate_data) override |
void | resetHierarchyConfiguration (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_level, int finest_level) override |
void | applyGradientDetector (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool uses_richardson_extrapolation_too) override |
void | putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override |
Public Member Functions inherited from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM > | |
virtual double | getLevelDt (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double dt_time, const bool initial_time) |
virtual double | advanceLevel (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< NDIM > > hierarchy, const double current_time, const double new_time, const bool first_step, const bool last_step, const bool regrid_advance=false) |
virtual void | resetTimeDependentData (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double new_time, const bool can_be_refined) |
virtual void | resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level) |
virtual void | applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< NDIM > > level, const double error_data_time, const int tag_index, const double deltat, const int error_coarsen_ratio, const bool initial_time, const bool uses_gradient_detector_too) |
virtual void | coarsenDataForRichardsonExtrapolation (const tbox::Pointer< hier::PatchHierarchy< NDIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< NDIM > > coarser_level, const double coarsen_data_time, const bool before_advance) |
Protected Attributes | |
IBHierarchyIntegrator * | d_ib_solver = nullptr |
bool | d_use_fixed_coupling_ops = false |
Class IBStrategy provides a generic interface for specifying the implementation details of a particular version of the IB method.
All classes inheriting from IBStrategy read the following values from the input database to determine which interpolation and spreading delta functions should be used:
interp_delta_fcn
: name of the interpolation kernel, provided as a string. Defaults to "IB_4"
. spread_delta_fcn
: name of the spreading kernel, provided as a string. Defaults to "IB_4"
. IB_delta_fcn
: overriding alias for the two previous entries - has the same default.
|
virtual |
Advance the positions of the Lagrangian structure using the standard 2nd-order Adams-Bashforth rule.
A default implementation is provided that emits an unrecoverable exception.
|
virtual |
Activate a previously inactivated structure or part to be used again in FSI calculations.
[in] | structure_number | Number of the structure/part. |
[in] | level_number | Level on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy. |
Reimplemented in IBAMR::IBMethod.
|
virtual |
Add the estimated computational work from the current object per cell into the specified workload_data_idx
.
An empty default implementation is provided.
Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBInterpolantMethod.
|
overridevirtual |
Set integer tags to "one" in cells where refinement of the given level should occur according to user-supplied feature detection criteria.
An empty default implementation is provided.
Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
Reimplemented in IBAMR::IIMethod, and IBAMR::IBStrategySet.
|
virtual |
Advance the positions of the Lagrangian structure using the (explicit) backward Euler method.
A default implementation is provided that emits an unrecoverable exception.
Reimplemented in IBAMR::IBMethod, IBAMR::IBInterpolantMethod, and IBAMR::CIBMethod.
|
virtual |
Begin redistributing Lagrangian data prior to regridding the patch hierarchy.
An empty default implementation is provided.
Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBInterpolantMethod.
|
virtual |
Compute the Lagrangian source/sink density at the specified time within the current time interval.
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.
|
pure virtual |
Compute the Lagrangian force at the specified time within the current time interval.
Implemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, and IBAMR::GeneralizedIBMethod.
|
virtual |
Complete redistributing Lagrangian data following regridding the patch hierarchy.
An empty default implementation is provided.
Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBInterpolantMethod.
|
pure virtual |
Advance the positions of the Lagrangian structure using the forward Euler method.
Implemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.
|
protected |
Get coarsen algorithm.
|
protected |
Get coarsen schedules.
|
protected |
Get ghost cell-filling refine algorithm.
|
protected |
Get ghost cell-filling refine schedules.
|
protected |
Return a pointer to a HierarchyMathOps object.
|
protected |
Return a pointer to the INSHierarchyIntegrator object being used with the IBHierarchyIntegrator class registered with this IBStrategy object.
|
virtual |
Determine whether or not the given structure or part is currently activated.
[in] | structure_number | Number of the structure/part. |
[in] | level_number | Level on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy. |
Reimplemented in IBAMR::IBMethod.
|
virtual |
Get the ratio of the maximum point displacement of all the structures owned by the current class to the cell width of the grid level on which the structure is assigned. This value is useful for determining if the Eulerian patch hierarchy needs to be regridded.
Reimplemented in IBAMR::IBStrategySet.
|
pure virtual |
Return the number of ghost cells required by the Lagrangian-Eulerian interaction routines.
Implemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBInterpolantMethod.
|
protected |
Return a pointer to the HierarchyDataOpsReal object associated with pressure-like variables.
|
protected |
Get data-prolonging refine algorithm.
|
protected |
Get data-prolonging refine schedules.
|
protected |
Return a pointer to the HierarchyDataOpsReal object associated with velocity-like variables.
|
virtual |
Indicate whether there are any internal fluid sources/sinks.
A default implementation is provided that returns false.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.
|
virtual |
Inactivate a structure or part. Such a structure will be ignored by FSI calculations: i.e., it will have its velocity set to zero and no forces will be spread from the structure to the Eulerian grid.
[in] | structure_number | Number of the structure/part. |
[in] | level_number | Level on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy. |
Reimplemented in IBAMR::IBMethod.
|
overridevirtual |
Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm.
An empty default implementation is provided.
Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, and IBAMR::IBStrategySet.
|
virtual |
Initialize Lagrangian data corresponding to the given AMR patch hierarchy at the start of a computation. If the computation is begun from a restart file, data may be read from the restart databases.
A patch data descriptor is provided for the Eulerian velocity in case initialization requires interpolating Eulerian data. Ghost cells for Eulerian data will be filled upon entry to this function.
An empty default implementation is provided.
Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.
|
virtual |
Compute the pressures at the positions of any distributed internal fluid sources or sinks.
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.
|
pure virtual |
Interpolate the Eulerian velocity to the curvilinear mesh at the specified time within the current time interval.
Implemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.
|
pure virtual |
Advance the positions of the Lagrangian structure using the (explicit) midpoint rule.
Implemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.
|
virtual |
Execute user-defined post-processing operations.
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.
|
virtual |
Method to clean up data following call(s) to integrateHierarchy().
An empty default implementation is provided.
Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.
|
virtual |
Execute user-defined routines just after solving the fluid equations.
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBLevelSetMethod, and IBAMR::ConstraintIBMethod.
|
virtual |
Method to prepare to advance data from current_time to new_time.
An empty default implementation is provided.
Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.
|
virtual |
Execute user-defined routines just before solving the fluid equations.
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBLevelSetMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.
|
overridevirtual |
Write out object state to the given database.
An empty default implementation is provided.
Implements SAMRAI::tbox::Serializable.
Reimplemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, and IBAMR::IBStrategySet.
|
protected |
Register a coarsen algorithm.
|
virtual |
Register Eulerian refinement or coarsening algorithms with the parent IBHierarchyIntegrator using the two versions of the protected methods IBStrategy::registerGhostfillRefineAlgorithm(), IBStrategy::registerProlongRefineAlgorithm(), and IBStrategy::registerCoarsenAlgorithm().
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.
|
virtual |
Register Eulerian variables with the parent IBHierarchyIntegrator with the VariableDatabase, or via the various versions of the protected method IBStrategy::registerVariable().
An empty default implementation is provided.
Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, IBAMR::ConstraintIBMethod, and IBAMR::CIBMethod.
|
protected |
Register a ghost cell-filling refine algorithm.
|
virtual |
Register the IBHierarchyIntegrator object that is using this strategy class.
Reimplemented in IBAMR::IBStrategySet, and IBAMR::IBLevelSetMethod.
|
virtual |
Register a load balancer and work load patch data index with the IB strategy object.
An empty default implementation is provided.
Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, and IBAMR::IBMethod.
|
protected |
Register a data-prolonging refine algorithm.
|
protected |
Register a state variable with the integrator. When a refine operator is specified, the data for the variable are automatically maintained as the patch hierarchy evolves.
All state variables are registered with three contexts: current, new, and scratch. The current context of a state variable is maintained from time step to time step and, if the necessary coarsen and refine operators are specified, as the patch hierarchy evolves.
|
protected |
Register a variable with the integrator that may not be maintained from time step to time step.
By default, variables are registered with the scratch context, which is deallocated after each time step.
|
overridevirtual |
Reset cached hierarchy dependent data.
An empty default implementation is provided.
Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.
Reimplemented in IBAMR::IIMethod, and IBAMR::IBStrategySet.
|
virtual |
Setup the tag buffer.
A default implementation is provided that sets the tag buffer to be at least the minimum ghost cell width.
Reimplemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBInterpolantMethod.
void IBAMR::IBStrategy::setUseFixedLEOperators | ( | bool | use_fixed_coupling_ops = true | ) |
Indicate whether "fixed" interpolation and spreading operators should be used during Lagrangian-Eulerian interaction.
|
virtual |
Indicate that multistep time stepping will be used.
A default implementation is provided that emits an unrecoverable exception.
[in] | n_previous_steps | Number of previous solution values that can be used by the multistep scheme. |
|
virtual |
Spread the Lagrangian source/sink density to the Cartesian grid at the specified time within the current time interval.
An empty default implementation is provided.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.
|
pure virtual |
Spread the Lagrangian force to the Cartesian grid at the specified time within the current time interval.
Implemented in IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.
|
pure virtual |
Advance the positions of the Lagrangian structure using the (explicit) trapezoidal rule.
Implemented in IBAMR::PenaltyIBMethod, IBAMR::IIMethod, IBAMR::IBStrategySet, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::GeneralizedIBMethod, and IBAMR::CIBMethod.
|
virtual |
Update the positions used for the "fixed" interpolation and spreading operators.
A default implementation is provided that emits an unrecoverable exception.
Reimplemented in IBAMR::IBStrategySet, IBAMR::IBMethod, and IBAMR::IBLevelSetMethod.
|
protected |
The IBHierarchyIntegrator object that is using this strategy class.
|
protected |
Whether to use "fixed" Lagrangian-Eulerian coupling operators.