Keywords: hyperelasticity, finite strain, large deformations, Newton's method, conjugate gradient, automatic differentiation
Figure 1: Cube loaded in torsion modeled with a hyperelastic material model and finite strain.
In this example we will solve a problem in a finite strain setting using an hyperelastic material model. In order to compute the stress we will use automatic differentiation, to solve the non-linear system we use Newton's method, and for solving the Newton increment we use conjugate gradients.
The weak form is expressed in terms of the first Piola-Kirchoff stress P as follows: Find u∈U such that
∫Ω[∇Xδu]:P(u) dΩ=∫Ωδu⋅b dΩ+∫ΓNδu⋅t dΓ∀δu∈U0,where u is the unknown displacement field, b is the body force acting on the reference domain, t is the traction acting on the Neumann part of the reference domain's boundary, and where U and U0 are suitable trial and test sets. Ω denotes the reference (sometimes also called initial or material) domain. Gradients are defined with respect to the reference domain, here denoted with an X. Formally this is expressed as (∇X∙)ij:=∂(∙)i∂Xj. Note that for large deformation problems it is also possible that gradients and integrals are defined on the deformed (sometimes also called current or spatial) domain, depending on the specific formulation.
The specific problem we will solve in this example is the cube from Figure 1: On one side we apply a rotation using Dirichlet boundary conditions, on the opposite side we fix the displacement with a homogeneous Dirichlet boundary condition, and on the remaining four sides we apply a traction in the normal direction of the surface. In addition, a body force is applied in one direction.
In addition to Ferrite.jl and Tensors.jl, this examples uses TimerOutputs.jl for timing the program and print a summary at the end, ProgressMeter.jl for showing a simple progress bar, and IterativeSolvers.jl for solving the linear system using conjugate gradients.
using Ferrite, Tensors, TimerOutputs, ProgressMeter, IterativeSolvers
The stress can be derived from an energy potential, defined in terms of the right Cauchy-Green tensor C=FT⋅F, where F=I+∇Xu is the deformation gradient. We shall use the compressible neo-Hookean model from Wikipedia with the potential
Ψ(C)=μ2(I1−3)⏟W(C)−μln(J)+λ2(J−1)2⏟U(J),where I1=tr(C) is the first invariant, J=√det(C) and μ and λ material parameters.
Extra details on compressible neo-Hookean formulations
The Neo-Hooke model is only a well defined terminology in the incompressible case. Thus, only W(C) specifies the neo-Hookean behavior, the volume penalty U(J) can vary in different formulations. In order to obtain a well-posed problem, it is crucial to choose a convex formulation of U(J). Other examples for U(J) can be found, e.g. in Hol:2000:nsm; Eq. (6.138) β−2(βlnJ+J−β−1)
where SimMie:1992:act; Eq. (2.37) published a non-generalized version with β=−2. This shows the possible variety of U(J) while all of them refer to compressible neo-Hookean models. Sometimes the modified first invariant ¯I1=I1I1/33 is used in W(C) instead of I1.
From the potential we obtain the second Piola-Kirchoff stress S as
S=2∂Ψ∂C,and the tangent of S as
∂S∂C=2∂2Ψ∂C2.Finally, for the finite element problem we need P and ∂P∂F, which can be obtained by using the following relations:
P=F⋅S,∂P∂F=Iˉ⊗S+2F⋅∂S∂C:FTˉ⊗I.Derivation of ∂P/∂F
Tip: See knutam.github.io/tensors for an explanation of the index notation used in this derivation. Using the product rule, the chain rule, and the relations P=F⋅S and C=FT⋅F, we obtain the following: ∂Pij∂Fkl=∂(FimSmj)∂Fkl=∂Fim∂FklSmj+Fim∂Smj∂Fkl=δikδmlSmj+Fim∂Smj∂Cno∂Cno∂Fkl=δikSlj+Fim∂Smj∂Cno∂(FTnpFpo)∂Fkl=δikSTjl+Fim∂Smj∂Cno(∂FTnp∂FklFpo+FTnp∂Fpo∂Fkl)=δikSjl+Fim∂Smj∂Cno(δnlδpkFpo+FTnpδpkδol)=δikSlj+Fim∂Smj∂Cno(FTokδnl+FTnkδol)=δikSjl+2Fim∂Smj∂CnoFTnkδol∂P∂F=Iˉ⊗S+2F⋅∂S∂C:FTˉ⊗I,
where we used the fact that S is symmetric (Slj=Sjl) and that ∂S∂C is minor symmetric (∂Smj∂Cno=∂Smj∂Con).
We can implement the material model as follows, where we utilize automatic differentiation for the stress and the tangent, and thus only define the potential:
struct NeoHooke
μ::Float64
λ::Float64
end
function Ψ(C, mp::NeoHooke)
μ = mp.μ
λ = mp.λ
Ic = tr(C)
J = sqrt(det(C))
return μ / 2 * (Ic - 3 - 2 * log(J)) + λ / 2 * (J - 1)^2
end
function constitutive_driver(C, mp::NeoHooke)
# Compute all derivatives in one function call
∂²Ψ∂C², ∂Ψ∂C = Tensors.hessian(y -> Ψ(y, mp), C, :all)
S = 2.0 * ∂Ψ∂C
∂S∂C = 2.0 * ∂²Ψ∂C²
return S, ∂S∂C
end;
As mentioned above, to deal with the non-linear weak form we first linearize the problem such that we can apply Newton's method, and then apply the FEM to discretize the problem. Skipping a detailed derivation, Newton's method can be expressed as: Given some initial guess for the degrees of freedom u_0, find a sequence u_k by iterating
u_k+1=u_k−Δu_kuntil some termination condition has been met. Therein we determine Δu_k from the linearized problem
K__(u_k)Δu_k=g_(u_k)where the global residual, g_, and the Jacobi matrix, K__=∂g_∂u_, are evaluated at the current guess u_k. The entries of g_ are given by
(g_)i=∫Ω[∇Xδui]:PdΩ−∫Ωδui⋅bdΩ−∫ΓNδui⋅t dΓ,and the entries of K__ are given by
(K__)ij=∫Ω[∇Xδui]:∂P∂F:[∇Xδuj]dΩ.A detailed derivation can be found in every continuum mechanics book, which has a chapter about finite elasticity theory. We used "Nonlinear solid mechanics: a continuum approach for engineering science." by Hol:2000:nsm; Chapter 8 as a reference.
The element routine for assembling the residual and tangent stiffness is implemented as usual, with loops over quadrature points and shape functions:
function assemble_element!(ke, ge, cell, cv, fv, mp, ue, ΓN)
# Reinitialize cell values, and reset output arrays
reinit!(cv, cell)
fill!(ke, 0.0)
fill!(ge, 0.0)
b = Vec{3}((0.0, -0.5, 0.0)) # Body force
tn = 0.1 # Traction (to be scaled with surface normal)
ndofs = getnbasefunctions(cv)
for qp in 1:getnquadpoints(cv)
dΩ = getdetJdV(cv, qp)
# Compute deformation gradient F and right Cauchy-Green tensor C
∇u = function_gradient(cv, qp, ue)
F = one(∇u) + ∇u
C = tdot(F) # F' ⋅ F
# Compute stress and tangent
S, ∂S∂C = constitutive_driver(C, mp)
P = F ⋅ S
I = one(S)
∂P∂F = otimesu(I, S) + 2 * F ⋅ ∂S∂C ⊡ otimesu(F', I)
# Loop over test functions
for i in 1:ndofs
# Test function and gradient
δui = shape_value(cv, qp, i)
∇δui = shape_gradient(cv, qp, i)
# Add contribution to the residual from this test function
ge[i] += (∇δui ⊡ P - δui ⋅ b) * dΩ
∇δui∂P∂F = ∇δui ⊡ ∂P∂F # Hoisted computation
for j in 1:ndofs
∇δuj = shape_gradient(cv, qp, j)
# Add contribution to the tangent
ke[i, j] += (∇δui∂P∂F ⊡ ∇δuj) * dΩ
end
end
end
# Surface integral for the traction
for facet in 1:nfacets(cell)
if (cellid(cell), facet) in ΓN
reinit!(fv, cell, facet)
for q_point in 1:getnquadpoints(fv)
t = tn * getnormal(fv, q_point)
dΓ = getdetJdV(fv, q_point)
for i in 1:ndofs
δui = shape_value(fv, q_point, i)
ge[i] -= (δui ⋅ t) * dΓ
end
end
end
end
return
end;
Assembling global residual and tangent is also done in the usual way, just looping over the elements, call the element routine and assemble in the the global matrix K and residual g.
function assemble_global!(K, g, dh, cv, fv, mp, u, ΓN)
n = ndofs_per_cell(dh)
ke = zeros(n, n)
ge = zeros(n)
# start_assemble resets K and g
assembler = start_assemble(K, g)
# Loop over all cells in the grid
@timeit "assemble" for cell in CellIterator(dh)
global_dofs = celldofs(cell)
ue = u[global_dofs] # element dofs
@timeit "element assemble" assemble_element!(ke, ge, cell, cv, fv, mp, ue, ΓN)
assemble!(assembler, global_dofs, ke, ge)
end
return
end;
Finally, we define a main function which sets up everything and then performs Newton iterations until convergence.
function solve()
reset_timer!()
# Generate a grid
N = 10
L = 1.0
left = zero(Vec{3})
right = L * ones(Vec{3})
grid = generate_grid(Tetrahedron, (N, N, N), left, right)
# Material parameters
E = 10.0
ν = 0.3
μ = E / (2(1 + ν))
λ = (E * ν) / ((1 + ν) * (1 - 2ν))
mp = NeoHooke(μ, λ)
# Finite element base
ip = Lagrange{RefTetrahedron, 1}()^3
qr = QuadratureRule{RefTetrahedron}(1)
qr_facet = FacetQuadratureRule{RefTetrahedron}(1)
cv = CellValues(qr, ip)
fv = FacetValues(qr_facet, ip)
# DofHandler
dh = DofHandler(grid)
add!(dh, :u, ip) # Add a displacement field
close!(dh)
function rotation(X, t)
θ = pi / 3 # 60°
x, y, z = X
return t * Vec{3}(
(
0.0,
L / 2 - y + (y - L / 2) * cos(θ) - (z - L / 2) * sin(θ),
L / 2 - z + (y - L / 2) * sin(θ) + (z - L / 2) * cos(θ),
)
)
end
dbcs = ConstraintHandler(dh)
# Add a homogeneous boundary condition on the "clamped" edge
dbc = Dirichlet(:u, getfacetset(grid, "right"), (x, t) -> [0.0, 0.0, 0.0], [1, 2, 3])
add!(dbcs, dbc)
dbc = Dirichlet(:u, getfacetset(grid, "left"), (x, t) -> rotation(x, t), [1, 2, 3])
add!(dbcs, dbc)
close!(dbcs)
t = 0.5
Ferrite.update!(dbcs, t)
# Neumann part of the boundary
ΓN = union(
getfacetset(grid, "top"),
getfacetset(grid, "bottom"),
getfacetset(grid, "front"),
getfacetset(grid, "back"),
)
# Pre-allocation of vectors for the solution and Newton increments
_ndofs = ndofs(dh)
un = zeros(_ndofs) # previous solution vector
u = zeros(_ndofs)
Δu = zeros(_ndofs)
ΔΔu = zeros(_ndofs)
apply!(un, dbcs)
# Create sparse matrix and residual vector
K = allocate_matrix(dh)
g = zeros(_ndofs)
# Perform Newton iterations
newton_itr = -1
NEWTON_TOL = 1.0e-8
NEWTON_MAXITER = 30
prog = ProgressMeter.ProgressThresh(NEWTON_TOL; desc = "Solving:")
while true
newton_itr += 1
# Construct the current guess
u .= un .+ Δu
# Compute residual and tangent for current guess
assemble_global!(K, g, dh, cv, fv, mp, u, ΓN)
# Apply boundary conditions
apply_zero!(K, g, dbcs)
# Compute the residual norm and compare with tolerance
normg = norm(g)
ProgressMeter.update!(prog, normg; showvalues = [(:iter, newton_itr)])
if normg < NEWTON_TOL
break
elseif newton_itr > NEWTON_MAXITER
error("Reached maximum Newton iterations, aborting")
end
# Compute increment using conjugate gradients
@timeit "linear solve" IterativeSolvers.cg!(ΔΔu, K, g; maxiter = 1000)
apply_zero!(ΔΔu, dbcs)
Δu .-= ΔΔu
end
# Save the solution
@timeit "export" begin
VTKGridFile("hyperelasticity", dh) do vtk
write_solution(vtk, dh, u)
end
end
print_timer(title = "Analysis with $(getncells(grid)) elements", linechars = :ascii)
return u
end
solve (generic function with 1 method)
Run the simulation
u = solve();
Solving: (thresh = 1e-08, value = 0.000259699) Solving: Time: 0:00:00 (6 iterations) iter: 5 ------------------------------------------------------------------------------- Analysis with 6000 elements Time Allocations ----------------------- ------------------------ Tot / % measured: 2.12s / 52.2% 147MiB / 44.5% Section ncalls time %tot avg alloc %tot avg ------------------------------------------------------------------------------- export 1 923ms 83.3% 923ms 59.3MiB 90.9% 59.3MiB assemble 6 110ms 9.9% 18.3ms 5.50MiB 8.4% 938KiB element assemble 36.0k 64.5ms 5.8% 1.79μs 0.00B 0.0% 0.00B linear solve 5 74.6ms 6.7% 14.9ms 473KiB 0.7% 94.6KiB -------------------------------------------------------------------------------
This notebook was generated using Literate.jl.