IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
IBAMR::IBMethod Class Reference

Class IBMethod is an implementation of the abstract base class IBImplicitStrategy that provides functionality required by the standard IB method. More...

#include </home/runner/work/IBAMR/IBAMR/include/ibamr/IBMethod.h>

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

Public Member Functions

 IBMethod (std::string object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, bool register_for_restart=true)
 Constructor.
 
 ~IBMethod ()
 Destructor.
 
void registerIBLagrangianForceFunction (SAMRAI::tbox::Pointer< IBLagrangianForceStrategy > ib_force_fcn)
 
void registerIBLagrangianSourceFunction (SAMRAI::tbox::Pointer< IBLagrangianSourceStrategy > ib_source_fcn)
 
void registerLInitStrategy (SAMRAI::tbox::Pointer< IBTK::LInitStrategy > l_initializer)
 
void freeLInitStrategy ()
 
void registerIBMethodPostProcessor (SAMRAI::tbox::Pointer< IBMethodPostProcessStrategy > post_processor)
 
IBTK::LDataManagergetLDataManager () const
 
SAMRAI::tbox::Pointer< IBInstrumentPanelgetIBInstrumentPanel () const
 
void registerLSiloDataWriter (SAMRAI::tbox::Pointer< IBTK::LSiloDataWriter > silo_writer)
 
const SAMRAI::hier::IntVector< NDIM > & getMinimumGhostCellWidth () const override
 
void setupTagBuffer (SAMRAI::tbox::Array< int > &tag_buffer, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) const override
 
virtual void inactivateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) override
 
virtual void activateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) override
 
virtual bool getLagrangianStructureIsActivated (int structure_number=0, int level_number=std::numeric_limits< int >::max()) const override
 
void preprocessIntegrateData (double current_time, double new_time, int num_cycles) override
 
void postprocessIntegrateData (double current_time, double new_time, int num_cycles) override
 
void createSolverVecs (Vec *X_vec, Vec *F_vec) override
 
void setupSolverVecs (Vec *X_vec, Vec *F_vec) override
 
void setUpdatedPosition (Vec &X_new_vec) override
 
void setLinearizedPosition (Vec &X_vec, double data_time) override
 
void computeResidual (Vec &R_vec) override
 
void computeLinearizedResidual (Vec &X_vec, Vec &R_vec) override
 
void updateFixedLEOperators () override
 
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) override
 
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) override
 
void forwardEulerStep (double current_time, double new_time) override
 
void backwardEulerStep (double current_time, double new_time) override
 
void midpointStep (double current_time, double new_time) override
 
void trapezoidalStep (double current_time, double new_time) override
 
void computeLagrangianForce (double data_time) override
 
void computeLinearizedLagrangianForce (Vec &X_vec, double data_time) override
 
void constructLagrangianForceJacobian (Mat &A, MatType mat_type, double data_time) override
 
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) override
 
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) override
 
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) override
 
bool hasFluidSources () const override
 
void computeLagrangianFluidSource (double data_time) override
 
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) override
 
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) override
 
void postprocessData () override
 
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) override
 
void registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx) override
 
void addWorkloadEstimate (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, const int workload_data_idx) override
 
void beginDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
void endDataRedistribution (SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > hierarchy, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) override
 
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 IBAMR::IBImplicitStrategy
 IBImplicitStrategy ()=default
 Constructor.
 
virtual ~IBImplicitStrategy ()=default
 Virtual destructor.
 
- Public Member Functions inherited from IBAMR::IBStrategy
 IBStrategy ()=default
 Constructor.
 
virtual void registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver)
 
virtual void registerEulerianVariables ()
 
virtual void registerEulerianCommunicationAlgorithms ()
 
virtual double getMaxPointDisplacement () const
 
void setUseFixedLEOperators (bool use_fixed_coupling_ops=true)
 
virtual void preprocessSolveFluidEquations (double current_time, double new_time, int cycle_num)
 
virtual void postprocessSolveFluidEquations (double current_time, double new_time, int cycle_num)
 
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 Member Functions

void getPositionData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **X_data, bool **X_needs_ghost_fill, double data_time)
 
void getLinearizedPositionData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **X_data, bool **X_needs_ghost_fill)
 
void getLECouplingPositionData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **X_LE_data, bool **X_LE_needs_ghost_fill, double data_time)
 
void getVelocityData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **U_data, double data_time)
 
void getLinearizedVelocityData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **U_data)
 
void getForceData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **F_data, bool **F_needs_ghost_fill, double data_time)
 
