Materials

Relevant module: edelweissfe.config.materiallibrary

Materials are defined in EdelweissFE using the *material or *advancedmaterial keyword. Arguments for the *material keyword are

*material : definition of a material

Option

Type

Description

id

string

A unique ID, which is used for referencing the material in EdelweissFE.

name

string

The name of of the material.

datalines

numpy float array

The material properties as a float vector, multiline possible.

provider

string

(optional) The material provider.

The material properties for materials using the *material keyword are assigned as seen here:

Use different non-advanced materials. Example:
*material, name=neohookewaplastic, id=Mat1, provider=edelweissmaterial
91304.34783, 100000., 260, 70, 320, 9

Arguments for the *advancedmaterial keyword are

*advancedmaterial : definition of a material

Option

Type

Description

id

string

A unique ID, which is used for referencing the material in EdelweissFE.

name

string

The name of of the material.

datalines

string

The material properties as a vector of strings, multiline possible.

provider

string

(optional) The material provider.

The material properties for materials using the *advancedmaterial keyword are assigned as seen here:

Use different advanced materials. Example:
*advancedmaterial, name=hyperplasticadvanced, id=hyperplasticadvanced, provider=edelweissmaterial
psi_e='mu/2 * (I1/J**(2/3) - 3) + K/8 * (J**2 + 1/J**2 - 2)'
mu=91304.34783, K=100000., fy0=260, HLin=70, dfy=320, delta=9
a='2,3,4'

If the material provider is not given, marmotmaterial is assumed. Materials are assigned to elements by means of Sections.

Provider marmotmaterial

Relevant module: Marmot.

Provider edelweiss

Relevant module: edelweissfe.materials

Currently available materials with this provider

Name

Description

Keyword

linearelastic

Linear elastic material.

*material

vonmises

Von Mises material.

*material

neohookewa

Neo-Hookean Pence-Gou formulation ‘a’ material.

*material

neohookewb

Neo-Hookean Pence-Gou formulation ‘b’ material.

*material

neohookewc

Neo-Hookean Pence-Gou formulation ‘c’ material.

*material

hyperelasticadvanced

Hyperelastic material with advanced defined energy density function using I1 and J.

*advancedmaterial

hyperelasticadvancedi2extended

Hyperelastic material with advanced defined energy density function using I1, I2, J and C itself.

*advancedmaterial

neohookewaplastic

Neo-Hookean Pence-Gou formulation ‘a’ material with J2 plasticity.

*material

neohookewbplastic

Neo-Hookean Pence-Gou formulation ‘b’ material with J2 plasticity.

*material

neohookewcplastic

Neo-Hookean Pence-Gou formulation ‘c’ material with J2 plasticity.

*material

hyperplasticadvanced

Hyperelastic-plastic material with advanced defined energy density function using I1 and J.

*advancedmaterial

Linear elastic material

This material uses a linear elastic relation between stress and strain. The linear elastic material can be used as 2D plane stress and plane strain material as well as 3D material. This material needs the following material properties in the correct order

  1. E - Elasticity module (Young’s modulus).

  2. \(\mathbf{\nu}\) - Poisson’s ratio.

For the 2D plane strain and 3D material the following law in voigt notation is used

\[\begin{split}\begin{bmatrix}\sigma_{11}\\ \sigma_{22} \\ \sigma_{33}\\ \sigma_{12}\\ \sigma_{23}\\ \sigma_{13}\end{bmatrix} = \frac{E}{(1+\nu)(1-2\nu)} \begin{bmatrix}(1-\nu) & \nu & \nu & 0 & 0 & 0 \\ & (1-\nu) & \nu & 0 & 0 & 0\\ & & (1-\nu) & 0 & 0 & 0\\ & & & \frac{1-2\nu}{2} & 0 & 0\\ & & & & \frac{1-2\nu}{2} & 0\\ \text{symm.}& & & & & \frac{1-2\nu}{2}\end{bmatrix} \begin{bmatrix}\varepsilon_{11}\\ \varepsilon_{22} \\ \varepsilon_{33}\\ \gamma_{12}\\ \gamma_{23}\\ \gamma_{13}\end{bmatrix}\end{split}\]

and for the 2D plane stress material the following law is used

