Hyperelasticity

Hyperelasticity#

Author: Jørgen S. Dokken and Garth N. Wells

This section shows how to solve the hyperelasticity problem for deformation of a beam.

We will also show how to create a constant boundary condition for a vector function space.

We start by importing DOLFINx and some additional dependencies. Then, we create a slender cantilever consisting of hexahedral elements and create the function space V for our unknown.

from dolfinx import log, default_scalar_type
from dolfinx.fem.petsc import NonlinearProblem
from dolfinx.nls.petsc import NewtonSolver
import pyvista
import numpy as np
import ufl

from mpi4py import MPI
from dolfinx import fem, mesh, plot
L = 20.0
domain = mesh.create_box(MPI.COMM_WORLD, [[0.0, 0.0, 0.0], [L, 1, 1]], [20, 5, 5], mesh.CellType.hexahedron)
V = fem.functionspace(domain, ("Lagrange", 2, (domain.geometry.dim, )))

We create two python functions for determining the facets to apply boundary conditions to

def left(x):
    return np.isclose(x[0], 0)


def right(x):
    return np.isclose(x[0], L)


fdim = domain.topology.dim - 1
left_facets = mesh.locate_entities_boundary(domain, fdim, left)
right_facets = mesh.locate_entities_boundary(domain, fdim, right)

Next, we create a marker based on these two functions

# Concatenate and sort the arrays based on facet indices. Left facets marked with 1, right facets with two
marked_facets = np.hstack([left_facets, right_facets])
marked_values = np.hstack([np.full_like(left_facets, 1), np.full_like(right_facets, 2)])
sorted_facets = np.argsort(marked_facets)
facet_tag = mesh.meshtags(domain, fdim, marked_facets[sorted_facets], marked_values[sorted_facets])

We then create a function for supplying the boundary condition on the left side, which is fixed.

u_bc = np.array((0,) * domain.geometry.dim, dtype=default_scalar_type)

To apply the boundary condition, we identity the dofs located on the facets marked by the MeshTag.

left_dofs = fem.locate_dofs_topological(V, facet_tag.dim, facet_tag.find(1))
bcs = [fem.dirichletbc(u_bc, left_dofs, V)]

Next, we define the body force on the reference configuration (B), and nominal (first Piola-Kirchhoff) traction (T).

B = fem.Constant(domain, default_scalar_type((0, 0, 0)))
T = fem.Constant(domain, default_scalar_type((0, 0, 0)))

Define the test and solution functions on the space \(V\)

v = ufl.TestFunction(V)
u = fem.Function(V)

Define kinematic quantities used in the problem

# Spatial dimension
d = len(u)

# Identity tensor
I = ufl.variable(ufl.Identity(d))

# Deformation gradient
F = ufl.variable(I + ufl.grad(u))

# Right Cauchy-Green tensor
C = ufl.variable(F.T * F)

# Invariants of deformation tensors
Ic = ufl.variable(ufl.tr(C))
J = ufl.variable(ufl.det(F))

Define the elasticity model via a stored strain energy density function \(\psi\), and create the expression for the first Piola-Kirchhoff stress:

# Elasticity parameters
E = default_scalar_type(1.0e4)
nu = default_scalar_type(0.3)
mu = fem.Constant(domain, E / (2 * (1 + nu)))
lmbda = fem.Constant(domain, E * nu / ((1 + nu) * (1 - 2 * nu)))
# Stored strain energy density (compressible neo-Hookean model)
psi = (mu / 2) * (Ic - 3) - mu * ufl.ln(J) + (lmbda / 2) * (ufl.ln(J))**2
# Stress
# Hyper-elasticity
P = ufl.diff(psi, F)

Comparison to linear elasticity

To illustrate the difference between linear and hyperelasticity, the following lines can be uncommented to solve the linear elasticity problem.

# P = 2.0 * mu * ufl.sym(ufl.grad(u)) + lmbda * ufl.tr(ufl.sym(ufl.grad(u))) * I

Define the variational form with traction integral over all facets with value 2. We set the quadrature degree for the integrals to 4.

metadata = {"quadrature_degree": 4}
ds = ufl.Measure('ds', domain=domain, subdomain_data=facet_tag, metadata=metadata)
dx = ufl.Measure("dx", domain=domain, metadata=metadata)
# Define form F (we want to find u such that F(u) = 0)
F = ufl.inner(ufl.grad(v), P) * dx - ufl.inner(v, B) * dx - ufl.inner(v, T) * ds(2)