void getLinearizedForceData (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **F_data, bool **F_needs_ghost_fill)
 
void reinitMidpointData (const std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &current_data, const std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &new_data, const std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &half_data)
 
void resetAnchorPointValues (std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > U_data, int coarsest_ln, int finest_ln)
 
PetscErrorCode computeForce (Vec X, Vec F)
 
- Protected Member Functions inherited from IBAMR::IBStrategy
INSHierarchyIntegratorgetINSHierarchyIntegrator () const
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > getVelocityHierarchyDataOps () const
 
SAMRAI::tbox::Pointer< SAMRAI::math::HierarchyDataOpsReal< NDIM, double > > getPressureHierarchyDataOps () const
 
SAMRAI::tbox::Pointer< IBTK::HierarchyMathOpsgetHierarchyMathOps () 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< IBTK::CartGridFunction > init_fcn=SAMRAI::tbox::Pointer< IBTK::CartGridFunction >(NULL))
 
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 >(NULL))
 
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
 

Static Protected Member Functions

static PetscErrorCode computeForce_SAMRAI (void *ctx, Vec X, Vec F)
 

Protected Attributes

bool d_do_log = false
 
SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > > d_hierarchy
 
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > d_gridding_alg
 
double d_current_time = std::numeric_limits<double>::quiet_NaN()
 
double d_new_time = std::numeric_limits<double>::quiet_NaN()
 
double d_half_time = std::numeric_limits<double>::quiet_NaN()
 
bool d_X_current_needs_ghost_fill = true
 
bool d_X_new_needs_ghost_fill = true
 
bool d_X_half_needs_ghost_fill = true
 
bool d_X_jac_needs_ghost_fill = true
 
bool d_X_LE_new_needs_ghost_fill = true
 
bool d_X_LE_half_needs_ghost_fill = true
 
bool d_F_current_needs_ghost_fill = true
 
bool d_F_new_needs_ghost_fill = true
 
bool d_F_half_needs_ghost_fill = true
 
bool d_F_jac_needs_ghost_fill = true
 
IBTK::LDataManagerd_l_data_manager
 
std::string d_interp_kernel_fcn = "IB_4"
 
std::string d_spread_kernel_fcn = "IB_4"
 
bool d_error_if_points_leave_domain = false
 
SAMRAI::hier::IntVector< NDIM > d_ghosts
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_current_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_new_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_half_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_jac_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_LE_new_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_X_LE_half_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_U_current_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_U_new_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_U_half_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_U_jac_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_F_current_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_F_new_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_F_half_data
 
std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > d_F_jac_data
 
std::vector< std::set< int > > d_anchor_point_local_idxs
 
SAMRAI::tbox::Pointer< IBInstrumentPaneld_instrument_panel
 
std::vector< double > d_total_flow_volume
 
SAMRAI::tbox::Pointer< IBTK::LInitStrategyd_l_initializer
 
SAMRAI::tbox::Pointer< IBLagrangianForceStrategyd_ib_force_fcn
 
bool d_ib_force_fcn_needs_init = true
 
SAMRAI::tbox::Pointer< IBLagrangianSourceStrategyd_ib_source_fcn
 
bool d_ib_source_fcn_needs_init = true
 
std::vector< std::vector< IBTK::Point > > d_X_src
 
std::vector< std::vector< double > > d_r_src
 
std::vector< std::vector< double > > d_P_src
 
std::vector< std::vector< double > > d_Q_src
 
std::vector< int > d_n_src
 
bool d_normalize_source_strength = false
 
SAMRAI::tbox::Pointer< IBMethodPostProcessStrategyd_post_processor
 
SAMRAI::tbox::Pointer< IBTK::LSiloDataWriterd_silo_writer
 
SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > d_load_balancer
 
int d_workload_idx = IBTK::invalid_index
 
std::string d_object_name
 
bool d_registered_for_restart
 
- Protected Attributes inherited from IBAMR::IBStrategy
IBHierarchyIntegratord_ib_solver = nullptr
 
bool d_use_fixed_coupling_ops = false
 

Detailed Description

Class IBMethod is an implementation of the abstract base class IBImplicitStrategy that provides functionality required by the standard IB method.

Member Function Documentation

◆ activateLagrangianStructure()

void IBAMR::IBMethod::activateLagrangianStructure ( int  structure_number = 0,
int  level_number = std::numeric_limits<int>::max() 
)
overridevirtual

Activate a previously inactivated structure/part to be used again in FSI calculations. See IBAMR::IBStrategy::activateLagrangianStructure().

