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
import matplotlib.pyplot as plt
pyvista.start_xvfb()
plotter = pyvista.Plotter()
plotter.open_gif("deformation.gif", fps=3)
topology, cells, geometry = plot.create_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 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}")
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)
plotter.update_coordinates(warped_n.points.copy(), render=False)
plotter.update_scalar_bar_range([0, 10])
plotter.update_scalars(magnitude.x.array)
plotter.write_frame()
plotter.close()
2023-03-27 08:06:40.978 ( 7.057s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:41.681 ( 7.760s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.702939, 0.640000, 0.070000 (PETSc Krylov solver)
2023-03-27 08:06:41.931 ( 8.010s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:42.490 ( 8.569s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559220, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:42.502 ( 8.582s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 22.2455 (tol = 1e-08) r (rel) = 0.134278(tol = 1e-08)
2023-03-27 08:06:42.736 ( 8.815s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:43.279 ( 9.358s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.543530, 0.510000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:43.293 ( 9.372s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 2.43261 (tol = 1e-08) r (rel) = 0.0146837(tol = 1e-08)
2023-03-27 08:06:43.527 ( 9.607s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:44.071 ( 10.150s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.543265, 0.520000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:44.083 ( 10.162s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 4.43158 (tol = 1e-08) r (rel) = 0.0267498(tol = 1e-08)
2023-03-27 08:06:44.313 ( 10.393s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:44.864 ( 10.944s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.551121, 0.510000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:44.879 ( 10.958s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 0.144189 (tol = 1e-08) r (rel) = 0.000870353(tol = 1e-08)
2023-03-27 08:06:45.120 ( 11.199s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:45.679 ( 11.759s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559441, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:45.692 ( 11.771s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 0.021424 (tol = 1e-08) r (rel) = 0.000129319(tol = 1e-08)
2023-03-27 08:06:45.922 ( 12.001s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:46.475 ( 12.555s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.553090, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:46.491 ( 12.570s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 7: r (abs) = 4.80067e-06 (tol = 1e-08) r (rel) = 2.89777e-08(tol = 1e-08)
2023-03-27 08:06:46.721 ( 12.800s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:47.265 ( 13.344s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.544057, 0.520000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:06:47.278 ( 13.357s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 8: r (abs) = 2.65191e-11 (tol = 1e-08) r (rel) = 1.60074e-13(tol = 1e-08)
2023-03-27 08:06:47.278 ( 13.357s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
Time step 1, Number of iterations 8, Load [ 0. 0. -1.5]
2023-03-27 08:06:47.851 ( 13.930s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:48.414 ( 14.493s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.563313, 0.520000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:06:48.662 ( 14.741s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:49.228 ( 15.307s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.565964, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:49.241 ( 15.320s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 17.3254 (tol = 1e-08) r (rel) = 0.117842(tol = 1e-08)
2023-03-27 08:06:49.471 ( 15.550s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:50.018 ( 16.098s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.547243, 0.520000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:50.031 ( 16.110s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 5.14882 (tol = 1e-08) r (rel) = 0.0350207(tol = 1e-08)
2023-03-27 08:06:50.271 ( 16.350s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:50.827 ( 16.906s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.556267, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:50.840 ( 16.919s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 7.24003 (tol = 1e-08) r (rel) = 0.0492445(tol = 1e-08)
2023-03-27 08:06:51.071 ( 17.150s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:51.630 ( 17.710s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559406, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:51.642 ( 17.721s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 0.777889 (tol = 1e-08) r (rel) = 0.00529096(tol = 1e-08)
2023-03-27 08:06:51.878 ( 17.957s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:52.431 ( 18.511s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.553662, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:52.445 ( 18.524s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 1.25525 (tol = 1e-08) r (rel) = 0.00853785(tol = 1e-08)
2023-03-27 08:06:52.683 ( 18.762s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:53.245 ( 19.324s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.561590, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:53.258 ( 19.337s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 7: r (abs) = 0.00849512 (tol = 1e-08) r (rel) = 5.77813e-05(tol = 1e-08)
2023-03-27 08:06:53.493 ( 19.572s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:54.047 ( 20.126s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.553690, 0.510000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:54.061 ( 20.140s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 8: r (abs) = 0.000192107 (tol = 1e-08) r (rel) = 1.30665e-06(tol = 1e-08)
2023-03-27 08:06:54.294 ( 20.374s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:54.845 ( 20.925s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.551013, 0.510000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:54.859 ( 20.938s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 9: r (abs) = 1.70369e-10 (tol = 1e-08) r (rel) = 1.1588e-12(tol = 1e-08)
2023-03-27 08:06:54.859 ( 20.938s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
Time step 2, Number of iterations 9, Load [ 0. 0. -3.]
2023-03-27 08:06:55.267 ( 21.347s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:55.824 ( 21.904s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.556921, 0.530000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:56.071 ( 22.151s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:56.628 ( 22.707s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.556160, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:06:56.641 ( 22.720s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 10.0011 (tol = 1e-08) r (rel) = 0.0887471(tol = 1e-08)
2023-03-27 08:06:56.885 ( 22.964s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:57.452 ( 23.531s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.567000, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:06:57.468 ( 23.548s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 5.33026 (tol = 1e-08) r (rel) = 0.0472992(tol = 1e-08)
2023-03-27 08:06:57.701 ( 23.781s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:58.257 ( 24.336s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.555298, 0.530000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:06:58.270 ( 24.350s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 11.9901 (tol = 1e-08) r (rel) = 0.106397(tol = 1e-08)
2023-03-27 08:06:58.502 ( 24.581s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:59.057 ( 25.137s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.555062, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:59.072 ( 25.152s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 2.29702 (tol = 1e-08) r (rel) = 0.0203831(tol = 1e-08)
2023-03-27 08:06:59.323 ( 25.403s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:06:59.883 ( 25.962s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559401, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:06:59.896 ( 25.975s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 3.90234 (tol = 1e-08) r (rel) = 0.0346282(tol = 1e-08)
2023-03-27 08:07:00.132 ( 26.211s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:00.694 ( 26.773s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.561995, 0.520000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:00.708 ( 26.787s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 7: r (abs) = 0.236535 (tol = 1e-08) r (rel) = 0.00209895(tol = 1e-08)
2023-03-27 08:07:00.939 ( 27.018s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:01.488 ( 27.567s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.548627, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:01.501 ( 27.581s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 8: r (abs) = 0.0427142 (tol = 1e-08) r (rel) = 0.000379034(tol = 1e-08)
2023-03-27 08:07:01.737 ( 27.817s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:02.306 ( 28.386s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.569044, 0.540000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:02.319 ( 28.399s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 9: r (abs) = 2.87798e-05 (tol = 1e-08) r (rel) = 2.55384e-07(tol = 1e-08)
2023-03-27 08:07:02.553 ( 28.633s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:03.104 ( 29.184s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.550988, 0.510000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:03.117 ( 29.197s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 10: r (abs) = 6.08768e-10 (tol = 1e-08) r (rel) = 5.40203e-12(tol = 1e-08)
2023-03-27 08:07:03.117 ( 29.197s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 10 iterations and 10 linear solver iterations.
Time step 3, Number of iterations 10, Load [ 0. 0. -4.5]
2023-03-27 08:07:03.515 ( 29.595s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:04.066 ( 30.145s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.550650, 0.510000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:04.311 ( 30.390s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:04.866 ( 30.945s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.555463, 0.500000, 0.060000 (PETSc Krylov solver)
2023-03-27 08:07:04.879 ( 30.959s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 5.50693 (tol = 1e-08) r (rel) = 0.0653918(tol = 1e-08)
2023-03-27 08:07:05.121 ( 31.200s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:05.679 ( 31.759s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.558693, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:05.695 ( 31.774s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 26.2489 (tol = 1e-08) r (rel) = 0.311692(tol = 1e-08)
2023-03-27 08:07:05.931 ( 32.010s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:06.489 ( 32.568s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557636, 0.510000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:06.502 ( 32.581s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 2.30927 (tol = 1e-08) r (rel) = 0.0274213(tol = 1e-08)
2023-03-27 08:07:06.734 ( 32.814s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:07.299 ( 33.378s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.564852, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:07.311 ( 33.390s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 14.0562 (tol = 1e-08) r (rel) = 0.16691(tol = 1e-08)
2023-03-27 08:07:07.550 ( 33.629s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:08.104 ( 34.183s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.554201, 0.510000, 0.060000 (PETSc Krylov solver)
2023-03-27 08:07:08.118 ( 34.197s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 0.222774 (tol = 1e-08) r (rel) = 0.00264532(tol = 1e-08)
2023-03-27 08:07:08.359 ( 34.438s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:08.909 ( 34.988s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.550086, 0.520000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:08.922 ( 35.001s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 7: r (abs) = 0.286671 (tol = 1e-08) r (rel) = 0.00340406(tol = 1e-08)
2023-03-27 08:07:09.158 ( 35.238s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:09.715 ( 35.794s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.556450, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:09.729 ( 35.809s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 8: r (abs) = 0.000321869 (tol = 1e-08) r (rel) = 3.82203e-06(tol = 1e-08)
2023-03-27 08:07:09.995 ( 36.074s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:10.567 ( 36.646s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.572021, 0.540000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:10.582 ( 36.661s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 9: r (abs) = 2.63797e-07 (tol = 1e-08) r (rel) = 3.13244e-09(tol = 1e-08)
2023-03-27 08:07:10.582 ( 36.661s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 9 iterations and 9 linear solver iterations.
Time step 4, Number of iterations 9, Load [ 0. 0. -6.]
2023-03-27 08:07:11.006 ( 37.085s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:11.564 ( 37.643s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557730, 0.510000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:11.809 ( 37.889s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:12.385 ( 38.464s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.575718, 0.530000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:12.398 ( 38.478s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 3.19462 (tol = 1e-08) r (rel) = 0.0496479(tol = 1e-08)
2023-03-27 08:07:12.632 ( 38.712s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:13.191 ( 39.270s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.558528, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:13.207 ( 39.286s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 7.71429 (tol = 1e-08) r (rel) = 0.119888(tol = 1e-08)
2023-03-27 08:07:13.443 ( 39.523s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:14.019 ( 40.098s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.575735, 0.520000, 0.070000 (PETSc Krylov solver)
2023-03-27 08:07:14.033 ( 40.112s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 0.850873 (tol = 1e-08) r (rel) = 0.0132235(tol = 1e-08)
2023-03-27 08:07:14.262 ( 40.341s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:14.811 ( 40.891s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.549244, 0.500000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:14.824 ( 40.903s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 0.371434 (tol = 1e-08) r (rel) = 0.0057725(tol = 1e-08)
2023-03-27 08:07:15.060 ( 41.139s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:15.618 ( 41.697s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557890, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:15.632 ( 41.712s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 0.00215066 (tol = 1e-08) r (rel) = 3.34236e-05(tol = 1e-08)
2023-03-27 08:07:15.868 ( 41.947s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:16.428 ( 42.508s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.560060, 0.540000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:16.441 ( 42.520s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 7: r (abs) = 2.54607e-06 (tol = 1e-08) r (rel) = 3.95687e-08(tol = 1e-08)
2023-03-27 08:07:16.678 ( 42.758s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:17.234 ( 43.313s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.555422, 0.530000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:17.247 ( 43.327s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 8: r (abs) = 2.62506e-13 (tol = 1e-08) r (rel) = 4.07963e-15(tol = 1e-08)
2023-03-27 08:07:17.247 ( 43.327s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 8 iterations and 8 linear solver iterations.
Time step 5, Number of iterations 8, Load [ 0. 0. -7.5]
2023-03-27 08:07:17.651 ( 43.730s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:18.204 ( 44.283s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.552721, 0.510000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:18.464 ( 44.544s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:19.019 ( 45.099s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.554733, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:19.033 ( 45.113s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 2.00649 (tol = 1e-08) r (rel) = 0.0395622(tol = 1e-08)
2023-03-27 08:07:19.270 ( 45.350s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:19.823 ( 45.903s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.552990, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:19.840 ( 45.919s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 4.60977 (tol = 1e-08) r (rel) = 0.0908914(tol = 1e-08)
2023-03-27 08:07:20.079 ( 46.158s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:20.639 ( 46.718s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559546, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:20.653 ( 46.733s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 0.185372 (tol = 1e-08) r (rel) = 0.00365501(tol = 1e-08)
2023-03-27 08:07:20.894 ( 46.974s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:21.485 ( 47.564s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.590371, 0.560000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:21.498 ( 47.577s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 0.024688 (tol = 1e-08) r (rel) = 0.000486777(tol = 1e-08)
2023-03-27 08:07:21.734 ( 47.813s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:22.322 ( 48.401s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.587834, 0.560000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:22.335 ( 48.415s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 5.69255e-06 (tol = 1e-08) r (rel) = 1.12241e-07(tol = 1e-08)
2023-03-27 08:07:22.579 ( 48.658s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:23.134 ( 49.213s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.554443, 0.520000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:23.146 ( 49.225s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 7: r (abs) = 2.61602e-11 (tol = 1e-08) r (rel) = 5.15804e-13(tol = 1e-08)
2023-03-27 08:07:23.146 ( 49.225s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 7 iterations and 7 linear solver iterations.
Time step 6, Number of iterations 7, Load [ 0. 0. -9.]
2023-03-27 08:07:23.545 ( 49.625s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:24.105 ( 50.184s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559107, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:24.357 ( 50.437s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:24.918 ( 50.998s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.561225, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:24.932 ( 51.012s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 1.38506 (tol = 1e-08) r (rel) = 0.0336622(tol = 1e-08)
2023-03-27 08:07:25.161 ( 51.240s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:25.724 ( 51.803s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.563132, 0.510000, 0.060000 (PETSc Krylov solver)
2023-03-27 08:07:25.737 ( 51.816s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 3.03739 (tol = 1e-08) r (rel) = 0.07382(tol = 1e-08)
2023-03-27 08:07:25.966 ( 52.046s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:26.521 ( 52.600s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.554639, 0.500000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:26.535 ( 52.614s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 0.0412386 (tol = 1e-08) r (rel) = 0.00100225(tol = 1e-08)
2023-03-27 08:07:26.770 ( 52.849s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:27.327 ( 53.406s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557202, 0.520000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:27.340 ( 53.419s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 0.00205057 (tol = 1e-08) r (rel) = 4.98364e-05(tol = 1e-08)
2023-03-27 08:07:27.578 ( 53.657s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:28.134 ( 54.214s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.556106, 0.530000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:28.147 ( 54.226s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 1.78861e-08 (tol = 1e-08) r (rel) = 4.347e-10(tol = 1e-08)
2023-03-27 08:07:28.147 ( 54.226s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
Time step 7, Number of iterations 6, Load [ 0. 0. -10.5]
2023-03-27 08:07:28.534 ( 54.613s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:29.090 ( 55.169s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.555948, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:29.334 ( 55.413s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:29.893 ( 55.973s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559528, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:29.906 ( 55.985s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 1.06336 (tol = 1e-08) r (rel) = 0.031085(tol = 1e-08)
2023-03-27 08:07:30.140 ( 56.220s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:30.695 ( 56.774s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.554641, 0.510000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:30.709 ( 56.788s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 2.0477 (tol = 1e-08) r (rel) = 0.0598598(tol = 1e-08)
2023-03-27 08:07:30.939 ( 57.018s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:31.499 ( 57.578s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559744, 0.510000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:31.512 ( 57.591s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 0.00897719 (tol = 1e-08) r (rel) = 0.000262427(tol = 1e-08)
2023-03-27 08:07:31.745 ( 57.825s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:32.305 ( 58.384s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.559744, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:32.318 ( 58.397s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 0.000167422 (tol = 1e-08) r (rel) = 4.89419e-06(tol = 1e-08)
2023-03-27 08:07:32.562 ( 58.641s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:33.119 ( 59.199s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557736, 0.550000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:33.135 ( 59.214s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 3.25075e-11 (tol = 1e-08) r (rel) = 9.50282e-13(tol = 1e-08)
2023-03-27 08:07:33.135 ( 59.214s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
Time step 8, Number of iterations 6, Load [ 0. 0. -12.]
2023-03-27 08:07:33.529 ( 59.608s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:34.087 ( 60.166s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557683, 0.510000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:34.326 ( 60.405s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:34.883 ( 60.963s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557747, 0.520000, 0.040000 (PETSc Krylov solver)
2023-03-27 08:07:34.896 ( 60.976s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 2: r (abs) = 0.898789 (tol = 1e-08) r (rel) = 0.0309666(tol = 1e-08)
2023-03-27 08:07:35.138 ( 61.217s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:35.708 ( 61.787s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.569686, 0.540000, 0.030000 (PETSc Krylov solver)
2023-03-27 08:07:35.723 ( 61.802s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 3: r (abs) = 1.38354 (tol = 1e-08) r (rel) = 0.0476679(tol = 1e-08)
2023-03-27 08:07:35.958 ( 62.037s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:36.529 ( 62.608s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.571211, 0.520000, 0.050000 (PETSc Krylov solver)
2023-03-27 08:07:36.542 ( 62.622s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 4: r (abs) = 0.00185096 (tol = 1e-08) r (rel) = 6.37724e-05(tol = 1e-08)
2023-03-27 08:07:36.773 ( 62.852s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
Time step 9, Number of iterations 6, Load [ 0. 0. -13.5]
2023-03-27 08:07:37.325 ( 63.404s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.552269, 0.540000, 0.020000 (PETSc Krylov solver)
2023-03-27 08:07:37.338 ( 63.417s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 5: r (abs) = 7.87183e-06 (tol = 1e-08) r (rel) = 2.71213e-07(tol = 1e-08)
2023-03-27 08:07:37.576 ( 63.655s) [main ] petsc.cpp:675 INFO| PETSc Krylov solver starting to solve system.
2023-03-27 08:07:38.133 ( 64.213s) [main ] TimeLogger.cpp:28 INFO| Elapsed wall, usr, sys time: 0.557644, 0.500000, 0.060000 (PETSc Krylov solver)
2023-03-27 08:07:38.147 ( 64.227s) [main ] NewtonSolver.cpp:36 INFO| Newton iteration 6: r (abs) = 2.81885e-13 (tol = 1e-08) r (rel) = 9.71198e-15(tol = 1e-08)
2023-03-27 08:07:38.147 ( 64.227s) [main ] NewtonSolver.cpp:255 INFO| Newton solver finished in 6 iterations and 6 linear solver iterations.
