MaterialModelsBase

MaterialModelsBase provides an implementation-agnostic interface for mechanical (stress-strain) material models. This facilitates interchanging material models between researchers by having a light-weight common interface.

Key features of this package are

  • Implement a material model for 3d, use stress states to get, e.g., uniaxial stress, plane stress, plane strain, and more without further implementation.
  • Conversion routine interface to convert e.g. parameters to a vector and back, allowing interfacing with optimization libraries for parameter identification.
  • Differentiation routine interface to allow taking derivatives wrt. material parameters etc. to enable gradient-based optimization for parameter identification.

In general, there are two types of user roles for this package (often the same person takes both roles):

  1. Someone writing code that use material models (e.g. finite element code)
  2. Someone implementing material models

For nr 1, a brief introduction on how to use material models following the MaterialModelsBase.jl interface is given below. For nr 2, a first introduction is provided as a tutorial.

Using material models

This section describes how to use a material model defined according to the described interface. As an example, we use the Zener material defined in the implementation tutorial. For normal usage different materials are defined in a material models package, such as MechanicalMaterialModels.jl.

If we would like to write a function that simulates the uniaxial response of a material model, we can write a function to do that as

using MaterialModelsBase, Tensors
function simulate_uniaxial(m::AbstractMaterial, ϵ11_history, time_history)
    state = initial_material_state(m)
    cache = allocate_material_cache(m)
    stress_state = UniaxialStress()
    σ11_history = zero(ϵ11_history)
    for i in eachindex(ϵ11_history, time_history)[2:end]
        Δt = time_history[i] - time_history[i-1]
        ϵ = SymmetricTensor{2,1}((ϵ11_history[i],))
        σ, dσdϵ, state = material_response(stress_state, m, ϵ, state, Δt, cache)
        σ11_history[i] = σ[1,1]
    end
    return σ11_history
end

To use this function, we define the material parameters for the Zener viscoelastic material model whose implementation is demonstrated here.

material = Zener(;K = 100.0, G0 = 20.0, G1 = 80.0, η1 = 10.0)

And then we define the strain and time history, before running the simulation,

ϵ11_history  = collect(range(0, 0.01; length=100))  # Ramp to 1 %
time_history = collect(range(0, 1; length=100))   # Constant time step
σ11_history  = simulate_uniaxial(material, ϵ11_history, time_history)

We can also plot the stress-strain result,

import CairoMakie as Plt
fig = Plt.Figure()
ax = Plt.Axis(fig[1,1]; xlabel = "strain [%]", ylabel = "stress [MPa]")
Plt.lines!(ax, ϵ11_history * 100, σ11_history)
fig
Example block output

This example used the stress iterations implemented in MaterialModelsBase.jl, see Stress States.

If these do not converge, a NoStressConvergence exception is thrown.

The package also contains the exception NoLocalConvergence, which should be thrown from inside implemented material routines to signal that something didn't converge and that the caller should consider to e.g. reduce the time step or handle the issue in some other way.

API

material_response

The main function is the material_response function that primarily dispatches on the AbstractMaterial input type. Two variants can be called, where the latter allows a reduced stress state, see Stress States for further details.

MaterialModelsBase.material_responseMethod
material_response(
    m::AbstractMaterial, 
    strain::Union{SecondOrderTensor,Vec}, 
    old::AbstractMaterialState, 
    Δt::Union{Number,Nothing}=nothing, 
    cache::AbstractMaterialCache=allocate_material_cache(m), 
    extras::AbstractExtraOutput=NoExtraOutput()
    )

Calculate the stress/traction, stiffness and updated state variables for the material m, given the strain input strain.

Mandatory arguments

  • m: Defines the material and its parameters
  • The second strain input can be, e.g.
    • ϵ::SymmetricTensor{2}: The total small strain tensor at the end of the increment
    • F::Tensor{2}: The total deformation gradient at the end of the increment
    • u::Vec: The deformation at the end of the increment (for cohesive elements)
  • old::AbstractMaterialState: The material state variables at the end of the last converged increment. The material initial material state can be obtained by the initial_material_state function.

Optional positional arguments

