IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class INSHierarchyIntegrator provides an abstract interface for a time integrator for the incompressible Navier-Stokes equations 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/INSHierarchyIntegrator.h>
Protected Attributes | |
bool | d_integrator_is_initialized = false |
TimeSteppingType | d_viscous_time_stepping_type = TRAPEZOIDAL_RULE |
TimeSteppingType | d_convective_time_stepping_type = ADAMS_BASHFORTH |
TimeSteppingType | d_init_convective_time_stepping_type = MIDPOINT_RULE |
StokesSpecifications | d_problem_coefs |
std::vector< SAMRAI::tbox::Pointer< AdvDiffHierarchyIntegrator > > | d_adv_diff_hier_integrators |
SAMRAI::tbox::Pointer< SAMRAI::pdat::FaceVariable< NDIM, double > > | d_U_adv_diff_var |
double | d_cfl_current = std::numeric_limits<double>::quiet_NaN() |
double | d_cfl_max = 1.0 |
bool | d_using_vorticity_tagging = false |
SAMRAI::tbox::Array< double > | d_Omega_rel_thresh |
SAMRAI::tbox::Array< double > | d_Omega_abs_thresh |
bool | d_normalize_pressure = false |
bool | d_normalize_velocity = false |
bool | d_creeping_flow = false |
double | d_regrid_max_div_growth_factor = 1.1 |
double | d_U_scale = 1.0 |
double | d_P_scale = 1.0 |
double | d_F_scale = 1.0 |
double | d_Q_scale = 1.0 |
double | d_Omega_scale = 1.0 |
double | d_Div_U_scale = 1.0 |
double | d_EE_scale = 1.0 |
bool | d_output_U = true |
bool | d_output_P = true |
bool | d_output_F = false |
bool | d_output_Q = false |
bool | d_output_Omega = true |
bool | d_output_Div_U = true |
bool | d_output_EE = false |
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > | d_U_var |
std::string | d_U_coarsen_type = "CONSERVATIVE_COARSEN" |
std::string | d_U_refine_type = "CONSERVATIVE_LINEAR_REFINE" |
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > | d_P_var |
std::string | d_P_coarsen_type = "CONSERVATIVE_COARSEN" |
std::string | d_P_refine_type = "LINEAR_REFINE" |
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > | d_F_var |
std::string | d_F_coarsen_type = "CONSERVATIVE_COARSEN" |
std::string | d_F_refine_type = "CONSERVATIVE_LINEAR_REFINE" |
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > | d_Q_var |
std::string | d_Q_coarsen_type = "CONSERVATIVE_COARSEN" |
std::string | d_Q_refine_type = "CONSTANT_REFINE" |
SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | d_U_init |
SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | d_P_init |
SAMRAI::solv::LocationIndexRobinBcCoefs< NDIM > | d_default_bc_coefs |
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > | d_bc_coefs |
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > | d_U_bc_coefs |
std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > | d_U_star_bc_coefs |
TractionBcType | d_traction_bc_type = TRACTION |
SAMRAI::solv::RobinBcCoefStrategy< NDIM > * | d_P_bc_coef |
std::unique_ptr< SAMRAI::solv::RobinBcCoefStrategy< NDIM > > | d_Phi_bc_coef |
SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | d_F_fcn |
SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | d_Q_fcn |
SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolation > | d_U_bdry_bc_fill_op |
SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolation > | d_P_bdry_bc_fill_op |
SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolation > | d_Q_bdry_bc_fill_op |
SAMRAI::tbox::Pointer< IBTK::HierarchyGhostCellInterpolation > | d_no_fill_op |
bool | d_use_div_sink_drag_term = false |
int | d_coarsest_reset_ln |
int | d_finest_reset_ln |
std::string | d_convective_op_type = "DEFAULT" |
ConvectiveDifferencingType | d_convective_difference_form = ADVECTIVE |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_convective_op_input_db |
SAMRAI::tbox::Pointer< ConvectiveOperator > | d_convective_op |
bool | d_convective_op_needs_init |
std::string | d_velocity_solver_type |
std::string | d_velocity_precond_type |
std::string | d_velocity_sub_precond_type |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_velocity_solver_db |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_velocity_precond_db |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_velocity_sub_precond_db |
SAMRAI::tbox::Pointer< IBTK::PoissonSolver > | d_velocity_solver |
bool | d_velocity_solver_needs_init |
std::string | d_pressure_solver_type |
std::string | d_pressure_precond_type |
std::string | d_pressure_sub_precond_type |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_pressure_solver_db |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_pressure_precond_db |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_pressure_sub_precond_db |
SAMRAI::tbox::Pointer< IBTK::PoissonSolver > | d_pressure_solver |
bool | d_pressure_solver_needs_init |
std::string | d_regrid_projection_solver_type |
std::string | d_regrid_projection_precond_type |
std::string | d_regrid_projection_sub_precond_type |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_regrid_projection_solver_db |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_regrid_projection_precond_db |
SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_regrid_projection_sub_precond_db |
![]() | |
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 |
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 INSHierarchyIntegrator provides an abstract interface for a time integrator for the incompressible Navier-Stokes equations on an AMR grid hierarchy, along with basic data management for variables defined on that hierarchy.
IBAMR::INSHierarchyIntegrator::~INSHierarchyIntegrator | ( | ) |
The destructor for class INSHierarchyIntegrator unregisters the integrator object with the restart manager when the object is so registered.
|
protected |
The constructor for class INSHierarchyIntegrator sets some default values, reads in configuration information from input and restart databases, and registers the integrator object with the restart manager when requested.
This constructor sets the default coarsen and refine operator types to:
The other constructor allows these default values to be overridden.
|
protected |
The constructor for class INSHierarchyIntegrator sets some default values, reads in configuration information from input and restart databases, and registers the integrator object with the restart manager when requested.
Pointer< FaceVariable< NDIM, double > > IBAMR::INSHierarchyIntegrator::getAdvectionVelocityVariable | ( | ) | const |
Return a pointer to a fluid velocity variable that can be used to advect quantities registered with an advection-diffusion solver.
Data are allocated for this variable using the current context. Patch data for this variable are allocated only when an advection-diffusion solver is registered with the Navier-Stokes solver.
Return a pointer to the body force variable.
|
pure virtual |
Get the convective operator being used by this solver class.
If the time integrator is configured to solve the time-dependent (creeping) Stokes equations, then the returned pointer will be nullptr.
If the convective operator has not already been constructed, then this function will initialize a default convective operator.
Implemented in IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, and IBAMR::INSVCStaggeredHierarchyIntegrator.
TimeSteppingType IBAMR::INSHierarchyIntegrator::getConvectiveTimeSteppingType | ( | ) | const |
Get the type of convective time integration scheme being employed by the incompressible flow solver.
Different implementations may support different time stepping schemes.
|
virtual |
Return the current CFL number (i.e., the CFL number computed from the current velocity field). This is typically computed at the end of each time step.
Return a pointer to the source strength variable.
TimeSteppingType IBAMR::INSHierarchyIntegrator::getInitialConvectiveTimeSteppingType | ( | ) | const |
Get the type of convective time integration scheme being employed by the incompressible flow solver during the initial time step.
Different implementations may support different time stepping schemes.
std::vector< RobinBcCoefStrategy< NDIM > * > IBAMR::INSHierarchyIntegrator::getIntermediateVelocityBoundaryConditions | ( | ) | const |
Get a vector of pointers to the intermediate velocity boundary condition specification objects.
|
overrideprotectedvirtual |
Return the maximum stable time step size.
Reimplemented from IBTK::HierarchyIntegrator.
|
protected |
Compute the maximum vorticity magnitude at any given point.
Omega_idx
.
|
overridevirtual |
Returns the number of cycles to perform for the present time step.
Reimplemented from IBTK::HierarchyIntegrator.
Reimplemented in IBAMR::INSVCStaggeredConservativeHierarchyIntegrator.
|
virtual |
Get a pointer to the pressure boundary condition specification object.
|
pure virtual |
Get the subdomain solver for the pressure subsystem. Such solvers can be useful in constructing block preconditioners.
If the pressure subdomain solver has not already been constructed, then this function will initialize a default solver.
Implemented in IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, and IBAMR::INSVCStaggeredHierarchyIntegrator.
Return a pointer to the fluid pressure state variable.
RobinBcCoefStrategy< NDIM > * IBAMR::INSHierarchyIntegrator::getProjectionBoundaryConditions | ( | ) | const |
Get a pointer to the projection Poisson problem boundary condition specification object.
|
protectedpure virtual |
Determine the largest stable timestep on an individual patch.
Implemented in IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, and IBAMR::INSVCStaggeredHierarchyIntegrator.
|
protected |
Determine the largest stable timestep on an individual patch level.
const StokesSpecifications * IBAMR::INSHierarchyIntegrator::getStokesSpecifications | ( | ) | const |
Get a const pointer to the problem coefficients object used by the solver.
|
virtual |
Get a vector of pointers to the velocity boundary condition specification objects.
Return a pointer to the variable that specifies the divergence of the velocity.
|
pure virtual |
Get the subdomain solver for the velocity subsystem. Such solvers can be useful in constructing block preconditioners.
If the velocity subdomain solver has not already been constructed, then this function will initialize a default solver.
Implemented in IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, and IBAMR::INSVCStaggeredHierarchyIntegrator.
Return a pointer to the fluid velocity variable.
TimeSteppingType IBAMR::INSHierarchyIntegrator::getViscousTimeSteppingType | ( | ) | const |
Get the type of viscous time integration scheme being employed by the incompressible flow solver.
Different implementations may support different time stepping schemes.
|
overridevirtual |
Finish postprocessing the hierarchy by computing the current CFL number.
Reimplemented from IBTK::HierarchyIntegrator.
Reimplemented in IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, IBAMR::INSVCStaggeredHierarchyIntegrator, and IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator.
|
overrideprotectedvirtual |
Write out specialized object state to the given database.
Reimplemented from IBTK::HierarchyIntegrator.
Reimplemented in IBAMR::INSStaggeredHierarchyIntegrator.
void IBAMR::INSHierarchyIntegrator::registerAdvDiffHierarchyIntegrator | ( | SAMRAI::tbox::Pointer< AdvDiffHierarchyIntegrator > | adv_diff_hier_integrator | ) |
Register an advection-diffusion solver with this incompressible Navier-Stokes solver.
void IBAMR::INSHierarchyIntegrator::registerBodyForceFunction | ( | SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | F_fcn | ) |
Supply a body force.
void IBAMR::INSHierarchyIntegrator::registerFluidSourceFunction | ( | SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | Q_fcn | ) |
Supply a fluid source/sink distribution.
void IBAMR::INSHierarchyIntegrator::registerPhysicalBoundaryConditions | ( | const std::vector< SAMRAI::solv::RobinBcCoefStrategy< NDIM > * > & | bc_coefs | ) |
Supply a physical boundary conditions specificaion for the velocity field. Boundary conditions take the form of
CartesianGeometry
database in the input file, specifically with the flag periodic_dimension
.void IBAMR::INSHierarchyIntegrator::registerPressureInitialConditions | ( | SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | P_init | ) |
Supply initial conditions for the pressure.
void IBAMR::INSHierarchyIntegrator::registerVelocityDivergenceFunction | ( | SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | Q_fcn | ) |
Supply a CartGridFunction that specifies
Q_fcn | is not specified, then ![]() |
void IBAMR::INSHierarchyIntegrator::registerVelocityInitialConditions | ( | SAMRAI::tbox::Pointer< IBTK::CartGridFunction > | U_init | ) |
Supply initial conditions for the velocity field.
|
protectedpure virtual |
Pure virtual method to project the velocity field following a regridding operation.
Implemented in IBAMR::INSCollocatedHierarchyIntegrator, IBAMR::INSStaggeredHierarchyIntegrator, IBAMR::INSVCStaggeredConservativeHierarchyIntegrator, and IBAMR::INSVCStaggeredNonConservativeHierarchyIntegrator.
void IBAMR::INSHierarchyIntegrator::setConvectiveOperator | ( | SAMRAI::tbox::Pointer< ConvectiveOperator > | convective_op | ) |
Register an operator to compute the convective acceleration term u*grad u.
If the supplied operator is nullptr, then the integrator will solve the time-dependent (creeping) Stokes equations instead of the Navier-Stokes equations.
void IBAMR::INSHierarchyIntegrator::setConvectiveOperatorNeedsInit | ( | ) |
Indicate that the convective operator should be (re-)initialized before the next time step.
void IBAMR::INSHierarchyIntegrator::setConvectiveTimeSteppingType | ( | TimeSteppingType | convective_time_stepping_type | ) |
Set the type of convective time integration scheme being employed by the incompressible flow solver.
Different implementations may support different time stepping schemes.
void IBAMR::INSHierarchyIntegrator::setInitialConvectiveTimeSteppingType | ( | TimeSteppingType | init_convective_time_stepping_type | ) | const |
Set the type of convective time integration scheme being employed by the incompressible flow solver during the initial time step.
Different implementations may support different time stepping schemes.
void IBAMR::INSHierarchyIntegrator::setPressureSubdomainSolver | ( | SAMRAI::tbox::Pointer< IBTK::PoissonSolver > | pressure_solver | ) |
Register a solver for the pressure subsystem.
void IBAMR::INSHierarchyIntegrator::setPressureSubdomainSolverNeedsInit | ( | ) |
Indicate that the velocity subdomain solver should be (re-)initialized before the next time step.
void IBAMR::INSHierarchyIntegrator::setStokesSpecifications | ( | StokesSpecifications | problem_coefs | ) |
Set the problem coefficients used by the solver.
void IBAMR::INSHierarchyIntegrator::setVelocitySubdomainSolver | ( | SAMRAI::tbox::Pointer< IBTK::PoissonSolver > | velocity_solver | ) |
Register a solver for the velocity subsystem.
void IBAMR::INSHierarchyIntegrator::setVelocitySubdomainSolverNeedsInit | ( | ) |
Indicate that the velocity subdomain solver should be (re-)initialized before the next time step.
void IBAMR::INSHierarchyIntegrator::setViscousTimeSteppingType | ( | TimeSteppingType | viscous_time_stepping_type | ) |
Set the type of viscous time integration scheme being employed by the incompressible flow solver.
Different implementations may support different time stepping schemes.
|
protected |
Tag cells based on the vorticity magnitude.
|
protectedvirtual |
Update the current CFL number (i.e., at the end of a timestep).
|
protected |
Current CFL number.
|
protected |
The maximum CFL number.
|
protected |
Hierarchy operators and solvers and related configuration data.
|
protected |
Enum indicating the time integration employed for the explicit discretization of the convective terms.
|
protected |
This boolean value determines whether the convective acceleration term is included in the momentum equation. (If it is not, this solver effectively solves the so-called creeping Stokes equations.)
|
protected |
Enum indicating the time integration employed for the explicit discretization of the convective terms during the initial time step.
|
protected |
This boolean value determines whether the pressure is normalized to have zero mean (i.e., discrete integral) at the end of each timestep.
|
protected |
This boolean value determines whether the velocity is normalized to have zero mean (i.e., discrete integral) at the end of each timestep.
This parameter only affects the case in which rho=0 (i.e. the steady Stokes equations).
|
protected |
Problem coeficients.
|
protected |
Threshold that determines whether the velocity field needs to be reprojected following adaptive regridding.
|
protected |
Objects to set initial conditions, boundary conditions, body forces, and fluid source/sink distributions.
|
protected |
Double precision values are (optional) factors used to rescale the velocity, pressure, and force for plotting.
Boolean values indicates whether to output various quantities for visualization.
|
protected |
Fluid solver variables.
|
protected |
Cell tagging criteria based on the relative and absolute magnitudes of the local vorticity.
|
protected |
Enum indicating the time integration employed for the implicit discretization of the viscous terms.