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-10-02 17:13:27.920 ( 6.536s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:30.442 ( 9.057s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:32.444 ( 11.059s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 22.2455 (tol = 1e-08) r (rel) = 0.134278(tol = 1e-08)
2024-10-02 17:13:32.775 ( 11.391s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:34.775 ( 13.390s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 2.43261 (tol = 1e-08) r (rel) = 0.0146837(tol = 1e-08)
2024-10-02 17:13:35.106 ( 13.721s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:37.105 ( 15.721s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 4.43158 (tol = 1e-08) r (rel) = 0.0267498(tol = 1e-08)
2024-10-02 17:13:37.441 ( 16.057s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:39.444 ( 18.059s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.144189 (tol = 1e-08) r (rel) = 0.000870353(tol = 1e-08)
2024-10-02 17:13:39.776 ( 18.391s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:41.782 ( 20.397s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 0.021424 (tol = 1e-08) r (rel) = 0.000129319(tol = 1e-08)
2024-10-02 17:13:42.119 ( 20.735s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:44.115 ( 22.730s) [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-10-02 17:13:44.457 ( 23.072s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:46.465 ( 25.080s) [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-10-02 17:13:46.465 ( 25.080s) [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-10-02 17:13:46.971 ( 25.587s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:49.318 ( 27.933s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:51.320 ( 29.935s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 17.3254 (tol = 1e-08) r (rel) = 0.117842(tol = 1e-08)
2024-10-02 17:13:51.659 ( 30.274s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:53.659 ( 32.275s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 5.14882 (tol = 1e-08) r (rel) = 0.0350207(tol = 1e-08)
2024-10-02 17:13:54.002 ( 32.617s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:56.004 ( 34.620s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 7.24003 (tol = 1e-08) r (rel) = 0.0492445(tol = 1e-08)
2024-10-02 17:13:56.346 ( 34.961s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:13:58.384 ( 36.999s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.777889 (tol = 1e-08) r (rel) = 0.00529096(tol = 1e-08)
2024-10-02 17:13:58.722 ( 37.337s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:00.724 ( 39.340s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 1.25525 (tol = 1e-08) r (rel) = 0.00853785(tol = 1e-08)
2024-10-02 17:14:01.056 ( 39.672s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:03.080 ( 41.695s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 7: r (abs) = 0.00849512 (tol = 1e-08) r (rel) = 5.77812e-05(tol = 1e-08)
2024-10-02 17:14:03.422 ( 42.038s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:05.421 ( 44.036s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 8: r (abs) = 0.000192107 (tol = 1e-08) r (rel) = 1.30665e-06(tol = 1e-08)
2024-10-02 17:14:05.758 ( 44.373s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
Time step 2, Number of iterations 9, Load [ 0. 0. -3.]
2024-10-02 17:14:07.768 ( 46.384s) [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-10-02 17:14:07.768 ( 46.384s) [main ] NewtonSolver.cpp:252 INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
2024-10-02 17:14:08.158 ( 46.773s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:10.499 ( 49.115s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:12.510 ( 51.126s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 10.0011 (tol = 1e-08) r (rel) = 0.0887471(tol = 1e-08)
2024-10-02 17:14:12.838 ( 51.453s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:14.842 ( 53.457s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 5.33026 (tol = 1e-08) r (rel) = 0.0472992(tol = 1e-08)
2024-10-02 17:14:15.181 ( 53.797s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:17.189 ( 55.804s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 11.9901 (tol = 1e-08) r (rel) = 0.106397(tol = 1e-08)
2024-10-02 17:14:17.522 ( 56.137s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:19.530 ( 58.145s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 2.29702 (tol = 1e-08) r (rel) = 0.0203831(tol = 1e-08)
2024-10-02 17:14:19.858 ( 58.473s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:21.872 ( 60.488s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 3.90234 (tol = 1e-08) r (rel) = 0.0346282(tol = 1e-08)
2024-10-02 17:14:22.211 ( 60.826s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:24.216 ( 62.831s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 7: r (abs) = 0.236535 (tol = 1e-08) r (rel) = 0.00209895(tol = 1e-08)
2024-10-02 17:14:24.552 ( 63.167s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:26.551 ( 65.167s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 8: r (abs) = 0.0427142 (tol = 1e-08) r (rel) = 0.000379034(tol = 1e-08)
2024-10-02 17:14:26.885 ( 65.501s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:28.891 ( 67.507s) [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-10-02 17:14:29.231 ( 67.846s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:31.237 ( 69.852s) [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-10-02 17:14:31.237 ( 69.852s) [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-10-02 17:14:31.626 ( 70.242s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:33.962 ( 72.578s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:35.960 ( 74.576s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 5.50693 (tol = 1e-08) r (rel) = 0.0653918(tol = 1e-08)
2024-10-02 17:14:36.299 ( 74.914s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:38.296 ( 76.912s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 26.2489 (tol = 1e-08) r (rel) = 0.311692(tol = 1e-08)
2024-10-02 17:14:38.633 ( 77.249s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:40.651 ( 79.267s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 2.30927 (tol = 1e-08) r (rel) = 0.0274213(tol = 1e-08)
2024-10-02 17:14:40.981 ( 79.597s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:42.981 ( 81.597s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 14.0562 (tol = 1e-08) r (rel) = 0.16691(tol = 1e-08)
2024-10-02 17:14:43.320 ( 81.936s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:45.330 ( 83.945s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 0.222774 (tol = 1e-08) r (rel) = 0.00264532(tol = 1e-08)
2024-10-02 17:14:45.668 ( 84.283s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:47.666 ( 86.281s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 7: r (abs) = 0.286671 (tol = 1e-08) r (rel) = 0.00340406(tol = 1e-08)
2024-10-02 17:14:48.000 ( 86.615s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:50.000 ( 88.615s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 8: r (abs) = 0.000321869 (tol = 1e-08) r (rel) = 3.82203e-06(tol = 1e-08)
2024-10-02 17:14:50.338 ( 88.954s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:52.344 ( 90.959s) [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-10-02 17:14:52.344 ( 90.959s) [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-10-02 17:14:52.732 ( 91.347s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:55.076 ( 93.692s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:57.079 ( 95.694s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 3.19462 (tol = 1e-08) r (rel) = 0.0496479(tol = 1e-08)
2024-10-02 17:14:57.411 ( 96.027s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:14:59.414 ( 98.029s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 7.71429 (tol = 1e-08) r (rel) = 0.119888(tol = 1e-08)
2024-10-02 17:14:59.741 ( 98.356s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:01.753 ( 100.368s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.850873 (tol = 1e-08) r (rel) = 0.0132235(tol = 1e-08)
2024-10-02 17:15:02.083 ( 100.698s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:04.082 ( 102.698s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.371434 (tol = 1e-08) r (rel) = 0.0057725(tol = 1e-08)
2024-10-02 17:15:04.418 ( 103.033s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:06.410 ( 105.026s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 0.00215066 (tol = 1e-08) r (rel) = 3.34236e-05(tol = 1e-08)
2024-10-02 17:15:06.743 ( 105.358s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:08.760 ( 107.375s) [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-10-02 17:15:09.088 ( 107.703s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:11.120 ( 109.735s) [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-10-02 17:15:11.120 ( 109.735s) [main ] NewtonSolver.cpp:252 INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
Time step 5, Number of iterations 8, Load [ 0. 0. -7.5]
2024-10-02 17:15:11.527 ( 110.142s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:13.869 ( 112.484s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:15.868 ( 114.484s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 2.00649 (tol = 1e-08) r (rel) = 0.0395622(tol = 1e-08)
2024-10-02 17:15:16.202 ( 114.817s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:18.209 ( 116.824s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 4.60977 (tol = 1e-08) r (rel) = 0.0908914(tol = 1e-08)
2024-10-02 17:15:18.546 ( 117.162s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:20.564 ( 119.179s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.185372 (tol = 1e-08) r (rel) = 0.00365501(tol = 1e-08)
2024-10-02 17:15:20.903 ( 119.518s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:22.922 ( 121.537s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.024688 (tol = 1e-08) r (rel) = 0.000486777(tol = 1e-08)
2024-10-02 17:15:23.261 ( 121.876s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:25.259 ( 123.874s) [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-10-02 17:15:25.598 ( 124.213s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
Time step 6, Number of iterations 7, Load [ 0. 0. -9.]
2024-10-02 17:15:27.597 ( 126.212s) [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-10-02 17:15:27.597 ( 126.212s) [main ] NewtonSolver.cpp:252 INFO| Newton solver finished in 7 iterations and 7 linear solver iterations.
2024-10-02 17:15:27.999 ( 126.614s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:30.326 ( 128.941s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:32.353 ( 130.968s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 1.38506 (tol = 1e-08) r (rel) = 0.0336622(tol = 1e-08)
2024-10-02 17:15:32.692 ( 131.308s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:34.696 ( 133.311s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 3.03739 (tol = 1e-08) r (rel) = 0.07382(tol = 1e-08)
2024-10-02 17:15:35.036 ( 133.651s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:37.036 ( 135.652s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.0412386 (tol = 1e-08) r (rel) = 0.00100225(tol = 1e-08)
2024-10-02 17:15:37.364 ( 135.979s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:39.366 ( 137.981s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.00205057 (tol = 1e-08) r (rel) = 4.98364e-05(tol = 1e-08)
2024-10-02 17:15:39.705 ( 138.320s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:41.715 ( 140.330s) [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-10-02 17:15:41.715 ( 140.330s) [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-10-02 17:15:42.117 ( 140.732s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:44.463 ( 143.078s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:46.463 ( 145.078s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 1.06336 (tol = 1e-08) r (rel) = 0.031085(tol = 1e-08)
2024-10-02 17:15:46.801 ( 145.416s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:48.799 ( 147.414s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 2.0477 (tol = 1e-08) r (rel) = 0.0598598(tol = 1e-08)
2024-10-02 17:15:49.126 ( 147.742s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:51.132 ( 149.747s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.00897719 (tol = 1e-08) r (rel) = 0.000262427(tol = 1e-08)
2024-10-02 17:15:51.458 ( 150.074s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:53.461 ( 152.076s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.000167422 (tol = 1e-08) r (rel) = 4.89419e-06(tol = 1e-08)
2024-10-02 17:15:53.800 ( 152.416s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:55.819 ( 154.434s) [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-10-02 17:15:55.819 ( 154.434s) [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-10-02 17:15:56.221 ( 154.836s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:15:58.558 ( 157.173s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:16:00.570 ( 159.186s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 0.898789 (tol = 1e-08) r (rel) = 0.0309666(tol = 1e-08)
2024-10-02 17:16:00.912 ( 159.528s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:16:02.911 ( 161.526s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 1.38354 (tol = 1e-08) r (rel) = 0.0476679(tol = 1e-08)
2024-10-02 17:16:03.239 ( 161.854s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:16:05.242 ( 163.858s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.00185096 (tol = 1e-08) r (rel) = 6.37724e-05(tol = 1e-08)
2024-10-02 17:16:05.576 ( 164.191s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:16:07.578 ( 166.194s) [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-10-02 17:16:07.919 ( 166.535s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-10-02 17:16:09.918 ( 168.533s) [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-10-02 17:16:09.918 ( 168.534s) [main ] NewtonSolver.cpp:252 INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
Time step 9, Number of iterations 6, Load [ 0. 0. -13.5]