\[\begin{split}\begin{bmatrix}\sigma_{11}\\ \sigma_{22} \\ \sigma_{12}\end{bmatrix} = \frac{E}{1-\nu^2} \begin{bmatrix}1 & \nu & 0\\ & 1 & 0\\ \text{symm.}& & \frac{1-\nu}{2}\end{bmatrix} \begin{bmatrix}\varepsilon_{11}\\ \varepsilon_{22} \\ \gamma_{12}\end{bmatrix}.\end{split}\]

For the second case the third strain component gets calculated by using

\[\varepsilon_{33} = -\frac{\nu}{1-\nu} (\varepsilon_{11}+\varepsilon_{22}).\]
class edelweissfe.materials.linearelastic.linearelastic.LinearElasticMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computePlaneStress(stress, dStress_dStrain, ...)

Computes the stresses for a 2D material with plane stress.

computeStress(stress, dStress_dStrain, ...)

Computes the stresses for a 3D material.

computeStress2D(stress, dStress_dStrain, ...)

Computes the stresses for a 2D material with plane strain.

computeUniaxialStress(stress, ...)

Computes the stresses for a uniaxial stress state.

elasticityMatrix()

Initalize a 3D material elasticity matrix.

elasticityMatrix2D()

Initalize a 2D plane strain material elasticity matrix.

elasticityMatrixPlaneStress()