As the varitional form is non-linear and written on residual form, we use the non-linear problem class from DOLFINx to set up required structures to use a Newton solver.

problem = NonlinearProblem(F, u, bcs)

and then create and customize the Newton solver

solver = NewtonSolver(domain.comm, problem)

# Set Newton solver options
solver.atol = 1e-8
solver.rtol = 1e-8
solver.convergence_criterion = "incremental"

We create a function to plot the solution at each time step.

pyvista.start_xvfb()
plotter = pyvista.Plotter()
plotter.open_gif("deformation.gif", fps=3)

topology, cells, geometry = plot.vtk_mesh(u.function_space)
function_grid = pyvista.UnstructuredGrid(topology, cells, geometry)

values = np.zeros((geometry.shape[0], 3))
values[:, :len(u)] = u.x.array.reshape(geometry.shape[0], len(u))
function_grid["u"] = values
function_grid.set_active_vectors("u")

# Warp mesh by deformation
warped = function_grid.warp_by_vector("u", factor=1)
warped.set_active_vectors("u")

# Add mesh to plotter and visualize
actor = plotter.add_mesh(warped, show_edges=True, lighting=False, clim=[0, 10])

# Compute magnitude of displacement to visualize in GIF
Vs = fem.functionspace(domain, ("Lagrange", 2))
magnitude = fem.Function(Vs)
us = fem.Expression(ufl.sqrt(sum([u[i]**2 for i in range(len(u))])), Vs.element.interpolation_points())
magnitude.interpolate(us)
warped["mag"] = magnitude.x.array

Finally, we solve the problem over several time steps, updating the z-component of the traction

log.set_log_level(log.LogLevel.INFO)
tval0 = -1.5
for n in range(1, 10):
    T.value[2] = n * tval0
    num_its, converged = solver.solve(u)
    assert (converged)
    u.x.scatter_forward()
    print(f"Time step {n}, Number of iterations {num_its}, Load {T.value}")
    function_grid["u"][:, :len(u)] = u.x.array.reshape(geometry.shape[0], len(u))
    magnitude.interpolate(us)
    warped.set_active_scalars("mag")
    warped_n = function_grid.warp_by_vector(factor=1)
    warped.points[:, :] = warped_n.points
    warped.point_data["mag"][:] = magnitude.x.array
    plotter.update_scalar_bar_range([0, 10])
    plotter.write_frame()
