|
IBAMR
IBAMR version 0.19.
|
Class IBFEMethod is an implementation of the abstract base class IBStrategy that provides functionality required by the IB method with finite element elasticity. Much of the setup for finite element computations is done by the IBAMR::FEMechanicsBase: see the documentation of that class for additional information on input parameters. More...
#include <ibamr/IBFEMethod.h>

Classes | |
| struct | LagBodySourceFcnData |
Public Types | |
| using | LagBodySourceFcnPtr = IBTK::ScalarMeshFcnPtr |
| 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 | PK1StressFcnPtr = IBTK::TensorMeshFcnPtr |
| using | LagBodyForceFcnPtr = IBTK::VectorMeshFcnPtr |
| using | LagSurfacePressureFcnPtr = IBTK::ScalarSurfaceFcnPtr |
| using | LagSurfaceForceFcnPtr = IBTK::VectorSurfaceFcnPtr |
| using | VolumetricEnergyDerivativeFcn = double(*)(double) |
Public Member Functions | |
| IBFEMethod (const std::string &object_name, const 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 for a single-part model. More... | |
| IBFEMethod (const std::string &object_name, const 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 for a multi-part model. More... | |
| IBFEMethod ()=delete | |
| Deleted default constructor. More... | |
| IBFEMethod (const IBFEMethod &from)=delete | |
| Deleted copy constructor. More... | |
| IBFEMethod & | operator= (const IBFEMethod &that)=delete |
| Deleted assignment operator. More... | |
| ~IBFEMethod () override=default | |
| Defaulted destructor. More... | |
| const std::string & | getSourceSystemName () const |
| IBTK::FEDataManager * | getFEDataManager (unsigned int part=0) const |
| virtual void | registerStaticPressurePart (PressureProjectionType projection_type=CONSISTENT_PROJECTION, FEMechanicsBase::VolumetricEnergyDerivativeFcn U_prime_fcn=nullptr, unsigned int part=0) override |
| virtual void | registerStressNormalizationPart (unsigned int part=0) |
| virtual void | registerLagBodySourceFunction (const LagBodySourceFcnData &data, unsigned int part=0) |
| LagBodySourceFcnData | getLagBodySourceFunction (unsigned int part=0) const |
| virtual void | registerDirectForcingKinematics (const SAMRAI::tbox::Pointer< IBAMR::IBFEDirectForcingKinematics > &data, unsigned int part=0) |
| const SAMRAI::hier::IntVector< NDIM > & | getMinimumGhostCellWidth () const override |
| void | setupTagBuffer (SAMRAI::tbox::Array< int > &tag_buffer, SAMRAI::tbox::Pointer< SAMRAI::mesh::GriddingAlgorithm< NDIM > > gridding_alg) const override |
| virtual double | getMaxPointDisplacement () const override |
| virtual void | inactivateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) override |
| virtual void | activateLagrangianStructure (int structure_number=0, int level_number=std::numeric_limits< int >::max()) override |
| virtual bool | getLagrangianStructureIsActivated (int structure_number=0, int level_number=std::numeric_limits< int >::max()) const override |
| void | preprocessIntegrateData (double current_time, double new_time, int num_cycles) override |
| void | postprocessIntegrateData (double current_time, double new_time, int num_cycles) override |
| void | 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 | setUseMultistepTimeStepping (unsigned int n_previous_steps=1) override |
| void | forwardEulerStep (double current_time, double new_time) override |
| void | backwardEulerStep (double current_time, double new_time) override |
| void | midpointStep (double current_time, double new_time) override |
| void | trapezoidalStep (double current_time, double new_time) override |
| void | AB2Step (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 |
| bool | hasFluidSources () const override |
| void | computeLagrangianFluidSource (double data_time) override |
| void | spreadFluidSource (int q_data_idx, IBTK::RobinPhysBdryPatchStrategy *q_phys_bdry_op, const std::vector< SAMRAI::tbox::Pointer< SAMRAI::xfer::RefineSchedule< NDIM > > > &q_prolongation_scheds, double data_time) override |
| IBTK::FEDataManager::InterpSpec | getDefaultInterpSpec () const |
| IBTK::FEDataManager::SpreadSpec | getDefaultSpreadSpec () const |
| void | setWorkloadSpec (const IBTK::FEDataManager::WorkloadSpec &workload_spec, unsigned int part=0) |
| 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 | 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, 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 |
| SAMRAI::tbox::Pointer< SAMRAI::hier::BasePatchHierarchy< NDIM > > | getScratchHierarchy () |
| libMesh::EquationSystems * | getEquationSystems (unsigned int part=0) const |
| const std::string & | getCurrentCoordinatesSystemName () const |
| const std::string & | getDisplacementSystemName () const |
| const std::string & | getForceSystemName () const |
| const std::string & | getPressureSystemName () const |
| const std::string & | getVelocitySystemName () const |
| std::shared_ptr< IBTK::FEData > | getFEData (unsigned int part=0) const |
| virtual void | registerInitialCoordinateMappingFunction (const CoordinateMappingFcnData &data, unsigned int part=0) |
| CoordinateMappingFcnData | getInitialCoordinateMappingFunction (unsigned int part=0) const |
| virtual void | registerInitialVelocityFunction (const InitialVelocityFcnData &data, unsigned int part=0) |
| InitialVelocityFcnData | getInitialVelocityFunction (unsigned int part=0) const |
| virtual void | registerPK1StressFunction (const PK1StressFcnData &data, unsigned int part=0) |
| std::vector< PK1StressFcnData > | getPK1StressFunction (unsigned int part=0) const |
| virtual void | registerLagBodyForceFunction (const LagBodyForceFcnData &data, unsigned int part=0) |
| LagBodyForceFcnData | getLagBodyForceFunction (unsigned int part=0) const |
| virtual void | registerLagSurfacePressureFunction (const LagSurfacePressureFcnData &data, unsigned int part=0) |
| LagSurfacePressureFcnData | getLagSurfacePressureFunction (unsigned int part=0) const |
| virtual void | registerLagSurfaceForceFunction (const LagSurfaceForceFcnData &data, unsigned int part=0) |
| LagSurfaceForceFcnData | getLagSurfaceForceFunction (unsigned int part=0) const |
| virtual void | registerDynamicPressurePart (PressureProjectionType projection_type=CONSISTENT_PROJECTION, VolumetricEnergyDerivativeFcn d2U_dJ2_fcn=nullptr, unsigned int part=0) |
| bool | partHasPressure (unsigned int part=0) |
| virtual void | initializeFEEquationSystems () |
| virtual void | initializeFEData () |
| virtual void | reinitializeFEData () |
| virtual void | writeFEDataToRestartFile (const std::string &restart_dump_dirname, unsigned int time_step_number) |
| virtual void | registerIBHierarchyIntegrator (IBHierarchyIntegrator *ib_solver) |
| virtual void | registerEulerianCommunicationAlgorithms () |
| void | setUseFixedLEOperators (bool use_fixed_coupling_ops=true) |
| virtual void | updateFixedLEOperators () |
| 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 | SOURCE_SYSTEM_NAME |
| static const std::string | COORDS_SYSTEM_NAME |
| static const std::string | COORD_MAPPING_SYSTEM_NAME |
| static const std::string | FORCE_SYSTEM_NAME |
| static const std::string | PRESSURE_SYSTEM_NAME |
| static const std::string | VELOCITY_SYSTEM_NAME |
Protected Member Functions | |
| virtual void | doInitializeFEEquationSystems () override |
| void | doInitializeFESystemVectors () override |
| virtual void | doInitializeFEData (bool use_present_data) override |
| void | computeStressNormalization (libMesh::PetscVector< double > &P_vec, libMesh::PetscVector< double > &X_vec, double data_time, unsigned int part) |
| Compute the stress normalization field. More... | |
| void | spreadTransmissionForceDensity (int f_data_idx, libMesh::PetscVector< double > &X_ghost_vec, double data_time, unsigned int part) |
| Spread the transmission force density along the physical boundary of the Lagrangian structure. More... | |
| void | imposeJumpConditions (int f_data_idx, libMesh::PetscVector< double > &F_ghost_vec, libMesh::PetscVector< double > &X_ghost_vec, double data_time, unsigned int part) |
| Impose jump conditions determined from the interior and transmission force densities along the physical boundary of the Lagrangian structure. More... | |
| int | getCoarsestPatchLevelNumber () const |
| int | getFinestPatchLevelNumber () const |
| SAMRAI::xfer::RefineSchedule< NDIM > & | getProlongationSchedule (int level_number, int coarse_data_idx, int fine_data_idx) |
| void | computeStaticPressure (libMesh::PetscVector< double > &P_vec, libMesh::PetscVector< double > &X_vec, double data_time, unsigned int part) |
| Compute the static pressure field. More... | |
| void | computeDynamicPressureRateOfChange (libMesh::PetscVector< double > &dP_dt_vec, libMesh::PetscVector< double > &X_vec, libMesh::PetscVector< double > &U_vec, double data_time, unsigned int part) |
| Compute the dynamic pressure time rate of change. More... | |
| virtual void | assembleInteriorForceDensityRHS (libMesh::PetscVector< double > &F_rhs_vec, libMesh::PetscVector< double > &X_vec, libMesh::PetscVector< double > *P_vec, double data_time, unsigned int part) |
| virtual void | initializeCoordinates (unsigned int part) |
| virtual void | updateCoordinateMapping (unsigned int part) |
| virtual void | initializeVelocity (unsigned int part) |
| virtual std::string | getLibMeshRestartFileName (const std::string &restart_dump_dirname, unsigned int time_step_number, unsigned int part, const std::string &extension) const |
| 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 |
Static Protected Member Functions | |
| static void | setup_system_vectors (libMesh::EquationSystems *equation_systems, const std::vector< std::string > &system_names, const std::vector< std::string > &vector_names) |
| static void | setup_system_vector (libMesh::System &system, const std::string &vector_name) |
Private Member Functions | |
| void | commonConstructor (const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, int max_levels) |
| void | getFromInput (const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &db, bool is_from_restart) |
| void | getFromRestart () |
| void | assertStructureOnFinestLevel () const |
| void | reinitElementMappings () |
By default, the libMesh data is partitioned once at the beginning of the computation by libMesh's default partitioner.
Like other classes inheriting from IBStrategy, most options regarding the actual IB method implementation can be specified with the provided input database. Parameters starting with IB_ set and override those with the same name starting with interp_ or spread_: e.g., IB_delta_fcn overrides both interp_delta_fcn and spread_delta_fcn.
interp_quad_type: Quadrature type for interpolation, provided as a string. Can be any quadrature type known to libMesh. Defaults to "QGAUSS". spread_quad_type: Quadrature type for spreading, provided as a string. Parsed in the same was as interp_quad_type. IB_quad_type: overriding alias for the two previous entries - has the same default. interp_use_adaptive_quadrature: Whether or not the current deformation of each element should be considered when determining which quadrature rule to use. Defaults to TRUE. interp_point_density: Multiplier on the number of points to use for doing IB calculations. Defaults to 2. spread_point_density: Same as above, but for spreading. IB_point_density: overriding alias for the previous entry - has the same default. interp_point_density: Parameter for adaptively computing the number of quadrature points in a quadrature rule. Defaults to 2.0. See IBTK::getQuadratureKey() for a detailed description. spread_point_density: Same as above, but for spreading. IB_point_density: overriding alias for the two previous entries - has the same default. interp_use_consistent_mass_matrix: Whether or not mass lumping should be applied when solving the L2 projection for computing the velocity of the structure. Defaults to FALSE. Note that no linear system is solved when computing forces so this parameter does not have a spreading equivalent. use_consistent_mass_matrix: Overriding alias of the previous entry. IB_use_consistent_mass_matrix: Overriding alias of the previous entry. interp_use_nodal_quadrature: Whether or not nodal quadrature should be used, which is essentially interpolation instead of projection. This is an experimental feature. Defaults to FALSE. spread_use_nodal_quadrature: Same as above, but for spreading. IB_use_nodal_quadrature: overriding alias for the two previous entries - has the same default. This feature is experimental: at the present time the default settings have the best performance and are the correct choice.
This class can repartition libMesh data in a way that matches SAMRAI's distribution of patches; put another way, if a certain region of space on the finest level is assigned to processor N, then all libMesh nodes and elements within that region will also be assigned to processor N. The actual partitioning here is done by the IBTK::BoxPartitioner class. See the discussion in IBTK::HierarchyIntegrator and IBTK::FEDataManager for descriptions on how this partitioning is performed.
The choice of libMesh partitioner depends on the libmesh_partitioner_type parameter in the input database:
If libmesh_partitioner_type is LIBMESH_DEFAULT then this class will never repartition libMesh data, since the default libMesh partitioner is already used at the beginning of the computation and, since no degrees of freedom are added or removed, any subsequent partitioning would have no effect.
libmesh_partitioner_type is SAMRAI_BOX then this class will always repartition the libMesh data with IBTK::BoxPartitioner every time the Eulerian data is regridded. The default value for libmesh_partitioner_type is LIBMESH_DEFAULT. The intent of these choices is to automatically use the fairest (that is, partitioning based on equal work when computing force densities and L2 projections) partitioner.
The main computational expenses of this class are IBFEMethod::interpolateVelocity() and IBFEMethod::spreadForce(). These two methods compute at IB points placed inside the patches owned on the current processor: i.e., they use the Eulerian partitioning of the domain. This partitioning scales very poorly at higher processor counts with some Lagrangian geometries since the Eulerian partitioning places equal number of cells, which do not necessarily coincide with IB points, on different processors: i.e., some processors will have a large number of IB points and some may have zero.
To get around this, this class can optionally work with a different partitioning of the Eulerian data that is partitioned so that each processor has roughly the same number of IB points, or some more elaborate partitioning scheme that takes into account the number of mesh nodes as well. This class will set up this scratch hierarchy and manage its state (see IBFEMethod::d_scratch_hierarchy). The scratch hierarchy can be set up by adding the following parameters to the input database:
i.e., providing use_scratch_hierarchy = TRUE (the default is FALSE) turns on the scratch hierarchy and the remaining parameters determine how patches are generated and load balanced. The extra argument type to LoadBalancer specifies whether an IBTK::MergingLoadBalancer (chosen by "MERGING") or the default SAMRAI LoadBalancer (chosen by "DEFAULT") is used. Since IBTK::MergingLoadBalancer is usually what one wants "MERGING" is the default. The merging option is better since it reduces the total number of elements which end up in patch ghost regions since some patches will be merged together.
The parameter workload_quad_point_weight is the multiplier assigned to an IB point when calculating the work per processor. Similarly, the parameter workload_duplicated_node_weight is the multiplier assigned to each node of every element (i.e., each node is counted more than once): see IBTK::FEDataManager::WorkloadSpec for more information.
For efficiency reasons this class only associates elements with patches if they can interact with the patches (e.g., the points used for velocity interpolation can lie inside the patch). Hence, over time, this association needs to change as the structure moves. The frequency of reassociation is controlled by two things:
patch_association_cfl to the input database. If you set IBHierarchyIntegrator::d_regrid_cfl_interval to a value larger than 1 and use AMR this class may not work correctly since, in that case, elements may move outside the patch level they are associated with before new cells are tagged to be on that level.
The logging options set by this class are propagated to the owned IBTK::FEDataManager objects.
enable_logging: set to TRUE to enable logging.. Defaults to false. skip_initial_workload_log: For testing purposes (see d_skip_initial_workload_log) it is necessary to disable some output: this option disables logging of workload data (quadrature point counts, etc.) before the first time step if set to TRUE. Defaults to false. It is often useful to add your own libMesh data to the EquationSystems object used by IBAMR objects. One such example would be defining fields of fibers to give direction to add anisotropy to solid models. Since a libMesh Mesh object also stores all the degrees of freedom of all the systems, there can only be one EquationSystems object for each Mesh object. If you set up your own systems in this way then they will be automatically saved and loaded from restart data.
Since initialization of this class occurs in multiple stages and IBFEMethod assumes that it is ultimately responsible for setting up the EquationSystems object (by calling EquationSystems::init()) other libMesh Systems must be added in a specific order. Additionally, some parts must only be set up if restart information is not available.
This code is based on an IBFE example and assumes that the IBAMR objects objects are already set up in the usual way.
IBFEMethod uses IBTK::FEDataManager to actually perform IB calculations with the finite element mesh. IBTK::FEDataManager objects are configured by setting spreading and interpolation parameters in the usual way (i.e., by providing the parameters described above in the input database). Options specific to the behavior of FEDataManager can be set by defining a database named FEDataManager inside the database provided to this class - see the documentation of FEDataManager for more information.
The caching of the IBFE restart data is not managed by SAMRAI's SAMRAI::tbox::RestartManager. It is instead handled by FEMechanicsBase::writeFEDataToRestartFile() given a restart_dump_dirname and time_step_number. Each instance of IBFEMethod is registered for restart by default, but the this option can be turned off. During a restart, the data is handled by the SAMRAI::tbox::RestartManager automatically to reinitiate the IBFEMethod.
Typedef specifying interface for Lagrangian mass source/sink distribution function.
|
inherited |
Typedef specifying interface for coordinate mapping function.
|
inherited |
Typedef specifying interface for initial velocity specification function.
|
inherited |
Typedef specifying interface for PK1 stress tensor function.
|
inherited |
Typedef specifying interface for Lagrangian body force distribution function.
|
inherited |
Typedef specifying interface for Lagrangian pressure force distribution function.
|
inherited |
Typedef specifying interface for Lagrangian surface force distribution function.
|
inherited |
Function signature for specifying the energy functional that determines the pressure.
| IBAMR::IBFEMethod::IBFEMethod | ( | const std::string & | object_name, |
| const 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::IBFEMethod::IBFEMethod | ( | const std::string & | object_name, |
| const 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 |
||
| ) |
|
delete |
|
delete |
|
overridedefault |
|
delete |
|
inline |
Return the source system name.
| IBTK::FEDataManager* IBAMR::IBFEMethod::getFEDataManager | ( | unsigned int | part = 0 | ) | const |
Return a pointer to the finite element data manager object for the specified part.
|
overridevirtual |
Indicate that a part should include a static pressure.
Reimplemented from IBAMR::FEMechanicsBase.
|
virtual |
Indicate that a part should use stress normalization.
|
virtual |
Register the (optional) function to compute a mass source/sink distribution on the Lagrangian finite element mesh.
| LagBodySourceFcnData IBAMR::IBFEMethod::getLagBodySourceFunction | ( | unsigned int | part = 0 | ) | const |
Get the Lagrangian body source function data.
|
virtual |
Register the (optional) direct forcing kinematics object with the finite element mesh.
|
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 |
Same as the base class.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Inactivate a structure/part. See IBAMR::IBStrategy::inactivateLagrangianStructure().
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Activate a previously inactivated structure/part to be used again in FSI calculations. See IBAMR::IBStrategy::activateLagrangianStructure().
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Determine whether or not the given structure or part is currently activated. See IBAMR::IBStrategy::getLagrangianStructureIsActivated().
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Method to prepare to advance data from current_time to new_time.
Reimplemented from IBAMR::FEMechanicsBase.
|
overridevirtual |
Method to clean up data following call(s) to integrateHierarchy().
Reimplemented from IBAMR::FEMechanicsBase.
|
overridevirtual |
Interpolate the Eulerian velocity to the curvilinear mesh at the specified time within the current time interval.
Implements IBAMR::IBStrategy.
|
overridevirtual |
Indicate that multistep time stepping will be used.
| [in] | n_previous_steps | Number of previous solution values that can be used by the multistep scheme. |
Reimplemented from IBAMR::IBStrategy.
Advance the positions of the Lagrangian structure using the forward Euler method.
Implements IBAMR::IBStrategy.
Advance the positions of the Lagrangian structure using the (explicit) backward Euler method.
Reimplemented from IBAMR::IBStrategy.
Advance the positions of the Lagrangian structure using the (explicit) midpoint rule.
Implements IBAMR::IBStrategy.
Advance the positions of the Lagrangian structure using the (explicit) trapezoidal rule.
Implements IBAMR::IBStrategy.
Advance the positions of the Lagrangian structure using the standard 2nd-order Adams-Bashforth rule.
Reimplemented from 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.
|
overridevirtual |
Indicate whether there are any internal fluid sources/sinks.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Compute the Lagrangian source/sink density at the specified time within the current time interval.
Reimplemented from IBAMR::IBStrategy.
|
overridevirtual |
Spread the Lagrangian source/sink density to the Cartesian grid at the specified time within the current time interval.
Reimplemented from IBAMR::IBStrategy.
| IBTK::FEDataManager::InterpSpec IBAMR::IBFEMethod::getDefaultInterpSpec | ( | ) | const |
Get the default interpolation spec object used by the class.
| IBTK::FEDataManager::SpreadSpec IBAMR::IBFEMethod::getDefaultSpreadSpec | ( | ) | const |
Get the default spread spec object used by the class.
| void IBAMR::IBFEMethod::setWorkloadSpec | ( | const IBTK::FEDataManager::WorkloadSpec & | workload_spec, |
| unsigned int | part = 0 |
||
| ) |
Set the workload spec object used with a particular mesh part.
| void IBAMR::IBFEMethod::setInterpSpec | ( | const IBTK::FEDataManager::InterpSpec & | interp_spec, |
| unsigned int | part = 0 |
||
| ) |
Set the interpolation spec object used with a particular mesh part.
| void IBAMR::IBFEMethod::setSpreadSpec | ( | const IBTK::FEDataManager::SpreadSpec & | spread_spec, |
| unsigned int | part = 0 |
||
| ) |
Set the spread spec object used with a particular mesh part.
|
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 |
this function only exists for compatibility with the base class and does nothing: data reinitialization is handled by endDataRedistribution() instead.
The reasoning is this: since this class stores data only on particular levels (at the present time, the structure is always on the finest level) setting up level data is nontrivial when generating the initial grid (i.e., when tagging cells that contain interaction points for refinement). In a sense there is no level data to compute until we are done regridding.
|
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.
The name here is misleading, but SAMRAI expects us to use one of two tagging methods to refine the grid, and IBAMR consistently uses gradient detection: hence this function has the same name but tags cells in a different way.
|
overridevirtual |
Write out object state to the given database.
Reimplemented from IBAMR::FEMechanicsBase.
| SAMRAI::tbox::Pointer<SAMRAI::hier::BasePatchHierarchy<NDIM> > IBAMR::IBFEMethod::getScratchHierarchy | ( | ) |
Return a pointer to the scratch hierarchy used by this object. See the main documentation of this class for more information.
|
overrideprotectedvirtual |
Do the actual work in initializeFEEquationSystems.
Reimplemented from IBAMR::FEMechanicsBase.
|
overrideprotectedvirtual |
Do the actual work of setting up libMesh system vectors.
Reimplemented from IBAMR::FEMechanicsBase.
|
overrideprotectedvirtual |
Do the actual work in reinitializeFEData and initializeFEData. if use_present_data is true then the current content of the solution vectors is used: more exactly, the coordinates and velocities (computed by initializeCoordinates and initializeVelocity) are considered as being up to date, as is the direct forcing kinematic data.
Reimplemented from IBAMR::FEMechanicsBase.
|
protected |
|
protected |
|
protected |
|
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 |
Get the schedule used to prolong force values. Data is read from coarse_data_idx on level level_number and written into fine_data_idx on level level_number + 1.
|
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.
|
private |
At the present time this class and FEDataManager assume that the finite element mesh is always on the finest grid level. This function explicitly asserts that this condition is met.
|
private |
Convenience function that reinitializes the patch-to-element mappings on all relevant FEDataManagers (i.e., for all parts and, if enabled, on the scratch hierarchy).
|
inherited |
Return a pointer to the equations system object for the specified part.
|
inlineinherited |
Return the current coordinates system name.
|
inlineinherited |
Return the displacement system name.
|
inlineinherited |
Return the force system name.
|
inlineinherited |
Return the pressure system name.
|
inlineinherited |
Return the velocity system name.
|
inherited |
Return a pointer to the FEData object for the specified part.
|
virtualinherited |
Register the (optional) function used to initialize the physical coordinates from the Lagrangian coordinates.
|
inherited |
Get the initial coordinate mapping function data.
|
virtualinherited |
Register the (optional) function used to initialize the velocity of the solid mesh.
|
inherited |
Get the initial velocity function data.
|
virtualinherited |
Register the (optional) function to compute the first Piola-Kirchhoff stress tensor, used to compute the forces on the Lagrangian finite element mesh.
|
inherited |
Get the PK1 stress function data.
|
virtualinherited |
Register the (optional) function to compute body force distributions on the Lagrangian finite element mesh.
|
inherited |
Get the Lagrangian body force function data.
|
virtualinherited |
Register the (optional) function to compute surface pressure distributions on the Lagrangian finite element mesh.
|
inherited |
Get the Lagrangian surface pressure function data.
|
virtualinherited |
Register the (optional) function to compute surface force distributions on the Lagrangian finite element mesh.
|
inherited |
Get the Lagrangian surface force function data.
|
virtualinherited |
Indicate that a part should include a dynamic pressure.
The pressure is determined via (P-dot, Q) = (J U''(J) (FF : Grad U), Q), using either a consistent or lumped mass matrix, or via a locally stabilized projection of the form (P, Q) + epsilon (P - Pi P, Q - Pi Q) = (U'(J), Q), in which P is the pressure and Q is an arbitrary test function.
Users can provide a function to evaluate U''(J). If no function is provided, we default to using U(J) = -kappa (J ln(J) − J + 1), so that U'(J) = -kappa ln J and U''(J) = -kappa J^{-1}. (Ref: C.H. Liu, G. Hofstetter, H.A. Mang, 3D finite element analysis of rubber-like materials at finite strains, Eng. Comput. 11 (2) (1994) 111–128.)
The sign convention used in the implementation generates a PK1 stress of the form PP = -J P FF^{-T}.
Indicate whether there is a static or dynamic pressure associated with the specified FE mesh part.
|
virtualinherited |
Initialize the FE equation systems objects. This method must be called prior to calling initializeFEData().
|
virtualinherited |
Initialize FE data. This method must be called prior to calling IBHierarchyIntegrator::initializePatchHierarchy().
|
virtualinherited |
Reinitialize FE data by calling reinit on each part's EquationSystem, reassembling the system matrices, and setting boundary conditions.
|
virtualinherited |
For technical reasons this class does not use SAMRAI's RestartManager, so restart files must be separately written for the FE objects. This function saves the solutions to the defined EquationSystems in an xdr file in restart_dump_dirname for each FE part. An example snippet is included below to show the distinct FE restart data saving step. The data will then be automatically read back into the system along with the RestartManager data during restart.
|
protectedinherited |
|
protectedinherited |
|
protectedvirtualinherited |
Assemble the RHS for the interior elastic density, possibly splitting off the normal component of the transmission force along the physical boundary of the Lagrangian structure.
|
protectedvirtualinherited |
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.
|
protectedvirtualinherited |
Compute dX = x - X, useful mainly for visualization purposes.
|
protectedvirtualinherited |
Initialize the velocity field using the supplied initial velocity specification function. If no function is provided, the initial velocity is taken to be zero.
|
protectedvirtualinherited |
Get the libMesh restart file name.
|
staticprotectedinherited |
Convenience function to setup system vectors and, if necessary, convert PARALLEL vectors into GHOSTED vectors for a collection of Systems.
|
staticprotectedinherited |
Convenience function to setup a system vector and, if necessary, convert a PARALLEL vector into a GHOSTED vector.
|
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.
|
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 |
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 |
|
protected |
Cached input databases.
|
protected |
|
protected |
Whether or not the initial (i.e., before the regrid prior to timestepping) workload calculations should be logged. This output is generally not stable between machines and so this is usually disabled in tests.
|
protected |
Whether or not we have started time integration. This is only used to determine whether or not we print some initial logging output: see d_skip_initial_workload_log for more information.
|
protected |
Boolean controlling whether or not the scratch hierarchy should be used.
|
protected |
Pointers to the patch hierarchy and gridding algorithm objects associated with this object.
|
protected |
|
protected |
|
protected |
Scratch data caching object. Used both by this class and the FEDataManager objects responsible for working on the primary partitioning.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Maximum level number in the patch hierarchy.
|
protected |
CFL-like number used to determine when we should call reinitElementMappings() based on maximum structure point displacement. More exactly: this class will call that function if the maximum displacement of the structure (calculated by comparing the position vector as of the last reassociation to the current position vector) exceeds dx * d_patch_association_cfl, where dx is the smallest Eulerian cell width.
Note that this is not a regridding, in the sense that the grid changes: instead only the association between patches and elements changes.
@seealso IBHierarchyIntegrator::d_regrid_cfl_interval
|
protected |
Indexing information determining whether a given part is active or not. The default state for each part is to be active. Parts are active unless inactivated via inactivateLagrangianStructure().
|
protected |
FEDataManager objects associated with the primary hierarchy (i.e., d_hierarchy). These are used by some other objects (such as IBFEPostProcessor); IBFEMethod keeps them up to date (i.e., reinitializing data after regrids).
|
protected |
FEDataManager objects that use the scratch hierarchy instead of d_hierarchy. These are only used internally by IBFEMethod and are not intended to be accessed by any other object.
|
protected |
The FEDataManager objects that are actually used in computations. This vector will be equal to either d_primary_fe_data_managers or d_scratch_fe_data_managers, dependent on which is actually used in IB calculations.
|
protected |
|
protected |
Schedules for prolonging data during spreading. The keys are the level number, the patch data index for the coarse level data, and the patch data index which will be filled with fine level data.
|
protected |
|
protected |
Vectors of pointers to the systems for each part (for sources; other systems are handled by the base class).
|
protected |
Object managing access to libMesh system vectors for the source/sink strength.
|
protected |
Object managing access to libMesh system vectors for the position that include IB ghosting information.
|
protected |
Object managing access to libMesh system vectors for the velocity that include IB ghosting information.
|
protected |
Object managing access to libMesh system vectors for the force that include IB ghosting information.
|
protected |
Object managing access to libMesh system vectors for the source/sink strength that include IB ghosting information.
|
protected |
IBFE method parameters.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Data related to handling stress normalization.
|
protected |
|
protected |
|
protected |
Objects used to impose direct forcing kinematics.
|
protected |
Functions used to compute source/sink strength on the Lagrangian mesh.
|
protected |
|
protected |
|
protected |
The optional second (i.e., scratch) hierarchy. This hierarchy is only used for evaluation of IB terms, i.e., in IBFEMethod::interpolateVelocity(), IBFEMethod::spreadForce(), and IBFEMethod::spreadFluidSource().
|
protected |
Data related to multistep time stepping.
|
protected |
|
staticinherited |
|
staticinherited |
|
staticinherited |
|
staticinherited |
|
staticinherited |
|
protectedinherited |
Cached input databases.
|
protectedinherited |
Indicates whether the integrator should output logging messages.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Meshes provided to this object. These are set up and managed outside this class. These meshes are modified by FEMechanicsBase since this class creates several libMesh Systems (and hence stores DoF information in these meshes).
|
protectedinherited |
The libMesh Systems set up by this system (for example, for velocity projection) consist of one variable per spatial component. By default, libMesh assumes that all variables in a given System couple to each other which, since we only ever solve projection problems in this class, is not the case. Hence we can save some memory by explicitly informing libMesh that the variables in a system only couple to themselves by providing a diagonal coupling matrix to each System.
|
protectedinherited |
EquationSystems objects, one per part. These contain the actual matrices and solution vectors for each relevant libMesh system.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Vectors of pointers to the systems for each part (for position, velocity, force, and pressure).
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Object managing access to libMesh system vectors for the position.
|
protectedinherited |
Object managing access to libMesh system vectors for the velocity.
|
protectedinherited |
Object managing access to libMesh system vectors for the force.
|
protectedinherited |
Object managing access to libMesh system vectors for the pressure.
|
protectedinherited |
Whether or not the libMesh equation systems objects have been initialized (i.e., whether or not initializeFEEquationSystems has been called).
|
protectedinherited |
Whether or not all finite element data (including that initialized by initializeFEEquationSystems), such system matrices, is available.
|
protectedinherited |
Type of partitioner to use. See the main documentation of this class for more information.
|
protectedinherited |
Whether or not to use AMR in the finite element discretization. This feature is not yet implemented and currently defaults to false.
|
protectedinherited |
Method parameters.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Functions used to compute the initial coordinates of the Lagrangian mesh.
|
protectedinherited |
Functions used to compute the initial coordinates of the Lagrangian mesh.
|
protectedinherited |
Functions used to compute the first Piola-Kirchhoff stress tensor.
|
protectedinherited |
Functions used to compute additional body and surface forces on the Lagrangian mesh.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Data related to handling static pressures.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
Data related to handling dynamic pressures.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
The object name is used as a handle to databases stored in restart files and for error reporting purposes.
|
protectedinherited |
A boolean value indicating whether the class is registered with the restart database.
|
protectedinherited |
Directory and time step number to use when restarting.
|
protectedinherited |
|
protectedinherited |
Restart file type for libMesh equation systems (e.g. xda or xdr).
|
protectedinherited |
The IBHierarchyIntegrator object that is using this strategy class.
|
protectedinherited |
Whether to use "fixed" Lagrangian-Eulerian coupling operators.
1.8.17