Reimplemented from IBAMR::IBStrategy.

◆ addWorkloadEstimate()

void IBAMR::IBMethod::addWorkloadEstimate ( SAMRAI::tbox::Pointer< SAMRAI::hier::PatchHierarchy< NDIM > >  hierarchy,
const int  workload_data_idx 
)
overridevirtual

Add the estimated computational work from the current object per cell into the specified workload_data_idx.

Reimplemented from IBAMR::IBStrategy.

◆ applyGradientDetector()

void IBAMR::IBMethod::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 
)
overridevirtual

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

See also
SAMRAI::mesh::StandardTagAndInitStrategy::applyGradientDetector

Reimplemented from SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.

◆ backwardEulerStep()

void IBAMR::IBMethod::backwardEulerStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the backward Euler method.

Reimplemented from IBAMR::IBStrategy.

◆ beginDataRedistribution()

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

Begin redistributing Lagrangian data prior to regridding the patch hierarchy.

Reimplemented from IBAMR::IBStrategy.

◆ computeLagrangianFluidSource()

void IBAMR::IBMethod::computeLagrangianFluidSource ( double  data_time)
overridevirtual

Compute the Lagrangian source/sink density at the specified time within the current time interval.

Reimplemented from IBAMR::IBStrategy.

◆ computeLagrangianForce()

void IBAMR::IBMethod::computeLagrangianForce ( double  data_time)
overridevirtual

Compute the Lagrangian force at the specified time within the current time interval.

Implements IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ computeLinearizedLagrangianForce()

void IBAMR::IBMethod::computeLinearizedLagrangianForce ( Vec &  X_vec,
double  data_time 
)
overridevirtual

Compute the Lagrangian force of the linearized problem for the specified configuration of the updated position vector.

Implements IBAMR::IBImplicitStrategy.

◆ computeLinearizedResidual()

void IBAMR::IBMethod::computeLinearizedResidual ( Vec &  X_vec,
Vec &  R_vec 
)
overridevirtual

Compute the linearized residual for the given intermediate position vector.

Implements IBAMR::IBImplicitStrategy.

◆ computeResidual()

void IBAMR::IBMethod::computeResidual ( Vec &  R_vec)
overridevirtual

Compute the residual on the specified level of the patch hierarchy.

Implements IBAMR::IBImplicitStrategy.

◆ constructInterpOp()

void IBAMR::IBMethod::constructInterpOp ( Mat &  J,
void(*)(const double, double *)  spread_fnc,
int  stencil_width,
const std::vector< int > &  num_dofs_per_proc,
int  dof_index_idx,
double  data_time 
)
overridevirtual

Construct the IB interpolation operator.

Implements IBAMR::IBImplicitStrategy.

◆ constructLagrangianForceJacobian()

void IBAMR::IBMethod::constructLagrangianForceJacobian ( Mat &  A,
MatType  mat_type,
double  data_time 
)
overridevirtual

Construct the linearized Lagrangian force Jacobian.

Implements IBAMR::IBImplicitStrategy.

◆ createSolverVecs()

void IBAMR::IBMethod::createSolverVecs ( Vec *  X_vec,
Vec *  F_vec 
)
overridevirtual

Create solution and rhs data on the specified level of the patch hierarchy.

Implements IBAMR::IBImplicitStrategy.

◆ endDataRedistribution()

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

Complete redistributing Lagrangian data following regridding the patch hierarchy.

Reimplemented from IBAMR::IBStrategy.

◆ forwardEulerStep()

void IBAMR::IBMethod::forwardEulerStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the forward Euler method.

Implements IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ freeLInitStrategy()

void IBAMR::IBMethod::freeLInitStrategy ( )

Free references to Lagrangian initialization objects.

◆ getForceData()

void IBAMR::IBMethod::getForceData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  F_data,
bool **  F_needs_ghost_fill,
double  data_time 
)
protected

Get the current structure force data.

◆ getIBInstrumentPanel()

Pointer< IBInstrumentPanel > IBAMR::IBMethod::getIBInstrumentPanel ( ) const

Return a pointer to the instrumentation manager object.

◆ getLagrangianStructureIsActivated()

bool IBAMR::IBMethod::getLagrangianStructureIsActivated ( int  structure_number = 0,
int  level_number = std::numeric_limits<int>::max() 
) const
overridevirtual

Determine whether or not the given structure or part is currently activated. See IBAMR::IBStrategy::getLagrangianStructureIsActivated().

Reimplemented from IBAMR::IBStrategy.

