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()
2024-05-13 18:35:32.893 (   6.575s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:35.376 (   9.058s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:37.343 (  11.025s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 22.2455 (tol = 1e-08) r (rel) = 0.134278(tol = 1e-08)
2024-05-13 18:35:37.675 (  11.356s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:39.648 (  13.330s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 2.43261 (tol = 1e-08) r (rel) = 0.0146837(tol = 1e-08)
2024-05-13 18:35:39.982 (  13.664s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:41.952 (  15.634s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 4.43158 (tol = 1e-08) r (rel) = 0.0267498(tol = 1e-08)
2024-05-13 18:35:42.276 (  15.957s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:44.249 (  17.931s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 0.144189 (tol = 1e-08) r (rel) = 0.000870353(tol = 1e-08)
2024-05-13 18:35:44.580 (  18.262s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:46.551 (  20.233s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 0.021424 (tol = 1e-08) r (rel) = 0.000129319(tol = 1e-08)
2024-05-13 18:35:46.875 (  20.556s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:48.843 (  22.525s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 7: r (abs) = 4.80065e-06 (tol = 1e-08) r (rel) = 2.89776e-08(tol = 1e-08)
2024-05-13 18:35:49.172 (  22.854s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:51.137 (  24.818s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 8: r (abs) = 2.65767e-11 (tol = 1e-08) r (rel) = 1.60422e-13(tol = 1e-08)
2024-05-13 18:35:51.137 (  24.818s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
Time step 1, Number of iterations 8, Load [ 0.   0.  -1.5]
2024-05-13 18:35:51.652 (  25.333s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:53.943 (  27.625s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:55.909 (  29.591s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 17.3254 (tol = 1e-08) r (rel) = 0.117842(tol = 1e-08)
2024-05-13 18:35:56.235 (  29.917s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:35:58.205 (  31.887s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 5.14882 (tol = 1e-08) r (rel) = 0.0350207(tol = 1e-08)
2024-05-13 18:35:58.542 (  32.224s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:00.520 (  34.201s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 7.24003 (tol = 1e-08) r (rel) = 0.0492445(tol = 1e-08)
2024-05-13 18:36:00.848 (  34.530s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:02.817 (  36.499s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 0.777889 (tol = 1e-08) r (rel) = 0.00529096(tol = 1e-08)
2024-05-13 18:36:03.154 (  36.836s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:05.127 (  38.809s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 1.25525 (tol = 1e-08) r (rel) = 0.00853785(tol = 1e-08)
2024-05-13 18:36:05.462 (  39.144s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:07.429 (  41.111s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 7: r (abs) = 0.00849512 (tol = 1e-08) r (rel) = 5.77812e-05(tol = 1e-08)
2024-05-13 18:36:07.765 (  41.447s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:09.737 (  43.419s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 8: r (abs) = 0.000192107 (tol = 1e-08) r (rel) = 1.30665e-06(tol = 1e-08)
2024-05-13 18:36:10.073 (  43.755s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:12.054 (  45.736s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 9: r (abs) = 1.70902e-10 (tol = 1e-08) r (rel) = 1.16242e-12(tol = 1e-08)
2024-05-13 18:36:12.054 (  45.736s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
Time step 2, Number of iterations 9, Load [ 0.  0. -3.]
2024-05-13 18:36:12.442 (  46.124s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:14.732 (  48.414s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:16.695 (  50.376s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 10.0011 (tol = 1e-08) r (rel) = 0.0887471(tol = 1e-08)
2024-05-13 18:36:17.030 (  50.712s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:18.988 (  52.670s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 5.33026 (tol = 1e-08) r (rel) = 0.0472992(tol = 1e-08)
2024-05-13 18:36:19.321 (  53.003s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:21.292 (  54.974s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 11.9901 (tol = 1e-08) r (rel) = 0.106397(tol = 1e-08)
2024-05-13 18:36:21.626 (  55.308s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:23.598 (  57.280s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 2.29702 (tol = 1e-08) r (rel) = 0.0203831(tol = 1e-08)
2024-05-13 18:36:23.932 (  57.614s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:25.903 (  59.585s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 3.90234 (tol = 1e-08) r (rel) = 0.0346282(tol = 1e-08)
2024-05-13 18:36:26.222 (  59.904s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:28.194 (  61.875s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 7: r (abs) = 0.236535 (tol = 1e-08) r (rel) = 0.00209895(tol = 1e-08)
2024-05-13 18:36:28.518 (  62.199s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:30.484 (  64.165s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 8: r (abs) = 0.0427142 (tol = 1e-08) r (rel) = 0.000379034(tol = 1e-08)
2024-05-13 18:36:30.807 (  64.489s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:32.776 (  66.458s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 9: r (abs) = 2.87798e-05 (tol = 1e-08) r (rel) = 2.55384e-07(tol = 1e-08)
2024-05-13 18:36:33.100 (  66.782s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:35.073 (  68.755s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 10: r (abs) = 6.09146e-10 (tol = 1e-08) r (rel) = 5.40539e-12(tol = 1e-08)
2024-05-13 18:36:35.073 (  68.755s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 10 iterations and 10 linear solver iterations.
Time step 3, Number of iterations 10, Load [ 0.   0.  -4.5]
2024-05-13 18:36:35.457 (  69.139s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:37.734 (  71.416s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:39.689 (  73.371s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 5.50693 (tol = 1e-08) r (rel) = 0.0653918(tol = 1e-08)
2024-05-13 18:36:40.007 (  73.689s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:41.962 (  75.643s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 26.2489 (tol = 1e-08) r (rel) = 0.311692(tol = 1e-08)
2024-05-13 18:36:42.289 (  75.970s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:44.250 (  77.932s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 2.30927 (tol = 1e-08) r (rel) = 0.0274213(tol = 1e-08)
2024-05-13 18:36:44.586 (  78.267s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:46.538 (  80.220s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 14.0562 (tol = 1e-08) r (rel) = 0.16691(tol = 1e-08)
2024-05-13 18:36:46.861 (  80.543s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:48.824 (  82.506s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 0.222774 (tol = 1e-08) r (rel) = 0.00264532(tol = 1e-08)
2024-05-13 18:36:49.158 (  82.840s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:51.127 (  84.809s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 7: r (abs) = 0.286671 (tol = 1e-08) r (rel) = 0.00340406(tol = 1e-08)
2024-05-13 18:36:51.460 (  85.142s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:53.437 (  87.118s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 8: r (abs) = 0.000321869 (tol = 1e-08) r (rel) = 3.82203e-06(tol = 1e-08)
2024-05-13 18:36:53.770 (  87.452s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:55.735 (  89.416s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 9: r (abs) = 2.63796e-07 (tol = 1e-08) r (rel) = 3.13244e-09(tol = 1e-08)
2024-05-13 18:36:55.735 (  89.416s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
Time step 4, Number of iterations 9, Load [ 0.  0. -6.]
2024-05-13 18:36:56.122 (  89.804s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:36:58.407 (  92.089s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:00.414 (  94.096s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 3.19462 (tol = 1e-08) r (rel) = 0.0496479(tol = 1e-08)
2024-05-13 18:37:00.748 (  94.429s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:02.712 (  96.393s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 7.71429 (tol = 1e-08) r (rel) = 0.119888(tol = 1e-08)
2024-05-13 18:37:03.045 (  96.726s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:05.008 (  98.690s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 0.850873 (tol = 1e-08) r (rel) = 0.0132235(tol = 1e-08)
2024-05-13 18:37:05.337 (  99.018s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:07.307 ( 100.989s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 0.371434 (tol = 1e-08) r (rel) = 0.0057725(tol = 1e-08)
2024-05-13 18:37:07.634 ( 101.316s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:09.603 ( 103.284s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 0.00215066 (tol = 1e-08) r (rel) = 3.34236e-05(tol = 1e-08)
2024-05-13 18:37:09.938 ( 103.620s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:11.907 ( 105.589s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 7: r (abs) = 2.54607e-06 (tol = 1e-08) r (rel) = 3.95687e-08(tol = 1e-08)
2024-05-13 18:37:12.238 ( 105.920s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
Time step 5, Number of iterations 8, Load [ 0.   0.  -7.5]
2024-05-13 18:37:14.208 ( 107.890s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 8: r (abs) = 1.73699e-13 (tol = 1e-08) r (rel) = 2.69948e-15(tol = 1e-08)
2024-05-13 18:37:14.209 ( 107.890s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
2024-05-13 18:37:14.603 ( 108.285s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:16.910 ( 110.592s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:18.882 ( 112.564s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 2.00649 (tol = 1e-08) r (rel) = 0.0395622(tol = 1e-08)
2024-05-13 18:37:19.210 ( 112.891s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:21.180 ( 114.862s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 4.60977 (tol = 1e-08) r (rel) = 0.0908914(tol = 1e-08)
2024-05-13 18:37:21.513 ( 115.195s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:23.482 ( 117.164s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 0.185372 (tol = 1e-08) r (rel) = 0.00365501(tol = 1e-08)
2024-05-13 18:37:23.815 ( 117.496s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:25.787 ( 119.469s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 0.024688 (tol = 1e-08) r (rel) = 0.000486777(tol = 1e-08)
2024-05-13 18:37:26.121 ( 119.803s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
Time step 6, Number of iterations 7, Load [ 0.  0. -9.]
2024-05-13 18:37:28.096 ( 121.777s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 5.69255e-06 (tol = 1e-08) r (rel) = 1.12241e-07(tol = 1e-08)
2024-05-13 18:37:28.429 ( 122.111s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:30.398 ( 124.080s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 7: r (abs) = 2.62006e-11 (tol = 1e-08) r (rel) = 5.166e-13(tol = 1e-08)
2024-05-13 18:37:30.398 ( 124.080s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 7 iterations and 7 linear solver iterations.
2024-05-13 18:37:30.798 ( 124.480s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:33.093 ( 126.775s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:35.061 ( 128.743s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 1.38506 (tol = 1e-08) r (rel) = 0.0336622(tol = 1e-08)
2024-05-13 18:37:35.394 ( 129.076s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:37.364 ( 131.046s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 3.03739 (tol = 1e-08) r (rel) = 0.07382(tol = 1e-08)
2024-05-13 18:37:37.698 ( 131.380s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:39.668 ( 133.350s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 0.0412386 (tol = 1e-08) r (rel) = 0.00100225(tol = 1e-08)
2024-05-13 18:37:39.989 ( 133.671s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:41.961 ( 135.643s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 0.00205057 (tol = 1e-08) r (rel) = 4.98364e-05(tol = 1e-08)
2024-05-13 18:37:42.289 ( 135.971s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:44.265 ( 137.947s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 1.78867e-08 (tol = 1e-08) r (rel) = 4.34714e-10(tol = 1e-08)
2024-05-13 18:37:44.265 ( 137.947s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
Time step 7, Number of iterations 6, Load [  0.    0.  -10.5]
2024-05-13 18:37:44.661 ( 138.343s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:46.966 ( 140.648s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:48.937 ( 142.619s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 1.06336 (tol = 1e-08) r (rel) = 0.031085(tol = 1e-08)
2024-05-13 18:37:49.272 ( 142.953s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:51.242 ( 144.924s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 2.0477 (tol = 1e-08) r (rel) = 0.0598598(tol = 1e-08)
2024-05-13 18:37:51.576 ( 145.258s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:53.547 ( 147.229s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 0.00897719 (tol = 1e-08) r (rel) = 0.000262427(tol = 1e-08)
2024-05-13 18:37:53.882 ( 147.564s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:55.852 ( 149.534s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 0.000167422 (tol = 1e-08) r (rel) = 4.89419e-06(tol = 1e-08)
2024-05-13 18:37:56.182 ( 149.864s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:37:58.155 ( 151.837s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 3.24488e-11 (tol = 1e-08) r (rel) = 9.48564e-13(tol = 1e-08)
2024-05-13 18:37:58.155 ( 151.837s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
Time step 8, Number of iterations 6, Load [  0.   0. -12.]
2024-05-13 18:37:58.552 ( 152.234s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:38:00.870 ( 154.551s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:38:02.844 ( 156.526s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 2: r (abs) = 0.898789 (tol = 1e-08) r (rel) = 0.0309666(tol = 1e-08)
2024-05-13 18:38:03.178 ( 156.859s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:38:05.146 ( 158.827s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 3: r (abs) = 1.38354 (tol = 1e-08) r (rel) = 0.0476679(tol = 1e-08)
2024-05-13 18:38:05.477 ( 159.158s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:38:07.445 ( 161.127s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 4: r (abs) = 0.00185096 (tol = 1e-08) r (rel) = 6.37724e-05(tol = 1e-08)
2024-05-13 18:38:07.776 ( 161.458s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
2024-05-13 18:38:09.744 ( 163.426s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 5: r (abs) = 7.87183e-06 (tol = 1e-08) r (rel) = 2.71213e-07(tol = 1e-08)
2024-05-13 18:38:10.068 ( 163.749s) [main            ]              petsc.cpp:700   INFO| PETSc Krylov solver starting to solve system.
Time step 9, Number of iterations 6, Load [  0.    0.  -13.5]
2024-05-13 18:38:12.037 ( 165.719s) [main            ]       NewtonSolver.cpp:38    INFO| Newton iteration 6: r (abs) = 2.05992e-13 (tol = 1e-08) r (rel) = 7.09718e-15(tol = 1e-08)
2024-05-13 18:38:12.037 ( 165.719s) [main            ]       NewtonSolver.cpp:252   INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
gif