# Hyperelasticity Solve a compressible neo-Hookean model in 3D. +++ ## UFL form file The UFL file is implemented in {download}`demo_hyperelasticity/hyperelasticity.py`. ````{admonition} UFL form implemented in python :class: dropdown ```python """Hyperelasticity variational forms.""" ``` The first step is to define the variational problem at hand. We are interested in solving for a discrete vector field in three dimensions, so first we need the appropriate finite element space and trial and test functions on this space: ```python from basix.ufl import element from ufl import ( Coefficient, Constant, FunctionSpace, Identity, Mesh, TestFunction, TrialFunction, derivative, det, diff, ds, dx, grad, inner, ln, tr, variable, ) ``` ```python # Function spaces e = element("Lagrange", "tetrahedron", 1, shape=(3,)) mesh = Mesh(e) V = FunctionSpace(mesh, e) ``` ```python # Trial and test functions du = TrialFunction(V) # Incremental displacement v = TestFunction(V) # Test function ``` Note that `element` with `shape=(3,)` creates a finite element space of vector fields. Next, we will be needing functions for the boundary source `B`, the traction `T` and the displacement solution itself `u`: ```python # Functions u = Coefficient(V) # Displacement from previous iteration B = Constant(mesh, shape=(3,)) # Body force per unit volume T = Constant(mesh, shape=(3,)) # Traction force on the boundary ``` Now, we can define the kinematic quantities involved in the model: ```python # Kinematics d = len(u) I = Identity(d) # Identity tensor # noqa: E741 F = variable(I + grad(u)) # Deformation gradient C = F.T * F # Right Cauchy-Green tensor ``` ```python # Invariants of deformation tensors Ic = tr(C) J = det(F) ``` Before defining the energy density and thus the total potential energy, it only remains to specify constants for the elasticity parameters: ```python # Elasticity parameters E = 10.0 nu = 0.3 mu = E / (2 * (1 + nu)) lmbda = E * nu / ((1 + nu) * (1 - 2 * nu)) ``` Both the first variation of the potential energy, and the Jacobian of the variation, can be automatically computed by a call to `derivative`: ```python # Stored strain energy density (compressible neo-Hookean model) psi = (mu / 2) * (Ic - 3) - mu * ln(J) + (lmbda / 2) * (ln(J)) ** 2 ``` ```python # Total potential energy Pi = psi * dx - inner(B, u) * dx - inner(T, u) * ds ``` ```python # First variation of Pi (directional derivative about u in the direction # of v) F_form = derivative(Pi, u, v) ``` ```python # Compute Jacobian of F J_form = derivative(F_form, u, du) ``` ```python # Compute Cauchy stress sigma = (1 / J) * diff(psi, F) * F.T ``` ```python forms = [F_form, J_form] elements = [e] expressions = [(sigma, [[0.25, 0.25, 0.25]])] ``` ```` +++ ## C++ program ```cpp #include "hyperelasticity.h" #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include ``` ```cpp using namespace dolfinx; using T = PetscScalar; using U = typename dolfinx::scalar_value_t; ``` ```cpp /// Hyperelastic problem class class HyperElasticProblem { public: /// Constructor HyperElasticProblem(fem::Form& L, fem::Form& J, const std::vector>& bcs) : _l(L), _j(J), _bcs(bcs.begin(), bcs.end()), _b_vec(L.function_spaces()[0]->dofmap()->index_map, L.function_spaces()[0]->dofmap()->index_map_bs()), _matJ(la::petsc::Matrix(fem::petsc::create_matrix(J, "aij"), false)), _solver(L.function_spaces()[0]->dofmap()->index_map->comm()) { auto map = L.function_spaces()[0]->dofmap()->index_map; const int bs = L.function_spaces()[0]->dofmap()->index_map_bs(); std::int32_t size_local = bs * map->size_local(); std::vector ghosts(map->ghosts().begin(), map->ghosts().end()); std::int64_t size_global = bs * map->size_global(); VecCreateGhostBlockWithArray(map->comm(), bs, size_local, size_global, ghosts.size(), ghosts.data(), _b_vec.array().data(), &_b); // Create linear solver. Default to LU. _solver.set_options_prefix("nls_solve_"); la::petsc::options::set("nls_solve_ksp_type", "preonly"); la::petsc::options::set("nls_solve_pc_type", "lu"); _solver.set_from_options(); } /// Destructor virtual ~HyperElasticProblem() { assert(_b); VecDestroy(&_b); } /// @brief Newton Solver /// @param x Solution vector /// @return Iteration count and flag indicating convergence std::pair solve(Vec x) { int iteration = 0; PetscReal residual0 = 0; auto converged = [&iteration, &residual0, this](const Vec r) -> std::pair { PetscReal residual = 0; VecNorm(r, NORM_2, &residual); // Relative residual const double relative_residual = residual / residual0; // Output iteration number and residual spdlog::info("Newton iteration {}" ": r (abs) = {} (tol = {}), r (rel) = {} (tol = {})", iteration, residual, atol, relative_residual, rtol); // Return true if convergence criterion is met bool converged = relative_residual < rtol or residual < atol; return {residual, converged}; }; assert(_b); F(x); auto [residual, newton_converged] = converged(_b); _solver.set_operators(_matJ.mat(), _matJ.mat()); Vec dx; MatCreateVecs(_matJ.mat(), &dx, nullptr); int max_it = 50; int krylov_iterations = 0; // Start iterations while (!newton_converged and iteration < max_it) { // Compute Jacobian assert(_matJ.mat()); J(x, _matJ.mat()); // Perform linear solve and update total number of Krylov iterations krylov_iterations += _solver.solve(dx, _b); // Update solution double relaxation_parameter = 1.0; VecAXPY(x, -relaxation_parameter, dx); // Increment iteration count ++iteration; // Compute F F(x); // Initialize residual0 if (iteration == 1) VecNorm(dx, NORM_2, &residual0); // Test for convergence std::tie(residual, newton_converged) = converged(_b); } if (not newton_converged) throw std::runtime_error("Newton solver did not converge."); spdlog::info("Newton solver finished in {} iterations and {} linear " "solver iterations.", iteration, krylov_iterations); VecDestroy(&dx); return {iteration, newton_converged}; } /// Compute F at current point x void F(const Vec x) { VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD); VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD); // Assemble b and update ghosts std::span b(_b_vec.array()); std::ranges::fill(b, 0); fem::assemble_vector(b, _l); VecGhostUpdateBegin(_b, ADD_VALUES, SCATTER_REVERSE); VecGhostUpdateEnd(_b, ADD_VALUES, SCATTER_REVERSE); // Set bcs fem::petsc::set_bc(_b, _bcs, x, -1); } /// Compute J = F' at current point x void J(const Vec, Mat A) { MatZeroEntries(A); fem::assemble_matrix(la::petsc::Matrix::set_block_fn(A, ADD_VALUES), _j, _bcs); MatAssemblyBegin(A, MAT_FLUSH_ASSEMBLY); MatAssemblyEnd(A, MAT_FLUSH_ASSEMBLY); fem::set_diagonal(la::petsc::Matrix::set_fn(A, INSERT_VALUES), *_j.function_spaces()[0], _bcs); MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY); MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY); } /// @brief Relative convergence tolerance. double rtol = 1e-9; /// @brief Absolute convergence tolerance. double atol = 1e-10; ``` ```cpp private: fem::Form& _l; fem::Form& _j; std::vector>> _bcs; la::Vector _b_vec; Vec _b = nullptr; // Jacobian matrix la::petsc::Matrix _matJ; // Linear solver dolfinx::la::petsc::KrylovSolver _solver; }; ``` ```cpp int main(int argc, char* argv[]) { init_logging(argc, argv); PetscInitialize(&argc, &argv, nullptr, nullptr); // Set the logging thread name to show the process rank int mpi_rank = dolfinx::MPI::rank(MPI_COMM_WORLD); std::string fmt = "[%Y-%m-%d %H:%M:%S.%e] [RANK " + std::to_string(mpi_rank) + "] [%l] %v"; spdlog::set_pattern(fmt); { // Inside the `main` function, we begin by defining a tetrahedral // mesh of the domain and the function space on this mesh. Here, we // choose to create a unit cube mesh with 25 ( = 24 + 1) vertices in // one direction and 17 ( = 16 + 1) vertices in the other two // directions. With this mesh, we initialize the (finite element) // function space defined by the generated code. // Create mesh and define function space auto mesh = std::make_shared>(mesh::create_box( MPI_COMM_WORLD, {{{0.0, 0.0, 0.0}, {1.0, 1.0, 1.0}}}, {10, 10, 10}, mesh::CellType::tetrahedron, mesh::create_cell_partitioner(mesh::GhostMode::none, 2))); auto element = basix::create_element( basix::element::family::P, basix::cell::type::tetrahedron, 1, basix::element::lagrange_variant::unset, basix::element::dpc_variant::unset, false); auto V = std::make_shared>(fem::create_functionspace( mesh, std::make_shared>( element, std::vector{3}))); auto B = std::make_shared>(std::vector{0, 0, 0}); auto traction = std::make_shared>(std::vector{0, 0, 0}); // Define solution function auto u = std::make_shared>(V); fem::Form a = fem::create_form(*form_hyperelasticity_J_form, {V, V}, {{"u", u}}, {{"B", B}, {"T", traction}}, {}, {}); fem::Form L = fem::create_form(*form_hyperelasticity_F_form, {V}, {{"u", u}}, {{"B", B}, {"T", traction}}, {}, {}); auto u_rotation = std::make_shared>(V); u_rotation->interpolate( [](auto x) -> std::pair, std::vector> { constexpr U scale = 0.005; // Center of rotation constexpr U x1_c = 0.5; constexpr U x2_c = 0.5; // Large angle of rotation (60 degrees) constexpr U theta = std::numbers::pi / 3; // New coordinates std::vector fdata(3 * x.extent(1), 0); md::mdspan> f( fdata.data(), 3, x.extent(1)); for (std::size_t p = 0; p < x.extent(1); ++p) { U x1 = x(1, p); U x2 = x(2, p); f(1, p) = scale * (x1_c + (x1 - x1_c) * std::cos(theta) - (x2 - x2_c) * std::sin(theta) - x1); f(2, p) = scale * (x2_c + (x1 - x1_c) * std::sin(theta) - (x2 - x2_c) * std::cos(theta) - x2); } return {std::move(fdata), {3, x.extent(1)}}; }); // Create Dirichlet boundary conditions auto bdofs_left = fem::locate_dofs_geometrical( *V, [](auto x) { constexpr U eps = 1.0e-6; std::vector marker(x.extent(1), false); for (std::size_t p = 0; p < x.extent(1); ++p) { if (std::abs(x(0, p)) < eps) marker[p] = true; } return marker; }); auto bdofs_right = fem::locate_dofs_geometrical( *V, [](auto x) { constexpr U eps = 1.0e-6; std::vector marker(x.extent(1), false); for (std::size_t p = 0; p < x.extent(1); ++p) { if (std::abs(x(0, p) - 1) < eps) marker[p] = true; } return marker; }); std::vector bcs = {fem::DirichletBC(std::vector{0, 0, 0}, bdofs_left, V), fem::DirichletBC(u_rotation, bdofs_right)}; HyperElasticProblem problem(L, a, bcs); problem.rtol = 10 * std::numeric_limits::epsilon(); problem.atol = 10 * std::numeric_limits::epsilon(); la::petsc::Vector _u(la::petsc::create_vector_wrap(*u->x()), false); auto [niter, success] = problem.solve(_u.vec()); std::cout << "Number of Newton iterations: " << niter << std::endl; // Compute Cauchy stress. Construct appropriate Basix element for // stress. fem::Expression sigma_expression = fem::create_expression( *expression_hyperelasticity_sigma, {{"u", u}}, {}, {}); constexpr auto family = basix::element::family::P; auto cell_type = mesh::cell_type_to_basix_type(mesh->topology()->cell_type()); constexpr int k = 0; constexpr bool discontinuous = true; basix::FiniteElement S_element = basix::create_element( family, cell_type, k, basix::element::lagrange_variant::unset, basix::element::dpc_variant::unset, discontinuous); auto S = std::make_shared>(fem::create_functionspace( mesh, std::make_shared>( S_element, std::vector{3, 3}))); fem::Function sigma(S); sigma.name = "cauchy_stress"; sigma.interpolate(sigma_expression); // Save solution in VTK format io::VTKFile file_u(mesh->comm(), "u.pvd", "w"); file_u.write({*u}, 0); // Save Cauchy stress in XDMF format io::XDMFFile file_sigma(mesh->comm(), "sigma.xdmf", "w"); file_sigma.write_mesh(*mesh); file_sigma.write_function(sigma, 0); } PetscFinalize(); return 0; } ```