◆ getLDataManager()

LDataManager * IBAMR::IBMethod::getLDataManager ( ) const

Return a pointer to the Lagrangian data manager object.

◆ getLECouplingPositionData()

void IBAMR::IBMethod::getLECouplingPositionData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  X_LE_data,
bool **  X_LE_needs_ghost_fill,
double  data_time 
)
protected

Get the current interpolation/spreading position data.

◆ getLinearizedForceData()

void IBAMR::IBMethod::getLinearizedForceData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  F_data,
bool **  F_needs_ghost_fill 
)
protected

Get the linearized structure force data.

◆ getLinearizedPositionData()

void IBAMR::IBMethod::getLinearizedPositionData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  X_data,
bool **  X_needs_ghost_fill 
)
protected

Get the linearized structure position data.

◆ getLinearizedVelocityData()

void IBAMR::IBMethod::getLinearizedVelocityData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  U_data)
protected

Get the linearized structure velocity data.

◆ getMinimumGhostCellWidth()

const IntVector< NDIM > & IBAMR::IBMethod::getMinimumGhostCellWidth ( ) const
overridevirtual

Return the number of ghost cells required by the Lagrangian-Eulerian interaction routines.

Implements IBAMR::IBStrategy.

◆ getPositionData()

void IBAMR::IBMethod::getPositionData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  X_data,
bool **  X_needs_ghost_fill,
double  data_time 
)
protected

Get the current structure position data.

◆ getVelocityData()

void IBAMR::IBMethod::getVelocityData ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > **  U_data,
double  data_time 
)
protected

Get the current structure velocity data.

◆ hasFluidSources()

bool IBAMR::IBMethod::hasFluidSources ( ) const
overridevirtual

Indicate whether there are any internal fluid sources/sinks.

Reimplemented from IBAMR::IBStrategy.

◆ inactivateLagrangianStructure()

void IBAMR::IBMethod::inactivateLagrangianStructure ( int  structure_number = 0,
int  level_number = std::numeric_limits<int>::max() 
)
overridevirtual

Inactivate a structure/part. See IBAMR::IBStrategy::inactivateLagrangianStructure().

Reimplemented from IBAMR::IBStrategy.

◆ initializeLevelData()

void IBAMR::IBMethod::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 
)
overridevirtual

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

See also
SAMRAI::mesh::StandardTagAndInitStrategy::initializeLevelData

Implements SAMRAI::mesh::StandardTagAndInitStrategy< NDIM >.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ initializePatchHierarchy()

void IBAMR::IBMethod::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 
)
overridevirtual

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.

Reimplemented from IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ interpolateLinearizedVelocity()

void IBAMR::IBMethod::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 
)
overridevirtual

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.

Implements IBAMR::IBImplicitStrategy.

◆ interpolatePressure()

void IBAMR::IBMethod::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 
)
overridevirtual

Compute the pressures at the positions of any distributed internal fluid sources or sinks.

Reimplemented from IBAMR::IBStrategy.

◆ interpolateVelocity()

void IBAMR::IBMethod::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 
)
overridevirtual

Interpolate the Eulerian velocity to the curvilinear mesh at the specified time within the current time interval.

Implements IBAMR::IBStrategy.

◆ midpointStep()

void IBAMR::IBMethod::midpointStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the midpoint rule.

Implements IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ postprocessData()

void IBAMR::IBMethod::postprocessData ( )
overridevirtual

Execute user-defined post-processing operations.

Reimplemented from IBAMR::IBStrategy.

◆ postprocessIntegrateData()

void IBAMR::IBMethod::postprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
overridevirtual

Method to clean up data following call(s) to integrateHierarchy().

Reimplemented from IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ preprocessIntegrateData()

void IBAMR::IBMethod::preprocessIntegrateData ( double  current_time,
double  new_time,
int  num_cycles 
)
overridevirtual

Method to prepare to advance data from current_time to new_time.

Reimplemented from IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ putToDatabase()

void IBAMR::IBMethod::putToDatabase ( SAMRAI::tbox::Pointer< SAMRAI::tbox::Database db)
overridevirtual

Write out object state to the given database.

Implements SAMRAI::tbox::Serializable.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ registerIBLagrangianForceFunction()

void IBAMR::IBMethod::registerIBLagrangianForceFunction ( SAMRAI::tbox::Pointer< IBLagrangianForceStrategy ib_force_fcn)

Supply a Lagrangian force object.

◆ registerIBLagrangianSourceFunction()

