|
IBAMR
IBAMR version 0.19.
|
Class FEMechanicsBase provides core finite element mechanics functionality and data management. More...
#include <ibamr/FEMechanicsBase.h>

Classes | |
| struct | CoordinateMappingFcnData |
| struct | InitialVelocityFcnData |
| struct | LagBodyForceFcnData |
| struct | LagSurfaceForceFcnData |
| struct | LagSurfacePressureFcnData |
| struct | PK1StressFcnData |
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 | PK1StressFcnPtr = IBTK::TensorMeshFcnPtr |
| using | LagBodyForceFcnPtr = IBTK::VectorMeshFcnPtr |
| using | LagSurfacePressureFcnPtr = IBTK::ScalarSurfaceFcnPtr |
| using | LagSurfaceForceFcnPtr = IBTK::VectorSurfaceFcnPtr |
| using | VolumetricEnergyDerivativeFcn = double(*)(double) |
Public Member Functions | |
| FEMechanicsBase (const std::string &object_name, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, libMesh::MeshBase *mesh, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0) | |
| FEMechanicsBase (const std::string &object_name, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, const std::vector< libMesh::MeshBase * > &meshes, bool register_for_restart=true, const std::string &restart_read_dirname="", unsigned int restart_restore_number=0) | |
| FEMechanicsBase ()=delete | |
| FEMechanicsBase (const FEMechanicsBase &from)=delete | |
| FEMechanicsBase & | operator= (const FEMechanicsBase &that)=delete |
| virtual | ~FEMechanicsBase () override |
| 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 | registerStaticPressurePart (PressureProjectionType projection_type=CONSISTENT_PROJECTION, VolumetricEnergyDerivativeFcn dU_dJ_fcn=nullptr, unsigned int part=0) |
| 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 | preprocessIntegrateData (double current_time, double new_time, int num_cycles) |
| virtual void | postprocessIntegrateData (double current_time, double new_time, int num_cycles) |
| virtual void | initializeFEEquationSystems () |
| virtual void | initializeFEData () |
| virtual void | reinitializeFEData () |
| virtual void | putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) override |
| virtual void | writeFEDataToRestartFile (const std::string &restart_dump_dirname, unsigned int time_step_number) |
Static Public Attributes | |
| 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 () |
| virtual void | doInitializeFESystemVectors () |
| virtual void | doInitializeFEData (bool use_present_data) |
| 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 |
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 std::string &object_name, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db, const std::vector< libMesh::MeshBase * > &meshes, bool register_for_restart, const std::string &restart_read_dirname, unsigned int restart_restore_number) |
| void | getFromInput (const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &db, bool is_from_restart) |
| void | getFromRestart () |
| using IBAMR::FEMechanicsBase::CoordinateMappingFcnPtr = void (*)(libMesh::Point& x, const libMesh::Point& X, void* ctx) |
Typedef specifying interface for coordinate mapping function.
| using IBAMR::FEMechanicsBase::InitialVelocityFcnPtr = void (*)(libMesh::VectorValue<double>& U0, const libMesh::Point& X0, void* ctx) |
Typedef specifying interface for initial velocity specification function.
Typedef specifying interface for PK1 stress tensor function.
Typedef specifying interface for Lagrangian body force distribution function.
Typedef specifying interface for Lagrangian pressure force distribution function.
Typedef specifying interface for Lagrangian surface force distribution function.
Function signature for specifying the energy functional that determines the pressure.
| IBAMR::FEMechanicsBase::FEMechanicsBase | ( | const std::string & | object_name, |
| const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > & | input_db, | ||
| libMesh::MeshBase * | mesh, | ||
| bool | register_for_restart = true, |
||
| const std::string & | restart_read_dirname = "", |
||
| unsigned int | restart_restore_number = 0 |
||
| ) |
Constructor for a single-part model.
| IBAMR::FEMechanicsBase::FEMechanicsBase | ( | const std::string & | object_name, |
| const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > & | input_db, | ||
| const std::vector< libMesh::MeshBase * > & | meshes, | ||
| bool | register_for_restart = true, |
||
| const std::string & | restart_read_dirname = "", |
||
| unsigned int | restart_restore_number = 0 |
||
| ) |
Constructor for a multi-part model.
|
delete |
Deleted default constructor.
|
delete |
Deleted copy constructor.
|
overridevirtual |
Destructor.
|
delete |
Deleted assignment operator.
| libMesh::EquationSystems* IBAMR::FEMechanicsBase::getEquationSystems | ( | unsigned int | part = 0 | ) | const |
Return a pointer to the equations system object for the specified part.
|
inline |
Return the current coordinates system name.
|
inline |
Return the displacement system name.
|
inline |
Return the force system name.
|
inline |
Return the pressure system name.
|
inline |
Return the velocity system name.
| std::shared_ptr<IBTK::FEData> IBAMR::FEMechanicsBase::getFEData | ( | unsigned int | part = 0 | ) | const |
Return a pointer to the FEData object for the specified part.
|
virtual |
Register the (optional) function used to initialize the physical coordinates from the Lagrangian coordinates.
| CoordinateMappingFcnData IBAMR::FEMechanicsBase::getInitialCoordinateMappingFunction | ( | unsigned int | part = 0 | ) | const |
Get the initial coordinate mapping function data.
|
virtual |
Register the (optional) function used to initialize the velocity of the solid mesh.
| InitialVelocityFcnData IBAMR::FEMechanicsBase::getInitialVelocityFunction | ( | unsigned int | part = 0 | ) | const |
Get the initial velocity function data.
|
virtual |
Register the (optional) function to compute the first Piola-Kirchhoff stress tensor, used to compute the forces on the Lagrangian finite element mesh.
| std::vector<PK1StressFcnData> IBAMR::FEMechanicsBase::getPK1StressFunction | ( | unsigned int | part = 0 | ) | const |
Get the PK1 stress function data.
|
virtual |
Register the (optional) function to compute body force distributions on the Lagrangian finite element mesh.
| LagBodyForceFcnData IBAMR::FEMechanicsBase::getLagBodyForceFunction | ( | unsigned int | part = 0 | ) | const |
Get the Lagrangian body force function data.
|
virtual |
Register the (optional) function to compute surface pressure distributions on the Lagrangian finite element mesh.
| LagSurfacePressureFcnData IBAMR::FEMechanicsBase::getLagSurfacePressureFunction | ( | unsigned int | part = 0 | ) | const |
Get the Lagrangian surface pressure function data.
|
virtual |
Register the (optional) function to compute surface force distributions on the Lagrangian finite element mesh.
| LagSurfaceForceFcnData IBAMR::FEMechanicsBase::getLagSurfaceForceFunction | ( | unsigned int | part = 0 | ) | const |
Get the Lagrangian surface force function data.
|
virtual |
Indicate that a part should include a static pressure.
The pressure is determined via (P, Q) = (U'(J), 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. (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}.
Reimplemented in IBAMR::IBFEMethod.
|
virtual |
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.
|
virtual |
Method to prepare to advance data from current_time to new_time.
Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.
|
virtual |
Method to clean up data following call(s) to integrateHierarchy().
Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.
|
virtual |
Initialize the FE equation systems objects. This method must be called prior to calling initializeFEData().
|
virtual |
Initialize FE data. This method must be called prior to calling IBHierarchyIntegrator::initializePatchHierarchy().
|
virtual |
Reinitialize FE data by calling reinit on each part's EquationSystem, reassembling the system matrices, and setting boundary conditions.
|
overridevirtual |
Write out object state to the given database.
Implements SAMRAI::tbox::Serializable.
Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.
|
virtual |
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.
|
protectedvirtual |
Do the actual work in initializeFEEquationSystems.
Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.
|
protectedvirtual |
Do the actual work of setting up libMesh system vectors.
Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.
|
protectedvirtual |
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.
Reimplemented in IBAMR::IBFEMethod, and IBAMR::FEMechanicsExplicitIntegrator.
|
protected |
|
protected |
|
protectedvirtual |
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.
|
protectedvirtual |
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.
|
protectedvirtual |
Compute dX = x - X, useful mainly for visualization purposes.
|
protectedvirtual |
Initialize the velocity field using the supplied initial velocity specification function. If no function is provided, the initial velocity is taken to be zero.
|
protectedvirtual |
Get the libMesh restart file name.
|
staticprotected |
Convenience function to setup system vectors and, if necessary, convert PARALLEL vectors into GHOSTED vectors for a collection of Systems.
|
staticprotected |
Convenience function to setup a system vector and, if necessary, convert a PARALLEL vector into a GHOSTED vector.
|
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.
|
static |
|
static |
|
static |
|
static |
|
static |
|
protected |
Cached input databases.
|
protected |
Indicates whether the integrator should output logging messages.
|
protected |
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).
|
protected |
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.
|
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 |
|
protected |
Vectors of pointers to the systems for each part (for position, velocity, force, and pressure).
|
protected |
|
protected |
|
protected |
|
protected |
Object managing access to libMesh system vectors for the position.
|
protected |
Object managing access to libMesh system vectors for the velocity.
|
protected |
Object managing access to libMesh system vectors for the force.
|
protected |
Object managing access to libMesh system vectors for the pressure.
|
protected |
Whether or not the libMesh equation systems objects have been initialized (i.e., whether or not initializeFEEquationSystems has been called).
|
protected |
Whether or not all finite element data (including that initialized by initializeFEEquationSystems), such system matrices, is available.
|
protected |
Type of partitioner to use. See the main documentation of this class for more information.
|
protected |
Whether or not to use AMR in the finite element discretization. This feature is not yet implemented and currently defaults to false.
|
protected |
Method parameters.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Functions used to compute the initial coordinates of the Lagrangian mesh.
|
protected |
Functions used to compute the initial coordinates of the Lagrangian mesh.
|
protected |
Functions used to compute the first Piola-Kirchhoff stress tensor.
|
protected |
Functions used to compute additional body and surface forces on the Lagrangian mesh.
|
protected |
|
protected |
|
protected |
Data related to handling static pressures.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Data related to handling dynamic pressures.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
The object name is used as a handle to databases stored in restart files and for error reporting purposes.
|
protected |
A boolean value indicating whether the class is registered with the restart database.
|
protected |
Directory and time step number to use when restarting.
|
protected |
|
protected |
Restart file type for libMesh equation systems (e.g. xda or xdr).
1.8.17