IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class IBImplicitStrategy provides a generic interface for specifying the implementation details of a particular implicit version of the IB method. More...
#include </home/runner/work/IBAMR/IBAMR/include/ibamr/IBImplicitStrategy.h>
Public Member Functions | |
IBImplicitStrategy ()=default | |
Constructor. | |
virtual | ~IBImplicitStrategy ()=default |
Virtual destructor. | |
virtual void | createSolverVecs (Vec *X_vec, Vec *F_vec)=0 |
virtual void | setupSolverVecs (Vec *X_vec, Vec *F_vec)=0 |
virtual void | setUpdatedPosition (Vec &X_new_vec)=0 |
virtual void | setLinearizedPosition (Vec &X_vec, double data_time)=0 |
virtual void | computeResidual (Vec &R_vec)=0 |
virtual void | computeLinearizedResidual (Vec &X_vec, Vec &R_vec)=0 |
virtual void | interpolateLinearizedVelocity (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 | computeLinearizedLagrangianForce (Vec &X_vec, double data_time)=0 |
virtual void | constructLagrangianForceJacobian (Mat &A, MatType mat_type, double data_time)=0 |
virtual void | spreadLinearizedForce (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 void | constructInterpOp (Mat &J, void(*spread_fnc)(const double, double *), int stencil_width, const std::vector< int > &num_dofs_per_proc, int dof_index_idx, double data_time)=0 |
Public Member Functions inherited from IBAMR::IBStrategy | |
IBStrategy ()=default | |
Constructor. | |
virtual void | registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver) |
virtual void | registerEulerianVariables () |
virtual void | registerEulerianCommunicationAlgorithms () |
virtual const SAMRAI::hier::IntVector< NDIM > & | getMinimumGhostCellWidth () const =0 |
virtual void | setupTagBuffer (SAMRAI::tbox::Array< int > &tag_buffer, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) const |
virtual void | inactivateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) |
virtual void | activateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) |
virtual bool | getLagrangianStructureIsActivated (int structure_number=0, int level_number=std::numeric_limits< int >::max()) const |
virtual double | getMaxPointDisplacement () const |
virtual void | preprocessIntegrateData (double current_time, double new_time, int num_cycles) |
virtual void | postprocessIntegrateData (double current_time, double new_time, int num_cycles) |
void | setUseFixedLEOperators (bool use_fixed_coupling_ops=true) |
virtual void | updateFixedLEOperators () |
virtual void | interpolateVelocity (int u_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &u_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &u_ghost_fill_scheds, double data_time)=0 |
virtual void | forwardEulerStep (double current_time, double new_time)=0 |
virtual void | backwardEulerStep (double current_time, double new_time) |
virtual void | midpointStep (double current_time, double new_time)=0 |
virtual void | trapezoidalStep (double current_time, double new_time)=0 |
virtual void | computeLagrangianForce (double data_time)=0 |
virtual void | spreadForce (int f_data_idx, IBTK::RobinPhysBdryPatchStrategy *f_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &f_prolongation_scheds, double data_time)=0 |
virtual bool | hasFluidSources () const |
virtual void | computeLagrangianFluidSource (double data_time) |
virtual void | spreadFluidSource (int q_data_idx, IBTK::RobinPhysBdryPatchStrategy *q_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &q_prolongation_scheds, double data_time) |
virtual void | interpolatePressure (int p_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &p_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &p_ghost_fill_scheds, double data_time) |
virtual void | preprocessSolveFluidEquations (double current_time, double new_time, int cycle_num) |
virtual void | postprocessSolveFluidEquations (double current_time, double new_time, int cycle_num) |
virtual void | postprocessData () |
virtual void | initializePatchHierarchy (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg, int u_data_idx, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::CoarsenSchedule< NDIM > > > &u_synch_scheds, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &u_ghost_fill_scheds, int integrator_step, double init_data_time, bool initial_time) |
virtual void | registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx) |
virtual void | addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx) |
virtual void | beginDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) |
virtual void | endDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) |
void | initializeLevelData (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double init_data_time, bool can_be_refined, bool initial_time, SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchLevel< NDIM > > old_level, bool allocate_data) override |
void | resetHierarchyConfiguration (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int coarsest_level, int finest_level) override |
void | applyGradientDetector (SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > hierarchy, int level_number, double error_data_time, int tag_index, bool initial_time, bool uses_richardson_extrapolation_too) override |
void | putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override |
Public Member Functions inherited from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM > | |
virtual double | getLevelDt (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double dt_time, const bool initial_time) |
virtual double | advanceLevel (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< NDIM > > hierarchy, const double current_time, const double new_time, const bool first_step, const bool last_step, const bool regrid_advance=false) |
virtual void | resetTimeDependentData (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level, const double new_time, const bool can_be_refined) |
virtual void | resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< NDIM > > level) |
virtual void | applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< NDIM > > level, const double error_data_time, const int tag_index, const double deltat, const int error_coarsen_ratio, const bool initial_time, const bool uses_gradient_detector_too) |
virtual void | coarsenDataForRichardsonExtrapolation (const tbox::Pointer< hier::PatchHierarchy< NDIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< NDIM > > coarser_level, const double coarsen_data_time, const bool before_advance) |
Class IBImplicitStrategy provides a generic interface for specifying the implementation details of a particular implicit version of the IB method.
|
pure virtual |
Compute the Lagrangian force of the linearized problem for the specified configuration of the updated position vector.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Compute the linearized residual for the given intermediate position vector.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Compute the nonlinear residual.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Construct the IB interpolation operator.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Construct the linearized Lagrangian force Jacobian.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Create solution and rhs data.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Interpolate the Eulerian velocity to the curvilinear mesh at the specified time within the current time interval for use in evaluating the residual of the linearized problem.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Set the value of the intermediate position vector used in evaluating the linearized problem.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Set the value of the updated position vector.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Setup solution and rhs data.
Implemented in IBAMR::IBMethod.
|
pure virtual |
Spread the Lagrangian force of the linearized problem to the Cartesian grid at the specified time within the current time interval.
Implemented in IBAMR::IBMethod.