Initalize a 2D plane stress material elasticity matrix.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computePlaneStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a 3D material.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeStress2D(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a 2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

elasticityMatrix()[source]

Initalize a 3D material elasticity matrix.

Returns

The elasticity matrix.

Return type

np.ndarray

elasticityMatrix2D()[source]

Initalize a 2D plane strain material elasticity matrix.

Returns

The elasticity matrix.

Return type

np.ndarray

elasticityMatrixPlaneStress()[source]

Initalize a 2D plane stress material elasticity matrix.

Returns

The elasticity matrix.

Return type

np.ndarray

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

Von Mises material

This material uses the same linear elastic law as the linear elastic material until plasticity is reached. The von Mises material can be used as 2D plane strain and 3D material. This material needs the following material properties in the correct order

  1. E - Elasticity module (Young’s modulus).

  2. \(\mathbf{\nu}\) - Poisson’s ratio.

  3. \(\mathbf{f_{y0}}\) - Yield stress.

  4. \(\mathbf{H_{lin}}\) - Linear plastic hardening parameter.

  5. \(\mathbf{\Delta f_y}\) - Multiplicator for nonlinear isotropic hardening.

  6. \(\mathbf{\delta}\) - Exponent for nonlinear isotropic hardening.

The law used for nonlinear isotropic hardening in this material is \(f_y(\kappa)=f_{y0}+H_{lin}\kappa+\Delta f_{y} e^{-\delta\kappa}\) with \(\kappa\) as the hardening parameter. Plasticity is reached once the yield function f is \(f(\kappa)=||\mathbf{s}||-\sqrt{\frac{2}{3}}f_y(\kappa) > 0\) with \(\mathbf{s}\) as the deviatoric stress. With plasticity reached the material calculates a new \(\Delta\kappa\) using Newton’s method which is afterwards added to the hardening parameter \(\kappa\). In the end the full material tangent and the back projected stress get calculated. For the back projected stress

\[\mathbf{\sigma} = \mathbf{\sigma}^{trial} - 2G\sqrt{\frac{3}{2}}\Delta\kappa\frac{\mathbf{s}}{||\mathbf{s}||}\]

is used with

\[G = \frac{E}{2(1 + \nu)}\]

and \(\mathbf{s}\) as the deviatoric stress.

class edelweissfe.materials.vonmises.vonmises.VonMisesMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computePlaneStress(stress, dStress_dStrain, ...)

Computes the stresses for a plane stress material.

computeStress(stress, dStress_dStrain, ...)

Computes the stresses for a 3D material.

computeUniaxialStress(stress, ...)

Computes the stresses for a uniaxial stress state.

elasticityMatrix()

Initalize a 3D material elasticity matrix.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computePlaneStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a plane stress material.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a 3D material.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

elasticityMatrix()[source]

Initalize a 3D material elasticity matrix.

Returns

The elasticity matrix.

Return type

np.ndarray

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

Elastic Neo-Hookean W(I1, J) Pence-Gou [a] materials

These materials need the following parameters as input data in the correct order:

  1. \(\mu\) - Shear modulus.

  2. \(\kappa\) - Bulk modulus.

These materials base on a energy density function \(W(I_1, J)\), where \(I_1=\text{trace}(\mathbf{b})\) is the first and \(J=\det(\mathbf{F})\) is the second invariant of the left Cauchy-Green tensor \(\mathbf{b}=\mathbf{FF}^\text{T}\) with the deformation gradient \(\mathbf{F}\), the Kirchhoff stress is then:

\[\tau = 2 \frac{\partial W}{\partial I_1} \mathbf{b} + J\frac{\partial W}{\partial J} \mathbf{I}.\]

The tangent modulus is then given by

\[C^{\tau \text{F}}_{ijkL} = 2b_{ij}\frac{\partial^2 W}{\partial I_1\partial F_{kL}} + 2\frac{\partial W}{\partial I_1} \frac{\partial b_{ij}}{\partial F_{kL}} +\frac{\partial}{\partial F_{kL}}\left[J\frac{\partial W}{\partial J}\right] \delta_{ij}.\]

The \(\mathbf{W_a}\) material is defined by

\[W_a (I_1, J) = \frac{\mu}{2} (I_1 - 3) + \left (\frac{\kappa}{2} - \frac{\mu}{3} \right) (J - 1)^2 - \mu \ln (J).\]

The stress can then be calculated with

\[\tau = \mu \mathbf{b} + \left[\left(\kappa - \frac{2\mu}{3}\right) (J^2-J) - \mu\right] \mathbf{I}\]

and the tangent modulus is

\[C^{\tau \text{F}}_{ijkL} = \mu\frac{\partial b_{ij}}{\partial F_{kL}} + \left(\kappa - \frac{2\mu}{3}\right) (2J^2-J) \delta_{ij} (F^{-1})_{Lk}.\]
class edelweissfe.materials.neohooke.neohookepencegouformulationa.NeoHookeanWaMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

The \(\mathbf{W_b}\) material is defined by

\[W_b (I_1, J) = \frac{\mu}{2} \left(\frac{I_1}{J^{2/3}} - 3\right) + \frac{\kappa}{8} \left (J^2 + \frac{1}{J^2} - 2 \right).\]

The stress can then be calculated with

\[\tau = \frac{\mu}{J^{2/3}} \mathbf{b} + \left[ \frac{\kappa}{4} \left( J^2 - \frac{1}{J^2} \right) -\frac{\mu I_1}{3 J^{2/3}}\right] \mathbf{I}\]

and the tangent modulus is

\[C^{\tau \text{F}}_{ijkL} = \frac{\mu}{J^{2/3}}\frac{\partial b_{ij}}{\partial F_{kL}} - \frac{2\mu}{3J^{2/3}} b_{ij} (F^{-1})_{Lk} - \frac{2\mu}{3J^{2/3}}\delta_{ij} F_{kL} + \left[ \frac{2\mu I_1}{9 J^{2/3}} + \frac{\kappa}{2} \left( J^2 + \frac{1}{J^2} \right) \right] \delta_{ij}\left(F^{-1}\right)_{Lk}.\]
class edelweissfe.materials.neohooke.neohookepencegouformulationb.NeoHookeanWbMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

The \(\mathbf{W_c}\) material is defined by

\[W_c (I_1, J) = \frac{\mu}{2} \left(I_1 - 3\right) + \frac{3\mu^2}{3\kappa - 2\mu} \left(J^{\frac{2}{3} - \frac{\kappa}{\mu}} - 1\right).\]

The stress can then be calculated with

\[\tau = \mu \left(\mathbf{b}-J^{\frac{2}{3} - \frac{\kappa}{\mu}}\mathbf{I}\right)\]

and the tangent modulus is

\[C^{\tau \text{F}}_{ijkL} = \mu\frac{\partial b_{ij}}{\partial F_{kL}} + J^{\frac{2}{3}-\frac{\kappa}{\mu}}\left(\kappa - \frac{2\mu}{3}\right) \delta_{ij} (F^{-1})_{Lk}.\]
class edelweissfe.materials.neohooke.neohookepencegouformulationc.NeoHookeanWcMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

Advanced hyperelastic materials

These materials need the following parameters as input data using keywords:

  1. mu - Shear modulus.

  2. K - Bulk modulus.

  3. psi_e - Energy density function.

  4. a - (Optional) your own parameters in psi_e (assignment with commas and quotation mark: a = '1.,2.,3').

The right Cauchy-Green tensor is defined as \(\mathbf{C}=\mathbf{F}^\text{T}\mathbf{F}\). These materials use the Kirchhoff stress and its corresponding tangent modulus.

Hyperelastic advanced material

This material class allows the user to input their own energy density function, as seen for the Neo-Hookean materials in the last section. This can be achieved by using the extra input parameter W=’f(C)’ in the input file. This material needs \(\kappa\) and \(\mu\) as input data. The energy density function must be a function of the form

\[W(\mu, \kappa, \mathbf{C}, \mathbf{a}) = W(\mu, \kappa, I_1, J, \mathbf{a})\]

with \(I_1=\text{trace}(\mathbf{C})\), \(J=\det(\mathbf{F})\) and the advanced parameters \(\mathbf{a}\) referred to as a[i] for the i-th parameter in the function in the input file.

Note

This material class uses the num-dual-package for automatic differentiation. This package needs to be installed before usage:

pip install num_dual

Hyperelastic advanced I2 extended material

This material class allows the user to input their own energy density function. Compared to the HyperelasticAdvancedMaterial, this material allows the second invariant \(I_2\) and \(\mathbf{C}\) itself to be used in the energy density function. This can be achieved by using the extra input parameter W=’f(C)’ in the input file. This material needs \(\kappa\) and \(\mu\) as input data. The energy density function must be a function of the form

\[W(\mu, \kappa, \mathbf{C}, \mathbf{a}) = W(\mu, \kappa, \mathbf{C}, I_1, I_2, J, \mathbf{a})\]

with \(I_1=\text{trace}(\mathbf{C})\), \(I_2=\frac{1}{2}\left((\text{trace}(\mathbf{C}))^2 - \text{trace}(\mathbf{C}^2)\right)\), \(J=\det(\mathbf{F})\) and the advanced parameters \(\mathbf{a}\) referred to as a[i] for the i-th parameter in the function in the input file. The right Cauchy-Green tensor itself may also be used as an input.

Note

This material class uses the autograd and num-dual-package for automatic differentiation depending on how complicated the energy density function is. These two packages need to be installed before usage:

pip install num_dual && pip install autograd

class edelweissfe.materials.hyperelasticadvanced.hyperelasticadvancedi2extended.HyperelasticAdvancedI2ExtendedMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (dict) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

setEnergyFunction(materialEnergy)

Sets the energy density function for the custom material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

setEnergyFunction(materialEnergy)[source]

Sets the energy density function for the custom material.

Parameters

materialEnergy (str) – Energy density function for the material.

property materialProperties: dict

The properties the material has.

Elastic-plastic Neo-Hookean W(I1, J) Pence-Gou [a] materials

These materials use the same hyperelastic formulations as for the Elastic Neo-Hookean W(I1, J) Pence-Gou [a] materials with plasticity using the spatial Hencky-strain \(\mathbf{e}\) and the same yield and hardening functions as for the Von Mises material added onto them. This material needs the following material properties in the correct order:

  1. \(\mathbf{\mu}\) - Shear modulus.

  2. \(\mathbf{\kappa}\) - Bulk modulus.

  3. \(\mathbf{f_{y0}}\) - Yield stress.

  4. \(\mathbf{H_{lin}}\) - Linear plastic hardening parameter.

  5. \(\mathbf{\Delta f_y}\) - Multiplicator for nonlinear isotropic hardening.

  6. \(\mathbf{\delta}\) - Exponent for nonlinear isotropic hardening.

These three material formulations solve for the residual:

\[\begin{split}\mathbf{R}(\bar{\mathbf{x}}) = \begin{bmatrix}\left(\mathbf{e}_{n+1}^\text{e}\right)_\text{v} - \left(\mathbf{e}_{n+1}^\text{e,trial}\right)_\text{v} + \sqrt{\frac{3}{2}}\Delta\bar{\kappa}\frac{\left(\tau^\text{dev}\right)_\text{v}}{||\tau^\text{dev}||} \\ ||\tau^\text{dev}(\tau(\mathbf{e}^\text{e}_{n+1}))|| - \sqrt{\frac{2}{3}} f_\text{y}(\bar{\kappa}^\text{old}+\Delta\bar{\kappa}) \end{bmatrix}\end{split}\]

and give back the Kirchhoff stress and its corresponding tangent modulus. The full algorithm for these materials can be found in [b], chapter 14.

class edelweissfe.materials.neohookeplastic.neohookepencegouformulationaplastic.NeoHookeanWaPlasticMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

class edelweissfe.materials.neohookeplastic.neohookepencegouformulationbplastic.NeoHookeanWbPlasticMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

class edelweissfe.materials.neohookeplastic.neohookepencegouformulationcplastic.NeoHookeanWcPlasticMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

Advanced hyperelastic-plastic material

This material class allows the user to input their own energy density function. This can be achieved by using the extra input parameter W=’f(C)’ in the input file. This material needs the following material properties in the correct order The algorithm from the section Elastic-plastic Neo-Hookean W(I1, J) Pence-Gou [a] materials is used. For this material, only the elastic part is differentiated by using automatic differentiation, the plastic part still uses analytical derivatives. This material needs the following parameters as input data using keywords:

  1. mu - Shear modulus.

  2. K - Bulk modulus.

  3. fy0 - Yield stress.

  4. HLin - Linear plastic hardening parameter.

  5. dfy - Multiplicator for nonlinear isotropic hardening.

  6. delta - Exponent for nonlinear isotropic hardening.

  7. psi_e - Energy density function.

  8. a - (Optional) your own parameters in psi_e (assignment with commas and quotation mark: a = '1.,2.,3').

The energy density function must be a function of the form

\[W(\mu, \kappa, I_1, J, \mathbf{a})\]

with \(\mathbf{C}\) itself and \(I_2\) being not allowed.

Note

This material class uses the num-dual-package for automatic differentiation. This package needs to be installed before usage:

pip install num_dual

Implementing your own materials

Relevant module for hypoelastic materials: edelweissfe.materials.base.basehypoelasticmaterial

class edelweissfe.materials.base.basehypoelasticmaterial.BaseHypoElasticMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computePlaneStress(stress, dStress_dStrain, ...)

Computes the stresses for a 2D material with plane stress.

computeStress(stress, dStress_dStrain, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computeUniaxialStress(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

abstract assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

abstract computePlaneStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

abstract computeStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

abstract computeUniaxialStress(stress, dStress_dStrain, dStrain, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dStrain (ndarray) – Matrix containing dStress/dStrain.

  • dStrain (ndarray) – Strain vector increment at time step t to t+dTime.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

abstract getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

abstract getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

Relevant module for hyperelastic materials: edelweissfe.materials.base.basehyperelasticmaterial

class edelweissfe.materials.base.basehyperelasticmaterial.BaseHyperElasticMaterial(materialProperties)[source]

Initialize.

Attributes
materialProperties

The properties the material has.

Parameters

materialProperties (ndarray) –

Methods

assignCurrentStateVars(currentStateVars)

Assign new current state vars.

computeKirchhoff(stress, ...)

Computes the stresses for a 3D material/2D material with plane strain.

computePlaneKirchhoff(stress, ...)

Computes the stresses for a 2D material with plane stress.

computeUniaxialKirchhoff(stress, ...)

Computes the stresses for a uniaxial stress state.

getNumberOfRequiredStateVars()

Returns number of needed material state Variables per integration point in the material.

getResult(result)

Get the result, as a persistent view which is continiously updated by the material.

abstract assignCurrentStateVars(currentStateVars)[source]

Assign new current state vars.

Parameters

currentStateVars (ndarray) – Array containing the material state vars.

abstract computeKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 3D material/2D material with plane strain.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

abstract computePlaneKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a 2D material with plane stress.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

abstract computeUniaxialKirchhoff(stress, dStress_dDeformationGradient, deformationGradient, time, dTime)[source]

Computes the stresses for a uniaxial stress state.

Parameters
  • stress (ndarray) – Vector containing the stresses.

  • dStress_dDeformationGradient (ndarray) – Matrix containing dStress/dStrain.

  • deformationGradient (ndarray) – The deformation gradient at time step t.

  • time (float) – Array of step time and total time.

  • dTime (float) – Current time step size.

abstract getNumberOfRequiredStateVars()[source]

Returns number of needed material state Variables per integration point in the material.

Returns

Number of needed material state Vars.

Return type

int

abstract getResult(result)[source]

Get the result, as a persistent view which is continiously updated by the material.

Parameters

result (str) – The name of the result.

Returns

The result.

Return type

float

property materialProperties: ndarray

The properties the material has.

[a] Thomas J. Pence and Kun Gou, “On compressible versions of the incompressible neo-Hookean material”, Mathematics and Mechanics of Solids, 20(2): 157-182, 2015

[b] EA de Souza Neto, D Perić and DRJ Owen, ”Computational methods for plasticity - Theory and applications”, John Wiley & Sons Ltd, Engineering, Swansea University Bay Campus, Fabian Way, Swansea, SA1 8EN, 2008