Deflation techniques for computing multiple solutions of nonlinear problems¶
Nonlinear problems can have multiple solutions. Deflation is an approach for computing multiple solutions of nonlinear problems. In this demo we show how several solutions of the same nonlinear problem can be computed.
The demo was contributed by Patrick Farrell.
Deflation [FBF15] is a numerical technique for computing multiple solutions of nonlinear problems. Imagine we begin Newton’s method from some initial guess \(u_0\) and find a first solution \(u_1\). Under mild regularity conditions we can then deflate the solution \(u_1\), removing it from the nonlinear problem, while leaving all other solutions. We can then initialise Newton’s method from \(u_0\) and, if it converges, find a second solution \(u_2\). This process can then be repeated until no new solutions are found from the available initial guesses.
We demonstrate the use of deflation in Firedrake on the Liouville–Bratu–Gelfand equation, a classical problem in numerical bifurcation analysis:
If \(\Omega = (0, 1)\), then for \(\lambda \in (0, \lambda^\star)\) the equation has two solutions, for \(\lambda \in \{0, \lambda^\star\}\) it has one solution, and for \(\lambda > \lambda^\star\) it has no solutions. Here \(\lambda^\star \approx 3.51\) is a constant with a known analytical expression. We will fix \(\lambda = 2\) and prescribe \(6x(1-x)\) as our initial guess.
We implement the usual weak formulation of the equation in Firedrake as standard:
from firedrake import *
mesh = UnitIntervalMesh(10)
V = FunctionSpace(mesh, "CG", 3)
x = SpatialCoordinate(mesh)[0]
u = Function(V)
guess = Function(V).interpolate(6*x*(1-x))
v = TestFunction(V)
# For this value of lambda we expect two solutions
lmbda = Constant(2)
F = - inner(grad(u), grad(v))*dx + lmbda*inner(exp(u), v)*dx
bcs = DirichletBC(V, 0, "on_boundary")
problem = NonlinearVariationalProblem(F, u, bcs)
Applying deflation requires two ingredients: the DeflatedSNES
nonlinear solver, and a Deflation
object. The Deflation
object records the solutions to be deflated, and specifies the sense of distance to use in deflation. In this example we use the metric induced by the \(L^2(\Omega)\) inner product:
sp = {"snes_type": "python",
"snes_python_type": "firedrake.DeflatedSNES",
"deflated_snes_type": "newtonls",
"deflated_snes_monitor": None,
"deflated_snes_linesearch_type": "basic",
"deflated_ksp_type": "preonly",
"deflated_pc_type": "lu"}
deflation = Deflation(op=lambda x, y: inner(x-y, x-y)*dx)
appctx = {"deflation": deflation}
solver = NonlinearVariationalSolver(problem, solver_parameters=sp, appctx=appctx)
We now find the first solution:
u.assign(guess)
solver.solve()
The first solution has now been deflated automatically in the Deflation
object. If we reset our initial guess and solve again, we find the second solution:
u.assign(guess)
solver.solve()
We can check that the two solutions are distinct:
# Prints 'Norm of difference: 1.7514003250270025'
(first, second) = deflation.roots
print(f"Norm of difference: {norm(first - second)}")
assert norm(first - second) > 1
We can plot the two solutions:
import matplotlib.pyplot as plt
ax = plt.gca()
plot(first, linestyle='-', edgecolor='tab:blue', axes=ax)
plot(second, linestyle='--', edgecolor='tab:red', axes=ax)
plt.show()
This results in the plot below (the first in blue, the second in red):

A Python script version of this demo can be found here.
References
Patrick E. Farrell, Ásgeir Birkisson, and Simon W. Funke. Deflation techniques for finding distinct solutions of nonlinear partial differential equations. SIAM Journal on Scientific Computing, 37(4):A2026–A2045, 2015. doi:10.1137/140984798.