void IBAMR::IBMethod::registerIBLagrangianSourceFunction ( SAMRAI::tbox::Pointer< IBLagrangianSourceStrategy ib_source_fcn)

Supply a Lagrangian source object.

◆ registerIBMethodPostProcessor()

void IBAMR::IBMethod::registerIBMethodPostProcessor ( SAMRAI::tbox::Pointer< IBMethodPostProcessStrategy post_processor)

Supply a post processor object.

◆ registerLInitStrategy()

void IBAMR::IBMethod::registerLInitStrategy ( SAMRAI::tbox::Pointer< IBTK::LInitStrategy l_initializer)

Supply a Lagrangian initialization object.

◆ registerLoadBalancer()

void IBAMR::IBMethod::registerLoadBalancer ( SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > >  load_balancer,
int  workload_data_idx 
)
overridevirtual

Register a load balancer and work load patch data index with the IB strategy object.

Deprecated:
This method is no longer necessary with the current workload estimation scheme.

Reimplemented from IBAMR::IBStrategy.

◆ registerLSiloDataWriter()

void IBAMR::IBMethod::registerLSiloDataWriter ( SAMRAI::tbox::Pointer< IBTK::LSiloDataWriter silo_writer)

Register a Lagrangian Silo data writer so this class will write plot files that may be postprocessed with the VisIt visualization tool.

◆ reinitMidpointData()

void IBAMR::IBMethod::reinitMidpointData ( const std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &  current_data,
const std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &  new_data,
const std::vector< SAMRAI::tbox::Pointer< IBTK::LData > > &  half_data 
)
protected

Interpolate the current and new data to obtain values at the midpoint of the time interval.

◆ resetAnchorPointValues()

void IBAMR::IBMethod::resetAnchorPointValues ( std::vector< SAMRAI::tbox::Pointer< IBTK::LData > >  U_data,
int  coarsest_ln,
int  finest_ln 
)
protected

Set the elements of the Lagrangian vector to zero at anchored nodes of the curvilinear mesh.

◆ resetHierarchyConfiguration()

void IBAMR::IBMethod::resetHierarchyConfiguration ( SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > >  hierarchy,
int  coarsest_level,
int  finest_level 
)
overridevirtual

◆ setLinearizedPosition()

void IBAMR::IBMethod::setLinearizedPosition ( Vec &  X_vec,
double  data_time 
)
overridevirtual

Set the value of the intermediate position vector used in evaluating the linearized problem.

Implements IBAMR::IBImplicitStrategy.

◆ setUpdatedPosition()

void IBAMR::IBMethod::setUpdatedPosition ( Vec &  X_new_vec)
overridevirtual

Set the value of the updated position vector.

Implements IBAMR::IBImplicitStrategy.

◆ setupSolverVecs()

void IBAMR::IBMethod::setupSolverVecs ( Vec *  X_vec,
Vec *  F_vec 
)
overridevirtual

Setup solution and rhs data on the specified level of the patch hierarchy.

Implements IBAMR::IBImplicitStrategy.

◆ setupTagBuffer()

void IBAMR::IBMethod::setupTagBuffer ( SAMRAI::tbox::Array< int > &  tag_buffer,
SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > >  gridding_alg 
) const
overridevirtual

Setup the tag buffer.

Reimplemented from IBAMR::IBStrategy.

◆ spreadFluidSource()

void IBAMR::IBMethod::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 
)
overridevirtual

Spread the Lagrangian source/sink density to the Cartesian grid at the specified time within the current time interval.

Reimplemented from IBAMR::IBStrategy.

◆ spreadForce()

void IBAMR::IBMethod::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 
)
overridevirtual

Spread the Lagrangian force to the Cartesian grid at the specified time within the current time interval.

Implements IBAMR::IBStrategy.

◆ spreadLinearizedForce()

void IBAMR::IBMethod::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 
)
overridevirtual

Spread the Lagrangian force of the linearized problem to the Cartesian grid at the specified time within the current time interval.

Implements IBAMR::IBImplicitStrategy.

◆ trapezoidalStep()

void IBAMR::IBMethod::trapezoidalStep ( double  current_time,
double  new_time 
)
overridevirtual

Advance the positions of the Lagrangian structure using the trapezoidal rule.

Implements IBAMR::IBStrategy.

Reimplemented in IBAMR::PenaltyIBMethod.

◆ updateFixedLEOperators()

void IBAMR::IBMethod::updateFixedLEOperators ( )
overridevirtual

Update the positions used for the "fixed" interpolation and spreading operators.

Reimplemented from IBAMR::IBStrategy.


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