{ "cells": [ { "cell_type": "markdown", "source": [ "# Hyperelasticity\n", "\n", "**Keywords**: *hyperelasticity*, *finite strain*, *large deformations*, *Newton's method*,\n", "*conjugate gradient*, *automatic differentiation*\n", "\n", "\n", "\n", "*Figure 1*: Cube loaded in torsion modeled with a hyperelastic material model and\n", "finite strain." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Introduction\n", "\n", "In this example we will solve a problem in a finite strain setting using an\n", "hyperelastic material model. In order to compute the stress we will use automatic\n", "differentiation, to solve the non-linear system we use Newton's\n", "method, and for solving the Newton increment we use conjugate gradients.\n", "\n", "The weak form is expressed in terms of the first Piola-Kirchoff stress $\\mathbf{P}$\n", "as follows: Find $\\mathbf{u} \\in \\mathbb{U}$ such that\n", "\n", "$$\n", "\\int_{\\Omega} [\\nabla_{\\mathbf{X}} \\delta \\mathbf{u}] : \\mathbf{P}(\\mathbf{u})\\ \\mathrm{d}\\Omega =\n", "\\int_{\\Omega} \\delta \\mathbf{u} \\cdot \\mathbf{b}\\ \\mathrm{d}\\Omega + \\int_{\\Gamma_\\mathrm{N}}\n", "\\delta \\mathbf{u} \\cdot \\mathbf{t}\\ \\mathrm{d}\\Gamma\n", "\\quad \\forall \\delta \\mathbf{u} \\in \\mathbb{U}^0,\n", "$$\n", "\n", "where $\\mathbf{u}$ is the unknown displacement field, $\\mathbf{b}$ is the body force acting\n", "on the reference domain, $\\mathbf{t}$ is the traction acting on the Neumann part of the reference\n", "domain's boundary, and where $\\mathbb{U}$ and $\\mathbb{U}^0$ are suitable trial and test sets.\n", "$\\Omega$ denotes the reference (sometimes also called *initial* or *material*) domain.\n", "Gradients are defined with respect to the reference domain, here denoted with an $\\mathbf{X}$.\n", "Formally this is expressed as $(\\nabla_{\\mathbf{X}} \\bullet)_{ij} := \\frac{\\partial(\\bullet)_i}{\\partial X_j}$.\n", "Note that for large deformation problems it is also possible that gradients and integrals\n", "are defined on the deformed (sometimes also called *current* or *spatial*) domain, depending\n", "on the specific formulation.\n", "\n", "The specific problem we will solve in this example is the cube from Figure 1: On one side\n", "we apply a rotation using Dirichlet boundary conditions, on the opposite side we fix the\n", "displacement with a homogeneous Dirichlet boundary condition, and on the remaining four\n", "sides we apply a traction in the normal direction of the surface. In addition, a body\n", "force is applied in one direction.\n", "\n", "In addition to Ferrite.jl and Tensors.jl, this examples uses\n", "[TimerOutputs.jl](https://github.com/KristofferC/TimerOutputs.jl) for timing the program\n", "and print a summary at the end,\n", "[ProgressMeter.jl](https://github.com/timholy/ProgressMeter.jl) for showing a simple\n", "progress bar, and\n", "[IterativeSolvers.jl](https://github.com/JuliaLinearAlgebra/IterativeSolvers.jl) for solving\n", "the linear system using conjugate gradients." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Ferrite, Tensors, TimerOutputs, ProgressMeter, IterativeSolvers" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "## Hyperelastic material model\n", "\n", "The stress can be derived from an energy potential, defined in\n", "terms of the right Cauchy-Green tensor $\\mathbf{C} = \\mathbf{F}^{\\mathrm{T}} \\cdot \\mathbf{F}$,\n", "where $\\mathbf{F} = \\mathbf{I} + \\nabla_{\\mathbf{X}} \\mathbf{u}$ is the deformation gradient.\n", "We shall use the compressible neo-Hookean model from [Wikipedia](https://en.wikipedia.org/wiki/Neo-Hookean_solid) with the potential\n", "\n", "$$\n", "\\Psi(\\mathbf{C}) = \\underbrace{\\frac{\\mu}{2} (I_1 - 3)}_{W(\\mathbf{C})} \\underbrace{- {\\mu} \\ln(J) + \\frac{\\lambda}{2} (J - 1)^2}_{U(J)},\n", "$$\n", "\n", "where $I_1 = \\mathrm{tr}(\\mathbf{C})$ is the first invariant, $J = \\sqrt{\\det(\\mathbf{C})}$\n", "and $\\mu$ and $\\lambda$ material parameters." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "> **Extra details on compressible neo-Hookean formulations**\n", ">\n", "> The Neo-Hooke model is only a well defined terminology in the incompressible case.\n", "> Thus, only $W(\\mathbf{C})$ specifies the neo-Hookean behavior, the volume penalty $U(J)$ can vary in different formulations.\n", "> In order to obtain a well-posed problem, it is crucial to choose a convex formulation of $U(J)$.\n", "> Other examples for $U(J)$ can be found, e.g. in [Hol:2000:nsm; Eq. (6.138)](@cite)\n", "> $$\n", "> \\beta^{-2} (\\beta \\ln J + J^{-\\beta} -1)\n", "> $$\n", "> where [SimMie:1992:act; Eq. (2.37)](@cite) published a non-generalized version with $\\beta=-2$.\n", "> This shows the possible variety of $U(J)$ while all of them refer to compressible neo-Hookean models.\n", "> Sometimes the modified first invariant $\\overline{I}_1=\\frac{I_1}{I_3^{1/3}}$ is used in $W(\\mathbf{C})$ instead of $I_1$." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "From the potential we obtain the second Piola-Kirchoff stress $\\mathbf{S}$ as\n", "\n", "$$\n", "\\mathbf{S} = 2 \\frac{\\partial \\Psi}{\\partial \\mathbf{C}},\n", "$$\n", "\n", "and the tangent of $\\mathbf{S}$ as\n", "\n", "$$\n", "\\frac{\\partial \\mathbf{S}}{\\partial \\mathbf{C}} = 2 \\frac{\\partial^2 \\Psi}{\\partial \\mathbf{C}^2}.\n", "$$\n", "\n", "Finally, for the finite element problem we need $\\mathbf{P}$ and\n", "$\\frac{\\partial \\mathbf{P}}{\\partial \\mathbf{F}}$, which can be\n", "obtained by using the following relations:\n", "\n", "$$\n", "\\begin{align*}\n", "\\mathbf{P} &= \\mathbf{F} \\cdot \\mathbf{S},\\\\\n", "\\frac{\\partial \\mathbf{P}}{\\partial \\mathbf{F}} &= \\mathbf{I} \\bar{\\otimes} \\mathbf{S} + 2\\, \\mathbf{F} \\cdot\n", "\\frac{\\partial \\mathbf{S}}{\\partial \\mathbf{C}} : \\mathbf{F}^\\mathrm{T} \\bar{\\otimes} \\mathbf{I}.\n", "\\end{align*}\n", "$$" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "> **Derivation of $\\partial \\mathbf{P} / \\partial \\mathbf{F}$**\n", ">\n", "> *Tip:* See [knutam.github.io/tensors](https://knutam.github.io/tensors/Theory/IndexNotation/) for\n", "> an explanation of the index notation used in this derivation.\n", "> Using the product rule, the chain rule, and the relations $\\mathbf{P} = \\mathbf{F} \\cdot\n", "> \\mathbf{S}$ and $\\mathbf{C} = \\mathbf{F}^\\mathrm{T} \\cdot \\mathbf{F}$, we obtain the\n", "> following:\n", "> $$\n", "> \\begin{aligned}\n", "> \\frac{\\partial P_{ij}}{\\partial F_{kl}} &=\n", "> \\frac{\\partial (F_{im}S_{mj})}{\\partial F_{kl}} \\\\ &=\n", "> \\frac{\\partial F_{im}}{\\partial F_{kl}}S_{mj} +\n", "> F_{im}\\frac{\\partial S_{mj}}{\\partial F_{kl}} \\\\ &=\n", "> \\delta_{ik}\\delta_{ml} S_{mj} +\n", "> F_{im}\\frac{\\partial S_{mj}}{\\partial C_{no}}\\frac{\\partial C_{no}}{\\partial F_{kl}} \\\\ &=\n", "> \\delta_{ik}S_{lj} +\n", "> F_{im}\\frac{\\partial S_{mj}}{\\partial C_{no}}\n", "> \\frac{\\partial (F^\\mathrm{T}_{np}F_{po})}{\\partial F_{kl}} \\\\ &=\n", "> \\delta_{ik}S^\\mathrm{T}_{jl} +\n", "> F_{im}\\frac{\\partial S_{mj}}{\\partial C_{no}}\n", "> \\left(\n", "> \\frac{\\partial F^\\mathrm{T}_{np}}{\\partial F_{kl}}F_{po} +\n", "> F^\\mathrm{T}_{np}\\frac{\\partial F_{po}}{\\partial F_{kl}}\n", "> \\right) \\\\ &=\n", "> \\delta_{ik}S_{jl} +\n", "> F_{im}\\frac{\\partial S_{mj}}{\\partial C_{no}}\n", "> (\\delta_{nl} \\delta_{pk} F_{po} + F^\\mathrm{T}_{np}\\delta_{pk} \\delta_{ol}) \\\\ &=\n", "> \\delta_{ik}S_{lj} +\n", "> F_{im}\\frac{\\partial S_{mj}}{\\partial C_{no}}\n", "> (F^\\mathrm{T}_{ok} \\delta_{nl} + F^\\mathrm{T}_{nk} \\delta_{ol}) \\\\ &=\n", "> \\delta_{ik}S_{jl} +\n", "> 2\\, F_{im} \\frac{\\partial S_{mj}}{\\partial C_{no}}\n", "> F^\\mathrm{T}_{nk} \\delta_{ol} \\\\\n", "> \\frac{\\partial \\mathbf{P}}{\\partial \\mathbf{F}} &=\n", "> \\mathbf{I}\\bar{\\otimes}\\mathbf{S} +\n", "> 2\\, \\mathbf{F} \\cdot \\frac{\\partial \\mathbf{S}}{\\partial \\mathbf{C}}\n", "> : \\mathbf{F}^\\mathrm{T} \\bar{\\otimes} \\mathbf{I},\n", "> \\end{aligned}\n", "> $$\n", "> where we used the fact that $\\mathbf{S}$ is symmetric ($S_{lj} = S_{jl}$) and that\n", "> $\\frac{\\partial \\mathbf{S}}{\\partial \\mathbf{C}}$ is *minor* symmetric ($\\frac{\\partial\n", "> S_{mj}}{\\partial C_{no}} = \\frac{\\partial S_{mj}}{\\partial C_{on}}$)." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "### Implementation of material model using automatic differentiation\n", "We can implement the material model as follows, where we utilize automatic differentiation\n", "for the stress and the tangent, and thus only define the potential:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "struct NeoHooke\n", " μ::Float64\n", " λ::Float64\n", "end\n", "\n", "function Ψ(C, mp::NeoHooke)\n", " μ = mp.μ\n", " λ = mp.λ\n", " Ic = tr(C)\n", " J = sqrt(det(C))\n", " return μ / 2 * (Ic - 3 - 2 * log(J)) + λ / 2 * (J - 1)^2\n", "end\n", "\n", "function constitutive_driver(C, mp::NeoHooke)\n", " # Compute all derivatives in one function call\n", " ∂²Ψ∂C², ∂Ψ∂C = Tensors.hessian(y -> Ψ(y, mp), C, :all)\n", " S = 2.0 * ∂Ψ∂C\n", " ∂S∂C = 2.0 * ∂²Ψ∂C²\n", " return S, ∂S∂C\n", "end;" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "## Newton's method\n", "\n", "As mentioned above, to deal with the non-linear weak form we first linearize\n", "the problem such that we can apply Newton's method, and then apply the FEM to\n", "discretize the problem. Skipping a detailed derivation, Newton's method can\n", "be expressed as:\n", "Given some initial guess for the degrees of freedom $\\underline{u}^0$, find a sequence\n", "$\\underline{u}^{k}$ by iterating\n", "\n", "$$\n", "\\underline{u}^{k+1} = \\underline{u}^{k} - \\Delta \\underline{u}^{k}\n", "$$\n", "\n", "until some termination condition has been met. Therein we determine $\\Delta \\underline{u}^{k}$\n", "from the linearized problem\n", "\n", "$$\n", "\\underline{\\underline{K}}(\\underline{u}^{k}) \\Delta \\underline{u}^{k} = \\underline{g}(\\underline{u}^{k})\n", "$$\n", "\n", "where the global residual, $\\underline{g}$, and the Jacobi matrix,\n", "$\\underline{\\underline{K}} = \\frac{\\partial \\underline{g}}{\\partial \\underline{u}}$, are\n", "evaluated at the current guess $\\underline{u}^k$. The entries of $\\underline{g}$ are given\n", "by\n", "\n", "$$\n", "(\\underline{g})_{i} = \\int_{\\Omega} [\\nabla_{\\mathbf{X}} \\delta \\mathbf{u}_{i}] :\n", "\\mathbf{P} \\, \\mathrm{d} \\Omega - \\int_{\\Omega} \\delta \\mathbf{u}_{i} \\cdot \\mathbf{b} \\,\n", "\\mathrm{d} \\Omega - \\int_{\\Gamma_\\mathrm{N}} \\delta \\mathbf{u}_i \\cdot \\mathbf{t}\\\n", "\\mathrm{d}\\Gamma,\n", "$$\n", "\n", "and the entries of $\\underline{\\underline{K}}$ are given by\n", "\n", "$$\n", "(\\underline{\\underline{K}})_{ij} = \\int_{\\Omega} [\\nabla_{\\mathbf{X}} \\delta\n", "\\mathbf{u}_{i}] : \\frac{\\partial \\mathbf{P}}{\\partial \\mathbf{F}} : [\\nabla_{\\mathbf{X}}\n", "\\delta \\mathbf{u}_{j}] \\, \\mathrm{d} \\Omega.\n", "$$\n", "\n", "\n", "A detailed derivation can be found in every continuum mechanics book, which has a\n", "chapter about finite elasticity theory. We used \"Nonlinear solid mechanics: a continuum\n", "approach for engineering science.\" by [Hol:2000:nsm; Chapter 8](@citet) as a reference.\n", "\n", "## Finite element assembly\n", "\n", "The element routine for assembling the residual and tangent stiffness is implemented\n", "as usual, with loops over quadrature points and shape functions:" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "function assemble_element!(ke, ge, cell, cv, fv, mp, ue, ΓN)\n", " # Reinitialize cell values, and reset output arrays\n", " reinit!(cv, cell)\n", " fill!(ke, 0.0)\n", " fill!(ge, 0.0)\n", "\n", " b = Vec{3}((0.0, -0.5, 0.0)) # Body force\n", " tn = 0.1 # Traction (to be scaled with surface normal)\n", " ndofs = getnbasefunctions(cv)\n", "\n", " for qp in 1:getnquadpoints(cv)\n", " dΩ = getdetJdV(cv, qp)\n", " # Compute deformation gradient F and right Cauchy-Green tensor C\n", " ∇u = function_gradient(cv, qp, ue)\n", " F = one(∇u) + ∇u\n", " C = tdot(F) # F' ⋅ F\n", " # Compute stress and tangent\n", " S, ∂S∂C = constitutive_driver(C, mp)\n", " P = F ⋅ S\n", " I = one(S)\n", " ∂P∂F = otimesu(I, S) + 2 * F ⋅ ∂S∂C ⊡ otimesu(F', I)\n", "\n", " # Loop over test functions\n", " for i in 1:ndofs\n", " # Test function and gradient\n", " δui = shape_value(cv, qp, i)\n", " ∇δui = shape_gradient(cv, qp, i)\n", " # Add contribution to the residual from this test function\n", " ge[i] += (∇δui ⊡ P - δui ⋅ b) * dΩ\n", "\n", " ∇δui∂P∂F = ∇δui ⊡ ∂P∂F # Hoisted computation\n", " for j in 1:ndofs\n", " ∇δuj = shape_gradient(cv, qp, j)\n", " # Add contribution to the tangent\n", " ke[i, j] += (∇δui∂P∂F ⊡ ∇δuj) * dΩ\n", " end\n", " end\n", " end\n", "\n", " # Surface integral for the traction\n", " for facet in 1:nfacets(cell)\n", " if (cellid(cell), facet) in ΓN\n", " reinit!(fv, cell, facet)\n", " for q_point in 1:getnquadpoints(fv)\n", " t = tn * getnormal(fv, q_point)\n", " dΓ = getdetJdV(fv, q_point)\n", " for i in 1:ndofs\n", " δui = shape_value(fv, q_point, i)\n", " ge[i] -= (δui ⋅ t) * dΓ\n", " end\n", " end\n", " end\n", " end\n", " return\n", "end;" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "Assembling global residual and tangent is also done in the usual way, just looping over\n", "the elements, call the element routine and assemble in the the global matrix K and\n", "residual g." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "function assemble_global!(K, g, dh, cv, fv, mp, u, ΓN)\n", " n = ndofs_per_cell(dh)\n", " ke = zeros(n, n)\n", " ge = zeros(n)\n", "\n", " # start_assemble resets K and g\n", " assembler = start_assemble(K, g)\n", "\n", " # Loop over all cells in the grid\n", " @timeit \"assemble\" for cell in CellIterator(dh)\n", " global_dofs = celldofs(cell)\n", " ue = u[global_dofs] # element dofs\n", " @timeit \"element assemble\" assemble_element!(ke, ge, cell, cv, fv, mp, ue, ΓN)\n", " assemble!(assembler, global_dofs, ke, ge)\n", " end\n", " return\n", "end;" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "Finally, we define a main function which sets up everything and then performs Newton\n", "iterations until convergence." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "solve (generic function with 1 method)" }, "metadata": {}, "execution_count": 5 } ], "cell_type": "code", "source": [ "function solve()\n", " reset_timer!()\n", "\n", " # Generate a grid\n", " N = 10\n", " L = 1.0\n", " left = zero(Vec{3})\n", " right = L * ones(Vec{3})\n", " grid = generate_grid(Tetrahedron, (N, N, N), left, right)\n", "\n", " # Material parameters\n", " E = 10.0\n", " ν = 0.3\n", " μ = E / (2(1 + ν))\n", " λ = (E * ν) / ((1 + ν) * (1 - 2ν))\n", " mp = NeoHooke(μ, λ)\n", "\n", " # Finite element base\n", " ip = Lagrange{RefTetrahedron, 1}()^3\n", " qr = QuadratureRule{RefTetrahedron}(1)\n", " qr_facet = FacetQuadratureRule{RefTetrahedron}(1)\n", " cv = CellValues(qr, ip)\n", " fv = FacetValues(qr_facet, ip)\n", "\n", " # DofHandler\n", " dh = DofHandler(grid)\n", " add!(dh, :u, ip) # Add a displacement field\n", " close!(dh)\n", "\n", " function rotation(X, t)\n", " θ = pi / 3 # 60°\n", " x, y, z = X\n", " return t * Vec{3}(\n", " (\n", " 0.0,\n", " L / 2 - y + (y - L / 2) * cos(θ) - (z - L / 2) * sin(θ),\n", " L / 2 - z + (y - L / 2) * sin(θ) + (z - L / 2) * cos(θ),\n", " )\n", " )\n", " end\n", "\n", " dbcs = ConstraintHandler(dh)\n", " # Add a homogeneous boundary condition on the \"clamped\" edge\n", " dbc = Dirichlet(:u, getfacetset(grid, \"right\"), (x, t) -> [0.0, 0.0, 0.0], [1, 2, 3])\n", " add!(dbcs, dbc)\n", " dbc = Dirichlet(:u, getfacetset(grid, \"left\"), (x, t) -> rotation(x, t), [1, 2, 3])\n", " add!(dbcs, dbc)\n", " close!(dbcs)\n", " t = 0.5\n", " Ferrite.update!(dbcs, t)\n", "\n", " # Neumann part of the boundary\n", " ΓN = union(\n", " getfacetset(grid, \"top\"),\n", " getfacetset(grid, \"bottom\"),\n", " getfacetset(grid, \"front\"),\n", " getfacetset(grid, \"back\"),\n", " )\n", "\n", " # Pre-allocation of vectors for the solution and Newton increments\n", " _ndofs = ndofs(dh)\n", " un = zeros(_ndofs) # previous solution vector\n", " u = zeros(_ndofs)\n", " Δu = zeros(_ndofs)\n", " ΔΔu = zeros(_ndofs)\n", " apply!(un, dbcs)\n", "\n", " # Create sparse matrix and residual vector\n", " K = allocate_matrix(dh)\n", " g = zeros(_ndofs)\n", "\n", " # Perform Newton iterations\n", " newton_itr = -1\n", " NEWTON_TOL = 1.0e-8\n", " NEWTON_MAXITER = 30\n", " prog = ProgressMeter.ProgressThresh(NEWTON_TOL; desc = \"Solving:\")\n", "\n", " while true\n", " newton_itr += 1\n", " # Construct the current guess\n", " u .= un .+ Δu\n", " # Compute residual and tangent for current guess\n", " assemble_global!(K, g, dh, cv, fv, mp, u, ΓN)\n", " # Apply boundary conditions\n", " apply_zero!(K, g, dbcs)\n", " # Compute the residual norm and compare with tolerance\n", " normg = norm(g)\n", " ProgressMeter.update!(prog, normg; showvalues = [(:iter, newton_itr)])\n", " if normg < NEWTON_TOL\n", " break\n", " elseif newton_itr > NEWTON_MAXITER\n", " error(\"Reached maximum Newton iterations, aborting\")\n", " end\n", "\n", " # Compute increment using conjugate gradients\n", " @timeit \"linear solve\" IterativeSolvers.cg!(ΔΔu, K, g; maxiter = 1000)\n", "\n", " apply_zero!(ΔΔu, dbcs)\n", " Δu .-= ΔΔu\n", " end\n", "\n", " # Save the solution\n", " @timeit \"export\" begin\n", " VTKGridFile(\"hyperelasticity\", dh) do vtk\n", " write_solution(vtk, dh, u)\n", " end\n", " end\n", "\n", " print_timer(title = \"Analysis with $(getncells(grid)) elements\", linechars = :ascii)\n", " return u\n", "end" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "Run the simulation" ], "metadata": {} }, { "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "\rSolving: (thresh = 1e-08, value = 0.000259699)\u001b[K\r\n", " iter: 3\u001b[K\r\u001b[A\n", "\r\u001b[K\u001b[A\rSolving: Time: 0:00:00 (6 iterations)\u001b[K\r\n", " iter: 5\u001b[K\n", "-------------------------------------------------------------------------------\n", " Analysis with 6000 elements Time Allocations \n", " ----------------------- ------------------------\n", " Tot / % measured: 2.12s / 52.2% 147MiB / 44.5% \n", "\n", "Section ncalls time %tot avg alloc %tot avg\n", "-------------------------------------------------------------------------------\n", "export 1 923ms 83.3% 923ms 59.3MiB 90.9% 59.3MiB\n", "assemble 6 110ms 9.9% 18.3ms 5.50MiB 8.4% 938KiB\n", " element assemble 36.0k 64.5ms 5.8% 1.79μs 0.00B 0.0% 0.00B\n", "linear solve 5 74.6ms 6.7% 14.9ms 473KiB 0.7% 94.6KiB\n", "-------------------------------------------------------------------------------\n" ] } ], "cell_type": "code", "source": [ "u = solve();" ], "metadata": {}, "execution_count": 6 }, { "cell_type": "markdown", "source": [ "---\n", "\n", "*This notebook was generated using [Literate.jl](https://github.com/fredrikekre/Literate.jl).*" ], "metadata": {} } ], "nbformat_minor": 3, "metadata": { "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.11.3" }, "kernelspec": { "name": "julia-1.11", "display_name": "Julia 1.11.3", "language": "julia" } }, "nbformat": 4 }