IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
|
Class IIMethod is an implementation of the abstract base class IBStrategy that provides functionality required by the IB method with a finite element representation of a surface mesh. More...
#include </home/runner/work/IBAMR/IBAMR/include/ibamr/IIMethod.h>
Classes | |
struct | CoordinateMappingFcnData |
struct | InitialVelocityFcnData |
struct | LagSurfaceForceFcnData |
struct | LagSurfacePressureFcnData |
Public Types | |
using | CoordinateMappingFcnPtr = void(*)(libMesh::Point &x, const libMesh::Point &X, void *ctx) |
using | InitialVelocityFcnPtr = void(*)(libMesh::VectorValue< double > &U0, const libMesh::Point &X0, void *ctx) |
using | LagSurfacePressureFcnPtr = IBTK::ScalarSurfaceFcnPtr |
using | LagSurfaceForceFcnPtr = IBTK::VectorSurfaceFcnPtr |
Public Member Functions | |
IIMethod (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, libMesh::MeshBase *mesh, int max_level_number, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0) | |
Constructor. | |
IIMethod (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::vector< libMesh::MeshBase * > &meshes, int max_level_number, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0) | |
Constructor. | |
~IIMethod () | |
Destructor. | |
IBTK::FEDataManager * | getFEDataManager (unsigned int part=0) const |
void | registerDisconElemFamilyForPressureJump (unsigned int part, libMesh::FEFamily fe_family, libMesh::Order fe_order) |
void | registerDisconElemFamilyForViscousJump (unsigned int part, libMesh::FEFamily fe_family, libMesh::Order fe_order) |
void | registerDisconElemFamilyForTraction (unsigned int part, libMesh::FEFamily fe_family, libMesh::Order fe_order) |
void | registerTangentialVelocityMotion (unsigned int part=0) |
void | registerPressureJumpNormalization (unsigned int part=0) |
void | registerInitialCoordinateMappingFunction (const CoordinateMappingFcnData &data, unsigned int part=0) |
void | registerInitialVelocityFunction (const InitialVelocityFcnData &data, unsigned int part=0) |
void | registerLagSurfacePressureFunction (const LagSurfacePressureFcnData &data, unsigned int part=0) |
void | registerLagSurfaceForceFunction (const LagSurfaceForceFcnData &data, unsigned int part=0) |
const libMesh::VectorValue< double > & | getSurfaceForceIntegral (unsigned int part=0) const |
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 |
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 | 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 | computeFluidTraction (double current_time, unsigned int part=0) |
void | extrapolatePressureForTraction (int p_data_idx, double data_time, unsigned int part=0) |
void | calculateInterfacialFluidForces (int p_data_idx, double data_time) |
void | forwardEulerStep (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 | 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 |
IBTK::FEDataManager::InterpSpec | getDefaultInterpSpec () const |
IBTK::FEDataManager::SpreadSpec | getDefaultSpreadSpec () const |
void | setInterpSpec (const IBTK::FEDataManager::InterpSpec &interp_spec, unsigned int part=0) |
void | setSpreadSpec (const IBTK::FEDataManager::SpreadSpec &spread_spec, unsigned int part=0) |
void | initializeFEEquationSystems () |
void | initializeFEData () |
void | registerEulerianVariables () override |
Register Eulerian variables with the parent IBHierarchyIntegrator. | |
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 |
void | writeFEDataToRestartFile (const std::string &restart_dump_dirname, unsigned int time_step_number) |
Public Member Functions inherited from IBAMR::IBStrategy | |
IBStrategy ()=default | |
Constructor. | |
virtual void | registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver) |
virtual void | registerEulerianCommunicationAlgorithms () |
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 |
void | setUseFixedLEOperators (bool use_fixed_coupling_ops=true) |
virtual void | updateFixedLEOperators () |
virtual void | setUseMultistepTimeStepping (unsigned int n_previous_steps=1) |
virtual void | backwardEulerStep (double current_time, double new_time) |
virtual void | AB2Step (double current_time, double new_time) |
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 () |
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) |
Static Public Attributes | |
static const std::string | COORD_MAPPING_SYSTEM_NAME = "coordinate mapping system" |
static const std::string | COORDS_SYSTEM_NAME = "coordinates system" |
static const std::string | FORCE_SYSTEM_NAME = "IB force system" |
static const std::string | NORMAL_VELOCITY_SYSTEM_NAME = "normal velocity system" |
static const std::string | PRESSURE_IN_SYSTEM_NAME = "One sided interior pressure system" |
static const std::string | PRESSURE_JUMP_SYSTEM_NAME = "[[p]] system" |
static const std::string | PRESSURE_OUT_SYSTEM_NAME = "One sided exterior pressure system" |
static const std::string | TANGENTIAL_VELOCITY_SYSTEM_NAME = "tangential velocity system" |
static const std::string | TAU_IN_SYSTEM_NAME = "Interior traction system" |
static const std::string | TAU_OUT_SYSTEM_NAME = "Exterior traction system" |
static const std::string | VELOCITY_SYSTEM_NAME = "velocity system" |
static const std::string | WSS_IN_SYSTEM_NAME = "One sided interior wall shear stress system" |
static const std::string | WSS_OUT_SYSTEM_NAME = "One sided exterior wall shear stress system" |
static const std::array< std::string, NDIM > | VELOCITY_JUMP_SYSTEM_NAME |
Protected Member Functions | |
void | imposeJumpConditions (const int f_data_idx, libMesh::PetscVector< double > &P_jump_ghost_vec, std::array< libMesh::PetscVector< double > *, NDIM > &DU_jump_ghost_vec, libMesh::PetscVector< double > &X_ghost_vec, const double data_time, const unsigned int part) |
bool | checkDoubleCountingIntersection (int axis, const double *dx, const libMesh::VectorValue< double > &n, const libMesh::Point &x, const libMesh::Point &xi, const SAMRAI::pdat::SideIndex< NDIM > &i_s, const SAMRAI::pdat::SideIndex< NDIM > &i_s_prime, const std::vector< libMesh::Point > &candidate_coords, const std::vector< libMesh::Point > &candidate_ref_coords, const std::vector< libMesh::VectorValue< double > > &candidate_normals) |
Helper function for checking possible double-counting intesection points. | |
void | initializeCoordinates (unsigned int part) |
Initialize the physical coordinates using the supplied coordinate mapping function. If no function is provided, the initial coordinates are taken to be the Lagrangian coordinates. | |
void | updateCoordinateMapping (unsigned int part) |
Compute dX = x - X, useful mainly for visualization purposes. | |
void | initializeVelocity (unsigned int part) |
Initialize the velocity field using the supplied initial velocity specification function. If no function is provided, the initial velocity is taken to be zero. | |
Protected Member Functions inherited from IBAMR::IBStrategy | |
INSHierarchyIntegrator * | getINSHierarchyIntegrator () 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::HierarchyMathOps > | getHierarchyMathOps () const |
void | registerVariable (int ¤t_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=nullptr, const bool register_for_restart=true) |
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 >(nullptr), const bool register_for_restart=true) |
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 |
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 |
bool | d_is_initialized = false |
std::shared_ptr< IBTK::SAMRAIDataCache > | d_eulerian_data_cache |
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() |
std::vector< libMesh::MeshBase * > | d_meshes |
int | d_max_level_number |
std::vector< libMesh::EquationSystems * > | d_equation_systems |
std::vector< bool > | d_use_discon_elem_for_jumps = { false } |
std::vector< bool > | d_use_tangential_velocity = { false } |
std::vector< bool > | d_normalize_pressure_jump = { false } |
const unsigned int | d_num_parts = 1 |
std::vector< IBTK::FEDataManager * > | d_fe_data_managers |
SAMRAI::hier::IntVector< NDIM > | d_ghosts = 0 |
std::vector< libMesh::System * > | d_X_systems |
std::vector< libMesh::System * > | d_U_systems |
std::vector< libMesh::System * > | d_U_n_systems |
std::vector< libMesh::System * > | d_U_t_systems |
std::vector< libMesh::System * > | d_F_systems |
std::vector< libMesh::System * > | d_P_jump_systems |
std::vector< libMesh::System * > | d_WSS_in_systems |
std::vector< libMesh::System * > | d_WSS_out_systems |
std::vector< libMesh::System * > | d_P_in_systems |
std::vector< libMesh::System * > | d_P_out_systems |
std::vector< libMesh::System * > | d_TAU_in_systems |
std::vector< libMesh::System * > | d_TAU_out_systems |
std::vector< std::array< libMesh::System *, NDIM > > | d_DU_jump_systems |
std::vector< libMesh::PetscVector< double > * > | d_F_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_F_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_X_current_vecs |
std::vector< libMesh::PetscVector< double > * > | d_X_new_vecs |
std::vector< libMesh::PetscVector< double > * > | d_X_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_X0_vecs |
std::vector< libMesh::PetscVector< double > * > | d_X_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_current_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_new_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_n_current_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_n_new_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_n_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_t_current_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_t_new_vecs |
std::vector< libMesh::PetscVector< double > * > | d_U_t_half_vecs |
std::vector< std::array< libMesh::PetscVector< double > *, NDIM > > | d_DU_jump_half_vecs |
std::vector< std::array< libMesh::PetscVector< double > *, NDIM > > | d_DU_jump_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_P_jump_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_P_jump_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_P_in_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_P_in_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_P_out_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_P_out_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_WSS_in_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_WSS_in_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_WSS_out_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_WSS_out_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_TAU_in_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_TAU_in_IB_ghost_vecs |
std::vector< libMesh::PetscVector< double > * > | d_TAU_out_half_vecs |
std::vector< libMesh::PetscVector< double > * > | d_TAU_out_IB_ghost_vecs |
bool | d_fe_equation_systems_initialized = false |
bool | d_fe_data_initialized = false |
IBTK::FEDataManager::InterpSpec | d_default_interp_spec |
IBTK::FEDataManager::SpreadSpec | d_default_spread_spec |
IBTK::FEDataManager::WorkloadSpec | d_default_workload_spec |
std::vector< IBTK::FEDataManager::InterpSpec > | d_interp_spec |
std::vector< IBTK::FEDataManager::SpreadSpec > | d_spread_spec |
bool | d_use_pressure_jump_conditions = false |
bool | d_use_velocity_jump_conditions = false |
bool | d_use_u_interp_correction = false |
bool | d_compute_fluid_traction = false |
bool | d_perturb_fe_mesh_nodes = true |
std::vector< libMesh::FEFamily > | d_fe_family |
std::vector< libMesh::FEFamily > | d_viscous_jump_fe_family |
std::vector< libMesh::FEFamily > | d_pressure_jump_fe_family |
std::vector< libMesh::FEFamily > | d_traction_fe_family |
std::vector< libMesh::Order > | d_fe_order |
std::vector< libMesh::Order > | d_viscous_jump_fe_order |
std::vector< libMesh::Order > | d_pressure_jump_fe_order |
std::vector< libMesh::Order > | d_traction_fe_order |
std::vector< libMesh::QuadratureType > | d_default_quad_type |
std::vector< libMesh::Order > | d_default_quad_order |
bool | d_use_consistent_mass_matrix = true |
bool | d_use_direct_forcing = false |
double | d_exterior_calc_coef = 1.0 |
double | d_wss_calc_width = 1.05 |
double | d_p_calc_width = 1.3 |
double | d_fuzzy_tol = 1e-5 |
double | d_mesh_perturb_tol = 1e-4 |
std::vector< CoordinateMappingFcnData > | d_coordinate_mapping_fcn_data |
std::vector< InitialVelocityFcnData > | d_initial_velocity_fcn_data |
std::vector< LagSurfacePressureFcnData > | d_lag_surface_pressure_fcn_data |
std::vector< LagSurfaceForceFcnData > | d_lag_surface_force_fcn_data |
std::vector< libMesh::VectorValue< double > > | d_lag_surface_force_integral |
SAMRAI::tbox::Pointer< SAMRAI::hier::Variable< NDIM > > | d_p_var |
int | d_p_scratch_idx = IBTK::invalid_index |
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 |
std::string | d_libmesh_restart_read_dir |
int | d_libmesh_restart_restore_number |
std::string | d_libmesh_restart_file_extension = "xdr" |
Protected Attributes inherited from IBAMR::IBStrategy | |
IBHierarchyIntegrator * | d_ib_solver = nullptr |
bool | d_use_fixed_coupling_ops = false |
Class IIMethod is an implementation of the abstract base class IBStrategy that provides functionality required by the IB method with a finite element representation of a surface mesh.
Coupling schemes include both IB formulations (integral operations with regularized delta function kernels) and an immersed interface method (IIM) scheme (E. M. Kolahdouz, A. P. S. Bhalla, B. A. Craven, and B. E. Griffith. An immersed interface method for discrete surfaces. J Comput Phys, 400:108854 (37 pages), 2020).
using IBAMR::IIMethod::CoordinateMappingFcnPtr = void (*)(libMesh::Point& x, const libMesh::Point& X, void* ctx) |
Typedef specifying interface for coordinate mapping function.
using IBAMR::IIMethod::InitialVelocityFcnPtr = void (*)(libMesh::VectorValue<double>& U0, const libMesh::Point& X0, void* ctx) |
Typedef specifying interface for initial velocity specification function.
using IBAMR::IIMethod::LagSurfaceForceFcnPtr = IBTK::VectorSurfaceFcnPtr |
Typedef specifying interface for Lagrangian surface force distribution function.
using IBAMR::IIMethod::LagSurfacePressureFcnPtr = IBTK::ScalarSurfaceFcnPtr |
Typedef specifying interface for Lagrangian pressure force distribution function.
|
overridevirtual |
Add the estimated computational work from the current object (i.e., the work required by the owned Lagrangian objects) per cell into the specified workload_data_idx
.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Set integer tags to "one" in cells where refinement of the given level should occur according to user-supplied feature detection criteria.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Begin redistributing Lagrangian data prior to regridding the patch hierarchy.
Reimplemented from IBAMR::IBStrategy.
void IBAMR::IIMethod::calculateInterfacialFluidForces | ( | int | p_data_idx, |
double | data_time | ||
) |
A wrapper to compute the interfacial pressure and fluid traction from the jumps.
void IBAMR::IIMethod::computeFluidTraction | ( | double | current_time, |
unsigned int | part = 0 |
||
) |
Compute the fluid traction if all jump conditions are applied.
|
overridevirtual |
Compute the Lagrangian force at the specified time within the current time interval.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Complete redistributing Lagrangian data following regridding the patch hierarchy.
Reimplemented from IBAMR::IBStrategy.
void IBAMR::IIMethod::extrapolatePressureForTraction | ( | int | p_data_idx, |
double | data_time, | ||
unsigned int | part = 0 |
||
) |
Compute the interior/exterior pressure used in the calculation of the fluid traction (pressure jump condition needs to be applied).
|
overridevirtual |
Advance the positions of the Lagrangian structure using the forward Euler method.
Implements IBAMR::IBStrategy.
FEDataManager::InterpSpec IBAMR::IIMethod::getDefaultInterpSpec | ( | ) | const |
Get the default interpolation spec object used by the class.
FEDataManager::SpreadSpec IBAMR::IIMethod::getDefaultSpreadSpec | ( | ) | const |
Get the default spread spec object used by the class.
FEDataManager * IBAMR::IIMethod::getFEDataManager | ( | unsigned int | part = 0 | ) | const |
Return a pointer to the finite element data manager object for the specified part.
|
overridevirtual |
Return the number of ghost cells required by the Lagrangian-Eulerian interaction routines.
Implements IBAMR::IBStrategy.
const VectorValue< double > & IBAMR::IIMethod::getSurfaceForceIntegral | ( | unsigned int | part = 0 | ) | const |
The current value of the integrated surface force.
|
protected |
Impose the jump conditions.
void IBAMR::IIMethod::initializeFEData | ( | ) |
Initialize FE data. This method must be called prior to calling IBHierarchyIntegrator::initializePatchHierarchy().
void IBAMR::IIMethod::initializeFEEquationSystems | ( | ) |
Initialize the FE equation systems objects. This method must be called prior to calling initializeFEData().
|
overridevirtual |
Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm.
Reimplemented from IBAMR::IBStrategy.
|
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.
|
overridevirtual |
Interpolate the Eulerian velocity to the curvilinear mesh at the specified time within the current time interval.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Advance the positions of the Lagrangian structure using the (explicit) midpoint rule.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Method to clean up data following call(s) to integrateHierarchy().
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Method to prepare to advance data from current_time to new_time.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Write out object state to the given database.
Reimplemented from IBAMR::IBStrategy.
void IBAMR::IIMethod::registerDisconElemFamilyForPressureJump | ( | unsigned int | part, |
libMesh::FEFamily | fe_family, | ||
libMesh::Order | fe_order | ||
) |
Register relevant part to use discontinuous element type family for the calculation of jumps plus traction quantities. This option should be used for geometries with sharp corners. The acceptable options are CONSTANT MONOMIAL, FIRST order MONOMIAL, or FIRST order L2_LAGRANGE.
void IBAMR::IIMethod::registerDisconElemFamilyForTraction | ( | unsigned int | part, |
libMesh::FEFamily | fe_family, | ||
libMesh::Order | fe_order | ||
) |
Register relevant part to use discontinuous element type family for the calculation of jumps plus traction quantities. This option should be used for geometries with sharp corners. The acceptable options are CONSTANT MONOMIAL, FIRST order MONOMIAL, or FIRST order L2_LAGRANGE.
void IBAMR::IIMethod::registerDisconElemFamilyForViscousJump | ( | unsigned int | part, |
libMesh::FEFamily | fe_family, | ||
libMesh::Order | fe_order | ||
) |
Register relevant part to use discontinuous element type family for the calculation of jumps plus traction quantities. This option should be used for geometries with sharp corners. The acceptable options are CONSTANT MONOMIAL, FIRST order MONOMIAL, or FIRST order L2_LAGRANGE.
void IBAMR::IIMethod::registerInitialCoordinateMappingFunction | ( | const CoordinateMappingFcnData & | data, |
unsigned int | part = 0 |
||
) |
Register the (optional) function used to initialize the physical coordinates from the Lagrangian coordinates.
void IBAMR::IIMethod::registerInitialVelocityFunction | ( | const InitialVelocityFcnData & | data, |
unsigned int | part = 0 |
||
) |
Register the (optional) function used to initialize the velocity of the solid mesh.
void IBAMR::IIMethod::registerLagSurfaceForceFunction | ( | const LagSurfaceForceFcnData & | data, |
unsigned int | part = 0 |
||
) |
Register the (optional) function to compute surface force distributions on the Lagrangian finite element mesh.
void IBAMR::IIMethod::registerLagSurfacePressureFunction | ( | const LagSurfacePressureFcnData & | data, |
unsigned int | part = 0 |
||
) |
Register the (optional) function to compute surface pressure distributions on the Lagrangian finite element mesh.
|
overridevirtual |
Register a load balancer and work load patch data index with the IB strategy object.
Reimplemented from IBAMR::IBStrategy.
void IBAMR::IIMethod::registerPressureJumpNormalization | ( | unsigned int | part = 0 | ) |
Register relevant part for which we would like to normalize the pressure jump.
void IBAMR::IIMethod::registerTangentialVelocityMotion | ( | unsigned int | part = 0 | ) |
Register relevant part to only use the tangential component of the Lagrangian velocity for displacement. The lack of normal velocity needs to be compensated by a large frictional penalty force (direct forcing for the normal velocity)
|
overridevirtual |
Reset cached hierarchy dependent data.
Reimplemented from IBAMR::IBStrategy.
void IBAMR::IIMethod::setInterpSpec | ( | const IBTK::FEDataManager::InterpSpec & | interp_spec, |
unsigned int | part = 0 |
||
) |
Set the interpolation spec object used with a particular mesh part.
void IBAMR::IIMethod::setSpreadSpec | ( | const IBTK::FEDataManager::SpreadSpec & | spread_spec, |
unsigned int | part = 0 |
||
) |
Set the spread spec object used with a particular mesh part.
|
overridevirtual |
Setup the tag buffer.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Spread the Lagrangian force to the Cartesian grid at the specified time within the current time interval.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Advance the positions of the Lagrangian structure using the (explicit) trapezoidal rule.
Implements IBAMR::IBStrategy.
void IBAMR::IIMethod::writeFEDataToRestartFile | ( | const std::string & | restart_dump_dirname, |
unsigned int | time_step_number | ||
) |
Write the equation_systems data to a restart file in the specified directory.
|
static |