|
| 1 | +#= |
| 2 | +
|
| 3 | +# 227 : Obstacle Problem LVPP |
| 4 | +([source code](@__SOURCE_URL__)) |
| 5 | +
|
| 6 | +This example computes the solution ``u`` of the nonlinear obstacle problem that seeks |
| 7 | +the minimiser of the energy functional |
| 8 | +```math |
| 9 | +\begin{aligned} |
| 10 | + E(u) = \frac{1}{2} \int_\Omega \lvert \nabla u \rvert^2 dx - \int_\Omega f u dx |
| 11 | +\end{aligned} |
| 12 | +``` |
| 13 | +with some right-hand side ``f`` within the set of admissible functions that lie above an obstacle ``\chi`` |
| 14 | +```math |
| 15 | +\begin{aligned} |
| 16 | + \mathcal{K} := \lbrace u \in H^1_0(\Omega) : u \geq \chi \rbrace. |
| 17 | +\end{aligned} |
| 18 | +``` |
| 19 | +
|
| 20 | +Opposite to Example225 the solution is computed by the latent variable proximal point (LVPP) method |
| 21 | +that solves the problem via a series of nonlinear mixed problems that guarantee the decay of the |
| 22 | +energy. Given ``\alpha_k`` and initial guess ``u_0`` and ``\psi_0``, the subproblem for ``k \geq 1`` |
| 23 | +seeks a solution ``u_{k} \in V := H^1_0(\Omega)`` and ``\psi_{k} \in W := L^\infty(\Omega)`` such that |
| 24 | +```math |
| 25 | +\begin{aligned} |
| 26 | +\alpha_k (\nabla u_k, \nabla v) + (\psi_k, v) = (\alpha_k f + \psi_{k-1},v) |
| 27 | +&& \text{for all } v \in V\\ |
| 28 | +(u_k, w) - (\varchi + \exp(\psi_k), w) & = 0 && |
| 29 | +\text{for all } w \in W |
| 30 | +\end{aligned} |
| 31 | +``` |
| 32 | +The parameter ``\alpha_k`` is initialized with ``\alpha_0 = 1`` and updated according to |
| 33 | +``\alpha_k = min(max(r^(q^k) - α), 10^3)`` with ``r = q = 1.5``. The problem for each ``k`` |
| 34 | +is solved by the Newton method. This implements Algorithm 3 in the reference below. |
| 35 | +
|
| 36 | +
|
| 37 | +!!! reference |
| 38 | +
|
| 39 | + ''Proximal Galerkin: A Structure-Preserving Finite Element Method for Pointwise Bound Constraints'' |
| 40 | + Brendan Keith, Thomas M. Surowiec, Found Comput Math (2024) |
| 41 | + [>Link<](https://doi.org/10.1007/s10208-024-09681-8) |
| 42 | +
|
| 43 | +
|
| 44 | + |
| 45 | +=# |
| 46 | + |
| 47 | +module Example227_ObstacleProblemLVPP |
| 48 | + |
| 49 | +using ExtendableFEM |
| 50 | +using ExtendableFEMBase |
| 51 | +using ExtendableGrids |
| 52 | +using LinearAlgebra |
| 53 | +using Metis |
| 54 | +using Test #hide |
| 55 | + |
| 56 | +## define obstacle |
| 57 | +const b = 9 // 20 |
| 58 | +const d = sqrt(1 // 4 - b^2) |
| 59 | +function χ(x) |
| 60 | + r = sqrt(x[1]^2 + x[2]^2) |
| 61 | + if r <= b |
| 62 | + return sqrt(1 // 4 - r^2) |
| 63 | + else |
| 64 | + return d + b^2 / d - b * r / d |
| 65 | + end |
| 66 | +end |
| 67 | + |
| 68 | +## transformation of latent variable ψ to constrained variable u |
| 69 | +function ∇R!(result, input, qpinfo) |
| 70 | + return result[1] = χ(qpinfo.x) + exp(input[1]) |
| 71 | +end |
| 72 | + |
| 73 | +## boundary data for latent variable ψ (such that ̃u := χ + ∇R(ψ) satisfies Dirichlet boundary conditions) |
| 74 | +function bnd_ψ!(result, qpinfo) |
| 75 | + return result[1] = log(-χ(qpinfo.x)) |
| 76 | +end |
| 77 | + |
| 78 | +function main(; |
| 79 | + nrefs = 5, |
| 80 | + α0 = 1.0, |
| 81 | + order = 1, |
| 82 | + parallel = false, |
| 83 | + npart = 8, |
| 84 | + tol = 1.0e-12, |
| 85 | + Plotter = nothing, |
| 86 | + kwargs... |
| 87 | + ) |
| 88 | + |
| 89 | + ## choose initial mesh |
| 90 | + xgrid = uniform_refine(grid_unitsquare(Triangle2D; scale = (2, 2), shift = (-0.5, -0.5)), nrefs) |
| 91 | + if parallel |
| 92 | + xgrid = partition(xgrid, RecursiveMetisPartitioning(npart = npart)) |
| 93 | + end |
| 94 | + |
| 95 | + ## create finite element space |
| 96 | + FES = [FESpace{H1Pk{1, 2, order}}(xgrid), FESpace{H1Pk{1, 2, order}}(xgrid)] |
| 97 | + |
| 98 | + ## init proximal parameter |
| 99 | + α = α0 |
| 100 | + |
| 101 | + ## prepare Laplacian and mass matrix |
| 102 | + L = FEMatrix(FES[1], FES[1]) |
| 103 | + assemble!(L, BilinearOperator([grad(1)], [grad(1)]; parallel = parallel)) |
| 104 | + function scaled_laplacian!(A, b, args; assemble_matrix = true, assemble_rhs = true, kwargs...) |
| 105 | + return if assemble_matrix |
| 106 | + ## add Laplacian scaled by α |
| 107 | + ExtendableFEMBase.add!(A, L.entries; factor = α) |
| 108 | + end |
| 109 | + end |
| 110 | + M = FEMatrix(FES[1], FES[2]) |
| 111 | + b = FEVector(FES[1]) |
| 112 | + assemble!(M, BilinearOperator([id(1)], [id(1)]; parallel = parallel)) |
| 113 | + |
| 114 | + ## problem description |
| 115 | + PD = ProblemDescription() |
| 116 | + u = Unknown("u"; name = "solution") |
| 117 | + ψ = Unknown("ψ"; name = "latent variable") |
| 118 | + assign_unknown!(PD, u) |
| 119 | + assign_unknown!(PD, ψ) |
| 120 | + assign_operator!(PD, CallbackOperator(scaled_laplacian!, [u]; kwargs...)) |
| 121 | + assign_operator!(PD, BilinearOperator([id(u)], [(id(ψ))]; transposed_copy = 1, store = true, parallel = parallel, kwargs...)) |
| 122 | + assign_operator!(PD, NonlinearOperator(∇R!, [id(ψ)], [id(ψ)]; parallel = parallel, factor = -1, kwargs...)) |
| 123 | + assign_operator!(PD, HomogeneousBoundaryData(u; regions = 1:4, kwargs...)) |
| 124 | + assign_operator!(PD, InterpolateBoundaryData(ψ, bnd_ψ!; regions = 1:4, kwargs...)) |
| 125 | + assign_operator!(PD, LinearOperator(b, [u]; kwargs...)) |
| 126 | + |
| 127 | + ## solve |
| 128 | + sol = FEVector(FES; tags = PD.unknowns) |
| 129 | + sol_prev = FEVector(FES; tags = PD.unknowns) |
| 130 | + SC = nothing |
| 131 | + r, q = 3 // 2, 3 // 2 |
| 132 | + converged = false |
| 133 | + k = 0 |
| 134 | + while !converged |
| 135 | + k += 1 |
| 136 | + @info "Step $k: α = $(α)" |
| 137 | + |
| 138 | + ## save previous solution and update right-hand side |
| 139 | + b.entries .= M.entries * view(sol[ψ]) |
| 140 | + sol_prev.entries .= sol.entries |
| 141 | + |
| 142 | + ## solve nonlinear problem |
| 143 | + sol, SC = solve(PD, FES, SC; init = sol, maxiterations = 20, verbosity = -1, timeroutputs = :hide, return_config = true, kwargs...) |
| 144 | + niterations = length(ExtendableFEM.residuals(SC)) |
| 145 | + |
| 146 | + ## compute distance |
| 147 | + dist = norm(view(sol[u]) .- view(sol_prev[u])) |
| 148 | + @info "dist = $dist, niterations = $(niterations - 1)" |
| 149 | + if dist < tol |
| 150 | + converged = true |
| 151 | + else |
| 152 | + # increase proximal parameter |
| 153 | + α = min(max(r^(q^k) - α), 10^3) |
| 154 | + end |
| 155 | + end |
| 156 | + |
| 157 | + ## postprocess latent variable v = ϕ + exp ψ = ∇R!(ψ) |
| 158 | + u2 = Unknown("̃u"; name = "solution 2") |
| 159 | + append!(sol, FES[1]; tag = u2) |
| 160 | + lazy_interpolate!(sol[u2], sol, [id(ψ)]; postprocess = ∇R!) |
| 161 | + |
| 162 | + ## plot |
| 163 | + plt = plot([id(u), id(ψ), id(u2)], sol; Plotter = Plotter, ncols = 3) |
| 164 | + |
| 165 | + return sol, plt |
| 166 | +end |
| 167 | + |
| 168 | +generateplots = ExtendableFEM.default_generateplots(Example227_ObstacleProblemLVPP, "example227.png") #hide |
| 169 | +end # module |
0 commit comments