{ "cells": [ { "cell_type": "markdown", "source": [ "# Stokes flow\n", "\n", "**Keywords**: *periodic boundary conditions, multiple fields, mean value constraint*" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "![](stokes-flow.png)\n", "*Figure 1*: Left: Computational domain $\\Omega$ with boundaries $\\Gamma_1$,\n", "$\\Gamma_3$ (periodic boundary conditions) and $\\Gamma_2$, $\\Gamma_4$ (homogeneous\n", "Dirichlet boundary conditions). Right: Magnitude of the resulting velocity field." ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Introduction and problem formulation\n", "\n", "This example is a translation of the [step-45 example from\n", "deal.ii](https://www.dealii.org/current/doxygen/deal.II/step_45.html) which solves Stokes\n", "flow on a quarter circle. In particular it shows how to use periodic boundary conditions,\n", "how to solve a problem with multiple unknown fields, and how to enforce a specific mean\n", "value of the solution. For the mesh generation we use\n", "[`Gmsh.jl`](https://github.com/JuliaFEM/Gmsh.jl) and then use\n", "[`FerriteGmsh.jl`](https://github.com/Ferrite-FEM/FerriteGmsh.jl) to import the mesh into\n", "Ferrite's format.\n", "\n", "The strong form of Stokes flow with velocity $\\boldsymbol{u}$ and pressure $p$ can be\n", "written as follows:\n", "$$\n", "\\begin{align*}\n", "-\\Delta \\boldsymbol{u} + \\boldsymbol{\\nabla} p &= \\bigl(\\exp(-100||\\boldsymbol{x} - (0.75, 0.1)||^2), 0\\bigr) =:\n", "\\boldsymbol{b} \\quad \\forall \\boldsymbol{x} \\in \\Omega,\\\\\n", "-\\boldsymbol{\\nabla} \\cdot \\boldsymbol{u} &= 0 \\quad \\forall \\boldsymbol{x} \\in \\Omega,\n", "\\end{align*}\n", "$$\n", "where the domain is defined as $\\Omega = \\{\\boldsymbol{x} \\in (0, 1)^2:\n", "\\ ||\\boldsymbol{x}|| \\in (0.5, 1)\\}$, see *Figure 1*. For the velocity we use periodic\n", "boundary conditions on the inlet $\\Gamma_1$ and outlet $\\Gamma_3$:\n", "$$\n", "\\begin{align*}\n", "u_x(0,\\nu) &= -u_y(\\nu, 0) \\quad & \\nu\\ \\in\\ [0.5, 1],\\\\\n", "u_y(0,\\nu) &= u_x(\\nu, 0) \\quad & \\nu\\ \\in\\ [0.5, 1],\n", "\\end{align*}\n", "$$\n", "and homogeneous Dirichlet boundary conditions for $\\Gamma_2$ and $\\Gamma_4$:\n", "$$\n", "\\boldsymbol{u} = \\boldsymbol{0} \\quad \\forall \\boldsymbol{x}\\ \\in\\\n", "\\Gamma_2 \\cup \\Gamma_4 := \\{ \\boldsymbol{x}:\\ ||\\boldsymbol{x}|| \\in \\{0.5, 1\\}\\}.\n", "$$\n", "\n", "The corresponding weak form reads as follows: Find $(\\boldsymbol{u}, p) \\in \\mathbb{U}\n", "\\times \\mathrm{L}_2$ s.t.\n", "$$\n", "\\begin{align*}\n", "\\int_\\Omega \\Bigl[[\\delta\\boldsymbol{u}\\otimes\\boldsymbol{\\nabla}]:[\\boldsymbol{u}\\otimes\\boldsymbol{\\nabla}] -\n", "(\\boldsymbol{\\nabla}\\cdot\\delta\\boldsymbol{u})\\ p\\ \\Bigr] \\mathrm{d}\\Omega &=\n", "\\int_\\Omega \\delta\\boldsymbol{u} \\cdot \\boldsymbol{b}\\ \\mathrm{d}\\Omega \\quad \\forall\n", "\\delta \\boldsymbol{u} \\in \\mathbb{U},\\\\\n", "\\int_\\Omega - (\\boldsymbol{\\nabla}\\cdot\\boldsymbol{u})\\ \\delta p\\ \\mathrm{d}\\Omega &= 0\n", "\\quad \\forall \\delta p \\in \\mathrm{L}_2,\n", "\\end{align*}\n", "$$\n", "where $\\mathbb{U}$ is a suitable function space, that, in particular, enforces the\n", "Dirichlet boundary conditions, and the periodicity constraints.\n", "This formulation is a saddle point problem, and, just like the example with\n", "Incompressible Elasticity, we need our formulation to fulfill the [LBB\n", "condition](https://en.wikipedia.org/wiki/Ladyzhenskaya%E2%80%93Babu%C5%A1ka%E2%80%93Brezzi_condition).\n", "We ensure this by using a quadratic approximation for the velocity field, and a linear\n", "approximation for the pressure.\n", "\n", "With this formulation and boundary conditions for $\\boldsymbol{u}$ the pressure will\n", "only be determined up to a constant. We will therefore add an additional constraint which\n", "fixes this constant (see [deal.ii\n", "step-11](https://www.dealii.org/current/doxygen/deal.II/step_11.html) for some more\n", "discussion around this). In particular, we will enforce the mean value of the pressure on\n", "the boundary to be 0, i.e. $\\int_{\\Gamma} p\\ \\mathrm{d}\\Gamma = 0$. One option is to\n", "enforce this using a Lagrange multipler. This would give a contribution $\\lambda\n", "\\int_{\\Gamma} \\delta p\\ \\mathrm{d}\\Gamma$ to the second equation in the weak form above,\n", "and a third equation $\\delta\\lambda \\int_{\\Gamma} p\\ \\mathrm{d}\\Gamma = 0$ so that we\n", "can solve for $\\lambda$. However, since we in this case are not interested in computing\n", "$\\lambda$, and since the constraint is linear, we can directly embed this constraint\n", "using an `AffineConstraint` in Ferrite.\n", "\n", "After FE discretization we obtain a linear system of the form\n", "$\\underline{\\underline{K}}\\ \\underline{a} = \\underline{f}$, where\n", "$$\n", "\\underline{\\underline{K}} =\n", "\\begin{bmatrix}\n", "\\underline{\\underline{K}}_{uu} & \\underline{\\underline{K}}_{pu}^\\textrm{T} \\\\\n", "\\underline{\\underline{K}}_{pu} & \\underline{\\underline{0}}\n", "\\end{bmatrix}, \\quad\n", "\\underline{a} = \\begin{bmatrix}\n", "\\underline{a}_{u} \\\\\n", "\\underline{a}_{p}\n", "\\end{bmatrix}, \\quad\n", "\\underline{f} = \\begin{bmatrix}\n", "\\underline{f}_{u} \\\\\n", "\\underline{0}\n", "\\end{bmatrix},\n", "$$\n", "and where\n", "$$\n", "\\begin{align*}\n", "(\\underline{\\underline{K}}_{uu})_{ij} &= \\int_\\Omega [\\boldsymbol{\\phi}^u_i\\otimes\\boldsymbol{\\nabla}]:[\\boldsymbol{\\phi}^u_j\\otimes\\boldsymbol{\\nabla}] \\mathrm{d}\\Omega, \\\\\n", "(\\underline{\\underline{K}}_{pu})_{ij} &= \\int_\\Omega - (\\boldsymbol{\\nabla}\\cdot\\boldsymbol{\\phi}^u_j)\\ \\phi^p_i\\ \\mathrm{d}\\Omega, \\\\\n", "(\\underline{f}_{u})_{i} &= \\int_\\Omega \\boldsymbol{\\phi}^u_i \\cdot \\boldsymbol{b}\\ \\mathrm{d}\\Omega.\n", "\\end{align*}\n", "$$\n", "\n", "The affine constraint to enforce zero mean pressure on the boundary is obtained from\n", "$\\underline{\\underline{C}}_p\\ \\underline{a}_p = \\underline{0}$, where\n", "$$\n", "(\\underline{\\underline{C}}_p)_{1j} = \\int_{\\Gamma} \\phi^p_j\\ \\mathrm{d}\\Gamma.\n", "$$\n", "\n", "!!! note\n", " The constraint matrix $\\underline{\\underline{C}}_p$ is the same matrix we would have\n", " obtained when assembling the system with the Lagrange multiplier. In that case the\n", " full system would be\n", " $$\n", " \\underline{\\underline{K}} =\n", " \\begin{bmatrix}\n", " \\underline{\\underline{K}}_{uu} & \\underline{\\underline{K}}_{pu}^\\textrm{T} &\n", " \\underline{\\underline{0}}\\\\\n", " \\underline{\\underline{K}}_{pu} & \\underline{\\underline{0}} & \\underline{\\underline{C}}_p^\\mathrm{T} \\\\\n", " \\underline{\\underline{0}} & \\underline{\\underline{C}}_p & 0 \\\\\n", " \\end{bmatrix}, \\quad\n", " \\underline{a} = \\begin{bmatrix}\n", " \\underline{a}_{u} \\\\\n", " \\underline{a}_{p} \\\\\n", " \\underline{a}_{\\lambda}\n", " \\end{bmatrix}, \\quad\n", " \\underline{f} = \\begin{bmatrix}\n", " \\underline{f}_{u} \\\\\n", " \\underline{0} \\\\\n", " \\underline{0}\n", " \\end{bmatrix}.\n", " $$" ], "metadata": {} }, { "cell_type": "markdown", "source": [ "## Commented program\n", "\n", "What follows is a program spliced with comments." ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "using Ferrite, FerriteGmsh, Gmsh, Tensors, LinearAlgebra, SparseArrays" ], "metadata": {}, "execution_count": 1 }, { "cell_type": "markdown", "source": [ "### Geometry and mesh generation with `Gmsh.jl`\n", "\n", "In the `setup_grid` function below we use the\n", "[`Gmsh.jl`](https://github.com/JuliaFEM/Gmsh.jl) package for setting up the geometry and\n", "performing the meshing. We will not discuss this part in much detail but refer to the\n", "[Gmsh API documentation](https://gmsh.info/doc/texinfo/gmsh.html#Gmsh-API) instead. The\n", "most important thing to note is the mesh periodicity constraint that is applied between\n", "the \"inlet\" and \"outlet\" parts using `gmsh.model.set_periodic`. This is necessary to later\n", "on apply a periodicity constraint for the approximated velocity field." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "setup_grid (generic function with 2 methods)" }, "metadata": {}, "execution_count": 2 } ], "cell_type": "code", "source": [ "function setup_grid(h=0.05)\n", " # Initialize gmsh\n", " Gmsh.initialize()\n", " gmsh.option.set_number(\"General.Verbosity\", 2)\n", "\n", " # Add the points\n", " o = gmsh.model.geo.add_point(0.0, 0.0, 0.0, h)\n", " p1 = gmsh.model.geo.add_point(0.5, 0.0, 0.0, h)\n", " p2 = gmsh.model.geo.add_point(1.0, 0.0, 0.0, h)\n", " p3 = gmsh.model.geo.add_point(0.0, 1.0, 0.0, h)\n", " p4 = gmsh.model.geo.add_point(0.0, 0.5, 0.0, h)\n", "\n", " # Add the lines\n", " l1 = gmsh.model.geo.add_line(p1, p2)\n", " l2 = gmsh.model.geo.add_circle_arc(p2, o, p3)\n", " l3 = gmsh.model.geo.add_line(p3, p4)\n", " l4 = gmsh.model.geo.add_circle_arc(p4, o, p1)\n", "\n", " # Create the closed curve loop and the surface\n", " loop = gmsh.model.geo.add_curve_loop([l1, l2, l3, l4])\n", " surf = gmsh.model.geo.add_plane_surface([loop])\n", "\n", " # Synchronize the model\n", " gmsh.model.geo.synchronize()\n", "\n", " # Create the physical domains\n", " gmsh.model.add_physical_group(1, [l1], -1, \"Γ1\")\n", " gmsh.model.add_physical_group(1, [l2], -1, \"Γ2\")\n", " gmsh.model.add_physical_group(1, [l3], -1, \"Γ3\")\n", " gmsh.model.add_physical_group(1, [l4], -1, \"Γ4\")\n", " gmsh.model.add_physical_group(2, [surf])\n", "\n", " # Add the periodicity constraint using 4x4 affine transformation matrix,\n", " # see https://en.wikipedia.org/wiki/Transformation_matrix#Affine_transformations\n", " transformation_matrix = zeros(4, 4)\n", " transformation_matrix[1, 2] = 1 # -sin(-pi/2)\n", " transformation_matrix[2, 1] = -1 # cos(-pi/2)\n", " transformation_matrix[3, 3] = 1\n", " transformation_matrix[4, 4] = 1\n", " transformation_matrix = vec(transformation_matrix')\n", " gmsh.model.mesh.set_periodic(1, [l1], [l3], transformation_matrix)\n", "\n", " # Generate a 2D mesh\n", " gmsh.model.mesh.generate(2)\n", "\n", " # Save the mesh, and read back in as a Ferrite Grid\n", " grid = mktempdir() do dir\n", " path = joinpath(dir, \"mesh.msh\")\n", " gmsh.write(path)\n", " togrid(path)\n", " end\n", "\n", " # Finalize the Gmsh library\n", " Gmsh.finalize()\n", "\n", " return grid\n", "end" ], "metadata": {}, "execution_count": 2 }, { "cell_type": "markdown", "source": [ "### Degrees of freedom\n", "\n", "As mentioned in the introduction we will use a quadratic approximation for the velocity\n", "field and a linear approximation for the pressure to ensure that we fulfill the LBB\n", "condition. We create the corresponding FE values with interpolations `ipu` for the\n", "velocity and `ipp` for the pressure. Note that we use linear geometric interpolation\n", "(`ipg`) for both the velocity and pressure, this is because our grid contains linear\n", "triangles. We also construct face-values for the pressure since we need to integrate along\n", "the boundary when assembling the constraint matrix $\\underline{\\underline{C}}$." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "setup_fevalues (generic function with 1 method)" }, "metadata": {}, "execution_count": 3 } ], "cell_type": "code", "source": [ "function setup_fevalues(ipu, ipp, ipg)\n", " qr = QuadratureRule{2,RefTetrahedron}(2)\n", " cvu = CellVectorValues(qr, ipu, ipg)\n", " cvp = CellScalarValues(qr, ipp, ipg)\n", " qr_face = QuadratureRule{1,RefTetrahedron}(2)\n", " fvp = FaceScalarValues(qr_face, ipp, ipg)\n", " return cvu, cvp, fvp\n", "end" ], "metadata": {}, "execution_count": 3 }, { "cell_type": "markdown", "source": [ "The `setup_dofs` function creates the `DofHandler`, and adds the two fields: a\n", "vector field `:u` with interpolation `ipu`, and a scalar field `:p` with interpolation\n", "`ipp`." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "setup_dofs (generic function with 1 method)" }, "metadata": {}, "execution_count": 4 } ], "cell_type": "code", "source": [ "function setup_dofs(grid, ipu, ipp)\n", " dh = DofHandler(grid)\n", " add!(dh, :u, 2, ipu)\n", " add!(dh, :p, 1, ipp)\n", " close!(dh)\n", " return dh\n", "end" ], "metadata": {}, "execution_count": 4 }, { "cell_type": "markdown", "source": [ "### Boundary conditions and constraints\n", "\n", "Now it is time to setup the `ConstraintHandler` and add our boundary conditions and the\n", "mean value constraint. This is perhaps the most interesting section in this example, and\n", "deserves some attention.\n", "\n", "Let's first discuss the assembly of the constraint matrix $\\underline{\\underline{C}}$\n", "and how to create an `AffineConstraint` from it. This is done in the\n", "`setup_mean_constraint` function below. Assembling this is not so different from standard\n", "assembly in Ferrite: we loop over all the faces, loop over the quadrature points, and loop\n", "over the shape functions. Note that since there is only one constraint the matrix will\n", "only have one row.\n", "After assembling `C` we construct an `AffineConstraint` from it. We select the constrained\n", "dof to be the one with the highest weight (just to avoid selecting one with 0 or a very\n", "small weight), then move the remaining to the right hand side. As an example, consider the\n", "case where the constraint equation $\\underline{\\underline{C}}_p\\ \\underline{a}_p$ is\n", "$$\n", "w_{10} p_{10} + w_{23} p_{23} + w_{154} p_{154} = 0\n", "$$\n", "i.e. dofs 10, 23, and 154, are the ones located on the boundary (all other dofs naturally\n", "gives 0 contribution). If $w_{23}$ is the largest weight, then we select $p_{23}$ to\n", "be the constrained one, and thus reorder the constraint to the form\n", "$$\n", "p_{23} = -\\frac{w_{10}}{w_{23}} p_{10} -\\frac{w_{154}}{w_{23}} p_{154} + 0,\n", "$$\n", "which is the form the `AffineConstraint` constructor expects.\n", "\n", "!!! note\n", " If all nodes along the boundary are equidistant all the weights would be the same. In\n", " this case we can construct the constraint without having to do any integration by\n", " simply finding all degrees of freedom that are located along the boundary (and using 1\n", " as the weight). This is what is done in the [deal.ii step-11\n", " example](https://www.dealii.org/current/doxygen/deal.II/step_11.html)." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "setup_mean_constraint (generic function with 1 method)" }, "metadata": {}, "execution_count": 5 } ], "cell_type": "code", "source": [ "function setup_mean_constraint(dh, fvp)\n", " assembler = start_assemble()\n", " # All external boundaries\n", " set = union(\n", " getfaceset(dh.grid, \"Γ1\"),\n", " getfaceset(dh.grid, \"Γ2\"),\n", " getfaceset(dh.grid, \"Γ3\"),\n", " getfaceset(dh.grid, \"Γ4\"),\n", " )\n", " # Allocate buffers\n", " range_p = dof_range(dh, :p)\n", " element_dofs = zeros(Int, ndofs_per_cell(dh))\n", " element_dofs_p = view(element_dofs, range_p)\n", " element_coords = zeros(Vec{2}, 3)\n", " Ce = zeros(1, length(range_p)) # Local constraint matrix (only 1 row)\n", " # Loop over all the boundaries\n", " for (ci, fi) in set\n", " Ce .= 0\n", " getcoordinates!(element_coords, dh.grid, ci)\n", " reinit!(fvp, element_coords, fi)\n", " celldofs!(element_dofs, dh, ci)\n", " for qp in 1:getnquadpoints(fvp)\n", " dΓ = getdetJdV(fvp, qp)\n", " for i in 1:getnbasefunctions(fvp)\n", " Ce[1, i] += shape_value(fvp, qp, i) * dΓ\n", " end\n", " end\n", " # Assemble to row 1\n", " assemble!(assembler, [1], element_dofs_p, Ce)\n", " end\n", " C = end_assemble(assembler)\n", " # Create an AffineConstraint from the C-matrix\n", " _, J, V = findnz(C)\n", " _, constrained_dof = findmax(abs2, V)\n", " V ./= V[constrained_dof]\n", " mean_value_constraint = AffineConstraint(\n", " constrained_dof,\n", " Pair{Int,Float64}[J[i] => -V[i] for i in 1:length(J) if J[i] != constrained_dof],\n", " 0.0,\n", " )\n", " return mean_value_constraint\n", "end" ], "metadata": {}, "execution_count": 5 }, { "cell_type": "markdown", "source": [ "We now setup all the boundary conditions in the `setup_constraints` function below.\n", "Since the periodicity constraint for this example is between two boundaries which are not\n", "parallel to each other we need to i) compute the mapping between each mirror face and the\n", "corresponding image face (on the element level) and ii) describe the dof relation between\n", "dofs on these two faces. In Ferrite this is done by defining a transformation of entities\n", "on the image boundary such that they line up with the matching entities on the mirror\n", "boundary. In this example we consider the inlet $\\Gamma_1$ to be the image, and the\n", "outlet $\\Gamma_3$ to be the mirror. The necessary transformation to apply then becomes a\n", "rotation of $\\pi/2$ radians around the out-of-plane axis. We set up the rotation matrix\n", "`R`, and then compute the mapping between mirror and image faces using\n", "`collect_periodic_faces` where the rotation is applied to the coordinates. In the\n", "next step we construct the constraint using the `PeriodicDirichlet` constructor.\n", "We pass the constructor the computed mapping, and also the rotation matrix. This matrix is\n", "used to rotate the dofs on the mirror surface such that we properly constrain\n", "$\\boldsymbol{u}_x$-dofs on the mirror to $-\\boldsymbol{u}_y$-dofs on the image, and\n", "$\\boldsymbol{u}_y$-dofs on the mirror to $\\boldsymbol{u}_x$-dofs on the image.\n", "\n", "For the remaining part of the boundary we add a homogeneous Dirichlet boundary condition\n", "on both components of the velocity field. This is done using the `Dirichlet`\n", "constructor, which we have discussed in other tutorials." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "setup_constraints (generic function with 1 method)" }, "metadata": {}, "execution_count": 6 } ], "cell_type": "code", "source": [ "function setup_constraints(dh, fvp)\n", " ch = ConstraintHandler(dh)\n", " # Periodic BC\n", " R = rotation_tensor(π / 2)\n", " periodic_faces = collect_periodic_faces(dh.grid, \"Γ3\", \"Γ1\", x -> R ⋅ x)\n", " periodic = PeriodicDirichlet(:u, periodic_faces, R, [1, 2])\n", " add!(ch, periodic)\n", " # Dirichlet BC\n", " Γ24 = union(getfaceset(dh.grid, \"Γ2\"), getfaceset(dh.grid, \"Γ4\"))\n", " dbc = Dirichlet(:u, Γ24, (x, t) -> [0, 0], [1, 2])\n", " add!(ch, dbc)\n", " # Compute mean value constraint and add it\n", " mean_value_constraint = setup_mean_constraint(dh, fvp)\n", " add!(ch, mean_value_constraint)\n", " # Finalize\n", " close!(ch)\n", " update!(ch, 0)\n", " return ch\n", "end" ], "metadata": {}, "execution_count": 6 }, { "cell_type": "markdown", "source": [ "### Global and local assembly\n", "\n", "Assembly of the global system is also something that we have seen in many previous\n", "tutorials. One interesting thing to note here is that, since we have two unknown fields,\n", "we use the `dof_range` function to make sure we assemble the element contributions\n", "to the correct block of the local stiffness matrix `ke`." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "assemble_system! (generic function with 1 method)" }, "metadata": {}, "execution_count": 7 } ], "cell_type": "code", "source": [ "function assemble_system!(K, f, dh, cvu, cvp)\n", " assembler = start_assemble(K, f)\n", " ke = zeros(ndofs_per_cell(dh), ndofs_per_cell(dh))\n", " fe = zeros(ndofs_per_cell(dh))\n", " range_u = dof_range(dh, :u)\n", " ndofs_u = length(range_u)\n", " range_p = dof_range(dh, :p)\n", " ndofs_p = length(range_p)\n", " ϕᵤ = Vector{Vec{2,Float64}}(undef, ndofs_u)\n", " ∇ϕᵤ = Vector{Tensor{2,2,Float64,4}}(undef, ndofs_u)\n", " divϕᵤ = Vector{Float64}(undef, ndofs_u)\n", " ϕₚ = Vector{Float64}(undef, ndofs_p)\n", " for cell in CellIterator(dh)\n", " reinit!(cvu, cell)\n", " reinit!(cvp, cell)\n", " ke .= 0\n", " fe .= 0\n", " for qp in 1:getnquadpoints(cvu)\n", " dΩ = getdetJdV(cvu, qp)\n", " for i in 1:ndofs_u\n", " ϕᵤ[i] = shape_value(cvu, qp, i)\n", " ∇ϕᵤ[i] = shape_gradient(cvu, qp, i)\n", " divϕᵤ[i] = shape_divergence(cvu, qp, i)\n", " end\n", " for i in 1:ndofs_p\n", " ϕₚ[i] = shape_value(cvp, qp, i)\n", " end\n", " # u-u\n", " for (i, I) in pairs(range_u), (j, J) in pairs(range_u)\n", " ke[I, J] += ( ∇ϕᵤ[i] ⊡ ∇ϕᵤ[j] ) * dΩ\n", " end\n", " # u-p\n", " for (i, I) in pairs(range_u), (j, J) in pairs(range_p)\n", " ke[I, J] += ( -divϕᵤ[i] * ϕₚ[j] ) * dΩ\n", " end\n", " # p-u\n", " for (i, I) in pairs(range_p), (j, J) in pairs(range_u)\n", " ke[I, J] += ( -divϕᵤ[j] * ϕₚ[i] ) * dΩ\n", " end\n", " # rhs\n", " for (i, I) in pairs(range_u)\n", " x = spatial_coordinate(cvu, qp, getcoordinates(cell))\n", " b = exp(-100 * norm(x - Vec{2}((0.75, 0.1)))^2)\n", " bv = Vec{2}((b, 0.0))\n", " fe[I] += (ϕᵤ[i] ⋅ bv) * dΩ\n", " end\n", " end\n", " assemble!(assembler, celldofs(cell), ke, fe)\n", " end\n", " return K, f\n", "end" ], "metadata": {}, "execution_count": 7 }, { "cell_type": "markdown", "source": [ "### Running the simulation\n", "\n", "We now have all the puzzle pieces, and just need to define the main function, which puts\n", "them all together." ], "metadata": {} }, { "outputs": [ { "output_type": "execute_result", "data": { "text/plain": "main (generic function with 1 method)" }, "metadata": {}, "execution_count": 8 } ], "cell_type": "code", "source": [ "function main()\n", " # Grid\n", " h = 0.05 # approximate element size\n", " grid = setup_grid(h)\n", " # Interpolations\n", " ipu = Lagrange{2,RefTetrahedron,2}() # quadratic\n", " ipp = Lagrange{2,RefTetrahedron,1}() # linear\n", " # Dofs\n", " dh = setup_dofs(grid, ipu, ipp)\n", " # FE values\n", " ipg = Lagrange{2,RefTetrahedron,1}() # linear geometric interpolation\n", " cvu, cvp, fvp = setup_fevalues(ipu, ipp, ipg)\n", " # Boundary conditions\n", " ch = setup_constraints(dh, fvp)\n", " # Global tangent matrix and rhs\n", " coupling = [true true; true false] # no coupling between pressure test/trial functions\n", " K = create_sparsity_pattern(dh, ch; coupling=coupling)\n", " f = zeros(ndofs(dh))\n", " # Assemble system\n", " assemble_system!(K, f, dh, cvu, cvp)\n", " # Apply boundary conditions and solve\n", " apply!(K, f, ch)\n", " u = K \\ f\n", " apply!(u, ch)\n", " # Export the solution\n", " vtk_grid(\"stokes-flow\", grid) do vtk\n", " vtk_point_data(vtk, dh, u)\n", " end\n", " return\n", "end" ], "metadata": {}, "execution_count": 8 }, { "cell_type": "markdown", "source": [ "Run it!" ], "metadata": {} }, { "outputs": [], "cell_type": "code", "source": [ "main()" ], "metadata": {}, "execution_count": 9 }, { "cell_type": "markdown", "source": [ "The resulting magnitude of the velocity field is visualized in *Figure 1*." ], "metadata": {} }, { "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.8.5" }, "kernelspec": { "name": "julia-1.8", "display_name": "Julia 1.8.5", "language": "julia" } }, "nbformat": 4 }