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.

import numpy as np
import ufl

from petsc4py import PETSc
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.VectorFunctionSpace(domain, ("CG", 2))

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=PETSc.ScalarType)

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, PETSc.ScalarType((0, 0, 0)))
T = fem.Constant(domain, PETSc.ScalarType((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 = PETSc.ScalarType(1.0e4)
nu = PETSc.ScalarType(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 = fem.petsc.NonlinearProblem(F, u, bcs)

and then create and customize the Newton solver

from dolfinx import nls
solver = nls.petsc.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.

import pyvista
pyvista.set_jupyter_backend("ipygany")

grid = pyvista.UnstructuredGrid(*plot.create_vtk_mesh(domain, domain.topology.dim))

def plot_function(t, uh):
    """
    Create a figure of the concentration uh warped visualized in 3D at timet step t.
    """
    p = pyvista.Plotter()
    # Create grid defined by the function space for visualization of the function
    topology, cells, geometry = plot.create_vtk_mesh(uh.function_space)
    function_grid = pyvista.UnstructuredGrid(topology, cells, geometry)
    var_name = f"u({t})"
    values = np.zeros((geometry.shape[0], 3))
    values[:, :len(uh)] = uh.x.array.reshape(geometry.shape[0], len(uh))
    function_grid[var_name] = values
    function_grid.set_active_vectors(var_name)
    # Warp mesh by deformation
    warped = function_grid.warp_by_vector(var_name, factor=1)
    
    # Add mesh to plotter and visualize
    actor = p.add_mesh(warped)
    p.show_axes()
    if not pyvista.OFF_SCREEN:
       p.show()
    else:
        pyvista.start_xvfb()
        figure_as_array = p.screenshot(f"diffusion_{t:.2f}.png")
        # Clear plotter for next plot
        p.remove_actor(actor)

plot_function(0, u)

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

from dolfinx import log
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}")
    plot_function(n, u)
2022-11-21 08:02:32.152 (   3.012s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:32.865 (   3.724s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:33.223 (   4.083s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 22.2455 (tol = 1e-08) r (rel) = 0.134278(tol = 1e-08)
2022-11-21 08:02:33.451 (   4.310s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:33.810 (   4.669s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 2.43261 (tol = 1e-08) r (rel) = 0.0146837(tol = 1e-08)
2022-11-21 08:02:34.036 (   4.895s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:34.394 (   5.254s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 4.43158 (tol = 1e-08) r (rel) = 0.0267498(tol = 1e-08)
2022-11-21 08:02:34.620 (   5.480s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:34.981 (   5.840s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 0.144189 (tol = 1e-08) r (rel) = 0.000870353(tol = 1e-08)
2022-11-21 08:02:35.207 (   6.067s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:35.567 (   6.426s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 0.021424 (tol = 1e-08) r (rel) = 0.000129319(tol = 1e-08)
2022-11-21 08:02:35.794 (   6.653s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:36.153 (   7.013s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 7: r (abs) = 4.80067e-06 (tol = 1e-08) r (rel) = 2.89777e-08(tol = 1e-08)
2022-11-21 08:02:36.380 (   7.239s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
Time step 1, Number of iterations 8, Load [ 0.   0.  -1.5]
2022-11-21 08:02:36.740 (   7.599s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 8: r (abs) = 2.64464e-11 (tol = 1e-08) r (rel) = 1.59635e-13(tol = 1e-08)
2022-11-21 08:02:36.740 (   7.599s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
2022-11-21 08:02:37.034 (   7.894s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:37.618 (   8.478s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:37.979 (   8.838s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 17.3254 (tol = 1e-08) r (rel) = 0.117842(tol = 1e-08)
2022-11-21 08:02:38.205 (   9.064s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:38.563 (   9.423s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 5.14882 (tol = 1e-08) r (rel) = 0.0350207(tol = 1e-08)
2022-11-21 08:02:38.790 (   9.649s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:39.147 (  10.006s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 7.24003 (tol = 1e-08) r (rel) = 0.0492445(tol = 1e-08)
2022-11-21 08:02:39.373 (  10.232s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:39.731 (  10.591s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 0.777889 (tol = 1e-08) r (rel) = 0.00529096(tol = 1e-08)
2022-11-21 08:02:39.957 (  10.817s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:40.315 (  11.174s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 1.25525 (tol = 1e-08) r (rel) = 0.00853785(tol = 1e-08)
2022-11-21 08:02:40.541 (  11.400s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:40.898 (  11.757s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 7: r (abs) = 0.00849512 (tol = 1e-08) r (rel) = 5.77813e-05(tol = 1e-08)
2022-11-21 08:02:41.124 (  11.984s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:41.481 (  12.341s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 8: r (abs) = 0.000192107 (tol = 1e-08) r (rel) = 1.30665e-06(tol = 1e-08)
2022-11-21 08:02:41.708 (  12.568s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
Time step 2, Number of iterations 9, Load [ 0.  0. -3.]
2022-11-21 08:02:42.067 (  12.926s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 9: r (abs) = 1.70757e-10 (tol = 1e-08) r (rel) = 1.16144e-12(tol = 1e-08)
2022-11-21 08:02:42.067 (  12.926s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
2022-11-21 08:02:42.362 (  13.222s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:42.946 (  13.805s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:43.303 (  14.162s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 10.0011 (tol = 1e-08) r (rel) = 0.0887471(tol = 1e-08)
2022-11-21 08:02:43.530 (  14.389s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:43.887 (  14.747s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 5.33026 (tol = 1e-08) r (rel) = 0.0472992(tol = 1e-08)
2022-11-21 08:02:44.115 (  14.975s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:44.473 (  15.333s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 11.9901 (tol = 1e-08) r (rel) = 0.106397(tol = 1e-08)
2022-11-21 08:02:44.699 (  15.559s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:45.058 (  15.918s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 2.29702 (tol = 1e-08) r (rel) = 0.0203831(tol = 1e-08)
2022-11-21 08:02:45.285 (  16.144s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:45.644 (  16.503s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 3.90234 (tol = 1e-08) r (rel) = 0.0346282(tol = 1e-08)
2022-11-21 08:02:45.869 (  16.729s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:46.229 (  17.088s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 7: r (abs) = 0.236535 (tol = 1e-08) r (rel) = 0.00209895(tol = 1e-08)
2022-11-21 08:02:46.455 (  17.315s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:46.814 (  17.673s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 8: r (abs) = 0.0427142 (tol = 1e-08) r (rel) = 0.000379034(tol = 1e-08)
2022-11-21 08:02:47.040 (  17.900s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:47.400 (  18.259s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 9: r (abs) = 2.87798e-05 (tol = 1e-08) r (rel) = 2.55384e-07(tol = 1e-08)
Time step 3, Number of iterations 10, Load [ 0.   0.  -4.5]
2022-11-21 08:02:47.627 (  18.487s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:47.986 (  18.846s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 10: r (abs) = 6.09481e-10 (tol = 1e-08) r (rel) = 5.40836e-12(tol = 1e-08)
2022-11-21 08:02:47.986 (  18.846s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 10 iterations and 10 linear solver iterations.
2022-11-21 08:02:48.280 (  19.139s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:48.865 (  19.725s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:49.224 (  20.083s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 5.50693 (tol = 1e-08) r (rel) = 0.0653918(tol = 1e-08)
2022-11-21 08:02:49.451 (  20.311s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:49.808 (  20.667s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 26.2489 (tol = 1e-08) r (rel) = 0.311692(tol = 1e-08)
2022-11-21 08:02:50.034 (  20.893s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:50.393 (  21.252s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 2.30927 (tol = 1e-08) r (rel) = 0.0274213(tol = 1e-08)
2022-11-21 08:02:50.619 (  21.479s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:50.978 (  21.837s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 14.0562 (tol = 1e-08) r (rel) = 0.16691(tol = 1e-08)
2022-11-21 08:02:51.204 (  22.063s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:51.562 (  22.422s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 0.222774 (tol = 1e-08) r (rel) = 0.00264532(tol = 1e-08)
2022-11-21 08:02:51.790 (  22.649s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:52.146 (  23.005s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 7: r (abs) = 0.286671 (tol = 1e-08) r (rel) = 0.00340406(tol = 1e-08)
2022-11-21 08:02:52.372 (  23.232s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:52.730 (  23.589s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 8: r (abs) = 0.000321869 (tol = 1e-08) r (rel) = 3.82203e-06(tol = 1e-08)
Time step 4, Number of iterations 9, Load [ 0.  0. -6.]
2022-11-21 08:02:52.956 (  23.815s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:53.312 (  24.172s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 9: r (abs) = 2.63796e-07 (tol = 1e-08) r (rel) = 3.13244e-09(tol = 1e-08)
2022-11-21 08:02:53.312 (  24.172s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
2022-11-21 08:02:53.606 (  24.465s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:54.190 (  25.050s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:54.547 (  25.406s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 3.19462 (tol = 1e-08) r (rel) = 0.0496479(tol = 1e-08)
2022-11-21 08:02:54.773 (  25.632s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:55.131 (  25.991s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 7.71429 (tol = 1e-08) r (rel) = 0.119888(tol = 1e-08)
2022-11-21 08:02:55.358 (  26.217s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:55.714 (  26.573s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 0.850873 (tol = 1e-08) r (rel) = 0.0132235(tol = 1e-08)
2022-11-21 08:02:55.940 (  26.799s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:56.296 (  27.155s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 0.371434 (tol = 1e-08) r (rel) = 0.0057725(tol = 1e-08)
2022-11-21 08:02:56.522 (  27.381s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:56.877 (  27.737s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 0.00215066 (tol = 1e-08) r (rel) = 3.34236e-05(tol = 1e-08)
2022-11-21 08:02:57.103 (  27.963s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:57.460 (  28.320s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 7: r (abs) = 2.54607e-06 (tol = 1e-08) r (rel) = 3.95687e-08(tol = 1e-08)
Time step 5, Number of iterations 8, Load [ 0.   0.  -7.5]
2022-11-21 08:02:57.687 (  28.546s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:58.045 (  28.904s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 8: r (abs) = 1.51485e-13 (tol = 1e-08) r (rel) = 2.35425e-15(tol = 1e-08)
2022-11-21 08:02:58.045 (  28.904s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
2022-11-21 08:02:58.339 (  29.198s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:58.920 (  29.780s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:59.277 (  30.136s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 2.00649 (tol = 1e-08) r (rel) = 0.0395622(tol = 1e-08)
2022-11-21 08:02:59.503 (  30.362s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:02:59.859 (  30.719s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 4.60977 (tol = 1e-08) r (rel) = 0.0908914(tol = 1e-08)
2022-11-21 08:03:00.085 (  30.945s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:00.441 (  31.301s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 0.185372 (tol = 1e-08) r (rel) = 0.00365501(tol = 1e-08)
2022-11-21 08:03:00.668 (  31.528s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:01.026 (  31.885s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 0.024688 (tol = 1e-08) r (rel) = 0.000486777(tol = 1e-08)
2022-11-21 08:03:01.251 (  32.111s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:01.607 (  32.466s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 5.69254e-06 (tol = 1e-08) r (rel) = 1.12241e-07(tol = 1e-08)
Time step 6, Number of iterations 7, Load [ 0.  0. -9.]
2022-11-21 08:03:01.833 (  32.693s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:02.190 (  33.050s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 7: r (abs) = 2.67626e-11 (tol = 1e-08) r (rel) = 5.27682e-13(tol = 1e-08)
2022-11-21 08:03:02.190 (  33.050s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 7 iterations and 7 linear solver iterations.
2022-11-21 08:03:02.482 (  33.341s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:03.065 (  33.925s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:03.421 (  34.280s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 1.38506 (tol = 1e-08) r (rel) = 0.0336622(tol = 1e-08)
2022-11-21 08:03:03.647 (  34.507s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:04.002 (  34.862s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 3.03739 (tol = 1e-08) r (rel) = 0.07382(tol = 1e-08)
2022-11-21 08:03:04.228 (  35.088s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:04.584 (  35.443s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 0.0412386 (tol = 1e-08) r (rel) = 0.00100225(tol = 1e-08)
2022-11-21 08:03:04.810 (  35.669s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
Time step 7, Number of iterations 6, Load [  0.    0.  -10.5]
2022-11-21 08:03:05.167 (  36.026s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 0.00205057 (tol = 1e-08) r (rel) = 4.98364e-05(tol = 1e-08)
2022-11-21 08:03:05.393 (  36.253s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:05.749 (  36.609s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 1.78865e-08 (tol = 1e-08) r (rel) = 4.34709e-10(tol = 1e-08)
2022-11-21 08:03:05.749 (  36.609s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
2022-11-21 08:03:06.041 (  36.901s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:06.623 (  37.483s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:06.980 (  37.839s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 1.06336 (tol = 1e-08) r (rel) = 0.031085(tol = 1e-08)
2022-11-21 08:03:07.206 (  38.066s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:07.564 (  38.423s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 2.0477 (tol = 1e-08) r (rel) = 0.0598598(tol = 1e-08)
2022-11-21 08:03:07.790 (  38.650s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:08.147 (  39.007s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 0.00897719 (tol = 1e-08) r (rel) = 0.000262427(tol = 1e-08)
2022-11-21 08:03:08.374 (  39.233s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
Time step 8, Number of iterations 6, Load [  0.   0. -12.]
2022-11-21 08:03:08.733 (  39.593s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 0.000167422 (tol = 1e-08) r (rel) = 4.89419e-06(tol = 1e-08)
2022-11-21 08:03:08.960 (  39.819s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:09.317 (  40.177s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 3.23994e-11 (tol = 1e-08) r (rel) = 9.4712e-13(tol = 1e-08)
2022-11-21 08:03:09.317 (  40.177s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
2022-11-21 08:03:09.610 (  40.470s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:10.194 (  41.053s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:10.551 (  41.411s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 2: r (abs) = 0.898789 (tol = 1e-08) r (rel) = 0.0309666(tol = 1e-08)
2022-11-21 08:03:10.778 (  41.637s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:11.135 (  41.994s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 3: r (abs) = 1.38354 (tol = 1e-08) r (rel) = 0.0476679(tol = 1e-08)
2022-11-21 08:03:11.361 (  42.220s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:11.720 (  42.579s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 4: r (abs) = 0.00185096 (tol = 1e-08) r (rel) = 6.37724e-05(tol = 1e-08)
2022-11-21 08:03:11.946 (  42.805s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
Time step 9, Number of iterations 6, Load [  0.    0.  -13.5]
2022-11-21 08:03:12.302 (  43.162s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 5: r (abs) = 7.87183e-06 (tol = 1e-08) r (rel) = 2.71213e-07(tol = 1e-08)
2022-11-21 08:03:12.528 (  43.388s) [main            ]              petsc.cpp:677   INFO| PETSc Krylov solver starting to solve system.
2022-11-21 08:03:12.889 (  43.749s) [main            ]       NewtonSolver.cpp:36    INFO| Newton iteration 6: r (abs) = 2.75782e-13 (tol = 1e-08) r (rel) = 9.50169e-15(tol = 1e-08)
2022-11-21 08:03:12.889 (  43.749s) [main            ]       NewtonSolver.cpp:255   INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.