Solid namespace
Hyperelastic solid mechanics module for large-deformation problems.
The Solid module provides a layered framework for finite-strain hyperelasticity. It is organized into four sub-layers:
Kinematics
KinematicState encapsulates all finite-strain kinematic quantities at a quadrature point. Set the displacement gradient with setDisplacementGradient(H) and query derived tensors:
| Getter | Returns | Definition |
|---|---|---|
getDisplacementGradient() | Input gradient | |
getDeformationGradient() | ||
getDeformationGradientInverse() | Inverse | |
getDeformationGradientInverseTranspose() | Inverse transpose | |
getRightCauchyGreenTensor() | ||
getLeftCauchyGreenTensor() | ||
getJacobian() | ||
getLogJacobian() | Natural log of Jacobian |
Invariant classes compute derived quantities from the kinematic state:
IsotropicInvariants— , ,FiberInvariants— , (for anisotropic materials; constructed with a fiber direction vector)
Constitutive Point and Input Injection
ConstitutivePoint bundles a KinematicState with an optional Geometry:: (for spatial location) and a tag-based auxiliary data store. Auxiliary data (fiber directions, activation values, etc.) is injected via typed tags:
cp.set<Solid::Tags::FiberDirection>(fiberVec); auto& f = cp.get<Solid::Tags::FiberDirection>(); bool has = cp.has<Solid::Tags::Activation>();
The Input<Derived> CRTP base and InputFunction callable allow the user to inject auxiliary data at each quadrature point during assembly. Integrators accept input via setInput().
Constitutive Laws
Each law derives from HyperElasticLaw<Derived> (CRTP base) and provides:
- A
Cachestruct for precomputed invariant data setCache(cache, cp)— populate the cache from aConstitutivePointgetFirstPiolaKirchhoffStress(P, cache, cp)—getMaterialTangent(dP, cache, cp, dF)— directional derivativegetStrainEnergyDensity(cache, cp)— returns
| Law | Constructor | Energy density |
|---|---|---|
Hooke | Hooke(lambda, mu) | (linear stress-strain; not hyperelastic) |
NeoHookean | NeoHookean(lambda, mu) | |
SaintVenantKirchhoff | SaintVenantKirchhoff(lambda, mu) | |
MooneyRivlin | MooneyRivlin(c1, c2, kappa) |
where is the Green-Lagrange strain, are modified invariants, and is the bulk modulus.
All Lamé-parameter constructors accept and computed from Young's modulus and Poisson's ratio :
Integrators
Both integrators use setDisplacement() to set the linearization point and setInput() for auxiliary data injection. Quadrature order is auto-selected (2× FE order) unless overridden with setQuadratureOrder().
| Integrator | Form type | Assembles | Mathematical expression |
|---|---|---|---|
InternalForce<Law, FES> | Linear (residual) | Internal force vector | |
MaterialTangent<Law, Solution, FES> | Bilinear (tangent) | Tangent stiffness matrix |
Constructor signatures (CTAD-deduced):
Solid::InternalForce residual(law, v); // (law, TestFunction) residual.setDisplacement(u); // set linearization point Solid::MaterialTangent tangent(law, du, v); // (law, TrialFunction, TestFunction) tangent.setDisplacement(u); // set linearization point
Derived Fields
FirstPiolaKirchhoffStress<Law>— Post-process fieldCauchyStress<Law>— Post-processGreenLagrangeStrain— Post-process
Typical Usage
A quasi-static hyperelastic block under gravity using NeoHookean (from examples/Solid/BlockGravity.cpp):
#include <Rodin/Solid.h> #include <Rodin/Solver/NewtonSolver.h> #include <Rodin/Solver/SparseLU.h> // Material from Young's modulus and Poisson's ratio const Real E = 200.0, nu = 0.3; const Real lambda = E * nu / ((1.0 + nu) * (1.0 - 2.0 * nu)); const Real mu = E / (2.0 * (1.0 + nu)); Solid::NeoHookean law(lambda, mu); P1 Vh(mesh, mesh.getSpaceDimension()); // vector P1 space TrialFunction du(Vh); TestFunction v(Vh); GridFunction u(Vh); // accumulated displacement // Build Newton system: K δu = -F_int(u) + F_body Solid::MaterialTangent tangent(law, du, v); tangent.setDisplacement(u); Solid::InternalForce residual(law, v); residual.setDisplacement(u); auto bodyForce = VectorFunction{ Zero(), RealFunction(-10.0) }; Problem newton(du, v); newton = tangent + residual - Integral(bodyForce, v) + DirichletBC(du, VectorFunction{ Zero(), Zero() }).on(bottomBC); SparseLU linearSolver(newton); NewtonSolver solver(linearSolver); solver.setMaxIterations(50) .setAbsoluteTolerance(1e-10) .setRelativeTolerance(1e-8); solver.solve(u);
Namespaces
- namespace Tags
- Standard tag types for auxiliary constitutive data.
Classes
- class BodyForce
- Describes a body force per unit reference volume.
-
template<class LawDerived>class CauchyStress
- Computes the Cauchy stress from a constitutive law.
- class ConstitutivePoint
- Central data bundle for constitutive evaluation at a quadrature point.
- class FiberInvariants
- Fiber invariants for anisotropic hyperelasticity.
-
template<class LawDerived>class FirstPiolaKirchhoffStress
- Computes the first Piola-Kirchhoff stress from a constitutive law.
- class GreenLagrangeStrain
- Computes the Green-Lagrange strain tensor from a kinematic state.
- class Hooke
- Isotropic Hooke's law for linear elasticity.
-
template<class Derived>class HyperElasticLaw
- CRTP base class for hyperelastic constitutive laws.
-
template<class Derived>class Input
- CRTP base class for inputs.
-
template<class LawDerived, class FES>class InternalForce
- Linear form integrator for the internal force vector in hyperelastic problems.
- class IsotropicInvariants
- Isotropic invariants of the right Cauchy-Green tensor.
- class KinematicState
- Kinematic state for finite-strain continuum mechanics.
-
template<class LawDerived, class Solution, class FES>class MaterialTangent
- Local bilinear form integrator for the material tangent stiffness in hyperelastic problems.
- class MooneyRivlin
- Compressible Mooney-Rivlin hyperelastic law.
- class NeoHookean
- Compressible Neo-Hookean hyperelastic law.
- class SaintVenantKirchhoff
- Saint-Venant-Kirchhoff hyperelastic law.
- class TractionForce
- Describes a traction (surface force) boundary condition.
Typedefs
- using InputFunction = std::function<void(ConstitutivePoint&)>
- Type-erased callable for input injection into ConstitutivePoint.
Functions
-
template<class LawDerived, class TestFES>InternalForce(const LawDerived&, const Variational::
TestFunction<TestFES>&) -> InternalForce< LawDerived, TestFES > - CTAD deduction guide for InternalForce.
-
template<class LawDerived, class S, class TrialFES, class TestFES>MaterialTangent(const LawDerived&, const Variational::
TrialFunction<S, TrialFES>&, const Variational:: TestFunction<TestFES>&) -> MaterialTangent< LawDerived, S, TrialFES > - CTAD deduction guide for MaterialTangent.
Typedef documentation
using Rodin:: Solid:: InputFunction = std::function<void(ConstitutivePoint&)>
#include <Rodin/Solid/Local/Input.h>
Type-erased callable for input injection into ConstitutivePoint.
Function documentation
template<class LawDerived, class TestFES>
Rodin:: Solid:: InternalForce(const LawDerived&,
const Variational:: TestFunction<TestFES>&) -> InternalForce< LawDerived, TestFES >
CTAD deduction guide for InternalForce.
template<class LawDerived, class S, class TrialFES, class TestFES>
Rodin:: Solid:: MaterialTangent(const LawDerived&,
const Variational:: TrialFunction<S, TrialFES>&,
const Variational:: TestFunction<TestFES>&) -> MaterialTangent< LawDerived, S, TrialFES >
CTAD deduction guide for MaterialTangent.