|
PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
|
Finite Element Elasto-Dynamics initial value problem with Dirichlet boundary conditions using BDF (backward differentiation formula) as the time discretization. More...
#include <FemElastoDynamics.h>
Public Types | |
| using | MeshType = fem::Mesh<TElement, Dims, TScalar, TIndex> |
| Mesh type. | |
| using | ElementType = TElement |
| Element type. | |
| using | ElasticEnergyType = THyperElasticEnergy |
| Elastic energy type. | |
| using | BdfType = integration::Bdf<TScalar, TIndex> |
| BDF time integrator type. | |
| using | ScalarType = TScalar |
| Floating point scalar type. | |
| using | IndexType = TIndex |
| Integer index type. | |
Public Member Functions | |
| FemElastoDynamics ()=default | |
| Construct an empty FemElastoDynamics problem. | |
| FemElastoDynamics (Eigen::Ref< Eigen::Matrix< ScalarType, Eigen::Dynamic, Eigen::Dynamic > const > const &V, Eigen::Ref< Eigen::Matrix< IndexType, Eigen::Dynamic, Eigen::Dynamic > const > const &C) | |
| Construct an FemElastoDynamics problem on the mesh domain (V,C). | |
| void | Construct (Eigen::Ref< Eigen::Matrix< ScalarType, Eigen::Dynamic, Eigen::Dynamic > const > const &V, Eigen::Ref< Eigen::Matrix< IndexType, Eigen::Dynamic, Eigen::Dynamic > const > const &C) |
| Construct an FemElastoDynamics problem on the mesh domain (V,C). | |
| template<class TDerivedX0, class TDerivedV0> | |
| void | SetInitialConditions (Eigen::DenseBase< TDerivedX0 > const &x0, Eigen::DenseBase< TDerivedV0 > const &v0) |
| Set the initial conditions for the initial value problem. | |
| void | SetMassMatrix (ScalarType rho) |
| Compute, lump and set the mass matrix with homogeneous density \( \rho \). | |
| void | SetElasticEnergy (ScalarType mu, ScalarType lambda) |
| Set the elastic energy quadrature for a homogeneous material with Lame coefficients \( \mu \) and \( \lambda \). | |
| void | SetExternalLoad (Eigen::Vector< ScalarType, kDims > const &b) |
| Compute and set the external load vector given by fixed body forces \( b \). | |
| void | SetTimeIntegrationScheme (ScalarType dt=ScalarType(1e-2), int s=1) |
| Set the BDF (backward differentiation formula) time integration scheme. | |
| template<typename TDerivedDirichletMask> | |
| void | Constrain (Eigen::DenseBase< TDerivedDirichletMask > const &D) |
| Set the Dirichlet boundary conditions. | |
| template<class TDerivedEg, class TDerivedWg, class TDerivedXg, class TDerivedRhog> | |
| void | SetMassMatrix (Eigen::DenseBase< TDerivedEg > const &eg, Eigen::MatrixBase< TDerivedWg > const &wg, Eigen::MatrixBase< TDerivedXg > const &Xg, Eigen::MatrixBase< TDerivedRhog > const &rhog) |
| Compute, lump and set the mass matrix with variable density \( \rho(X) \) at quadrature points \( X_g \) of the given quadrature rule \( (w_g, X_g) \). | |
| template<class TDerivedEg, class TDerivedWg, class TDerivedXg, class TDerivedMug, class TDerivedLambdag> | |
| void | SetElasticEnergy (Eigen::DenseBase< TDerivedEg > const &eg, Eigen::DenseBase< TDerivedWg > const &wg, Eigen::MatrixBase< TDerivedXg > const &Xg, Eigen::DenseBase< TDerivedMug > const &mug, Eigen::DenseBase< TDerivedLambdag > const &lambdag) |
| Compute and set the elastic energy quadrature for a heterogeneous material with variable Lame coefficients \( \mu(X) \) and \( \lambda(X) \) at quadrature points \( X_g
\) of the given quadrature rule \( (w_g, X_g) \). | |
| template<class TDerivedEg, class TDerivedWg, class TDerivedXg, class TDerivedBg> | |
| void | SetExternalLoad (Eigen::DenseBase< TDerivedEg > const &eg, Eigen::MatrixBase< TDerivedWg > const &wg, Eigen::MatrixBase< TDerivedXg > const &Xg, Eigen::MatrixBase< TDerivedBg > const &bg) |
| Compute and set the external load vector given by variable body forces \( b(X) \) at quadrature points \( X_g \) of the given quadrature rule \( (w_g, X_g) \). | |
| void | SetupTimeIntegrationOptimization () |
| Set the BDF inertial target for elasto dynamics. | |
| void | ComputeElasticEnergy (int eElasticComputationFlags, fem::EHyperElasticSpdCorrection eSpdCorrectionFlags) |
| Compute the quadrature point elastic energies of the current configuration into Ug, Gg, Hg. | |
| auto | M () const |
| k-dimensional mass matrix | |
| auto | aext () const |
| External acceleration. | |
| bool | IsDirichletNode (IndexType node) const |
| Check if a node is Dirichlet constrained. | |
| bool | IsDirichletDof (IndexType i) const |
| Check if a coordinate is Dirichlet constrained. | |
| auto | DirichletNodes () const |
| Array of Dirichlet constrained nodes. | |
| auto | DirichletCoordinates () const |
| Dirichlet nodal positions. | |
| auto | DirichletVelocities () const |
| Dirichlet nodal velocities. | |
| auto | DirichletCoordinates () |
| Dirichlet nodal positions. | |
| auto | DirichletVelocities () |
| Dirichlet nodal velocities. | |
| auto | FreeNodes () const |
| Array of unconstrained nodes. | |
| auto | FreeCoordinates () const |
| Free nodal positions. | |
| auto | FreeVelocities () const |
| Free nodal velocities. | |
| auto | FreeCoordinates () |
| Free nodal positions. | |
| auto | FreeVelocities () |
| Free nodal velocities. | |
| void | Serialize (io::Archive &archive) const |
| Serialize to HDF5 group. | |
| void | Deserialize (io::Archive const &archive) |
| Deserialize from HDF5 group. | |
Public Attributes | |
| MeshType | mesh |
| FEM mesh. | |
| BdfType | bdf |
| BDF time integration scheme. | |
| Eigen::Matrix< ScalarType, kDims, Eigen::Dynamic > | fext |
kDims x |# nodes| matrix of external forces at nodes | |
| Eigen::Vector< ScalarType, Eigen::Dynamic > | m |
|# nodes| x 1 lumped mass matrix | |
| Eigen::Matrix< ScalarType, kDims, Eigen::Dynamic > | xtilde |
kDims x |# nodes| inertial targets | |
| Eigen::Matrix< ScalarType, kDims, Eigen::Dynamic > | x |
kDims x |# nodes| positions | |
| Eigen::Matrix< ScalarType, kDims, Eigen::Dynamic > | v |
kDims x |# nodes| velocities | |
| Eigen::Vector< IndexType, Eigen::Dynamic > | egU |
| Eigen::Vector< ScalarType, Eigen::Dynamic > | wgU |
|# quad.pts.| x 1 vector of quadrature weights for elastic potential | |
| Eigen::Matrix< ScalarType, Eigen::Dynamic, Eigen::Dynamic > | GNegU |
| Eigen::Matrix< ScalarType, Eigen::Dynamic, Eigen::Dynamic > | lamegU |
2 x |# quad.pts.| matrix of Lame coefficients at quadrature points | |
| Eigen::Vector< ScalarType, Eigen::Dynamic > | UgU |
|# quad.pts.| x 1 vector of elastic energy density at quadrature points | |
| Eigen::Matrix< ScalarType, kDims *ElementType::kNodes, Eigen::Dynamic > | GgU |
| Eigen::Matrix< ScalarType, kDims *ElementType::kNodes, Eigen::Dynamic > | HgU |
| IndexType | ndbc |
| Number of Dirichlet constrained nodes. | |
| Eigen::Vector< IndexType, Eigen::Dynamic > | dbc |
| Eigen::Vector< bool, Eigen::Dynamic > | dmask |
Static Public Attributes | |
| static auto constexpr | kDims = Dims |
| Dimensionality of the mesh. | |
Finite Element Elasto-Dynamics initial value problem with Dirichlet boundary conditions using BDF (backward differentiation formula) as the time discretization.
Represents the problem
\[\min_x \frac{1}{2} || x - \tilde{x} ||_M^2 + \tilde{\beta}_\text{bdf-s}^2 U(x) , \]
where \( M \) is the mass matrix, \( \tilde{x} \) is the BDF inertial target, \( \tilde{\beta}_\text{bdf-s} \) is the forcing term's coefficient, and \( U(x) \) is the hyper-elastic potential.
| TElement | Element type |
| Dims | Dimensionality of the mesh |
| THyperElasticEnergy | Hyper elastic energy type |
| TScalar | Floating point scalar type |
| TIndex | Integer index type |
|
inline |
Construct an FemElastoDynamics problem on the mesh domain (V,C).
All FemElastoDynamics quantities are initialized to sensible defaults, i.e. rest pose positions, zero velocities, homogeneous rubber-like material properties, and gravitational load.
| V | kDims x |# verts| matrix of mesh vertex positions |
| C | |# cell nodes| x |# cells| matrix of mesh cells |
|
inline |
External acceleration.
kDims x |# nodes| external acceleration
|
inline |
Compute the quadrature point elastic energies of the current configuration into Ug, Gg, Hg.
| eElasticComputationFlags | Flags for computing elastic potential, gradient, and/or hessian |
| eSpdCorrectionFlags | Flags for SPD correction of element hessians |
|
inline |
Set the Dirichlet boundary conditions.
| D | |# nodes| x 1 mask of Dirichlet boundary conditions s.t. D(i) == true if node i is constrained |
D.size() == mesh.X.cols()
|
inline |
Construct an FemElastoDynamics problem on the mesh domain (V,C).
All FemElastoDynamics quantities are initialized to sensible defaults, i.e. rest pose positions, zero velocities, homogeneous rubber-like material properties, and gravitational load.
| V | kDims x |# verts| matrix of mesh vertex positions |
| C | |# cell nodes| x |# cells| matrix of mesh cells |
|
inline |
Deserialize from HDF5 group.
| archive | Archive to deserialize from |
|
inline |
Dirichlet nodal positions.
kDims x ndbc matrix of Dirichlet constrained nodal positions
|
inline |
Dirichlet nodal positions.
kDims x ndbc matrix of Dirichlet constrained nodal positions
|
inline |
Array of Dirichlet constrained nodes.
ndbc x 1 array of Dirichlet constrained nodes
|
inline |
Dirichlet nodal velocities.
kDims x ndbc matrix of Dirichlet constrained nodal velocities
|
inline |
Dirichlet nodal velocities.
kDims x ndbc matrix of Dirichlet constrained nodal velocities
|
inline |
Free nodal positions.
kDims x |# nodes - ndbc| matrix of unconstrained nodal positions
|
inline |
Free nodal positions.
kDims x |# nodes - ndbc| matrix of unconstrained nodal positions
|
inline |
Array of unconstrained nodes.
|# nodes - ndbc| x 1 array of unconstrained nodes
|
inline |
Free nodal velocities.
kDims x |# nodes - ndbc| matrix of unconstrained nodal velocities
|
inline |
Free nodal velocities.
kDims x |# nodes - ndbc| matrix of unconstrained nodal velocities
|
inline |
Check if a coordinate is Dirichlet constrained.
| i | Coordinate index |
true if the coordinate is Dirichlet constrained, false otherwise
|
inline |
Check if a node is Dirichlet constrained.
| node | Node index |
true if the node is Dirichlet constrained, false otherwise
|
inline |
k-dimensional mass matrix
kDims * |#nodes| x 1 vector of the kDims-dimensional lumped mass matrix diagonal coefficients
|
inline |
Serialize to HDF5 group.
| archive | Archive to serialize to |
|
inline |
Compute and set the elastic energy quadrature for a heterogeneous material with variable Lame coefficients \( \mu(X) \) and \( \lambda(X) \) at quadrature points \( X_g \) of the given quadrature rule \( (w_g, X_g) \).
| TDerivedEg | Eigen dense expression type for element indices |
| TDerivedWg | Eigen dense expression type for quadrature weights |
| TDerivedXg | Eigen dense expression type for quadrature points |
| TDerivedMug | Eigen dense expression type for 1st Lame coefficients |
| TDerivedLambdag | Eigen dense expression type for 2nd Lame coefficients |
| eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
| wg | |# quad.pts.| x 1 vector of quadrature weights |
| Xg | |# dims| x |# quad.pts.| matrix of quadrature points |
| mug | |# quad.pts.| x 1 vector of 1st Lame coefficients at quadrature points |
| lambdag | |# quad.pts.| x 1 vector of 2nd Lame coefficients at quadrature points |
|
inline |
Set the elastic energy quadrature for a homogeneous material with Lame coefficients \( \mu \) and \( \lambda \).
| mu | 1st Lame coefficient \( \mu \) |
| lambda | 2nd Lame coefficient \( \lambda \) |
|
inline |
Compute and set the external load vector given by variable body forces \( b(X) \) at quadrature points \( X_g \) of the given quadrature rule \( (w_g, X_g) \).
| TDerivedEg | Eigen dense expression type for element indices |
| TDerivedWg | Eigen dense expression type for quadrature weights |
| TDerivedXg | Eigen dense expression type for quadrature points |
| TDerivedBg | Eigen dense expression type for body forces |
| eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
| wg | |# quad.pts.| x 1 vector of quadrature weights |
| Xg | |# dims| x |# quad.pts.| matrix of quadrature points |
| bg | kDims x |# quad.pts.| matrix of body forces at quadrature points |
|
inline |
Compute and set the external load vector given by fixed body forces \( b \).
| b | kDims x 1 fixed body forces |
|
inline |
Set the initial conditions for the initial value problem.
| x0 | kDims x |# nodes| matrix of initial nodal coordinates |
| v0 | kDims x |# nodes| matrix of initial nodal velocities |
|
inline |
Compute, lump and set the mass matrix with variable density \( \rho(X) \) at quadrature points \( X_g \) of the given quadrature rule \( (w_g, X_g) \).
| TDerivedEg | Eigen dense expression type for element indices |
| TDerivedWg | Eigen dense expression type for quadrature weights |
| TDerivedXg | Eigen dense expression type for quadrature points |
| TDerivedRhog | Eigen dense expression type for mass density |
| eg | |# quad.pts.| x 1 vector of element indices at quadrature points |
| wg | |# quad.pts.| x 1 vector of quadrature weights |
| Xg | |# dims| x |# quad.pts.| matrix of quadrature points |
| rhog | |# quad.pts.| x 1 vector of mass density at quadrature points |
|
inline |
Compute, lump and set the mass matrix with homogeneous density \( \rho \).
| rho | Mass density of the material |
|
inline |
Set the BDF (backward differentiation formula) time integration scheme.
| dt | Time step size |
| s | Step of the BDF scheme |
| Eigen::Vector<IndexType, Eigen::Dynamic> pbat::sim::dynamics::FemElastoDynamics< TElement, Dims, THyperElasticEnergy, TScalar, TIndex >::dbc |
|# nodes| x 1 concatenated vector of Dirichlet unconstrained and constrained nodes, partitioned as [ dbc(0 : |#nodes|-ndbc), dbc(|# nodes|-ndbc : |# nodes|) ]
| Eigen::Vector<bool, Eigen::Dynamic> pbat::sim::dynamics::FemElastoDynamics< TElement, Dims, THyperElasticEnergy, TScalar, TIndex >::dmask |
|# nodes| x 1 mask of Dirichlet boundary conditions s.t. dmask(i) == true if node i is constrained
| Eigen::Vector<IndexType, Eigen::Dynamic> pbat::sim::dynamics::FemElastoDynamics< TElement, Dims, THyperElasticEnergy, TScalar, TIndex >::egU |
|# quad.pts.| x 1 vector of element indices for quadrature points of elastic potential
| Eigen::Matrix<ScalarType, kDims * ElementType::kNodes, Eigen::Dynamic> pbat::sim::dynamics::FemElastoDynamics< TElement, Dims, THyperElasticEnergy, TScalar, TIndex >::GgU |
|# dims * # elem nodes| x |# quad.pts.| matrix of element elastic gradient vectors at quadrature points
| Eigen::Matrix<ScalarType, Eigen::Dynamic, Eigen::Dynamic> pbat::sim::dynamics::FemElastoDynamics< TElement, Dims, THyperElasticEnergy, TScalar, TIndex >::GNegU |
|ElementType::kNodes| x |kDims * # quad.pts.| matrix of shape function gradients at quadrature points
| Eigen::Matrix<ScalarType, kDims * ElementType::kNodes, Eigen::Dynamic> pbat::sim::dynamics::FemElastoDynamics< TElement, Dims, THyperElasticEnergy, TScalar, TIndex >::HgU |
|# dims * # elem nodes| x |# dims * # elem nodes * # quad.pts.| matrix of element elastic hessian matrices at quadrature points