|
| 1 | +# In this tutorial, we will implement a Hybridizable Discontinuous Galerkin (HDG) method |
| 2 | +# for solving the Poisson equation. The HDG method is an efficient variant of DG methods |
| 3 | +# that introduces an auxiliary variable m on mesh interfaces to reduce the global system size. |
| 4 | +# |
| 5 | +# ## HDG Discretization |
| 6 | +# |
| 7 | +# We consider the Poisson equation with Dirichlet boundary conditions: |
| 8 | +# |
| 9 | +# ```math |
| 10 | +# \begin{aligned} |
| 11 | +# -\Delta u &= f \quad \text{in} \quad \Omega\\ |
| 12 | +# u &= g \quad \text{in} \quad \partial\Omega |
| 13 | +# \end{aligned} |
| 14 | +# ``` |
| 15 | +# |
| 16 | +# The HDG method first rewrites the problem as a first-order system: |
| 17 | +# |
| 18 | +# ```math |
| 19 | +# \begin{aligned} |
| 20 | +# \boldsymbol{q} + \nabla u &= \boldsymbol{0} \quad \text{in} \quad \Omega\\ |
| 21 | +# \nabla \cdot \boldsymbol{q} &= f \quad \text{in} \quad \Omega\\ |
| 22 | +# u &= g \quad \text{on} \quad \partial\Omega |
| 23 | +# \end{aligned} |
| 24 | +# ``` |
| 25 | +# |
| 26 | +# The HDG discretization introduces three variables: |
| 27 | +# - ``\boldsymbol{q}_h``: the approximation to the flux ``\boldsymbol{q}`` |
| 28 | +# - ``u_h``: the approximation to the solution ``u`` |
| 29 | +# - ``m_h``: the approximation to the trace of ``u`` on element faces |
| 30 | +# |
| 31 | +# Numerical fluxes are defindes as |
| 32 | +# |
| 33 | +# ```math |
| 34 | +# \widehat{\boldsymbol{q}}_h = \boldsymbol{q}_h + \tau(u_h - m_h)\boldsymbol{n} |
| 35 | +# ``` |
| 36 | +# |
| 37 | +# where $\tau$ is a stabilization parameter. |
| 38 | +# |
| 39 | +# First, let's load the required Gridap packages: |
| 40 | + |
| 41 | +using Gridap |
| 42 | +using Gridap.Geometry |
| 43 | +using Gridap.FESpaces |
| 44 | +using Gridap.MultiField |
| 45 | +using Gridap.CellData |
| 46 | + |
| 47 | +# ## Manufactured Solution |
| 48 | +# |
| 49 | +# We use the method of manufactured solutions to verify our implementation. |
| 50 | +# We choose a solution u and derive the corresponding source term f: |
| 51 | + |
| 52 | +u(x) = sin(2*π*x[1])*sin(2*π*x[2]) |
| 53 | +q(x) = -∇(u)(x) # Define the flux q = -∇u |
| 54 | +f(x) = (∇ ⋅ q)(x) # Source term f = -Δu = -∇⋅(∇u)$ |
| 55 | + |
| 56 | +# ## Geometry |
| 57 | +# |
| 58 | +# We generate a D-dimensional simplicial mesh from a Cartesian grid: |
| 59 | + |
| 60 | +D = 2 # Problem dimension |
| 61 | +nc = Tuple(fill(8, D)) # 4 cells in each direction |
| 62 | +domain = Tuple(repeat([0, 1], D)) # Unit cube domain |
| 63 | +model = simplexify(CartesianDiscreteModel(domain,nc)) |
| 64 | + |
| 65 | +# From this mesh, we will require two triangulations where to define our HDG spaces: |
| 66 | +# |
| 67 | +# 1. A cell triangulation $\Omega$, for the volume variables |
| 68 | +# 2. A face triangulation $\Gamma$, for the skeleton variables |
| 69 | +# |
| 70 | +# These are given by |
| 71 | + |
| 72 | +Ω = Triangulation(ReferenceFE{D}, model) # Volume triangulation |
| 73 | +Γ = Triangulation(ReferenceFE{D-1}, model) # Skeleton triangulation |
| 74 | + |
| 75 | +# ## PatchTopology and PatchTriangulation |
| 76 | +# |
| 77 | +# A key aspect of hybrid methods is the use of static condensation, which is the |
| 78 | +# elimination of cell unknowns to reduce the size of the global system. |
| 79 | +# To achieve this, we need to be able to assemble and solve local problems on each cell, that |
| 80 | +# involve |
| 81 | +# - contributions from the cell itself |
| 82 | +# - contributions from the cell faces |
| 83 | +# To this end, Gridap provides a general framework for patch-assembly and solves. The idea |
| 84 | +# is to define a patch decomposition of the mesh (in this case, a patch is a cell and its sourrounding |
| 85 | +# faces). We can then gather contributions for each patch, solve the local problems, and |
| 86 | +# assemble the results into the global system. |
| 87 | +# |
| 88 | +# The following code creates the required `PatchTopology` for the problem at hand. We then |
| 89 | +# take d-dimensional slices of it by the means of `PatchTriangulation` and `PatchBoundaryTriangulation`. |
| 90 | +# These are the `Triangulation`s we will integrate our weakform over. |
| 91 | + |
| 92 | +ptopo = Geometry.PatchTopology(model) |
| 93 | +Ωp = Geometry.PatchTriangulation(model,ptopo) # Patch volume triangulation |
| 94 | +Γp = Geometry.PatchBoundaryTriangulation(model,ptopo) # Patch skeleton triangulation |
| 95 | + |
| 96 | +# ## FESpaces |
| 97 | +# |
| 98 | +# HDG uses three different finite element spaces: |
| 99 | +# 1. A vector-valued space for the flux q (Q) |
| 100 | +# 2. A scalar space for the solution u (V) |
| 101 | +# 3. A scalar space for the interface variable m (M) |
| 102 | +# |
| 103 | +# We then define discontinuous finite element spaces of the approppriate order, locally $\mathbb{P}^k$. |
| 104 | +# Note that only the skeletal space has Dirichlet boundary conditions. |
| 105 | + |
| 106 | +order = 1 # Polynomial order |
| 107 | +reffe_Q = ReferenceFE(lagrangian, VectorValue{D, Float64}, order; space=:P) |
| 108 | +reffe_V = ReferenceFE(lagrangian, Float64, order; space=:P) |
| 109 | +reffe_M = ReferenceFE(lagrangian, Float64, order; space=:P) |
| 110 | + |
| 111 | +V = TestFESpace(Ω, reffe_V; conformity=:L2) # Discontinuous vector space |
| 112 | +Q = TestFESpace(Ω, reffe_Q; conformity=:L2) # Discontinuous scalar space |
| 113 | +M = TestFESpace(Γ, reffe_M; conformity=:L2, dirichlet_tags="boundary") # Interface space |
| 114 | +N = TrialFESpace(M, u) |
| 115 | + |
| 116 | +# ## MultiField Structure |
| 117 | +# |
| 118 | +# Since we are doing static condensation, we need assemble by blocks. In particular, the |
| 119 | +# `StaticCondensationOperator` expects the variables to be groupped in two blocks: |
| 120 | +# - The eliminated variables (in this case, the volume variables q and u) |
| 121 | +# - The retained variables (in this case, the interface variable m) |
| 122 | +# We will assemble by blocks using the `BlockMultiFieldStyle` API. |
| 123 | + |
| 124 | +mfs = BlockMultiFieldStyle(2,(2,1)) # Special blocking for efficient static condensation |
| 125 | +X = MultiFieldFESpace([V, Q, N]; style=mfs) |
| 126 | + |
| 127 | +# ## Weak Form and integration |
| 128 | +# |
| 129 | + |
| 130 | +degree = 2*(order+1) # Integration degree |
| 131 | +dΩp = Measure(Ωp,degree) # Volume measure, on the patch triangulation |
| 132 | +dΓp = Measure(Γp,degree) # Surface measure, on the patch boundary triangulation |
| 133 | + |
| 134 | +τ = 1.0 # HDG stabilization parameter |
| 135 | + |
| 136 | +n = get_normal_vector(Γp) # Face normal vector |
| 137 | +Πn(u) = u⋅n # Normal component |
| 138 | +Π(u) = change_domain(u,Γp,DomainStyle(u)) # Project to skeleton |
| 139 | + |
| 140 | +a((uh,qh,sh),(vh,wh,lh)) = ∫( qh⋅wh - uh*(∇⋅wh) - qh⋅∇(vh) )dΩp + ∫(sh*Πn(wh))dΓp + |
| 141 | + ∫((Πn(qh) + τ*(Π(uh) - sh))*(Π(vh) + lh))dΓp |
| 142 | +l((vh,wh,lh)) = ∫( f*vh )dΩp |
| 143 | + |
| 144 | +# ## Static Condensation and Solution |
| 145 | +# |
| 146 | +# With all these ingredients, we can now build our statically-condensed operator ans |
| 147 | +# solve the problem. Note that we are solving a scatically-condensed system. We can |
| 148 | +# retrieve the internal `AffineFEOperator` from `op.sc_op`. |
| 149 | + |
| 150 | +op = MultiField.StaticCondensationOperator(ptopo,X,a,l) |
| 151 | +uh, qh, sh = solve(op) |
| 152 | + |
| 153 | +dΩ = Measure(Ω,degree) |
| 154 | +eh = uh - u |
| 155 | +l2_uh = sqrt(sum(∫(eh⋅eh)*dΩ)) |
| 156 | + |
| 157 | +writevtk(Ω,"results",cellfields=["uh"=>uh,"qh"=>qh,"eh"=>eh]) |
| 158 | + |
| 159 | +# ## Going Further |
| 160 | +# |
| 161 | +# By modifying the stabilisation term, HDG can also work on polytopal meshes. An driver |
| 162 | +# solving the same problem on a polytopal mesh is available [in the Gridap repository](https://github.com/gridap/Gridap.jl/blob/75efc9a7a7e286c27e7ca3ddef5468e591845484/test/GridapTests/HDGPolytopalTests.jl). |
| 163 | +# A tutorial for HHO on polytopal meshes is also available. |
| 164 | +# |
0 commit comments