Skip to content

Commit b740dbf

Browse files
author
chmerdon
committed
add LVPP method example for obstacle problem
1 parent ce1b3a8 commit b740dbf

File tree

1 file changed

+161
-0
lines changed

1 file changed

+161
-0
lines changed
Lines changed: 161 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,161 @@
1+
#=
2+
3+
# 225 : 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+
![](example227.png)
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+
function ∇R!(result, input, qpinfo)
69+
return result[1] = χ(qpinfo.x) + exp(input[1])
70+
end
71+
72+
function main(;
73+
nrefs = 5,
74+
α0 = 1.0,
75+
order = 1,
76+
parallel = false,
77+
npart = 8,
78+
Plotter = nothing,
79+
kwargs...
80+
)
81+
82+
## choose initial mesh
83+
xgrid = uniform_refine(grid_unitsquare(Triangle2D; scale = (2, 2), shift = (-0.5, -0.5)), nrefs)
84+
if parallel
85+
xgrid = partition(xgrid, RecursiveMetisPartitioning(npart = npart))
86+
end
87+
88+
## create finite element space
89+
FES = [FESpace{H1Pk{1, 2, order}}(xgrid), FESpace{H1Pk{1, 2, order}}(xgrid)]
90+
91+
## init proximal parameter
92+
α = α0
93+
94+
## prepare Laplacian and mass matrix
95+
L = FEMatrix(FES[1], FES[1])
96+
assemble!(L, BilinearOperator([grad(1)], [grad(1)]; parallel = parallel))
97+
function scaled_laplacian!(A, b, args; assemble_matrix = true, assemble_rhs = true, kwargs...)
98+
return if assemble_matrix
99+
## add Laplacian scaled by α
100+
ExtendableFEMBase.add!(A, L.entries; factor = α)
101+
end
102+
end
103+
M = FEMatrix(FES[1], FES[2])
104+
b = FEVector(FES[1])
105+
assemble!(M, BilinearOperator([id(1)], [id(1)]; parallel = parallel))
106+
107+
## problem description
108+
PD = ProblemDescription()
109+
u = Unknown("u"; name = "solution")
110+
ψ = Unknown("ψ"; name = "latent variable")
111+
assign_unknown!(PD, u)
112+
assign_unknown!(PD, ψ)
113+
assign_operator!(PD, CallbackOperator(scaled_laplacian!, [u]; kwargs...))
114+
assign_operator!(PD, BilinearOperator([id(u)], [(id(ψ))]; transposed_copy = 1, store = true, parallel = parallel, kwargs...))
115+
assign_operator!(PD, NonlinearOperator(∇R!, [id(ψ)], [id(ψ)]; parallel = parallel, factor = -1, kwargs...))
116+
assign_operator!(PD, HomogeneousBoundaryData(u; regions = 1:4, kwargs...))
117+
assign_operator!(PD, LinearOperator(b, [u]; kwargs...))
118+
119+
## solve
120+
sol = FEVector(FES; tags = PD.unknowns)
121+
sol_prev = FEVector(FES; tags = PD.unknowns)
122+
SC = nothing
123+
r, q = 3 // 2, 3 // 2
124+
converged = false
125+
k = 0
126+
while !converged
127+
k += 1
128+
@info "Step $k: α = $(α)"
129+
130+
## save previous solution and update right-hand side
131+
b.entries .= M.entries * view(sol[ψ])
132+
sol_prev.entries .= sol.entries
133+
134+
## solve nonlinear problem
135+
sol, SC = solve(PD, FES, SC; init = sol, maxiterations = 20, verbosity = -1, timeroutputs = :hide, return_config = true, kwargs...)
136+
niterations = length(ExtendableFEM.residuals(SC))
137+
138+
## compute distance
139+
dist = norm(view(sol[u]) .- view(sol_prev[u]))
140+
@info "dist = $dist, niterations = $(niterations - 1)"
141+
if dist < 1.0e-10
142+
converged = true
143+
else
144+
# increase proximal parameter
145+
α = min(max(r^(q^k) - α), 10^3)
146+
end
147+
end
148+
149+
## postprocess latent variable v = ϕ + exp ψ = ∇R!(ψ)
150+
u2 = Unknown("̃u"; name = "solution 2")
151+
append!(sol, FES[1]; tag = u2)
152+
lazy_interpolate!(sol[u2], sol, [id(ψ)]; postprocess = ∇R!)
153+
154+
## plot
155+
plt = plot([id(u), id(ψ), id(u2)], sol; Plotter = Plotter, ncols = 3)
156+
157+
return sol, plt
158+
end
159+
160+
generateplots = ExtendableFEM.default_generateplots(Example227_ObstacleProblemLVPP, "example227.png") #hide
161+
end # module

0 commit comments

Comments
 (0)