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-09-16 08:05:25.054 ( 6.730s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:27.653 ( 9.329s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:29.746 ( 11.421s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 22.2455 (tol = 1e-08) r (rel) = 0.134278(tol = 1e-08)
2024-09-16 08:05:30.072 ( 11.748s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:32.125 ( 13.801s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 2.43261 (tol = 1e-08) r (rel) = 0.0146837(tol = 1e-08)
2024-09-16 08:05:32.453 ( 14.129s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:34.520 ( 16.196s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 4.43158 (tol = 1e-08) r (rel) = 0.0267498(tol = 1e-08)
2024-09-16 08:05:34.857 ( 16.532s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:36.927 ( 18.603s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.144189 (tol = 1e-08) r (rel) = 0.000870353(tol = 1e-08)
2024-09-16 08:05:37.261 ( 18.937s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:39.321 ( 20.996s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 0.021424 (tol = 1e-08) r (rel) = 0.000129319(tol = 1e-08)
2024-09-16 08:05:39.655 ( 21.331s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:41.754 ( 23.430s) [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-09-16 08:05:42.082 ( 23.758s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
Time step 1, Number of iterations 8, Load [ 0. 0. -1.5]
2024-09-16 08:05:44.143 ( 25.819s) [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-09-16 08:05:44.143 ( 25.819s) [main ] NewtonSolver.cpp:252 INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
2024-09-16 08:05:44.657 ( 26.333s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:47.082 ( 28.758s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:49.162 ( 30.838s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 17.3254 (tol = 1e-08) r (rel) = 0.117842(tol = 1e-08)
2024-09-16 08:05:49.494 ( 31.170s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:51.547 ( 33.222s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 5.14882 (tol = 1e-08) r (rel) = 0.0350207(tol = 1e-08)
2024-09-16 08:05:51.875 ( 33.551s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:53.955 ( 35.630s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 7.24003 (tol = 1e-08) r (rel) = 0.0492445(tol = 1e-08)
2024-09-16 08:05:54.281 ( 35.957s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:56.368 ( 38.044s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.777889 (tol = 1e-08) r (rel) = 0.00529096(tol = 1e-08)
2024-09-16 08:05:56.743 ( 38.419s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:05:58.801 ( 40.476s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 1.25525 (tol = 1e-08) r (rel) = 0.00853785(tol = 1e-08)
2024-09-16 08:05:59.130 ( 40.806s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:01.212 ( 42.887s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 7: r (abs) = 0.00849512 (tol = 1e-08) r (rel) = 5.77812e-05(tol = 1e-08)
2024-09-16 08:06:01.546 ( 43.222s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:03.627 ( 45.303s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 8: r (abs) = 0.000192107 (tol = 1e-08) r (rel) = 1.30665e-06(tol = 1e-08)
2024-09-16 08:06:03.957 ( 45.633s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:06.044 ( 47.720s) [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-09-16 08:06:06.044 ( 47.720s) [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-09-16 08:06:06.434 ( 48.109s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:08.863 ( 50.539s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:10.936 ( 52.612s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 10.0011 (tol = 1e-08) r (rel) = 0.0887471(tol = 1e-08)
2024-09-16 08:06:11.267 ( 52.943s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:13.338 ( 55.014s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 5.33026 (tol = 1e-08) r (rel) = 0.0472992(tol = 1e-08)
2024-09-16 08:06:13.676 ( 55.351s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:15.779 ( 57.454s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 11.9901 (tol = 1e-08) r (rel) = 0.106397(tol = 1e-08)
2024-09-16 08:06:16.117 ( 57.793s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:18.184 ( 59.859s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 2.29702 (tol = 1e-08) r (rel) = 0.0203831(tol = 1e-08)
2024-09-16 08:06:18.510 ( 60.185s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:20.584 ( 62.260s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 3.90234 (tol = 1e-08) r (rel) = 0.0346282(tol = 1e-08)
2024-09-16 08:06:20.914 ( 62.590s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:23.003 ( 64.679s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 7: r (abs) = 0.236535 (tol = 1e-08) r (rel) = 0.00209895(tol = 1e-08)
2024-09-16 08:06:23.339 ( 65.015s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:25.418 ( 67.094s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 8: r (abs) = 0.0427142 (tol = 1e-08) r (rel) = 0.000379034(tol = 1e-08)
2024-09-16 08:06:25.750 ( 67.425s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:27.845 ( 69.520s) [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-09-16 08:06:28.181 ( 69.856s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:30.237 ( 71.913s) [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-09-16 08:06:30.237 ( 71.913s) [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-09-16 08:06:30.624 ( 72.300s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:33.038 ( 74.714s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:35.120 ( 76.796s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 5.50693 (tol = 1e-08) r (rel) = 0.0653918(tol = 1e-08)
2024-09-16 08:06:35.456 ( 77.132s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:37.536 ( 79.212s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 26.2489 (tol = 1e-08) r (rel) = 0.311692(tol = 1e-08)
2024-09-16 08:06:37.870 ( 79.546s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:39.941 ( 81.617s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 2.30927 (tol = 1e-08) r (rel) = 0.0274213(tol = 1e-08)
2024-09-16 08:06:40.277 ( 81.953s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:42.368 ( 84.044s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 14.0562 (tol = 1e-08) r (rel) = 0.16691(tol = 1e-08)
2024-09-16 08:06:42.707 ( 84.383s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:44.780 ( 86.455s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 0.222774 (tol = 1e-08) r (rel) = 0.00264532(tol = 1e-08)
2024-09-16 08:06:45.114 ( 86.790s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:47.234 ( 88.910s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 7: r (abs) = 0.286671 (tol = 1e-08) r (rel) = 0.00340406(tol = 1e-08)
2024-09-16 08:06:47.561 ( 89.236s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:49.628 ( 91.304s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 8: r (abs) = 0.000321869 (tol = 1e-08) r (rel) = 3.82203e-06(tol = 1e-08)
2024-09-16 08:06:49.964 ( 91.640s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:52.047 ( 93.722s) [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-09-16 08:06:52.047 ( 93.722s) [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-09-16 08:06:52.446 ( 94.122s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:54.858 ( 96.533s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:56.944 ( 98.620s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 3.19462 (tol = 1e-08) r (rel) = 0.0496479(tol = 1e-08)
2024-09-16 08:06:57.284 ( 98.960s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:06:59.345 ( 101.021s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 7.71429 (tol = 1e-08) r (rel) = 0.119888(tol = 1e-08)
2024-09-16 08:06:59.684 ( 101.360s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:01.760 ( 103.436s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.850873 (tol = 1e-08) r (rel) = 0.0132235(tol = 1e-08)
2024-09-16 08:07:02.086 ( 103.762s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:04.176 ( 105.852s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.371434 (tol = 1e-08) r (rel) = 0.0057725(tol = 1e-08)
2024-09-16 08:07:04.513 ( 106.189s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:06.605 ( 108.281s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 6: r (abs) = 0.00215066 (tol = 1e-08) r (rel) = 3.34236e-05(tol = 1e-08)
2024-09-16 08:07:06.983 ( 108.658s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:09.059 ( 110.735s) [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-09-16 08:07:09.396 ( 111.072s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:11.484 ( 113.160s) [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-09-16 08:07:11.484 ( 113.160s) [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-09-16 08:07:11.872 ( 113.548s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:14.275 ( 115.951s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:16.340 ( 118.016s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 2.00649 (tol = 1e-08) r (rel) = 0.0395622(tol = 1e-08)
2024-09-16 08:07:16.676 ( 118.352s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:18.757 ( 120.433s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 4.60977 (tol = 1e-08) r (rel) = 0.0908914(tol = 1e-08)
2024-09-16 08:07:19.084 ( 120.760s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:21.172 ( 122.848s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.185372 (tol = 1e-08) r (rel) = 0.00365501(tol = 1e-08)
2024-09-16 08:07:21.500 ( 123.176s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:23.609 ( 125.285s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.024688 (tol = 1e-08) r (rel) = 0.000486777(tol = 1e-08)
2024-09-16 08:07:23.946 ( 125.621s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:26.033 ( 127.708s) [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-09-16 08:07:26.359 ( 128.034s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:28.443 ( 130.119s) [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-09-16 08:07:28.443 ( 130.119s) [main ] NewtonSolver.cpp:252 INFO| Newton solver finished in 7 iterations and 7 linear solver iterations.
Time step 6, Number of iterations 7, Load [ 0. 0. -9.]
2024-09-16 08:07:28.840 ( 130.516s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:31.268 ( 132.944s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:33.352 ( 135.028s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 1.38506 (tol = 1e-08) r (rel) = 0.0336622(tol = 1e-08)
2024-09-16 08:07:33.683 ( 135.359s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:35.768 ( 137.444s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 3.03739 (tol = 1e-08) r (rel) = 0.07382(tol = 1e-08)
2024-09-16 08:07:36.097 ( 137.773s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:38.178 ( 139.854s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.0412386 (tol = 1e-08) r (rel) = 0.00100225(tol = 1e-08)
2024-09-16 08:07:38.516 ( 140.192s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:40.605 ( 142.281s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.00205057 (tol = 1e-08) r (rel) = 4.98364e-05(tol = 1e-08)
2024-09-16 08:07:40.944 ( 142.619s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:43.026 ( 144.702s) [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-09-16 08:07:43.026 ( 144.702s) [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-09-16 08:07:43.412 ( 145.088s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:45.814 ( 147.490s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:47.872 ( 149.547s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 1.06336 (tol = 1e-08) r (rel) = 0.031085(tol = 1e-08)
2024-09-16 08:07:48.203 ( 149.879s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:50.227 ( 151.902s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 2.0477 (tol = 1e-08) r (rel) = 0.0598598(tol = 1e-08)
2024-09-16 08:07:50.553 ( 152.229s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:52.579 ( 154.254s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.00897719 (tol = 1e-08) r (rel) = 0.000262427(tol = 1e-08)
2024-09-16 08:07:52.914 ( 154.589s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:54.925 ( 156.600s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 5: r (abs) = 0.000167422 (tol = 1e-08) r (rel) = 4.89419e-06(tol = 1e-08)
2024-09-16 08:07:55.249 ( 156.924s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:07:57.267 ( 158.943s) [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-09-16 08:07:57.268 ( 158.943s) [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-09-16 08:07:57.664 ( 159.339s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:08:00.043 ( 161.719s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:08:02.120 ( 163.796s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 2: r (abs) = 0.898789 (tol = 1e-08) r (rel) = 0.0309666(tol = 1e-08)
2024-09-16 08:08:02.456 ( 164.132s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:08:04.490 ( 166.165s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 3: r (abs) = 1.38354 (tol = 1e-08) r (rel) = 0.0476679(tol = 1e-08)
2024-09-16 08:08:04.824 ( 166.499s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:08:06.861 ( 168.537s) [main ] NewtonSolver.cpp:38 INFO| Newton iteration 4: r (abs) = 0.00185096 (tol = 1e-08) r (rel) = 6.37724e-05(tol = 1e-08)
2024-09-16 08:08:07.188 ( 168.863s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:08:09.263 ( 170.939s) [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-09-16 08:08:09.602 ( 171.277s) [main ] petsc.cpp:700 INFO| PETSc Krylov solver starting to solve system.
2024-09-16 08:08:11.705 ( 173.381s) [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-09-16 08:08:11.705 ( 173.381s) [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]