MarmotFiniteStrainMechanicsCore
-
template<int nDim>
struct AlgorithmicModuli - #include <MarmotMaterialFiniteStrain.h>
Algorithmic tangent moduli of a material.
Contains the algorithmic tangent moduli \(\frac{\partial \boldsymbol{\tau}}{\partial \boldsymbol{F}}\) with respect to the deformation gradient \(\boldsymbol{F}\).
- Template Parameters:
nDim – Number of spatial dimensions (2 or 3).
-
template<int nDim>
struct ConstitutiveResponse - #include <MarmotMaterialFiniteStrain.h>
Constitutive response of a material at given state.
Contains stress, density and elastic energy density.
- Template Parameters:
nDim – Number of spatial dimensions (2 or 3).
-
template<int nDim>
struct Deformation - #include <MarmotMaterialFiniteStrain.h>
Represents the deformation state of a material.
This struct holds the deformation gradient \(\boldsymbol{F}\) which describes the local deformation of a material.
-
class MarmotMaterialFiniteStrain : public MarmotMaterial
Inheritance diagram for MarmotMaterialFiniteStrain:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"2" [label="MarmotMaterial" tooltip="MarmotMaterial"]
"1" [label="MarmotMaterialFiniteStrain" tooltip="MarmotMaterialFiniteStrain" fillcolor="#BFBFBF"]
"1" -> "2" [dir=forward tooltip="public-inheritance"]
}](../../_images/graphviz-bbb6c1377508c8f915f56c319c8a82b6daa3a1e0.png)
Collaboration diagram for MarmotMaterialFiniteStrain:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"2" [label="MarmotMaterial" tooltip="MarmotMaterial"]
"1" [label="MarmotMaterialFiniteStrain" tooltip="MarmotMaterialFiniteStrain" fillcolor="#BFBFBF"]
"1" -> "2" [dir=forward tooltip="public-inheritance"]
}](../../_images/graphviz-bbb6c1377508c8f915f56c319c8a82b6daa3a1e0.png)
Abstract basic class for mechanical materials in the finite strain regime.
Public Functions
-
virtual void computeStress(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3>&, const TimeIncrement&) = 0
Updates the material state.
Computes the Kirchhoff \(\boldsymbol{\tau}\) stress from the deformation gradient \(\boldsymbol{F}\) and the time increment \(\Delta t\). It further updates the mass density \(\rho\) and the elastic energy density. Additionally, computes the algorithmic tangent moduli \(\frac{\partial \boldsymbol{\tau}}{\partial \boldsymbol{F}}\).
- Parameters:
response – [inout] ConstitutiveResponse instance
tangents – [out] AlgorithmicModuli instance
deformation – [in] Deformation instance
timeIncrement – [in] TimeIncrement instance
-
virtual void computeStress(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3> &deformation, const TimeIncrement &timeIncrement, const std::tuple<double, double, double> &eigenDeformation)
Computes the Kirchhoff stress given the deformation, time increment, and eigen deformation.
- Parameters:
response – [inout] ConstitutiveResponse instance
tangents – [out] AlgorithmicModuli instance
deformation – [in] Deformation instance
timeIncrement – [in] TimeIncrement instance
eigenDeformation – [in] Tuple representing eigen deformation in each spatial direction.
-
virtual void computePlaneStrain(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &algorithmicModuli, const Deformation<3> &deformation, const TimeIncrement &timeIncrement)
Compute stress under plane strain conditions.
It uses the general 3D computeStress function for a plane strain Deformation. The algorithmic tangent is modified according to plane strain conditions.
- Parameters:
response – [inout] ConstitutiveResponse instance
algorithmicModuli – [out] AlgorithmicModuli instance
deformation – [in] Deformation instance
timeIncrement – [in] TimeIncrement instance
-
virtual void computePlaneStrain(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &algorithmicModuli, const Deformation<3> &deformation, const TimeIncrement &timeIncrement, const std::tuple<double, double, double> &eigenDeformation)
Compute stress under plane strain conditions with eigen deformation.
It uses the general 3D computeStress function for a plane strain Deformation. The algorithmic tangent is modified according to plane strain conditions.
- Parameters:
response – [inout] ConstitutiveResponse instance
algorithmicModuli – [out] AlgorithmicModuli instance
deformation – [in] Deformation instance
timeIncrement – [in] TimeIncrement instance
eigenDeformation – [in] Tuple representing eigen deformation in each spatial direction.
-
virtual void computePlaneStress(ConstitutiveResponse<2> &response, AlgorithmicModuli<2> &algorithmicModuli, const Deformation<2> &deformation, const TimeIncrement &timeIncrement)
Compute stress under plane stress conditions.
It uses the general 3D computeStress function and iteratively finds the out-of-plane deformation. The algorithmic tangent is modified according to plane stress conditions.
- Parameters:
response – [inout] ConstitutiveResponse instance
algorithmicModuli – [out] AlgorithmicModuli instance
deformation – [in] Deformation instance
timeIncrement – [in] TimeIncrement instance
-
std::tuple<double, double, double> findEigenDeformationForEigenStress(const std::tuple<double, double, double> &initialGuess, const std::tuple<double, double, double> &eigenStress)
Find the eigen deformation that corresponds to a given eigen stress.
This function iteratively finds the eigen deformation that corresponds to a given eigen stress. This is used e.g. for geostatic stress initialization.
- Parameters:
initialGuess – Initial guess for the eigen deformation.
eigenStress – Target eigen stress.
- Returns:
Eigen deformation that corresponds to the given eigen stress.
-
virtual void computeStress(ConstitutiveResponse<3> &response, AlgorithmicModuli<3> &tangents, const Deformation<3>&, const TimeIncrement&) = 0
-
struct TimeIncrement
- #include <MarmotMaterialFiniteStrain.h>
Represents a time increment in a simulation.
This struct holds information about the current time and the time step size (dT).
-
namespace Fastor
-
namespace FastorIndices
-
namespace FastorStandardTensors
-
namespace Marmot
-
namespace ContinuumMechanics
-
namespace DeformationMeasures
Functions
-
template<typename T>
Tensor33t<T> rightCauchyGreen(const Tensor33t<T> &F) Computes the right Cauchy-Green tensor \(\boldsymbol{C}\).
The right Cauchy-Green deformation tensor is defined as:
\[ \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \]or in index notation:\[ C_{IJ} = F_{iI} F_{iJ} \]- Template Parameters:
T – Scalar type (e.g., float, double, autodiff::dual)
- Parameters:
F – Deformation gradient tensor
- Returns:
The right Cauchy-Green tensor
-
template<typename T>
Tensor33t<T> leftCauchyGreen(const Tensor33t<T> &F) Computes the left Cauchy-Green tensor \(\boldsymbol{b}\).
The left Cauchy-Green deformation tensor is defined as:
\[ \boldsymbol{b} = \boldsymbol{F} \boldsymbol{F}^T \]or in index notation:\[ b_{ij} = F_{iI} F_{jI} \]- Template Parameters:
T – Scalar type (e.g., float, double, autodiff::dual)
- Parameters:
F – Deformation gradient tensor
- Returns:
The left Cauchy-Green tensor b
-
template<typename T>
-
namespace FirstOrderDerived
Functions
-
template<typename T>
std::pair<Tensor33t<T>, Tensor3333t<T>> rightCauchyGreen(const Tensor33t<T> &F) Computes the right Cauchy-Green tensor \(\boldsymbol{C}\) and its derivative with respect to \(\boldsymbol{F}\).
The right Cauchy-Green deformation tensor is defined as:
\[ \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \]or in index notation:\[ C_{IJ} = F_{iI} F_{iJ} \]The derivative of \(\boldsymbol{C}\) with respect to \(\boldsymbol{F}\) is given by:\[ \frac{\partial C_{IJ}}{\partial F_{kK}} = \delta_{IK} F_{kJ} + F_{kI} \delta_{JK} \]- Template Parameters:
T – Scalar type (e.g., float, double, autodiff::dual)
- Parameters:
F – Deformation gradient tensor
- Returns:
A pair containing the right Cauchy-Green tensor C and its derivative with respect to F
-
template<typename T>
std::pair<Tensor33t<T>, Tensor3333t<T>> leftCauchyGreen(const Tensor33t<T> &F) Computes the left Cauchy-Green tensor \(\boldsymbol{b}\) and its derivative with respect to \(\boldsymbol{F}\).
The left Cauchy-Green deformation tensor is defined as:
\[ \boldsymbol{b} = \boldsymbol{F} \boldsymbol{F}^T \]or in index notation:\[ b_{ij} = F_{iI} F_{jI} \]The derivative of \(\boldsymbol{b}\) with respect to \(\boldsymbol{F}\) is given by:\[ \frac{\partial b_{ij}}{\partial F_{kK}} = \delta_{ik} F_{jK} + F_{iK} \delta_{jk} \]- Template Parameters:
T – Scalar type (e.g., float, double, autodiff::dual)
- Parameters:
F – Deformation gradient tensor
- Returns:
A pair containing the left Cauchy-Green tensor b and its derivative with respect to F
-
template<typename T>
-
namespace EnergyDensityFunctions
Functions
-
template<typename T>
T PenceGouPotentialA(const Tensor33t<T> &C, const double K, const double G) Hyperelastic Energy Density Function Wa acc. Pence & Gou (2015), Eq. (2.11)
The energy density function \(W_a\) is given as
\[ W_a = \frac{G}{2} (I_1 - 3) + \left(\frac{K}{2} - \frac{G}{3}\right) (J - 1)^2 - G \ln(J) \]where \( I_1 = \text{tr}(\boldsymbol{C}) \) is the first invariant of the right Cauchy-Green tensor \( \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \), \( J = \sqrt{\det(\boldsymbol{C})} = \det(\boldsymbol{F}) \) is the determinant of the deformation gradient, and \( K, G \) are the bulk and shear modulus, respectively.- Template Parameters:
T – Scalar type, e.g. double, float, etc.
- Parameters:
C – Right Cauchy-Green tensor
K – Bulk modulus
G – Shear modulus
- Returns:
Energy density
-
template<typename T>
T PenceGouPotentialB(const Tensor33t<T> &C, const double K, const double G) Hyperelastic Energy Density Function Wb acc. Pence & Gou (2015), Eq. (2.12)
The energy density function \(W_b\) is given as
\[ W_b = \frac{K}{8} \left(J - \frac{1}{J}\right)^2 + \frac{G}{2} \left(I_1 J^{-\frac{2}{3}} - 3\right) \]where \( I_1 = \text{tr}(\boldsymbol{C}) \) is the first invariant of the right Cauchy-Green tensor \( \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \), \( J = \sqrt{\det(\boldsymbol{C})} = \det(\boldsymbol{F}) \) is the determinant of the deformation gradient, and \( K, G \) are the bulk and shear modulus, respectively.- Template Parameters:
T – Scalar type, e.g. double, float, etc.
- Parameters:
C – Right Cauchy-Green tensor
K – Bulk modulus
G – Shear modulus
- Returns:
Energy density
-
template<typename T>
T PenceGouPotentialC(const Tensor33t<T> &C, const double K, const double G) Hyperelastic Energy Density Function Wc acc. Pence & Gou (2015), Eq. (2.13)
The energy density function \(W_c\) is given as
\[ W_c = \frac{G}{2} (I_1 - 3) + \frac{3 G^2}{3 K - 2 G} \left(J^{\frac{2}{3} - \frac{K}{G}} - 1\right) \]where \( I_1 = \text{tr}(\boldsymbol{C}) \) is the first invariant of the right Cauchy-Green tensor \( \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \), \( J = \sqrt{\det(\boldsymbol{C})} = \det(\boldsymbol{F}) \) is the determinant of the deformation gradient, and \( K, G \) are the bulk and shear modulus, respectively.- Template Parameters:
T – Scalar type, e.g. double, float, etc.
- Parameters:
C – Right Cauchy-Green tensor
K – Bulk modulus
G – Shear modulus
- Returns:
Energy density
-
template<typename T>
-
namespace FirstOrderDerived
Functions
-
template<typename T>
std::tuple<T, Tensor33t<T>> PenceGouPotentialB(const Tensor33t<T> &C, const double K, const double G) Hyperelastic Energy Density Function Wb acc. Pence & Gou (2015), Eq. (2.12) and its first derivative w.r.t. C.
The energy density function \(W_b\) is given as
\[ W_b = \frac{K}{8} \left(J - \frac{1}{J}\right)^2 + \frac{G}{2} \left(I_1 J^{-\frac{2}{3}} - 3\right) \]where \( I_1 = \text{tr}(\boldsymbol{C}) \) is the first invariant of the right Cauchy-Green tensor \( \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \), \( J = \sqrt{\det(\boldsymbol{C})} = \det(\boldsymbol{F}) \) is the determinant of the deformation gradient, and \( K, G \) are the bulk and shear modulus, respectively.Additionally, the first derivative w.r.t. C is computed as
\[ \frac{\partial W_b}{\partial \boldsymbol{C}} = \frac{\partial W_b}{\partial J} \frac{\partial J}{\partial \boldsymbol{C}} + \frac{\partial W_b}{\partial I_1} \frac{\partial I_1}{\partial \boldsymbol{C}} \]where\[ \frac{\partial J}{\partial \boldsymbol{C}} = \frac{1}{2} J \boldsymbol{C}^{-1} \]and\[ \frac{\partial I_1}{\partial \boldsymbol{C}} = \boldsymbol{I} \]- Template Parameters:
T – Scalar type, e.g. double, float, etc.
- Parameters:
C – Right Cauchy-Green tensor
K – Bulk modulus
G – Shear modulus
- Returns:
A tuple containing energy density and its first derivative w.r.t. C
-
template<typename T>
-
namespace SecondOrderDerived
Functions
-
template<typename T>
std::tuple<T, Tensor33t<T>, Tensor3333t<T>> PenceGouPotentialB(const Tensor33t<T> &C, const double K, const double G) Hyperelastic Energy Density Function Wb acc. Pence & Gou (2015), Eq. (2.12) and its first and second derivative w.r.t. C.
The energy density function \(W_b\) is given as
\[ W_b = \frac{K}{8} \left(J - \frac{1}{J}\right)^2 + \frac{G}{2} \left(I_1 J^{-\frac{2}{3}} - 3\right) \]where \( I_1 = \text{tr}(\boldsymbol{C}) \) is the first invariant of the right Cauchy-Green tensor \( \boldsymbol{C} = \boldsymbol{F}^T \boldsymbol{F} \), \( J = \sqrt{\det(\boldsymbol{C})} = \det(\boldsymbol{F}) \) is the determinant of the deformation gradient, and \( K, G \) are the bulk and shear modulus, respectively.Additionally, the first and second derivative w.r.t. \(\boldsymbol{C}\), i.e.,
\[\frac{\partial W_b}{\partial \boldsymbol{C}} \quad \text{and} \quad \frac{\partial^2 W_b}{\partial \boldsymbol{C} \partial \boldsymbol{C}} \]are computed.- Template Parameters:
T – Scalar type, e.g. double, float, etc.
- Parameters:
C – Right Cauchy-Green tensor
K – Bulk modulus
G – Shear modulus
- Returns:
A tuple containing energy density, its first and second derivative w.r.t. C
-
template<typename T>
-
namespace FiniteStrain
-
namespace Plasticity
-
namespace FlowIntegration
Functions
-
template<typename T>
Tensor33t<T> exponentialMap(const Tensor33t<T> &dGp) Computes the incremental plastic deformation gradient from the plastic velocity gradient using the exponential map.
The exponential map is given by
\[ \Delta F_{\bar I J} = \exp\left(\Delta \lambda \frac{\partial g}{\partial M_{\bar K L}}\right)_{J \bar{I}} \]where \( \Delta \lambda \frac{\partial g}{\partial M_{\bar K L}} \) is the incremental plastic velocity gradient.- Template Parameters:
T – Scalar type.
- Parameters:
dGp – Incremental plastic velocity gradient.
- Returns:
Incremental plastic deformation gradient.
-
template<typename T>
-
namespace FirstOrderDerived
Functions
-
std::pair<Tensor33d, Tensor3333d> explicitIntegration(const Tensor33d &deltaGp)
Computes the incremental plastic deformation gradient from the plastic velocity gradient using a first order approximation. The first order approximation is given by
\[ \Delta F_{\bar I J} = \left( \delta_{\bar K L} + \Delta \lambda \frac{\partial f}{\partial M_{\bar K L}}\right)_{J \bar{I}} \]where \( \Delta \lambda \frac{\partial f}{\partial M_{\bar K L}} \) is the incremental plastic velocity gradient.- Parameters:
deltaGp – Incremental plastic velocity gradient.
- Returns:
A pair of the incremental plastic deformation gradient and its derivative w.r.t. the plastic velocity gradient.
-
std::pair<Tensor33d, Tensor3333d> exponentialMap(const Tensor33d &deltaGp)
Computes the incremental plastic deformation gradient from the plastic velocity gradient using the exponential map. The exponential map is given by
\[ \Delta F_{\bar I J} = \exp\left(\Delta \lambda \frac{\partial f}{\partial M_{\bar I J}}\right)_{J \bar{I}} \]where \( \Delta \lambda \frac{\partial f}{\partial M_{\bar I J}} \) is the incremental plastic velocity gradient.- Parameters:
deltaGp – Incremental plastic velocity gradient.
- Returns:
A pair of the incremental plastic deformation gradient and its derivative w.r.t. the plastic velocity gradient.
-
std::pair<Tensor33d, Tensor3333d> explicitIntegration(const Tensor33d &deltaGp)
-
namespace StressMeasures
Functions
-
template<typename T>
Tensor33t<T> KirchhoffStressFromPK2(const Tensor33t<T> &PK2, const Tensor33t<T> &F) Computes the Kirchhoff stress from the 2nd Piola-Kirchhoff stress and the deformation gradient.
The Kirchhoff stress is computed as
\[ \boldsymbol{\tau} = \boldsymbol{F} \boldsymbol{S} \boldsymbol{F}^T \]or in index notation\[\tau_{ij} = F_{iI} S_{IJ} F_{jJ}. \]- Template Parameters:
T – Scalar type, e.g. double, float
- Parameters:
PK2 – 2nd Piola-Kirchhoff stress
F – Deformation gradient
- Returns:
Kirchhoff stress
-
template<typename T>
-
namespace FirstOrderDerived
Functions
-
template<typename T>
std::tuple<Tensor33t<T>, Tensor3333t<T>, Tensor3333t<T>> KirchhoffStressFromPK2(const Tensor33t<T> &PK2, const Tensor33t<T> &F) Computes the Kirchhoff stress from the 2nd Piola-Kirchhoff stress and the deformation gradient and the respective partial derivatives.
The Kirchoff stress is computed as
\[ \boldsymbol{\tau} = \boldsymbol{F} \boldsymbol{S} \boldsymbol{F}^T \]or in index notation\[\tau_{ij} = F_{iI} S_{IJ} F_{jJ}. \]Additionally, the derivatives with respect to \( \boldsymbol{S} \) and \( \boldsymbol{F} \) are computed in index notation as\[ \frac{\partial \tau_{ij}}{\partial S_{IJ}} = F_{iI} F_{jJ} \]and\[ \frac{\partial \tau_{ij}}{\partial F_{kL}} = \delta_{ik} S_{LJ} F_{jJ} + F_{iI} S_{IL} \delta_{jk} \]- Template Parameters:
T – Scalar type, e.g. double, float
- Parameters:
PK2 – 2nd Piola-Kirchhoff stress
F – Deformation gradient
- Returns:
A tuple containing Kirchhoff stress and its derivatives w.r.t \(\boldsymbol{S}\) and \(\boldsymbol{F}\)
-
template<typename T>
- file MarmotDeformationMeasures.h
- #include “Marmot/MarmotFastorTensorBasics.h”
Include dependency graph for MarmotDeformationMeasures.h:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotDeformationMeasures.h" tooltip="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotDeformationMeasures.h" fillcolor="#BFBFBF"]
"2" [label="Marmot/MarmotFastorTensorBasics.h" tooltip="Marmot/MarmotFastorTensorBasics.h"]
"1" -> "2" [dir=forward tooltip="include"]
}](../../_images/graphviz-ac09e225d3d8319c33bdfc9aa7d6ea630c4ca8cd.png)
- file MarmotEnergyDensityFunctions.h
- #include “Fastor/Fastor.h”#include “Marmot/MarmotFastorTensorBasics.h”
Include dependency graph for MarmotEnergyDensityFunctions.h:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotEnergyDensityFunctions.h" tooltip="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotEnergyDensityFunctions.h" fillcolor="#BFBFBF"]
"2" [label="Fastor/Fastor.h" tooltip="Fastor/Fastor.h"]
"3" [label="Marmot/MarmotFastorTensorBasics.h" tooltip="Marmot/MarmotFastorTensorBasics.h"]
"1" -> "2" [dir=forward tooltip="include"]
"1" -> "3" [dir=forward tooltip="include"]
}](../../_images/graphviz-ec3973084b4b414d358c2d8e274b38a65cb60dd8.png)
- file MarmotFiniteStrainPlasticity.h
- #include “Marmot/MarmotFastorTensorBasics.h”#include “Marmot/MarmotTensorExponential.h”
Include dependency graph for MarmotFiniteStrainPlasticity.h:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotFiniteStrainPlasticity.h" tooltip="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotFiniteStrainPlasticity.h" fillcolor="#BFBFBF"]
"2" [label="Marmot/MarmotFastorTensorBasics.h" tooltip="Marmot/MarmotFastorTensorBasics.h"]
"3" [label="Marmot/MarmotTensorExponential.h" tooltip="Marmot/MarmotTensorExponential.h"]
"1" -> "2" [dir=forward tooltip="include"]
"1" -> "3" [dir=forward tooltip="include"]
}](../../_images/graphviz-99506fe367f3e4090fa489c09a60387d6c0ccda5.png)
- file MarmotMaterialFiniteStrain.h
- #include “Fastor/Fastor.h”#include “Marmot/MarmotMaterial.h”#include <Fastor/tensor/Tensor.h>
Include dependency graph for MarmotMaterialFiniteStrain.h:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotMaterialFiniteStrain.h" tooltip="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotMaterialFiniteStrain.h" fillcolor="#BFBFBF"]
"2" [label="Fastor/Fastor.h" tooltip="Fastor/Fastor.h"]
"4" [label="Fastor/tensor/Tensor.h" tooltip="Fastor/tensor/Tensor.h"]
"3" [label="Marmot/MarmotMaterial.h" tooltip="Marmot/MarmotMaterial.h"]
"1" -> "2" [dir=forward tooltip="include"]
"1" -> "3" [dir=forward tooltip="include"]
"1" -> "4" [dir=forward tooltip="include"]
}](../../_images/graphviz-4efef05669b7358751411f15049b80b6ae98c426.png)
- file MarmotStressMeasures.h
- #include “Marmot/MarmotFastorTensorBasics.h”#include <Fastor/tensor/Tensor.h>
Include dependency graph for MarmotStressMeasures.h:
![digraph {
graph [bgcolor="#00000000"]
node [shape=rectangle style=filled fillcolor="#FFFFFF" font=Helvetica padding=2]
edge [color="#1414CE"]
"1" [label="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotStressMeasures.h" tooltip="/home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot/MarmotStressMeasures.h" fillcolor="#BFBFBF"]
"3" [label="Fastor/tensor/Tensor.h" tooltip="Fastor/tensor/Tensor.h"]
"2" [label="Marmot/MarmotFastorTensorBasics.h" tooltip="Marmot/MarmotFastorTensorBasics.h"]
"1" -> "2" [dir=forward tooltip="include"]
"1" -> "3" [dir=forward tooltip="include"]
}](../../_images/graphviz-a62a139e9b93e2f5d327a753aa86d78408e5f24f.png)
- dir /home/runner/work/Marmot/Marmot/modules/core
- dir /home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include
- dir /home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore/include/Marmot
- dir /home/runner/work/Marmot/Marmot/modules/core/MarmotFiniteStrainMechanicsCore
- dir /home/runner/work/Marmot/Marmot/modules