template<class LawDerived, class Solution, class FES>
Rodin::Solid::MaterialTangent class

Local bilinear form integrator for the material tangent stiffness in hyperelastic problems.

Template parameters
LawDerived The hyperelastic constitutive law type
Solution The solution type of the trial function
FES The finite element space type

Computes the element-level tangent stiffness matrix:

\[ K^e_{ij} = \int_{K} D\mathbf{P}[\nabla \phi_j] : \nabla \phi_i \, dX \]

Obtains the finite element basis from the FE space via getFiniteElement(), supports configurable quadrature order, and builds a ConstitutivePoint (composed over Geometry::Point) at each quadrature point for constitutive evaluation. An optional Input can inject auxiliary data (fiber directions, activation, etc.) into the ConstitutivePoint at each quadrature point.

Constructors, destructors, conversion operators

template<class S, class TrialFES, class TestFES>
MaterialTangent(const LawDerived& law, const Variational::TrialFunction<S, TrialFES>& u, const Variational::TestFunction<TestFES>& v)
Constructs the material tangent integrator.

Public functions

auto setDisplacement(const GridFunctionType& gf) -> MaterialTangent&
Sets the linearization point (current displacement GridFunction).
auto setQuadratureOrder(size_t order) -> MaterialTangent&
Sets the quadrature order.
auto setInput(InputFunction input) -> MaterialTangent&
Sets an input for auxiliary constitutive data.
auto copy() const -> MaterialTangent* final noexcept
Creates a copy of this integrator.
auto getLaw() const -> const LawType&
Gets the constitutive law.

Function documentation

template<class LawDerived, class Solution, class FES> template<class S, class TrialFES, class TestFES>
Rodin::Solid::MaterialTangent<LawDerived, Solution, FES>::MaterialTangent(const LawDerived& law, const Variational::TrialFunction<S, TrialFES>& u, const Variational::TestFunction<TestFES>& v)

Constructs the material tangent integrator.

Parameters
law The constitutive law (stored by value)
u The trial function
v The test function

template<class LawDerived, class Solution, class FES>
MaterialTangent& Rodin::Solid::MaterialTangent<LawDerived, Solution, FES>::setDisplacement(const GridFunctionType& gf)

Sets the linearization point (current displacement GridFunction).

Parameters
gf Reference to the displacement GridFunction
Returns Reference to this object for chaining

template<class LawDerived, class Solution, class FES>
MaterialTangent& Rodin::Solid::MaterialTangent<LawDerived, Solution, FES>::setQuadratureOrder(size_t order)

Sets the quadrature order.

Parameters
order Polynomial order for exact integration (0 = auto)
Returns Reference to this object for chaining

When set to 0 (the default), the quadrature order is determined automatically from the finite element approximation order as 2 * fe.getOrder().

template<class LawDerived, class Solution, class FES>
MaterialTangent& Rodin::Solid::MaterialTangent<LawDerived, Solution, FES>::setInput(InputFunction input)

Sets an input for auxiliary constitutive data.

Parameters
input A callable with signature void(ConstitutivePoint&)
Returns Reference to this object for chaining

The input is called at each quadrature point after the ConstitutivePoint has been constructed with geometric context and kinematics, allowing injection of fiber directions, activation parameters, region-wise material properties, etc.

template<class LawDerived, class Solution, class FES>
MaterialTangent* Rodin::Solid::MaterialTangent<LawDerived, Solution, FES>::copy() const final noexcept

Creates a copy of this integrator.

Returns Pointer to newly allocated copy