Skip to content

Commit d7096d1

Browse files
authored
Merge pull request #275 from mfem/mfem_48_dev_ex40
Mfem 48 dev ex40
2 parents 127ec4f + fa28724 commit d7096d1

File tree

3 files changed

+723
-1
lines changed

3 files changed

+723
-1
lines changed

examples/ex39.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -249,7 +249,7 @@ def run(order=1,
249249
if __name__ == "__main__":
250250
from mfem.common.arg_parser import ArgParser
251251

252-
parser = ArgParser(description='Ex1 (Laplace Problem)')
252+
parser = ArgParser(description='Ex39 (Named attribute)')
253253
parser.add_argument('-m', '--mesh',
254254
default='compass.msh',
255255
action='store', type=str,

examples/ex40.py

Lines changed: 327 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,327 @@
1+
'''
2+
MFEM example 40 (converted from ex40.cpp)
3+
4+
See c++ version in the MFEM library for more detail
5+
6+
Sample runs: python ex40.py -step 10.0 -gr 2.0
7+
python ex40.py -step 10.0 -gr 2.0 -o 3 -r 1
8+
python ex40.py -step 10.0 -gr 2.0 -r 4 -m ../data/l-shape.mesh
9+
python ex40.py -step 10.0 -gr 2.0 -r 2 -m ../data/fichera.mesh
10+
11+
Description: This example code demonstrates how to use MFEM to solve the
12+
eikonal equation,
13+
14+
|∇𝑢| = 1 in Ω, 𝑢 = 0 on ∂Ω.
15+
16+
The viscosity solution of this problem coincides with the unique optimum
17+
of the nonlinear program
18+
19+
maximize ∫_Ω 𝑢 d𝑥 subject to |∇𝑢| ≤ 1 in Ω, 𝑢 = 0 on ∂Ω, (⋆)
20+
21+
which is the foundation for method implemented below.
22+
23+
Following the proximal Galerkin methodology [1,2] (see also Example
24+
36), we construct a Legendre function for the closed unit ball
25+
𝐵₁ := {𝑥 ∈ Rⁿ | |𝑥| ≤ 1}. Our choice is the Hellinger entropy,
26+
27+
R(𝑥) = −( 1 − |𝑥|² )^{1/2},
28+
29+
although other choices are possible, each leading to a slightly
30+
different algorithm. We then adaptively regularize the optimization
31+
problem (⋆) with the Bregman divergence of the Hellinger entropy,
32+
33+
maximize ∫_Ω 𝑢 d𝑥 - αₖ⁻¹ D(∇𝑢,∇𝑢ₖ₋₁) subject to 𝑢 = 0 on Ω.
34+
35+
This results in a sequence of functions ( 𝜓ₖ , 𝑢ₖ ),
36+
37+
𝑢ₖ → 𝑢, 𝜓ₖ/|𝜓ₖ| → ∇𝑢 as k → ∞,
38+
39+
defined by the nonlinear saddle-point problems
40+
41+
Find 𝜓ₖ ∈ H(div,Ω) and 𝑢ₖ ∈ L²(Ω) such that
42+
( (∇R)⁻¹(𝜓ₖ) , τ ) + ( 𝑢ₖ , ∇⋅τ ) = 0 ∀ τ ∈ H(div,Ω)
43+
( ∇⋅𝜓ₖ , v ) = ( ∇⋅𝜓ₖ₋₁ - αₖ , v ) ∀ v ∈ L²(Ω)
44+
45+
where (∇R)⁻¹(𝜓) = 𝜓 / ( 1 + |𝜓|² )^{1/2} and αₖ = α₀rᵏ, where r ≥ 1
46+
is a prescribed growth rate. (r = 1 is the most stable.) The
47+
saddle-point problems are solved using a damped quasi-Newton method
48+
with a tunable regularization parameter 0 ≤ ϵ << 1.
49+
50+
[1] Keith, B. and Surowiec, T. (2024) Proximal Galerkin: A structure-
51+
preserving finite element method for pointwise bound constraints.
52+
Foundations of Computational Mathematics, 1–97.
53+
[2] Dokken, J., Farrell, P., Keith, B., Papadopoulos, I., and
54+
Surowiec, T. (2025) The latent variable proximal point algorithm
55+
for variational problems with inequality constraints. (To appear.)
56+
57+
'''
58+
import os
59+
from os.path import expanduser, join
60+
import numpy as np
61+
62+
import mfem.ser as mfem
63+
64+
if hasattr(mfem, "UMFPackSolver"):
65+
use_umfpack = True if not args.no_use_umfpack else False
66+
else:
67+
use_umfpack = False
68+
69+
70+
def run(meshfile="",
71+
order=1,
72+
max_it=5,
73+
ref_levels=3,
74+
alpha=1.0,
75+
growth_rate=1.0,
76+
newton_scaling=0.8,
77+
eps=1e-6,
78+
tol=1e-4,
79+
visualization=True):
80+
81+
# 2. Read the mesh from the mesh file.
82+
mesh = mfem.Mesh(meshfile, 1, 1)
83+
dim = mesh.Dimension()
84+
sdim = mesh.SpaceDimension()
85+
86+
# 3. Postprocess the mesh.
87+
# 3A. Refine the mesh to increase the resolution.
88+
for i in range(ref_levels):
89+
mesh.UniformRefinement()
90+
91+
# 3B. Interpolate the geometry after refinement to control geometry error.
92+
# NOTE: Minimum second-order interpolation is used to improve the accuracy.
93+
curvature_order = max(order, 2)
94+
mesh.SetCurvature(curvature_order)
95+
96+
# 4. Define the necessary finite element spaces on the mesh.
97+
RTfec = mfem.RT_FECollection(order, dim)
98+
RTfes = mfem.FiniteElementSpace(mesh, RTfec)
99+
100+
L2fec = mfem.L2_FECollection(order, dim)
101+
L2fes = mfem.FiniteElementSpace(mesh, L2fec)
102+
103+
print("Number of H(div) dofs: " + str(RTfes.GetTrueVSize()))
104+
print("Number of L² dofs: " + str(L2fes.GetTrueVSize()))
105+
106+
# 5. Define the offsets for the block matrices
107+
offsets = mfem.intArray([0, RTfes.GetVSize(), L2fes.GetVSize()])
108+
offsets.PartialSum()
109+
110+
x = mfem.BlockVector(offsets)
111+
rhs = mfem.BlockVector(offsets)
112+
x.Assign(0.0)
113+
rhs.Assign(0.0)
114+
115+
# 6. Define the solution vectors as a finite element grid functions
116+
# corresponding to the fespaces.
117+
118+
u_gf = mfem.GridFunction()
119+
delta_psi_gf = mfem.GridFunction()
120+
delta_psi_gf.MakeRef(RTfes, x, offsets[0])
121+
u_gf.MakeRef(L2fes, x, offsets[1])
122+
123+
psi_old_gf = mfem.GridFunction(RTfes)
124+
psi_gf = mfem.GridFunction(RTfes)
125+
u_old_gf = mfem.GridFunction(L2fes)
126+
127+
# 7. Define initial guesses for the solution variables.
128+
delta_psi_gf.Assign(0.0)
129+
psi_gf.Assign(0.0)
130+
u_gf.Assign(0.0)
131+
psi_old_gf.Assign(0.0)
132+
u_old_gf.Assign(0.0)
133+
134+
# 8. Prepare for glvis output.
135+
# 14. Send the solution by socket to a GLVis server.
136+
if visualization:
137+
sol_sock = mfem.socketstream("localhost", 19916)
138+
sol_sock.precision(8)
139+
140+
# 9. Coefficients to be used later.
141+
neg_alpha_cf = mfem.ConstantCoefficient(-1.0*alpha)
142+
zero_cf = mfem.ConstantCoefficient(0.0)
143+
144+
psigf_cf = mfem.VectorGridFunctionCoefficient(psi_gf)
145+
146+
@mfem.jit.vector(shape=(sdim,), dependency=(psigf_cf,))
147+
def Z(_ptx, psi_vals):
148+
norm = np.linalg.norm(psi_vals)
149+
phi = 1.0 / np.sqrt(1.0 + norm*norm)
150+
vvec = psi_vals*phi
151+
return vvec
152+
153+
@mfem.jit.matrix(shape=(sdim, sdim), dependency=(psigf_cf,))
154+
def DZ(_ptx, psi_vals):
155+
norm = np.linalg.norm(psi_vals)
156+
phi = 1.0 / np.sqrt(1.0 + norm*norm)
157+
158+
kmat = np.zeros((sdim, sdim))
159+
for i in range(sdim):
160+
kmat[i, i] = phi + eps
161+
for j in range(sdim):
162+
kmat[i, j] -= psi_vals[i]*psi_vals[j] * phi**3
163+
return kmat
164+
165+
neg_Z = mfem.ScalarVectorProductCoefficient(-1.0, Z)
166+
div_psi_cf = mfem.DivergenceGridFunctionCoefficient(psi_gf)
167+
div_psi_old_cf = mfem.DivergenceGridFunctionCoefficient(psi_old_gf)
168+
psi_old_minus_psi = mfem.SumCoefficient(
169+
div_psi_old_cf, div_psi_cf, 1.0, -1.0)
170+
171+
# 10. Assemble constant matrices/vectors to avoid reassembly in the loop.
172+
b0 = mfem.LinearForm()
173+
b1 = mfem.LinearForm()
174+
b0.MakeRef(RTfes, rhs.GetBlock(0), 0)
175+
b1.MakeRef(L2fes, rhs.GetBlock(1), 0)
176+
177+
b0.AddDomainIntegrator(mfem.VectorFEDomainLFIntegrator(neg_Z))
178+
b1.AddDomainIntegrator(mfem.DomainLFIntegrator(neg_alpha_cf))
179+
b1.AddDomainIntegrator(mfem.DomainLFIntegrator(psi_old_minus_psi))
180+
181+
a00 = mfem.BilinearForm(RTfes)
182+
a00.AddDomainIntegrator(mfem.VectorFEMassIntegrator(DZ))
183+
184+
a10 = mfem.MixedBilinearForm(RTfes, L2fes)
185+
a10.AddDomainIntegrator(mfem.VectorFEDivergenceIntegrator())
186+
a10.Assemble()
187+
a10.Finalize()
188+
A10 = a10.SpMat()
189+
A01 = mfem.Transpose(A10)
190+
191+
# 11. Iterate.
192+
total_iterations = 0
193+
increment_u = 0.1
194+
u_tmp = mfem.GridFunction(L2fes)
195+
196+
for k in range(max_it):
197+
u_tmp.Assign(u_old_gf)
198+
199+
print("\nOUTER ITERATION " + str(k+1))
200+
201+
for j in range(5):
202+
total_iterations += 1
203+
204+
b0.Assemble()
205+
b1.Assemble()
206+
207+
a00.Assemble(0)
208+
a00.Finalize(0)
209+
A00 = a00.SpMat()
210+
211+
# Construct Schur-complement preconditioner
212+
A00_diag = mfem.Vector(a00.Height())
213+
A00.GetDiag(A00_diag)
214+
A00_diag.Reciprocal()
215+
S = mfem.Mult_AtDA(A01, A00_diag)
216+
217+
# Python note:
218+
# owns_blocks should be 0 in Python
219+
# because wrapper class will delete blocks
220+
prec = mfem.BlockDiagonalPreconditioner(offsets)
221+
prec.owns_blocks = 0
222+
223+
prec.SetDiagonalBlock(0, mfem.DSmoother(A00))
224+
225+
if not use_umfpack:
226+
prec.SetDiagonalBlock(1, mfem.GSSmoother(S))
227+
else:
228+
prec.SetDiagonalBlock(1, mfem.UMFPackSolver(S))
229+
230+
A = mfem.BlockOperator(offsets)
231+
A.SetBlock(0, 0, A00)
232+
A.SetBlock(1, 0, A10)
233+
A.SetBlock(0, 1, A01)
234+
235+
mfem.MINRES(A, prec, rhs, x, 0, 2000, 1e-12)
236+
237+
del S
238+
239+
u_tmp -= u_gf # u_tmp = u_tmp - u_gf
240+
241+
Newton_update_size = u_tmp.ComputeL2Error(zero_cf)
242+
u_tmp.Assign(u_gf)
243+
244+
# Damped Newton update
245+
psi_gf.Add(newton_scaling, delta_psi_gf)
246+
a00.Update()
247+
248+
if visualization:
249+
sol_sock << "solution\n" << mesh << u_gf
250+
sol_sock << "window_title 'Discrete solution'"
251+
sol_sock.flush()
252+
253+
print("Newton_update_size = " + "{:g}".format(Newton_update_size))
254+
255+
if Newton_update_size < increment_u:
256+
break
257+
258+
u_tmp.Assign(u_gf)
259+
u_tmp -= u_old_gf
260+
increment_u = u_tmp.ComputeL2Error(zero_cf)
261+
262+
print("Number of Newton iterations = " + str(j+1))
263+
print("Increment (|| uₕ - uₕ_prvs||) = " + "{:g}".format(increment_u))
264+
265+
u_old_gf.Assign(u_gf)
266+
psi_old_gf.Assign(psi_gf)
267+
268+
if increment_u < tol or k == max_it-1:
269+
break
270+
271+
alpha *= max(growth_rate, 1.0)
272+
neg_alpha_cf.constant = -alpha
273+
274+
print("\n Outer iterations: " + str(k+1) +
275+
"\n Total iterations: " + str(total_iterations) +
276+
"\n Total dofs: " + str(RTfes.GetTrueVSize() + L2fes.GetTrueVSize()))
277+
278+
279+
if __name__ == "__main__":
280+
from mfem.common.arg_parser import ArgParser
281+
282+
parser = ArgParser(description='Ex40 (Eikonal queation)')
283+
parser.add_argument('-m', '--mesh',
284+
default='star.mesh',
285+
action='store', type=str,
286+
help='Mesh file to use.')
287+
parser.add_argument('-o', '--order',
288+
action='store', default=1, type=int,
289+
help="Finite element order (polynomial degree).")
290+
parser.add_argument('-r', '--refs',
291+
action='store', default=3, type=int,
292+
help="Number of h-refinements.")
293+
parser.add_argument('-mi', '--max-it',
294+
action='store', default=5, type=int,
295+
help="Maximum number of iterations")
296+
parser.add_argument("-tol", "--tol",
297+
action='store', default=1e-4, type=float,
298+
help="Stopping criteria based on the difference between.")
299+
parser.add_argument('-step', '--step',
300+
action='store', default=1.0, type=float,
301+
help="Initial size alpha")
302+
parser.add_argument("-gr", "--growth-rate",
303+
action='store', default=1.0, type=float,
304+
help="Growth rate of the step size alpha")
305+
parser.add_argument('-no-vis', '--no-visualization',
306+
action='store_true',
307+
default=False,
308+
help='Disable or disable GLVis visualization')
309+
310+
args = parser.parse_args()
311+
parser.print_options(args)
312+
313+
meshfile = expanduser(
314+
join(os.path.dirname(__file__), '..', 'data', args.mesh))
315+
visualization = not args.no_visualization
316+
in_alpha = args.step
317+
318+
run(meshfile=meshfile,
319+
order=args.order,
320+
max_it=args.max_it,
321+
ref_levels=args.refs,
322+
alpha=in_alpha,
323+
growth_rate=args.growth_rate,
324+
newton_scaling=0.8,
325+
eps=1e-6,
326+
tol=args.tol,
327+
visualization=visualization)

0 commit comments

Comments
 (0)