Implementing viscoelasticity with state variables #666
-
Hi @adtzlr , I really like the def ψ(x, mu1, mu2, alph1, alph2, a1, a2, m1, m2, K1, K2, bta1, bta2, eta0, etaInf, λ, dt):
F, p, ζ = x[0], x[1], x[2]
J = det(F)
Js = (λ + p + sqrt((λ + p) ** 2 + 4 * λ * (mu1 + mu2 + m1 + m2))) / (2 * λ)
C = F.T @ F
Cvinv = inv(ζ)
C_Cvinv = C @ Cvinv
Iv1 = trace(ζ)
I1 = trace(C)
Ie1 = trace(C_Cvinv)
W_I1 = 3 ** (1 - alph1) / 2 / alph1 * mu1 * (I1**alph1 - 3**alph1) + 3 ** (
1 - alph2
) / 2 / alph2 * mu2 * (I1**alph2 - 3**alph2)
W_Ie1 = 3 ** (1 - a1) / 2 / a1* m1 * (Ie1**a1- 3**a1) + 3 ** (
1 - a2
) / 2 / a2* m2* (Ie1**a2 - 3**a2)
ψs = p * (J - Js) + W_I1 - (mu1 + mu2 + m1 + m2) * log(Js) + λ / 2 * (Js - 1) ** 2
ζ_new = ζ + rk_update() # this update of the state variable
return gradient(ψs, F), gradient(ψs, p), ζ_new where the update function for the internal state variables ζ (if I understand these correctly) can be an explicit (like Runge-Kutta) or implicit (backward euler). The latter would need a Newton solver embedded within the constitutive relation and hence I don't want to get into that yet (or maybe its easily implementable). For now, consider that it is fully explicit, i.e. knowing at time-step def rk_update(C, Cv, Cn, Cvn):
pass Concretely, I am planning to implement the following constitutive behavior with and the ODE G(t, Cv) can be discretized using a time-discretization of your choice (in this case, say RK4 or RK5) Any thoughts on how this could be achieved within felupe with/without using matadi. |
Beta Was this translation helpful? Give feedback.
Replies: 7 comments 13 replies
-
And the constitutive model is here: https://pamies.cee.illinois.edu/Publications_files/CRM_2016.pdf Specifically, equations (36)-(37) along with (40)-(41) |
Beta Was this translation helpful? Give feedback.
-
Hi @bhaveshshrimali, you have to split / concatenate the state variables Edit: I think you need only one state variable in this model, Code:import matadi as mat
from matadi.math import (
det,
trace,
inv,
sqrt,
gradient,
log,
astensor,
asvoigt,
vertsplit,
vertcat,
)
F = mat.Variable("F", 3, 3)
p = mat.Variable("p", 1, 1)
ζ = mat.Variable("ζ", 12, 1)
def ψ(
x, mu1, mu2, alph1, alph2, a1, a2, m1, m2, K1, K2, bta1, bta2, eta0, etaInf, λ, dt
):
F, p, ζ = x[0], x[1], x[2]
Cn, Cvn = vertsplit(ζ, [0, 6, 12])
Cn = astensor(Cn)
Cvn = astensor(Cvn)
J = det(F)
Js = (λ + p + sqrt((λ + p) ** 2 + 4 * λ * (mu1 + mu2 + m1 + m2))) / (2 * λ)
C = F.T @ F
Cvinv = inv(Cv)
C_Cvinv = C @ Cvinv
Iv1 = trace(Cv)
I1 = trace(C)
Ie1 = trace(C_Cvinv)
W_I1 = 3 ** (1 - alph1) / 2 / alph1 * mu1 * (I1**alph1 - 3**alph1) + 3 ** (
1 - alph2
) / 2 / alph2 * mu2 * (I1**alph2 - 3**alph2)
W_Ie1 = 3 ** (1 - a1) / 2 / a1 * m1 * (Ie1**a1 - 3**a1) + 3 ** (
1 - a2
) / 2 / a2 * m2 * (Ie1**a2 - 3**a2)
ψs = p * (J - Js) + W_I1 - (mu1 + mu2 + m1 + m2) * log(Js) + λ / 2 * (Js - 1) ** 2
Cv = Cvn + rk_update(C, Cn, Cvn) # this update of the state variable
ζ_new = vertcat(asvoigt(C), asvoigt(Cv))
return gradient(ψs, F), gradient(ψs, p), ζ_new (this is still an incomplete function, I just added the state variable handling) CasADi (the backend of matadi) also comes with an ode function. I can't help you here as I've never used it but hopefully their docs give you some idea. Otherwise proceed it as you described it - use your own update function of state variables for the evolution equation. For my understanding, you have to do it at the beginning (like in my code block). Edit: I see, due to the explicit Runge-Kutta update the update relies on the old state variables. |
Beta Was this translation helpful? Give feedback.
-
I tried to implement it with pure FElupe using If you're interested, here is the code. Syntax is a bit different compared to matadi.import numpy as np
import felupe as fem
import tensortrax.math as tm
def viscoelastic_two_potential(C, ζn, dt, μ, α, m, a, η0, ηinf, β, K):
# extract old state variables
Cvn = tm.special.from_triu_1d(ζn, like=C)
# third invariant of the right Cauchy-Green deformation tensor
# first invariant of distortional part of the right Cauchy-Green deformation tensor
I3 = tm.linalg.det(C)
I1 = tm.trace(C) * I3 ** (-1 / 3)
def invariants(Cv):
"Invariants of the elastic and viscous parts of the deformation."
Ce = C @ tm.linalg.inv(Cv)
I1e = I3 ** (-1 / 3) * tm.trace(Ce)
I1v = tm.trace(Cv)
I2e = I3 ** (-2 / 3) * (I1e**2 - tm.trace(Ce @ Ce)) / 2
return I1e, I1v, I2e
def G(Cv):
"Evolution equation."
I1e, I1v, I2e = invariants(Cv)
J2Neq = (I1e**2 / 3 - I2e) * tm.sum(
[3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
) ** 2
u = [3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
ηK = ηinf + (η0 - ηinf + K[0] * (I1v ** β[0] - 3 ** β[0])) / (
1 + (K[1] * J2Neq) ** β[1]
)
return tm.sum(u) / ηK * (C - I1e / 3 * tm.linalg.inv(Cv))
# update the state variables
k1 = G(Cvn)
k2 = G(Cvn + k1 * dt / 2)
k3 = G(Cvn + (3 * k1 + k2) * dt / 16)
k4 = G(Cvn + k3 * dt / 2)
k5 = G(Cvn + 3 * (-k2 + 2 * k3 + 3 * k4) * dt / 16)
k6 = G(Cvn + (k1 + 4 * k2 + 6 * k3 - 12 * k4 + 8 * k5) * dt / 7)
Cv = Cvn + dt / 90 * (7 * k1 + 32 * k3 + 12 * k4 + 32 * k5 + 7 * k6)
I1e, I1v, I2e = invariants(Cv)
# strain energy functions of the equilibrium and non-equilibrium parts
p = [0, 1]
ψEq = [3 ** (1 - α[r]) / (2 * α[r]) * μ[r] * (I1 ** α[r] - 3 ** α[r]) for r in p]
ψNEq = [3 ** (1 - a[r]) / (2 * a[r]) * m[r] * (I1e ** a[r] - 3 ** a[r]) for r in p]
return tm.sum(ψEq) + tm.sum(ψNEq), tm.special.triu_1d(Cv)
λ = fem.math.linsteps([1, 2, 1.5], num=[10, 5])
dλ = np.abs(np.diff(λ)).mean()
dλdt = 0.00023
umat = fem.Hyperelastic(
viscoelastic_two_potential,
nstatevars=6,
dt=dλ / dλdt,
μ=[1.08, 0.0170],
α=[0.26, 7.68],
m=[1.57, 0.59],
a=[-10, 7.53],
η0=2.11,
ηinf=0.1,
β=[3.0, 1.929],
K=[442, 1289.49],
)
mesh = fem.Cube(n=2)
region = fem.RegionHexahedron(mesh)
field = fem.FieldContainer([fem.Field(region, dim=3)])
boundaries, loadcase = fem.dof.uniaxial(field, clamped=False)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)
solid.results.statevars[[0, 3, 5]] += 1.0
step = fem.Step(items=[solid], ramp={boundaries["move"]: λ - 1}, boundaries=boundaries)
job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate()
fig, ax = job.plot(
xlabel="Displacement $u$ in mm $\longrightarrow$",
ylabel="Normal Force $F$ in N $\longrightarrow$",
)
img = solid.imshow("Principal Values of Cauchy Stress") |
Beta Was this translation helpful? Give feedback.
-
Hmm. I don't find |
Beta Was this translation helpful? Give feedback.
-
The def G(Cv):
# ....
ηK = ηinf + (η0 - ηinf + K[0] * (I1v ** β[0] - 3 ** β[0])) / (
1 + (K[1] * J2Neq + 1e-6) ** β[1]
)
# ... Due to the Runge-Kutta method, small increments are necessary. Newton-iterations would definitely improve this. Shouldn't be too hard to implement. Edit: See #741 for details on implementing Newton iterations. Code (with matadi)import numpy as np
import felupe as fem
from pypardiso import spsolve
import matadi as mat
from matadi.math import (
det,
trace,
inv,
gradient,
)
F = mat.Variable("F", 3, 3)
ζ = mat.Variable("ζ", 3, 3)
def ψ(x, dt, μ, α, m, a, η0, ηinf, β, K):
F, ζn = x[0], x[1]
# extract old state variables
Cvn = ζn
J = det(F)
C = (F.T @ F) * J ** (-2 / 3)
I1 = trace(C)
def invariants(Cv):
"Invariants of the elastic and viscous parts of the deformation."
Ce = C @ inv(Cv)
I1e = trace(Ce)
I1v = trace(Cv)
I2e = (I1e**2 - trace(Ce @ Ce)) / 2
return I1e, I1v, I2e
I1e, I1v, I2e = invariants(Cvn)
def G(Cv):
"Evolution equation."
I1e, I1v, I2e = invariants(Cv)
J2Neq = (I1e**2 / 3 - I2e) * sum(
[3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
) ** 2
u = [3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
ηK = ηinf + (η0 - ηinf + K[0] * (I1v ** β[0] - 3 ** β[0])) / (
1 + (K[1] * J2Neq + 1e-6) ** β[1]
)
return sum(u) / ηK * (C - I1e / 3 * Cv)
k1 = G(Cvn)
k2 = G(Cvn + k1 * dt / 2)
k3 = G(Cvn + (3 * k1 + k2) * dt / 16)
k4 = G(Cvn + k3 * dt / 2)
k5 = G(Cvn + 3 * (-k2 + 2 * k3 + 3 * k4) * dt / 16)
k6 = G(Cvn + (k1 + 4 * k2 + 6 * k3 - 12 * k4 + 8 * k5) * dt / 7)
Cv = Cvn + dt / 90 * (7 * k1 + 32 * k3 + 12 * k4 + 32 * k5 + 7 * k6)
I1e, I1v, I2e = invariants(Cv)
# strain energy functions of the equilibrium and non-equilibrium parts
p = [0, 1]
ψEq = [3 ** (1 - α[r]) / (2 * α[r]) * μ[r] * (I1 ** α[r] - 3 ** α[r]) for r in p]
ψNEq = [3 ** (1 - a[r]) / (2 * a[r]) * m[r] * (I1e ** a[r] - 3 ** a[r]) for r in p]
return gradient(sum(ψEq) + sum(ψNEq), F), Cv
umat = mat.MaterialTensor(
x=[F, ζ],
fun=ψ,
statevars=1,
kwargs=dict(
dt=0.1,
μ=[1.08, 0.0170],
α=[0.26, 7.68],
m=[1.57, 0.59],
a=[-10, 7.53],
η0=2.11,
ηinf=0.1,
β=[3.0, 1.929],
K=[442, 1289.49],
),
compress=True,
)
view = fem.ViewMaterialIncompressible(
umat,
ux=fem.math.linsteps([1, 1.5, 1], num=[200]),
bx=None,
ps=None,
statevars=np.eye(3).reshape(3, 3, 1, 1)
)
view.plot(
show_kwargs=False,
)
mesh = fem.Cube(n=2)
region = fem.RegionHexahedron(mesh)
field = fem.FieldContainer([fem.Field(region, dim=3)])
boundaries, loadcase = fem.dof.uniaxial(field, clamped=False)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)
solid.results.statevars.T[:] = np.eye(3)
move = fem.math.linsteps([0, 0.5, 0, 0.5], num=200)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)
job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(solver=spsolve, tol=1e-2)
fig, ax = job.plot(
xlabel="Displacement $u$ in mm $\longrightarrow$",
ylabel="Normal Force $F$ in N $\longrightarrow$",
)
solid.imshow("Principal Values of Cauchy Stress") |
Beta Was this translation helpful? Give feedback.
-
Hi @bhaveshshrimali, I have something for you 🎁 🎁 . I know we've already discussed this question in detail and we found a solution, but I managed to implement the Newton-Rhapson iterations - finally! Unfortunately not with matadi, but with tensortrax. It's possible to use quite large increments with this constitutive model formulation. The trick is to un-couple the (Unfortunately I still have no idea how to handle that with matadi.) Code (this requires FElupe >=8.4.0)from tensortrax import function, jacobian
from tensortrax.math import trace
from tensortrax.math.linalg import inv
from tensortrax.math.special import dev, from_triu_1d, triu_1d
import felupe as fem
@fem.constitution.isochoric_volumetric_split
def viscoelastic_two_potential(C, ζn, dt, μ, α, m, a, η0, ηinf, β, K):
# extract old state variables
Cvn = from_triu_1d(ζn, like=C)
# first invariant of the right Cauchy-Green deformation tensor
I1 = trace(C)
def invariants(Cv, C):
"Invariants of the elastic and viscous parts of the deformation."
Ce = C @ inv(Cv)
I1e = trace(Ce)
I1v = trace(Cv)
I2e = (I1e**2 - trace(Ce @ Ce)) / 2
return I1e, I1v, I2e
def evolution(Cv, Cvn, C):
"Evolution equation."
I1e, I1v, I2e = invariants(Cv, C)
J2Neq = (I1e**2 / 3 - I2e) * sum(
[3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
) ** 2
u = [3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
ηK = ηinf + (η0 - ηinf + K[0] * (I1v ** β[0] - 3 ** β[0])) / (
1 + (K[1] * J2Neq + 1e-10) ** β[1]
)
G = sum(u) / ηK * dev(C @ inv(Cv)) @ Cv
return G - (Cv - Cvn) / dt
# update the state variables
Cv = fem.newtonrhapson(
x0=Cvn,
fun=function(evolution, ntrax=C.ntrax),
jac=jacobian(evolution, ntrax=C.ntrax),
solve=fem.math.solve_2d,
args=(Cvn, C.x),
verbose=0,
tol=1e-1,
).x
I1e, I1v, I2e = invariants(Cv, C)
# strain energy functions of the equilibrium and non-equilibrium parts
p = [0, 1]
ψEq = [3 ** (1 - α[r]) / (2 * α[r]) * μ[r] * (I1 ** α[r] - 3 ** α[r]) for r in p]
ψNEq = [3 ** (1 - a[r]) / (2 * a[r]) * m[r] * (I1e ** a[r] - 3 ** a[r]) for r in p]
return sum(ψEq) + sum(ψNEq), triu_1d(Cv)
dt_vals = 4.0
umat = fem.Hyperelastic(
viscoelastic_two_potential,
nstatevars=6,
dt=dt_vals,
μ=[13.54, 1.08],
α=[1.0, -2.474],
m=[5.42, 20.78],
a=[-10, 1.948],
η0=7014.0,
ηinf=0.1,
β=[1.852, 0.26],
K=[3507.0, 1.0],
# parallel=True,
)
mesh = fem.Cube(n=2)
region = fem.RegionHexahedron(mesh)
field = fem.FieldContainer([fem.Field(region, dim=3)])
boundaries, loadcase = fem.dof.uniaxial(field, clamped=False)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)
solid.results.statevars[[0, 3, 5]] += 1.0
ldot = 0.05
lfinal = 3.0
tfinal = (lfinal - 1.0) / ldot
num_steps = tfinal / dt_vals
move = fem.math.linsteps([0.0, lfinal - 1.00, 0.358], num=int(num_steps) + 1)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)
job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate()
fig, ax = job.plot(
xlabel="Displacement $u$ in mm $\longrightarrow$",
ylabel="Normal Force $F$ in N $\longrightarrow$",
marker="x",
)
img = solid.imshow("Principal Values of Cauchy Stress") If you're interested, we could add this to the examples section in the docs!? It's running fast in that basic setting, so that should be fine for the read-the-docs CI/CD-pipeline. I'd be great if you could open a PR with a new example as |
Beta Was this translation helpful? Give feedback.
-
Hi @bhaveshshrimali! In addition to the above tensortrax-based Code, it is now possible to use JAX, even in parallel on multiple CPU cores! 🚀 This works also well for Code (JAX, code based on main-branch, API could change until the release of v9.1.0)r"""
Uniaxial loading/unloading of a viscoelastic material (VHB 4910)
----------------------------------------------------------------
This example shows how to implement a constitutive material model for
rubber viscoelastic materials using a strain energy density function coupled
with an ODE for an internal (state) variable [1]_. The ODE is discretized using a
backward-Euler scheme and the resulting nonlinear algebraic equations for the
internal variable are solved using Newton's method [2]_.
"""
import os
os.environ["XLA_FLAGS"] = "--xla_force_host_platform_device_count=8"
import felupe as fem
import jax
import jax.numpy as jnp
from jax.scipy.optimize import minimize
import felupe.constitution.jax as mat
@mat.isochoric_volumetric_split
def viscoelastic_two_potential(C, ζn, dt, μ, α, m, a, η0, ηinf, β, K):
# extract old state variables
I = jnp.eye(3)
Cvn = ζn.reshape(3, 3) + I
def invariants(Cv, C):
"Invariants of the elastic and viscous parts of the deformation."
Ce = C @ jnp.linalg.inv(Cv)
I1e = jnp.trace(Ce)
I1v = jnp.trace(Cv)
I2e = (I1e**2 - jnp.trace(Ce @ Ce)) / 2
return I1e, I1v, I2e
def evolution(Cv, Cvn, C, ξ=1e-10):
"Evolution equation for the internal (state) variable Cv."
Cv = Cv.reshape(3, 3)
I1e, I1v, I2e = invariants(Cv, C)
J2Neq = (I1e**2 / 3 - I2e) * sum(
[3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
) ** 2
u = [3 ** (1 - a[r]) * m[r] * I1e ** (a[r] - 1) for r in [0, 1]]
ηK = ηinf + (η0 - ηinf + K[0] * (I1v ** β[0] - 3 ** β[0])) / (
1 + (K[1] * J2Neq + ξ) ** β[1]
)
dev = lambda C: C - jnp.trace(C) / 3 * jnp.eye(3)
G = sum(u) / ηK * dev(C @ jnp.linalg.inv(Cv)) @ Cv
residuals = G - (Cv - Cvn) / dt
return residuals.ravel() @ residuals.ravel()
# update the state variables
Cv = minimize(
evolution,
x0=Cvn.ravel(),
args=(Cvn, jax.lax.stop_gradient(C)),
method="BFGS",
).x.reshape(3, 3)
# invariants of the (elastic and viscous) right Cauchy-Green deformation tensor
I1 = jnp.trace(C)
I1e, I1v, I2e = invariants(Cv, C)
# strain energy functions of the equilibrium and non-equilibrium parts
p = [0, 1]
ψEq = [3 ** (1 - α[r]) / (2 * α[r]) * μ[r] * (I1 ** α[r] - 3 ** α[r]) for r in p]
ψNEq = [3 ** (1 - a[r]) / (2 * a[r]) * m[r] * (I1e ** a[r] - 3 ** a[r]) for r in p]
return sum(ψEq) + sum(ψNEq), (Cv - I).ravel()
dt_vals = 4.0
umat = mat.Hyperelastic(
viscoelastic_two_potential,
nstatevars=9,
dt=dt_vals,
μ=[13.54, 1.08],
α=[1.0, -2.474],
m=[5.42, 20.78],
a=[-10, 1.948],
η0=7014.0,
ηinf=0.1,
β=[1.852, 0.26],
K=[3507.0, 1.0],
parallel=True,
)
mesh = fem.Cube(n=6)
region = fem.RegionHexahedron(mesh)
field = fem.FieldContainer([fem.Field(region, dim=3)])
boundaries, loadcase = fem.dof.uniaxial(field, clamped=False)
solid = fem.SolidBodyNearlyIncompressible(umat, field, bulk=5000)
ldot = 0.05
lfinal = 3.0
tfinal = (lfinal - 1.0) / ldot
num_steps = tfinal / dt_vals
move = fem.math.linsteps([0.0, lfinal - 1.00, 0.358], num=int(num_steps) + 1)
step = fem.Step(items=[solid], ramp={boundaries["move"]: move}, boundaries=boundaries)
job = fem.CharacteristicCurve(steps=[step], boundary=boundaries["move"])
job.evaluate(tol=1e-2)
fig, ax = job.plot(
xlabel=r"Displacement $u$ in mm $\longrightarrow$",
ylabel=r"Normal Force $F$ in N $\longrightarrow$",
marker="x",
)
ax.set_title("Uniaxial loading/unloading\n of a viscoelastic material (VHB 4910)")
solid.plot("Principal Values of Cauchy Stress").show() |
Beta Was this translation helpful? Give feedback.
The$\beta$ -parameter seems to make the model unstable. Or I made a typo somewhere. With $\beta_r=1$ it works. If I add a small number like $\eta_K$ , the code is much more stable, even with the original parameters.
1e-6
inDue to the Runge-Kutta method, small increments are necessary. Newton-iterations would definitely improve this. Shouldn't be too hard to implement.
Edit: See #741 for details on implementing Newton iterations.
Code (with matadi)