|
IBAMR
IBAMR version 0.19.
|
Class IBFESurfaceMethod 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 <ibamr/IBFESurfaceMethod.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 | |
| IBFESurfaceMethod (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, libMesh::MeshBase *mesh, int max_levels, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0) | |
| Constructor. More... | |
| IBFESurfaceMethod (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::vector< libMesh::MeshBase * > &meshes, int max_levels, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0) | |
| Constructor. More... | |
| ~IBFESurfaceMethod () | |
| Destructor. More... | |
| IBTK::FEDataManager * | getFEDataManager (unsigned int part=0) const |
| 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 | 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. More... | |
| 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 | 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) |
| 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 () |
| virtual void | registerLoadBalancer (SAMRAI::tbox::Pointer< SAMRAI::mesh::LoadBalancer< NDIM > > load_balancer, int workload_data_idx) |
| virtual void | initializeLevelData (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double init_data_time, const bool can_be_refined, const bool initial_time, const tbox::Pointer< hier::BasePatchLevel< DIM > > old_level=tbox::Pointer< hier::BasePatchLevel< DIM > >(NULL), const bool allocate_data=true)=0 |
| virtual void | resetHierarchyConfiguration (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int coarsest_level, const int finest_level)=0 |
| virtual void | applyGradientDetector (const tbox::Pointer< hier::BasePatchHierarchy< DIM > > hierarchy, const int level_number, const double error_data_time, const int tag_index, const bool initial_time, const bool uses_richardson_extrapolation_too) |
| virtual double | getLevelDt (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const double dt_time, const bool initial_time) |
| virtual double | advanceLevel (const tbox::Pointer< hier::BasePatchLevel< DIM > > level, const tbox::Pointer< hier::BasePatchHierarchy< DIM > > 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< DIM > > level, const double new_time, const bool can_be_refined) |
| virtual void | resetDataToPreadvanceState (const tbox::Pointer< hier::BasePatchLevel< DIM > > level) |
| virtual void | applyRichardsonExtrapolation (const tbox::Pointer< hier::PatchLevel< DIM > > 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< DIM > > hierarchy, const int level_number, const tbox::Pointer< hier::PatchLevel< DIM > > coarser_level, const double coarsen_data_time, const bool before_advance) |
Static Public Attributes | |
| static const std::string | COORD_MAPPING_SYSTEM_NAME |
| static const std::string | COORDS_SYSTEM_NAME |
| static const std::string | FORCE_SYSTEM_NAME |
| static const std::string | NORMAL_VELOCITY_SYSTEM_NAME |
| static const std::string | PRESSURE_IN_SYSTEM_NAME |
| static const std::string | PRESSURE_JUMP_SYSTEM_NAME |
| static const std::string | PRESSURE_OUT_SYSTEM_NAME |
| static const std::string | TANGENTIAL_VELOCITY_SYSTEM_NAME |
| static const std::string | TAU_IN_SYSTEM_NAME |
| static const std::string | TAU_OUT_SYSTEM_NAME |
| static const std::string | VELOCITY_SYSTEM_NAME |
| static const std::string | WSS_IN_SYSTEM_NAME |
| static const std::string | WSS_OUT_SYSTEM_NAME |
| static const std::array< std::string, NDIM > | VELOCITY_JUMP_SYSTEM_NAME |
Protected Member Functions | |
| void | imposeJumpConditions (const int f_data_idx, libMesh::PetscVector< double > &DP_ghost_vec, libMesh::PetscVector< double > &X_ghost_vec, const double data_time, const unsigned int part) |
| int | getCoarsestPatchLevelNumber () const |
| int | getFinestPatchLevelNumber () const |
| 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. More... | |
| 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. More... | |
| void | updateCoordinateMapping (unsigned int part) |
| Compute dX = x - X, useful mainly for visualization purposes. More... | |
| 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. More... | |
| 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 |
Private Member Functions | |
| IBFESurfaceMethod ()=delete | |
| Default constructor. More... | |
| IBFESurfaceMethod (const IBFESurfaceMethod &from)=delete | |
| Copy constructor. More... | |
| IBFESurfaceMethod & | operator= (const IBFESurfaceMethod &that)=delete |
| Assignment operator. More... | |
| void | commonConstructor (const std::string &object_name, SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > input_db, const std::vector< libMesh::MeshBase * > &meshes, int max_levels, bool register_for_restart, const std::string &restart_read_dirname, unsigned int restart_restore_number) |
| void | getFromInput (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db, bool is_from_restart) |
| void | getFromRestart () |
Private Attributes | |
| SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | d_input_db |
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::IBFESurfaceMethod::CoordinateMappingFcnPtr = void (*)(libMesh::Point& x, const libMesh::Point& X, void* ctx) |
Typedef specifying interface for coordinate mapping function.
| using IBAMR::IBFESurfaceMethod::InitialVelocityFcnPtr = void (*)(libMesh::VectorValue<double>& U0, const libMesh::Point& X0, void* ctx) |
Typedef specifying interface for initial velocity specification function.
Typedef specifying interface for Lagrangian pressure force distribution function.
Typedef specifying interface for Lagrangian surface force distribution function.
| IBAMR::IBFESurfaceMethod::IBFESurfaceMethod | ( | const std::string & | object_name, |
| SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | input_db, | ||
| libMesh::MeshBase * | mesh, | ||
| int | max_levels, | ||
| bool | register_for_restart = true, |
||
| const std::string & | restart_read_dirname = "", |
||
| unsigned int | restart_restore_number = 0 |
||
| ) |
| IBAMR::IBFESurfaceMethod::IBFESurfaceMethod | ( | const std::string & | object_name, |
| SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > | input_db, | ||
| const std::vector< libMesh::MeshBase * > & | meshes, | ||
| int | max_levels, | ||
| bool | register_for_restart = true, |
||
| const std::string & | restart_read_dirname = "", |
||
| unsigned int | restart_restore_number = 0 |
||
| ) |
| IBAMR::IBFESurfaceMethod::~IBFESurfaceMethod | ( | ) |
|
privatedelete |
|
privatedelete |
| from | The value to copy to this object. |
| IBTK::FEDataManager* IBAMR::IBFESurfaceMethod::getFEDataManager | ( | unsigned int | part = 0 | ) | const |
Return a pointer to the finite element data manager object for the specified part.
| void IBAMR::IBFESurfaceMethod::registerInitialCoordinateMappingFunction | ( | const CoordinateMappingFcnData & | data, |
| unsigned int | part = 0 |
||
| ) |
Register the (optional) function used to initialize the physical coordinates from the Lagrangian coordinates.
| void IBAMR::IBFESurfaceMethod::registerInitialVelocityFunction | ( | const InitialVelocityFcnData & | data, |
| unsigned int | part = 0 |
||
| ) |
Register the (optional) function used to initialize the velocity of the solid mesh.
| void IBAMR::IBFESurfaceMethod::registerLagSurfacePressureFunction | ( | const LagSurfacePressureFcnData & | data, |
| unsigned int | part = 0 |
||
| ) |
Register the (optional) function to compute surface pressure distributions on the Lagrangian finite element mesh.
| void IBAMR::IBFESurfaceMethod::registerLagSurfaceForceFunction | ( | const LagSurfaceForceFcnData & | data, |
| unsigned int | part = 0 |
||
| ) |
Register the (optional) function to compute surface force distributions on the Lagrangian finite element mesh.
| const libMesh::VectorValue<double>& IBAMR::IBFESurfaceMethod::getSurfaceForceIntegral | ( | unsigned int | part = 0 | ) | const |
The current value of the integrated surface force.
|
overridevirtual |
Return the number of ghost cells required by the Lagrangian-Eulerian interaction routines.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Setup the tag buffer.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Method to prepare to advance data from current_time to new_time.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Method to clean up data following call(s) to integrateHierarchy().
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 forward Euler method.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Advance the positions of the Lagrangian structure using the (explicit) midpoint rule.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Advance the positions of the Lagrangian structure using the (explicit) trapezoidal rule.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Compute the Lagrangian force at the specified time within the current time interval.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Spread the Lagrangian force to the Cartesian grid at the specified time within the current time interval.
Implements IBAMR::IBStrategy.
| IBTK::FEDataManager::InterpSpec IBAMR::IBFESurfaceMethod::getDefaultInterpSpec | ( | ) | const |
Get the default interpolation spec object used by the class.
| IBTK::FEDataManager::SpreadSpec IBAMR::IBFESurfaceMethod::getDefaultSpreadSpec | ( | ) | const |
Get the default spread spec object used by the class.
| void IBAMR::IBFESurfaceMethod::setInterpSpec | ( | const IBTK::FEDataManager::InterpSpec & | interp_spec, |
| unsigned int | part = 0 |
||
| ) |
Set the interpolation spec object used with a particular mesh part.
| void IBAMR::IBFESurfaceMethod::setSpreadSpec | ( | const IBTK::FEDataManager::SpreadSpec & | spread_spec, |
| unsigned int | part = 0 |
||
| ) |
Set the spread spec object used with a particular mesh part.
| void IBAMR::IBFESurfaceMethod::initializeFEEquationSystems | ( | ) |
Initialize the FE equation systems objects. This method must be called prior to calling initializeFEData().
| void IBAMR::IBFESurfaceMethod::initializeFEData | ( | ) |
Initialize FE data. This method must be called prior to calling IBHierarchyIntegrator::initializePatchHierarchy().
|
overridevirtual |
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 |
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 |
Begin redistributing Lagrangian data prior to regridding the patch hierarchy.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Complete redistributing Lagrangian data following regridding the patch hierarchy.
Reimplemented from IBAMR::IBStrategy.
|
override |
Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm.
|
override |
Reset cached hierarchy dependent data.
|
override |
Set integer tags to "one" in cells where refinement of the given level should occur according to user-supplied feature detection criteria.
|
overridevirtual |
Write out object state to the given database.
Implements SAMRAI::tbox::Serializable.
| void IBAMR::IBFESurfaceMethod::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.
|
protected |
Impose the jump conditions.
|
protected |
Get the coarsest patch level number on which elements (including all parts) are assigned.
|
protected |
Get the finest patch level number on which elements (including all parts) are assigned.
|
protected |
|
protected |
|
protected |
|
protected |
|
privatedelete |
| that | The value to assign to this object. |
|
private |
Implementation of class constructor.
|
private |
Read input values from a given database.
|
private |
Read object state from the restart file and initialize class data members.
|
virtualinherited |
Register the IBHierarchyIntegrator object that is using this strategy class.
Reimplemented in IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Register Eulerian refinement or coarsening algorithms with the parent IBHierarchyIntegrator using the two versions of the protected methods IBStrategy::registerGhostfillRefineAlgorithm(), IBStrategy::registerProlongRefineAlgorithm(), and IBStrategy::registerCoarsenAlgorithm().
An empty default implementation is provided.
Reimplemented in IBAMR::CIBMethod, IBAMR::IBLevelSetMethod, IBAMR::IBInterpolantMethod, IBAMR::IBStrategySet, and IBAMR::GeneralizedIBMethod.
|
virtualinherited |
Inactivate a structure or part. Such a structure will be ignored by FSI calculations: i.e., it will have its velocity set to zero and no forces will be spread from the structure to the Eulerian grid.
| [in] | structure_number | Number of the structure/part. |
| [in] | level_number | Level on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy. |
Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBMethod.
|
virtualinherited |
Activate a previously inactivated structure or part to be used again in FSI calculations.
| [in] | structure_number | Number of the structure/part. |
| [in] | level_number | Level on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy. |
Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBMethod.
|
virtualinherited |
Determine whether or not the given structure or part is currently activated.
| [in] | structure_number | Number of the structure/part. |
| [in] | level_number | Level on which the structure lives. For some inheriting classes (e.g., IBAMR::IBMethod) the structure number alone is not enough to establish uniqueness. The default value is interpreted as the finest level in the patch hierarchy. |
Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBMethod.
|
virtualinherited |
Get the ratio of the maximum point displacement of all the structures owned by the current class to the cell width of the grid level on which the structure is assigned. This value is useful for determining if the Eulerian patch hierarchy needs to be regridded.
Reimplemented in IBAMR::IBFEMethod, and IBAMR::IBStrategySet.
|
inherited |
Indicate whether "fixed" interpolation and spreading operators should be used during Lagrangian-Eulerian interaction.
|
virtualinherited |
Update the positions used for the "fixed" interpolation and spreading operators.
A default implementation is provided that emits an unrecoverable exception.
Reimplemented in IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Indicate that multistep time stepping will be used.
A default implementation is provided that emits an unrecoverable exception.
| [in] | n_previous_steps | Number of previous solution values that can be used by the multistep scheme. |
Reimplemented in IBAMR::IBFEMethod.
|
virtualinherited |
Advance the positions of the Lagrangian structure using the (explicit) backward Euler method.
A default implementation is provided that emits an unrecoverable exception.
Reimplemented in IBAMR::IBFEMethod, IBAMR::CIBMethod, IBAMR::IBMethod, and IBAMR::IBInterpolantMethod.
Advance the positions of the Lagrangian structure using the standard 2nd-order Adams-Bashforth rule.
A default implementation is provided that emits an unrecoverable exception.
Reimplemented in IBAMR::IBFEMethod.
|
virtualinherited |
Indicate whether there are any internal fluid sources/sinks.
A default implementation is provided that returns false.
Reimplemented in IBAMR::IBFEMethod, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Compute the Lagrangian source/sink density at the specified time within the current time interval.
An empty default implementation is provided.
Reimplemented in IBAMR::IBFEMethod, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Spread the Lagrangian source/sink density to the Cartesian grid at the specified time within the current time interval.
An empty default implementation is provided.
Reimplemented in IBAMR::IBFEMethod, IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Compute the pressures at the positions of any distributed internal fluid sources or sinks.
An empty default implementation is provided.
Reimplemented in IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Execute user-defined routines just before solving the fluid equations.
An empty default implementation is provided.
Reimplemented in IBAMR::IBLevelSetMethod, IBAMR::IBStrategySet, IBAMR::CIBMethod, and IBAMR::ConstraintIBMethod.
|
virtualinherited |
Execute user-defined routines just after solving the fluid equations.
An empty default implementation is provided.
Reimplemented in IBAMR::IBLevelSetMethod, IBAMR::IBStrategySet, and IBAMR::ConstraintIBMethod.
|
virtualinherited |
Execute user-defined post-processing operations.
An empty default implementation is provided.
Reimplemented in IBAMR::IBMethod, IBAMR::IBLevelSetMethod, and IBAMR::IBStrategySet.
|
virtualinherited |
Register a load balancer and work load patch data index with the IB strategy object.
An empty default implementation is provided.
Reimplemented in IBAMR::IBMethod, IBAMR::IBStrategySet, and IBAMR::IMPMethod.
|
pure virtualinherited |
Initialize data on a new level after it is inserted into an AMR patch hierarchy by the gridding algorithm. The level number indicates that of the new level.
Generally, when data is set, it is interpolated from coarser levels in the hierarchy. If the old level pointer in the argument list is non-null, then data is copied from the old level to the new level on regions of intersection between those levels before interpolation occurs. In this case, the level number must match that of the old level. The specific operations that occur when initializing level data are determined by the particular solution methods in use; i.e., in the subclass of this abstract base class.
The boolean argument initial_time indicates whether the level is being introduced for the first time (i.e., at initialization time), or after some regrid process during the calculation beyond the initial hierarchy construction. This information is provided since the initialization of the data may be different in each of those circumstances. The can_be_refined boolean argument indicates whether the level is the finest allowable level in the hierarchy.
|
pure virtualinherited |
After hierarchy levels have changed and data has been initialized on the new levels, this routine can be used to reset any information needed by the solution method that is particular to the hierarchy configuration. For example, the solution procedure may cache communication schedules to amortize the cost of data movement on the AMR patch hierarchy. This function will be called by the gridding algorithm after the initialization occurs so that the algorithm-specific subclass can reset such things. Also, if the solution method must make the solution consistent across multiple levels after the hierarchy is changed, this process may be invoked by this routine. Of course the details of these processes are determined by the particular solution methods in use.
The level number arguments indicate the coarsest and finest levels in the current hierarchy configuration that have changed. It should be assumed that all intermediate levels have changed as well.
|
virtualinherited |
Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied gradient criteria. The double time argument is the regrid time. The integer "tag_index" argument is the patch descriptor index of the cell-centered integer tag array on each patch in the hierarchy. The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. Typically, this information is passed to the user's patch tagging routines since the error estimator or gradient detector may be different in each case.
The boolean uses_richardson_extrapolation_too is true when Richardson extrapolation error estimation is used in addition to the gradient detector, and false otherwise. This argument helps the user to manage multiple regridding criteria.
This routine is only when gradient detector is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
protectedinherited |
Return a pointer to the INSHierarchyIntegrator object being used with the IBHierarchyIntegrator class registered with this IBStrategy object.
|
protectedinherited |
Return a pointer to the HierarchyDataOpsReal object associated with velocity-like variables.
|
protectedinherited |
Return a pointer to the HierarchyDataOpsReal object associated with pressure-like variables.
|
protectedinherited |
Return a pointer to a HierarchyMathOps object.
|
protectedinherited |
Register a state variable with the integrator. When a refine operator is specified, the data for the variable are automatically maintained as the patch hierarchy evolves.
All state variables are registered with three contexts: current, new, and scratch. The current context of a state variable is maintained from time step to time step and, if the necessary coarsen and refine operators are specified, as the patch hierarchy evolves.
|
protectedinherited |
Register a variable with the integrator that may not be maintained from time step to time step.
By default, variables are registered with the scratch context, which is deallocated after each time step.
|
protectedinherited |
Register a ghost cell-filling refine algorithm.
|
protectedinherited |
Register a data-prolonging refine algorithm.
|
protectedinherited |
Register a coarsen algorithm.
|
protectedinherited |
Get ghost cell-filling refine algorithm.
|
protectedinherited |
Get data-prolonging refine algorithm.
|
protectedinherited |
Get coarsen algorithm.
|
protectedinherited |
Get ghost cell-filling refine schedules.
|
protectedinherited |
Get data-prolonging refine schedules.
|
protectedinherited |
Get coarsen schedules.
|
virtualinherited |
Determine time increment to advance data on level. The recompute_dt option specifies whether to compute the timestep using the current level data or to return the value stored by the time integrator. The default true setting means the timestep will be computed if no value is supplied.
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Advance data on all patches on specified patch level from current time (current_time) to new time (new_time). This routine is called only during time-dependent regridding procedures, such as Richardson extrapolation. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed. The boolean arguments are used to determine the state of the algorithm and the data when the advance routine is called. Note that this advance function is also used during normal time integration steps.
When this function is called, the level data required to begin the advance must be allocated and be defined appropriately. Typically, this is equivalent to what is needed to initialize a new level after regridding. Upon exiting this routine, both current and new data may exist on the level. This data is needed until level synchronization occurs, in general. Current and new data may be reset by calling the member function resetTimeDependentData().
This routine is called from two different points within the Richardson exptrapolation process: to advance a temporary level that is coarser than the hierarchy level on which error estimation is performed, and to advance the hierarchy level itself. In the first case, the values of the boolean flags are:
In the second case, the values of the boolean flags are:
|
virtualinherited |
Reset time-dependent data storage for the specified patch level.
This routine only applies when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Reset data on the patch level by destroying all patch data other than that which is needed to initialize the solution on that level. In other words, this is the data needed to begin a time integration step on the level.
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Set integer tags to "one" in cells where refinement of the given level should occur according to some user-supplied Richardson extrapolation criteria. The "error_data_time" argument is the regrid time. The "deltat" argument is the time increment to advance the solution on the level to be refined. Note that that level is finer than the level in the argument list, in general. The ratio between the argument level and the actual hierarchy level is given by the integer "coarsen ratio".
The integer "tag_index" argument is the patch descriptor index of the cell-centered integer tag array on each patch in the hierarchy.
The boolean argument initial_time indicates whether the level is being subject to refinement at the initial simulation time. If it is false, then the error estimation process is being invoked at some later time after the AMR hierarchy was initially constructed. Typically, this information is passed to the user's patch tagging routines since the application of the Richardson extrapolation process may be different in each case.
The boolean uses_gradient_detector_too is true when a gradient detector procedure is used in addition to Richardson extrapolation, and false otherwise. This argument helps the user to manage multiple regridding criteria.
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
virtualinherited |
Coarsen solution data from level to coarse_level for Richardson extrapolation. Note that this routine will be called twice during the Richardson extrapolation error estimation process, once to set data on the coarser level and once to coarsen data from after advancing the fine level. The init_coarse_level boolean argument indicates whether data is set on the coarse level by coarsening the "old" time level solution or by coarsening the "new" solution on the fine level (i.e., after it has been advanced).
This routine is only when Richardson extrapolation is being used. It is virtual with an empty implementation here (rather than pure virtual) so that users are not required to provide an implementation when the function is not needed.
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
EquationSystems objects, one per part. These contain the actual matrices and solution vectors for each relevant libMesh system.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Object managing access to libMesh system vectors for the position system.
|
protected |
Object managing access to libMesh system vectors for the velocity system.
|
protected |
Object managing access to libMesh system vectors for the normal velocity system.
|
protected |
Object managing access to libMesh system vectors for the tangential velocity system.
|
protected |
Object managing access to libMesh system vectors for the force system.
|
protected |
Object managing access to libMesh system vectors for the pressure jump system.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
private |
The input database. This is explicitly stored (and used outside the constructor) since the FEDataManager instances created by this class will also read part of it.
|
protectedinherited |
The IBHierarchyIntegrator object that is using this strategy class.
|
protectedinherited |
Whether to use "fixed" Lagrangian-Eulerian coupling operators.
1.8.17