|
IBAMR
IBAMR version 0.19.
|
Class FEMechanicsExplicitIntegrator is an implementation of the abstract base class FEMechanicsBase that provides a simple explicit elastodynamics time integrator with an interface that is similar to IBFEMethod to facilitate model re-use. More...
#include <ibamr/FEMechanicsExplicitIntegrator.h>

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 | |
| FEMechanicsExplicitIntegrator (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. More... | |
| FEMechanicsExplicitIntegrator (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. More... | |
| FEMechanicsExplicitIntegrator ()=delete | |
| Deleted default constructor. More... | |
| FEMechanicsExplicitIntegrator (const FEMechanicsExplicitIntegrator &from)=delete | |
| Deleted copy constructor. More... | |
| FEMechanicsExplicitIntegrator & | operator= (const FEMechanicsExplicitIntegrator &that)=delete |
| Deleted assignment operator. More... | |
| ~FEMechanicsExplicitIntegrator () override=default | |
| Defaulted destructor. More... | |
| 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 | forwardEulerStep (double current_time, double new_time) |
| void | modifiedEulerStep (double current_time, double new_time) |
| void | backwardEulerStep (double current_time, double new_time) |
| void | midpointStep (double current_time, double new_time) |
| void | trapezoidalStep (double current_time, double new_time) |
| void | modifiedTrapezoidalStep (double current_time, double new_time) |
| void | SSPRK3Step (double current_time, double new_time, unsigned int n_stages=4) |
| void | computeLagrangianForce (double data_time) |
| void | computeLagrangianForce (libMesh::PetscVector< double > &F_vec, libMesh::PetscVector< double > &X_vec, libMesh::PetscVector< double > &U_vec, libMesh::PetscVector< double > *P_vec, double data_time, unsigned int part) |
| void | putToDatabase (SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > db) 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 | initializeFEEquationSystems () |
| virtual void | initializeFEData () |
| virtual void | reinitializeFEData () |
| 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 | |
| void | doInitializeFEEquationSystems () override |
| void | doInitializeFESystemVectors () override |
| void | doInitializeFEData (bool use_present_data) override |
| void | doForwardEulerStep (libMesh::PetscVector< double > &X_new_vec, libMesh::PetscVector< double > &U_new_vec, libMesh::PetscVector< double > *P_new_vec, libMesh::PetscVector< double > &X_current_vec, libMesh::PetscVector< double > &U_current_vec, libMesh::PetscVector< double > *P_current_vec, libMesh::PetscVector< double > &F_current_vec, libMesh::PetscVector< double > *dP_dt_current_vec, double current_time, double new_time, unsigned int part) |
| 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 SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db) |
| void | getFromInput (const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &db, bool is_from_restart) |
| void | getFromRestart () |
|
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::FEMechanicsExplicitIntegrator::FEMechanicsExplicitIntegrator | ( | 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 |
||
| ) |
| IBAMR::FEMechanicsExplicitIntegrator::FEMechanicsExplicitIntegrator | ( | 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 |
||
| ) |
|
delete |
|
delete |
|
overridedefault |
|
delete |
|
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.
| void IBAMR::FEMechanicsExplicitIntegrator::forwardEulerStep | ( | double | current_time, |
| double | new_time | ||
| ) |
Advance the structural velocities and positions using forward Euler.
| void IBAMR::FEMechanicsExplicitIntegrator::modifiedEulerStep | ( | double | current_time, |
| double | new_time | ||
| ) |
Advance the structural velocities and positions using "modified" Euler.
NOTE: This scheme first advances the velocity using forward Euler, and then uses the updated velocity to advance the position of the structure.
| void IBAMR::FEMechanicsExplicitIntegrator::backwardEulerStep | ( | double | current_time, |
| double | new_time | ||
| ) |
Advance the structural velocities and positions using backward Euler.
NOTE: This is not a true backward Euler method; we simply use the most recent approximations to the updated velocity and position to compute the force at the end of the time step, and then use that force to increment the velocity and position via backward Euler.
Advance the structural velocities and positions using the explicit midpoint rule.
Advance the structural velocities and positions using the explicit SSP RK2 scheme, which is the explicit trapezoidal rule.
| void IBAMR::FEMechanicsExplicitIntegrator::modifiedTrapezoidalStep | ( | double | current_time, |
| double | new_time | ||
| ) |
Advance the structural velocities and positions using a modified explicit trapezoidal rule.
NOTE: In this scheme, the first stage uses a modified Euler scheme (
| void IBAMR::FEMechanicsExplicitIntegrator::SSPRK3Step | ( | double | current_time, |
| double | new_time, | ||
| unsigned int | n_stages = 4 |
||
| ) |
Advance the structural velocities and positions using the explicit SSP RK3 scheme.
The implementation supports either 3- or 4-stage variants. The default is the 4-stage version, which has a larger stability region.
| void IBAMR::FEMechanicsExplicitIntegrator::computeLagrangianForce | ( | double | data_time | ) |
Compute the Lagrangian force at the specified time within the current time interval.
| void IBAMR::FEMechanicsExplicitIntegrator::computeLagrangianForce | ( | libMesh::PetscVector< double > & | F_vec, |
| libMesh::PetscVector< double > & | X_vec, | ||
| libMesh::PetscVector< double > & | U_vec, | ||
| libMesh::PetscVector< double > * | P_vec, | ||
| double | data_time, | ||
| unsigned int | part | ||
| ) |
Compute the Lagrangian force at the specified time within the current time interval using the provided data vectors.
|
overridevirtual |
Write out object state to the given database.
Reimplemented from IBAMR::FEMechanicsBase.
|
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 |
Perform a forward Euler step.
|
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.
|
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 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.
|
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.
|
protected |
|
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).
1.8.17