{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Solid Mechanics Example" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [ "# dependencies\n", "using LFAToolkit\n", "using LinearAlgebra\n", "using Pkg\n", "Pkg.activate(\"./\")\n", "Pkg.instantiate()\n", "using Plots" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "code_folding": [ 20, 49 ] }, "outputs": [], "source": [ "# setup\n", "mesh = Mesh3D(1.0, 1.0, 1.0)\n", "finep = 3\n", "coarsep = 2\n", "numbercomponents = 1\n", "dimension = 3\n", "finebasis = TensorH1LagrangeBasis(finep, finep, numbercomponents, dimension)\n", "coarsebasis = TensorH1LagrangeBasis(coarsep, finep, numbercomponents, dimension)\n", "ctofbasis = TensorH1LagrangeBasis(coarsep, finep, numbercomponents, dimension, lagrangequadrature = true)\n", "\n", "# constants\n", "e = 1E6 # Young's modulus\n", "ν = 0.3 # Poisson's ratio\n", "K = e/(3*(1 - 2*ν)) # bulk modulus\n", "λ = e*ν/((1 + ν)*(1 - 2*ν)) # Lamé parameters\n", "μ = e/(2*(1 + ν))\n", "\n", "# state\n", "gradu = [1; 2; 3]*ones(1, 3);\n", "\n", "function neohookeanweakform(deltadu::Array{Float64}, w::Array{Float64})\n", " # dP = dF S + F dS\n", "\n", " # deformation gradient\n", " F = gradu + I\n", " J = det(F)\n", " # Green-Lagrange strain tensor\n", " E = (gradu*gradu' + gradu'*gradu)/2\n", " # right Cauchy-Green tensor\n", " C = 2*E + I\n", " C_inv = C^-1\n", " # second Piola-Kirchhoff\n", " S = λ*log(J)*C_inv + 2*μ*C_inv*E\n", "\n", " # delta du\n", " deltadu = deltadu'\n", " # dF\n", " dF = deltadu + I\n", " # deltaE\n", " deltaE = (deltadu*deltadu' + deltadu'*deltadu)/2\n", " # dS\n", " dS = λ*sum(C_inv.*deltaE)*C_inv + 2*(μ - λ*log(J))*C_inv*deltaE*C_inv\n", " # dP\n", " dP = dF*S + F*dS\n", "\n", " return [dP']\n", "end\n", "\n", "# linearized Neo-Hookean operators\n", "function makeoperator(basis::TensorBasis)\n", " inputs = [\n", " OperatorField(basis, [EvaluationMode.gradient], \"gradent of deformation\"),\n", " OperatorField(basis, [EvaluationMode.quadratureweights], \"quadrature weights\"),\n", " ]\n", " outputs = [\n", " OperatorField(\n", " basis,\n", " [EvaluationMode.gradient],\n", " \"test function gradient of deformation\",\n", " ),\n", " ]\n", " return Operator(neohookeanweakform, mesh, inputs, outputs)\n", "end\n", "fineoperator = makeoperator(finebasis)\n", "coarseoperator = makeoperator(coarsebasis)\n", "\n", "# Chebyshev smoother\n", "chebyshev = Chebyshev(fineoperator)\n", "\n", "# p-multigrid preconditioner\n", "multigrid =\n", " PMultigrid(fineoperator, coarseoperator, chebyshev, [ctofbasis, ctofbasis, ctofbasis])" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "code_folding": [ 8, 41 ] }, "outputs": [], "source": [ "# full operator symbols\n", "numberruns = 250\n", "maxeigenvalue = 0\n", "θ_min = -π/2\n", "θ_max = 3π/2\n", "\n", "# compute and plot smoothing factor\n", "# -- 1D --\n", "if dimension == 1\n", " # setup\n", " ω = [0.636]\n", " v = [1, 1]\n", " maxeigenvalues = zeros(numberruns)\n", "\n", " # compute\n", " for i in 1:numberruns\n", " θ = [θ_min + (θ_max - θ_min)*i/numberruns]\n", " if abs(θ[1] % 2π) > π/128\n", " A = computesymbols(multigrid, ω, v, θ)\n", " eigenvalues = [abs(val) for val in eigvals(A)]\n", " maxeigenvalues[i] = max(eigenvalues...)\n", " maxeigenvalue = max(maxeigenvalue, maxeigenvalues[i])\n", " end\n", " end\n", "\n", " # plot\n", " println(\"max eigenvalue: \", maxeigenvalue)\n", " xrange = θ_min/π:(θ_max - θ_min)/π/(numberruns-1):θ_max/π\n", " plot(\n", " xrange,\n", " xlabel=\"θ/π\",\n", " xtickfont=font(12, \"Courier\"),\n", " maxeigenvalues,\n", " ytickfont=font(12, \"Courier\"),\n", " ylabel=\"spectral radius\",\n", " linewidth=3,\n", " legend=:none,\n", " title=\"P-Multigrid Error Symbol\"\n", " )\n", " ylims!(0.0, max(maxeigenvalues...) * 1.1)\n", "# -- 2D --\n", "elseif dimension == 2\n", " # setup\n", " ω = [0.636]\n", " v = [1, 1]\n", " maxeigenvalues = zeros(numberruns, numberruns)\n", "\n", " # compute\n", " for i in 1:numberruns, j in 1:numberruns\n", " θ = [\n", " θ_min + (θ_max - θ_min)*i/numberruns,\n", " θ_min + (θ_max - θ_min)*j/numberruns\n", " ]\n", " if sqrt(abs(θ[1] % 2π)^2 + abs(θ[2] % 2π)^2) > π/128\n", " A = computesymbols(multigrid, ω, v, θ)\n", " eigenvalues = [abs(val) for val in eigvals(A)]\n", " maxeigenvalues[i, j] = max(eigenvalues...)\n", " maxeigenvalue = max(maxeigenvalue, maxeigenvalues[i, j])\n", " end\n", " end\n", "\n", " # plot\n", " println(\"max eigenvalue: \", maxeigenvalue)\n", " xrange = θ_min/π:(θ_max - θ_min)/π/(numberruns-1):θ_max/π\n", " heatmap(\n", " xrange,\n", " xlabel=\"θ/π\",\n", " xtickfont=font(12, \"Courier\"),\n", " xrange,\n", " ylabel=\"θ/π\",\n", " ytickfont=font(12, \"Courier\"),\n", " maxeigenvalues,\n", " title=\"P-Multigrid Error Symbol\",\n", " transpose=true,\n", " aspect_ratio=:equal\n", " )\n", " xlims!(θ_min/π, θ_max/π)\n", " ylims!(θ_min/π, θ_max/π)\n", "end" ] } ], "metadata": { "kernelspec": { "display_name": "Julia 1.5.2", "language": "julia", "name": "julia-1.5" }, "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.5.2" } }, "nbformat": 4, "nbformat_minor": 4 }