When calling the function, the following arguments are optional. When implementing a material, it is not necessary to implement these defaults, but the method signature should allow all arguments to be compatible with libraries relying on the interface. Typically, this can done by using args..., e.g., material_response(m::MyMat, ϵ, state, args...)

  • Δt: The time step in the current increment. Defaults: nothing.
  • cache::AbstractMaterialCache: Cache variables that can be used to avoid allocations during each call to the material_response function. Default: allocate_material_cache(m)
  • extras: Updated with requested extra output. Default: NoExtraOutput (Empty struct)

Outputs

  1. stress, is the stress measure that is energy conjugated to the strain (2nd) input.
  2. stiffness, is the derivative of the stress output wrt. the strain input.
  3. new_state, are the updated state variables
Note

The state given as input should not be mutated. That is, someone calling material_response multiple times with the same input variables should get the same output each time.

Common strain and stress pairs are

  • If the second input is the small strain tensor, $\boldsymbol{\epsilon}$ (ϵ::SymmetricTensor{2}), the outputs are
    • σ::SymmetricTensor{2}: Cauchy stress tensor, $\boldsymbol{\sigma}$
    • dσdϵ::SymmetricTensor{4}: Algorithmic tangent stiffness tensor, $\mathrm{d}\boldsymbol{\sigma}/\mathrm{d}\boldsymbol{\epsilon}$
  • If the second input is the deformation gradient $\boldsymbol{F}$ (F::Tensor{2})`, the outputs are
    • P::Tensor{2}: First Piola-Kirchhoff stress, $\boldsymbol{P}$
    • dPdF::Tensor{4}: Algorithmic tangent stiffness tensor, $\mathrm{d}\boldsymbol{P}/\mathrm{d}\boldsymbol{F}$
  • If the second input is deformation vector, $\boldsymbol{u}$ (u::Vec), the outputs are
    • t::Vec: Traction vector, $\boldsymbol{t}$
    • dtdu::SecondOrderTensor: Algorithmic tangent stiffness tensor, $\mathrm{d}\boldsymbol{t}/\mathrm{d}\boldsymbol{u}$
source
MaterialModelsBase.material_responseMethod
material_response(stress_state::AbstractStressState, m::AbstractMaterial, args...)

To be able to use material models implemented for 3d stress and strain states in lower-dimensional simulations, such as 2d plane stress, MaterialModelsBase.jl provides a set of stress states. For some states, such as plane stress, iterations will be performed to find the correct state. For other states, such as plane strain, the input is only padded with zeros and the out-of-plane components are removed from the output.

For someone implementing a material model, it is also possible to use dispatch on both the stress state and the material to provide an efficient implementation of a reduced stress state. Note that the interface expects the full strain tensor to be given as a fourth output in this case, but it is optional to implement this but such a deviation should be documented as it could cause problems for users of the material implementation.

The arguments are the same as for material_response(::AbstractMaterial). However, both a full and reduced strain input is accepted. For a full strain input, the out-of-plane components are used as an initial guess. For all cases, the full strain tensor giving the desired reduced response is given as a 4th output.

See also ReducedStressState.

source

State variables

To support state variables define

Cache

It is possible to pre-allocate a cache to avoid allocations during the material model calculation, to do this define

MaterialModelsBase.allocate_material_cacheFunction
allocate_material_cache(m::AbstractMaterial)

Return a cache::AbstractMaterialCache that can be used when calling the material to reduce allocations

Defaults to the empty NoMaterialCache()

source

Extra outputs

In some cases, it is necessary to define additional outputs to provide slight variations in the material calculation. In this case, a mutable input can be given as

MaterialModelsBase.AbstractExtraOutputType
AbstractExtraOutput

By allocating an AbstractExtraOutput type, this type can be mutated to extract additional information from the internal calculations in material_response only in cases when this is desired. E.g., when calculating derivatives or for multiphysics simulations. The concrete NoExtraOutput<:AbstractExtraOutput exists for the case when no additional output should be calculated.

source

Exceptions

Finally, the following exceptions are included

MaterialModelsBase.NoLocalConvergenceType
NoLocalConvergence(msg::String)
NoLocalConvergence(args...)

Throw if the material_response routine doesn't converge internally. Other arguments than a single ::String, are converted to String with string

source