plotter.close()
[2025-01-13 11:25:45.669] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:25:48.152] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:25:50.176] [info] Newton iteration 2: r (abs) = 22.24548907411152 (tol = 1e-08), r (rel) = 0.1342779429427 (tol = 1e-08)
[2025-01-13 11:25:50.445] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:25:52.468] [info] Newton iteration 3: r (abs) = 2.432613387004531 (tol = 1e-08), r (rel) = 0.014683710503896342 (tol = 1e-08)
[2025-01-13 11:25:52.735] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:25:54.758] [info] Newton iteration 4: r (abs) = 4.431579458563597 (tol = 1e-08), r (rel) = 0.0267498445055792 (tol = 1e-08)
[2025-01-13 11:25:55.034] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:25:57.048] [info] Newton iteration 5: r (abs) = 0.1441891503961589 (tol = 1e-08), r (rel) = 0.000870352746363483 (tol = 1e-08)
[2025-01-13 11:25:57.316] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:25:59.374] [info] Newton iteration 6: r (abs) = 0.02142395087661431 (tol = 1e-08), r (rel) = 0.00012931898434928528 (tol = 1e-08)
[2025-01-13 11:25:59.645] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:01.740] [info] Newton iteration 7: r (abs) = 4.800650765296232e-06 (tol = 1e-08), r (rel) = 2.8977628111600535e-08 (tol = 1e-08)
[2025-01-13 11:26:02.018] [info] PETSc Krylov solver starting to solve system.
Time step 1, Number of iterations 8, Load [ 0.   0.  -1.5]
[2025-01-13 11:26:04.022] [info] Newton iteration 8: r (abs) = 2.6576687328769084e-11 (tol = 1e-08), r (rel) = 1.6042186768064864e-13 (tol = 1e-08)
[2025-01-13 11:26:04.022] [info] Newton solver finished in 8 iterations and 8 linear solver iterations.
[2025-01-13 11:26:04.452] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:06.733] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:08.756] [info] Newton iteration 2: r (abs) = 17.325417844303768 (tol = 1e-08), r (rel) = 0.11784222265438753 (tol = 1e-08)
[2025-01-13 11:26:09.036] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:11.065] [info] Newton iteration 3: r (abs) = 5.148824049644568 (tol = 1e-08), r (rel) = 0.03502073517181963 (tol = 1e-08)
[2025-01-13 11:26:11.340] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:13.370] [info] Newton iteration 4: r (abs) = 7.240032557943796 (tol = 1e-08), r (rel) = 0.049244499404597956 (tol = 1e-08)
[2025-01-13 11:26:13.639] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:15.641] [info] Newton iteration 5: r (abs) = 0.7778888038797866 (tol = 1e-08), r (rel) = 0.005290963049257456 (tol = 1e-08)
[2025-01-13 11:26:15.916] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:17.921] [info] Newton iteration 6: r (abs) = 1.2552533176963954 (tol = 1e-08), r (rel) = 0.008537851281911266 (tol = 1e-08)
[2025-01-13 11:26:18.193] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:20.206] [info] Newton iteration 7: r (abs) = 0.008495123844168218 (tol = 1e-08), r (rel) = 5.77812485977187e-05 (tol = 1e-08)
[2025-01-13 11:26:20.481] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:22.504] [info] Newton iteration 8: r (abs) = 0.00019210706900551987 (tol = 1e-08), r (rel) = 1.3066538540467734e-06 (tol = 1e-08)
Time step 2, Number of iterations 9, Load [ 0.  0. -3.]
[2025-01-13 11:26:22.772] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:24.789] [info] Newton iteration 9: r (abs) = 1.709021236534951e-10 (tol = 1e-08), r (rel) = 1.1624242652424269e-12 (tol = 1e-08)
[2025-01-13 11:26:24.789] [info] Newton solver finished in 9 iterations and 9 linear solver iterations.
[2025-01-13 11:26:25.121] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:27.409] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:29.427] [info] Newton iteration 2: r (abs) = 10.001117236873005 (tol = 1e-08), r (rel) = 0.08874707902942588 (tol = 1e-08)
[2025-01-13 11:26:29.695] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:31.742] [info] Newton iteration 3: r (abs) = 5.330258435993069 (tol = 1e-08), r (rel) = 0.0472992022253551 (tol = 1e-08)
[2025-01-13 11:26:32.018] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:34.068] [info] Newton iteration 4: r (abs) = 11.990116716277255 (tol = 1e-08), r (rel) = 0.10639689652555262 (tol = 1e-08)
[2025-01-13 11:26:34.344] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:36.362] [info] Newton iteration 5: r (abs) = 2.2970192409224013 (tol = 1e-08), r (rel) = 0.020383097535809913 (tol = 1e-08)
[2025-01-13 11:26:36.634] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:38.639] [info] Newton iteration 6: r (abs) = 3.9023388029843677 (tol = 1e-08), r (rel) = 0.03462824821922908 (tol = 1e-08)
[2025-01-13 11:26:38.909] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:40.928] [info] Newton iteration 7: r (abs) = 0.23653543084877515 (tol = 1e-08), r (rel) = 0.0020989483552298553 (tol = 1e-08)
[2025-01-13 11:26:41.196] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:43.210] [info] Newton iteration 8: r (abs) = 0.04271420898976795 (tol = 1e-08), r (rel) = 0.000379033781037806 (tol = 1e-08)
Time step 3, Number of iterations 10, Load [ 0.   0.  -4.5][2025-01-13 11:26:43.478] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:45.477] [info] Newton iteration 9: r (abs) = 2.877976847117939e-05 (tol = 1e-08), r (rel) = 2.553835063090331e-07 (tol = 1e-08)
[2025-01-13 11:26:45.745] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:47.748] [info] Newton iteration 10: r (abs) = 6.091459406196152e-10 (tol = 1e-08), r (rel) = 5.405388383340817e-12 (tol = 1e-08)
[2025-01-13 11:26:47.748] [info] Newton solver finished in 10 iterations and 10 linear solver iterations.
[2025-01-13 11:26:48.071] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:50.356] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:52.367] [info] Newton iteration 2: r (abs) = 5.506929937904466 (tol = 1e-08), r (rel) = 0.06539183684430351 (tol = 1e-08)
[2025-01-13 11:26:52.643] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:54.687] [info] Newton iteration 3: r (abs) = 26.24891969494894 (tol = 1e-08), r (rel) = 0.3116918307271015 (tol = 1e-08)
[2025-01-13 11:26:54.955] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:56.995] [info] Newton iteration 4: r (abs) = 2.309270517276505 (tol = 1e-08), r (rel) = 0.027421347755981797 (tol = 1e-08)
[2025-01-13 11:26:57.264] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:26:59.283] [info] Newton iteration 5: r (abs) = 14.056244587477249 (tol = 1e-08), r (rel) = 0.16691035896085923 (tol = 1e-08)
[2025-01-13 11:26:59.561] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:01.584] [info] Newton iteration 6: r (abs) = 0.22277413250167966 (tol = 1e-08), r (rel) = 0.0026453232363483556 (tol = 1e-08)
[2025-01-13 11:27:01.860] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:03.904] [info] Newton iteration 7: r (abs) = 0.2866705220675385 (tol = 1e-08), r (rel) = 0.0034040585623003423 (tol = 1e-08)
[2025-01-13 11:27:04.183] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:06.210] [info] Newton iteration 8: r (abs) = 0.0003218693539167365 (tol = 1e-08), r (rel) = 3.8220257954677744e-06 (tol = 1e-08)
Time step 4, Number of iterations 9, Load [ 0.  0. -6.]
[2025-01-13 11:27:06.483] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:08.502] [info] Newton iteration 9: r (abs) = 2.63796397740629e-07 (tol = 1e-08), r (rel) = 3.132440614948943e-09 (tol = 1e-08)
[2025-01-13 11:27:08.502] [info] Newton solver finished in 9 iterations and 9 linear solver iterations.
[2025-01-13 11:27:08.836] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:11.118] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:13.130] [info] Newton iteration 2: r (abs) = 3.194619708714672 (tol = 1e-08), r (rel) = 0.04964789515160881 (tol = 1e-08)
[2025-01-13 11:27:13.401] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:15.415] [info] Newton iteration 3: r (abs) = 7.714285842397931 (tol = 1e-08), r (rel) = 0.1198884654809222 (tol = 1e-08)
[2025-01-13 11:27:15.690] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:17.773] [info] Newton iteration 4: r (abs) = 0.850873121067315 (tol = 1e-08), r (rel) = 0.01322350181050775 (tol = 1e-08)
[2025-01-13 11:27:18.045] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:20.107] [info] Newton iteration 5: r (abs) = 0.3714343540930969 (tol = 1e-08), r (rel) = 0.00577249737031741 (tol = 1e-08)
[2025-01-13 11:27:20.378] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:22.427] [info] Newton iteration 6: r (abs) = 0.002150656449271473 (tol = 1e-08), r (rel) = 3.342356074786805e-05 (tol = 1e-08)
[2025-01-13 11:27:22.700] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:24.724] [info] Newton iteration 7: r (abs) = 2.546067958388237e-06 (tol = 1e-08), r (rel) = 3.9568689413047105e-08 (tol = 1e-08)
Time step 5, Number of iterations 8, Load [ 0.   0.  -7.5]
[2025-01-13 11:27:24.997] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:27.005] [info] Newton iteration 8: r (abs) = 1.7369915835599454e-13 (tol = 1e-08), r (rel) = 2.6994754895100873e-15 (tol = 1e-08)
[2025-01-13 11:27:27.005] [info] Newton solver finished in 8 iterations and 8 linear solver iterations.
[2025-01-13 11:27:27.333] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:29.608] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:31.624] [info] Newton iteration 2: r (abs) = 2.0064883010780563 (tol = 1e-08), r (rel) = 0.039562203134381475 (tol = 1e-08)
[2025-01-13 11:27:31.895] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:33.927] [info] Newton iteration 3: r (abs) = 4.6097682272934195 (tol = 1e-08), r (rel) = 0.09089142803006316 (tol = 1e-08)
[2025-01-13 11:27:34.195] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:36.209] [info] Newton iteration 4: r (abs) = 0.18537214353432482 (tol = 1e-08), r (rel) = 0.003655007803444642 (tol = 1e-08)
[2025-01-13 11:27:36.484] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:38.516] [info] Newton iteration 5: r (abs) = 0.024688001493842943 (tol = 1e-08), r (rel) = 0.00048677668818530124 (tol = 1e-08)
[2025-01-13 11:27:38.790] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:40.813] [info] Newton iteration 6: r (abs) = 5.692545102089998e-06 (tol = 1e-08), r (rel) = 1.1224068715452317e-07 (tol = 1e-08)
[2025-01-13 11:27:41.088] [info] PETSc Krylov solver starting to solve system.
Time step 6, Number of iterations 7, Load [ 0.  0. -9.]
[2025-01-13 11:27:43.091] [info] Newton iteration 7: r (abs) = 2.6200570099858678e-11 (tol = 1e-08), r (rel) = 5.166002094157671e-13 (tol = 1e-08)
[2025-01-13 11:27:43.091] [info] Newton solver finished in 7 iterations and 7 linear solver iterations.
[2025-01-13 11:27:43.419] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:45.695] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:47.701] [info] Newton iteration 2: r (abs) = 1.3850621340169573 (tol = 1e-08), r (rel) = 0.0336621757045322 (tol = 1e-08)
[2025-01-13 11:27:47.970] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:49.979] [info] Newton iteration 3: r (abs) = 3.0373936007636577 (tol = 1e-08), r (rel) = 0.0738199930252921 (tol = 1e-08)
[2025-01-13 11:27:50.247] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:52.259] [info] Newton iteration 4: r (abs) = 0.04123862114946107 (tol = 1e-08), r (rel) = 0.0010022523010717154 (tol = 1e-08)
[2025-01-13 11:27:52.526] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:54.534] [info] Newton iteration 5: r (abs) = 0.002050566445156464 (tol = 1e-08), r (rel) = 4.9836412587848416e-05 (tol = 1e-08)
[2025-01-13 11:27:54.802] [info] PETSc Krylov solver starting to solve system.
Time step 7, Number of iterations 6, Load [  0.    0.  -10.5]
[2025-01-13 11:27:56.795] [info] Newton iteration 6: r (abs) = 1.7886722652583207e-08 (tol = 1e-08), r (rel) = 4.3471407233065016e-10 (tol = 1e-08)
[2025-01-13 11:27:56.795] [info] Newton solver finished in 6 iterations and 6 linear solver iterations.
[2025-01-13 11:27:57.120] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:27:59.381] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:01.393] [info] Newton iteration 2: r (abs) = 1.0633648536770233 (tol = 1e-08), r (rel) = 0.031085002055400286 (tol = 1e-08)
[2025-01-13 11:28:01.669] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:03.666] [info] Newton iteration 3: r (abs) = 2.047703133550441 (tol = 1e-08), r (rel) = 0.059859845748295164 (tol = 1e-08)
[2025-01-13 11:28:03.934] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:05.921] [info] Newton iteration 4: r (abs) = 0.008977192756641317 (tol = 1e-08), r (rel) = 0.0002624273825930661 (tol = 1e-08)
[2025-01-13 11:28:06.193] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:08.195] [info] Newton iteration 5: r (abs) = 0.0001674217868497817 (tol = 1e-08), r (rel) = 4.894187136567673e-06 (tol = 1e-08)
[2025-01-13 11:28:08.466] [info] PETSc Krylov solver starting to solve system.
Time step 8, Number of iterations 6, Load [  0.   0. -12.]
[2025-01-13 11:28:10.487] [info] Newton iteration 6: r (abs) = 3.2448763364689164e-11 (tol = 1e-08), r (rel) = 9.485642415194253e-13 (tol = 1e-08)
[2025-01-13 11:28:10.487] [info] Newton solver finished in 6 iterations and 6 linear solver iterations.
[2025-01-13 11:28:10.813] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:13.095] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:15.122] [info] Newton iteration 2: r (abs) = 0.89878861295836 (tol = 1e-08), r (rel) = 0.030966555668592563 (tol = 1e-08)
[2025-01-13 11:28:15.392] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:17.413] [info] Newton iteration 3: r (abs) = 1.3835366464849637 (tol = 1e-08), r (rel) = 0.04766789873082136 (tol = 1e-08)
[2025-01-13 11:28:17.681] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:19.686] [info] Newton iteration 4: r (abs) = 0.0018509604182366244 (tol = 1e-08), r (rel) = 6.377235759921814e-05 (tol = 1e-08)
[2025-01-13 11:28:19.957] [info] PETSc Krylov solver starting to solve system.
[2025-01-13 11:28:21.993] [info] Newton iteration 5: r (abs) = 7.871826891356691e-06 (tol = 1e-08), r (rel) = 2.712132331565423e-07 (tol = 1e-08)
[2025-01-13 11:28:22.268] [info] PETSc Krylov solver starting to solve system.
Time step 9, Number of iterations 6, Load [  0.    0.  -13.5]
[2025-01-13 11:28:24.294] [info] Newton iteration 6: r (abs) = 2.0599218895847923e-13 (tol = 1e-08), r (rel) = 7.097184471087099e-15 (tol = 1e-08)
[2025-01-13 11:28:24.294] [info] Newton solver finished in 6 iterations and 6 linear solver iterations.
gif