KASKADE 7 development version
Classes | Functions

Building blocks for viscoplastic problems. More...

Classes

class  Kaskade::Elastomechanics::MaterialLaws::ViscoPlasticEnergy< HyperelasticEnergy >
 An adaptor for using hyperelastic stored energies for viscoplasticity. More...
 

Functions

double Kaskade::Elastomechanics::yieldStrength (std::string const &name)
 Yield strengths of several materials in Pa. More...
 
template<int dim, class Scalar >
Scalar Kaskade::Elastomechanics::MaterialLaws::vonMisesStress (Dune::FieldMatrix< Scalar, dim, dim > const &stress)
 Computes the von Mises equivalent stress \( \sigma_v \). More...
 
template<int dim, class Scalar = double>
Dune::FieldVector< Scalar, dim *(dim+1)/2 > Kaskade::Elastomechanics::MaterialLaws::duvautLionsJ2Flow (Dune::FieldMatrix< Scalar, dim, dim > const &cauchyStress, Dune::FieldMatrix< Scalar, dim *(dim+1)/2, dim *(dim+1)/2 > const &C, double tau, double sigmaY)
 A simple Duvaut-Lions flow rule with \( J_2 \) (von Mises) yield surface. More...
 

Detailed Description

Building blocks for viscoplastic problems.

Function Documentation

◆ duvautLionsJ2Flow()

template<int dim, class Scalar = double>
Dune::FieldVector< Scalar, dim *(dim+1)/2 > Kaskade::Elastomechanics::MaterialLaws::duvautLionsJ2Flow ( Dune::FieldMatrix< Scalar, dim, dim > const &  cauchyStress,
Dune::FieldMatrix< Scalar, dim *(dim+1)/2, dim *(dim+1)/2 > const &  C,
double  tau,
double  sigmaY 
)

A simple Duvaut-Lions flow rule with \( J_2 \) (von Mises) yield surface.

The Duvaut-Lions model of viscoplasticity defines the flow rate, i.e. the time derivative of the viscoplastic strain, as

\[ \dot \epsilon^{\rm vp} = \tau^{-1} C^{-1} (\sigma - P\sigma), \]

where \( P \) is the closest point projector onto the admissible stress state. In \( J_2 \) viscoplasticity, the admissible stresses are characterized by \( \|\sigma\| \le \sqrt{2/3}\,\sigma_Y \), and the boundary of that set is known as von Mises yield surface.

Parameters
cauchyStressthe true stress tensor
Cthe strain to stress matrix in Voigt notation
tauthe relaxation time (>0)
sigmaYthe yield strength (>=0)
Returns
the flow rate \( \dot \epsilon^{\rm vp} \) in Voigt notation

◆ vonMisesStress()

template<int dim, class Scalar >
Scalar Kaskade::Elastomechanics::MaterialLaws::vonMisesStress ( Dune::FieldMatrix< Scalar, dim, dim > const &  stress)

Computes the von Mises equivalent stress \( \sigma_v \).

The equivalent von Mises stress is defined as

\[ 2\sigma_v^2 = (\sigma_{11}-\sigma_{22})^2 + (\sigma_{22}-\sigma_{33})^2 + (\sigma_{33}-\sigma_{11})^2 + 6 (\sigma_{23}^2+\sigma_{31}^2 + \sigma_{12}^2) \]

in terms of the Cauchy stress \( \sigma \).

Template Parameters
dimthe spatial dimension. For d=2, plane stress is assumed.
Scalarthe scalar field type of the stress tensor, usually double

Assuming there are finite element functions lambda and mu giving the material parameters, and u the displacement, a scalar finite element function sv of von Mises equivalent stress values can be obtained as follows:

interpolateGlobally<PlainAverage>(sv,makeFunctionView(u.space(), [&] (auto const& evaluator)
{
ElasticModulus modulus(lambda.value(evaluator.cell(),evaluator.xloc()),
mu.value(evaluator.cell(),evaluator.xloc()));
HyperelasticVariationalFunctional<StVenantKirchhoff<dim>,GreenLagrangeTensor<double,dim>> material(modulus);
material.setLinearizationPoint(u.derivative(evaluator));
return Dune::FieldVector<double,1>(vonMisesStress(material.cauchyStress()));
}));

◆ yieldStrength()

double Kaskade::Elastomechanics::yieldStrength ( std::string const &  name)

Yield strengths of several materials in Pa.

Throws a LookupException if the material is not available in the data base.