Nonlinear solver control
Generally, nonlinear systems in this package are solved using Newton's method. In many cases, the default settings provided by this package work well. However, the convergence of Newton's method is only guaranteed with initial values s7ufficiently close to the exact solution. This notebook describes how change the default settings for the solution of nonlinear problems with VoronoiFVM.jl.
begin
import Pkg as _Pkg
haskey(ENV, "PLUTO_PROJECT") && _Pkg.activate(ENV["PLUTO_PROJECT"])
using Revise
using VoronoiFVM
using ExtendableGrids
using Test
using PlutoUI
using LinearAlgebra
using GridVisualize
using CairoMakie
CairoMakie.activate!(; type = "png", visible = false)
GridVisualize.default_plotter!(CairoMakie)
end;
Define a nonlinear Poisson equation to have an example. Let \(Ω=(0,10)\) and define
$$\begin{aligned} -Δ u + e^u-e^{-u} & = 0 & \text{in}\; Ω \\ u(0)&=100\\ u(10)&=0 \end{aligned}$$
X = 0:0.001:1
0.0:0.001:1.0
flux(y, u, edge, data) = y[1] = u[1, 1] - u[1, 2];
function reaction(y, u, node, data)
eplus = exp(u[1])
eminus = 1 / eplus
y[1] = eplus - eminus
end
reaction (generic function with 1 method)
function bc(y, u, node, data)
boundary_dirichlet!(y, u, node; region = 1, value = 100)
boundary_dirichlet!(y, u, node; region = 2, value = 0.0)
end;
system = VoronoiFVM.System(X; flux = flux, reaction = reaction, bcondition = bc, species = 1);
Solution using default settings
begin
sol = solve(system; log = true)
hist = history(sol)
end;
scalarplot(system,
sol;
resolution = (500, 200),
xlabel = "x",
ylable = "y",
title = "solution")
With log=true
, the solve
method in addition to the solution records the solution history which after finished solution can be obtatined as history(sol)
.
The history can be plotted:
function plothistory(h)
scalarplot(1:length(h),
h;
resolution = (500, 200),
yscale = :log,
xlabel = "step",
ylabel = "||δu||_∞",
title = "Maximum norm of Newton update")
end;
plothistory(hist)
History can be summarized:
summary(hist)
(seconds = 13.6, tasm = 11.2, tlinsolve = 0.713, iters = 93, absnorm = 1.6e-12, relnorm = 1.62e-14, roundoff = 1.69e-13, factorizations = 1, liniters = 0)
History can be explored in detail:
93-element Vector{Any}: (update = 98.8, contraction = 1.0, round = 426.0) (update = 1.0, contraction = 0.0101, round = 0.0222) (update = 1.0, contraction = 1.0, round = 0.0223) (update = 1.0, contraction = 1.0, round = 0.0225) (update = 1.0, contraction = 1.0, round = 0.0227) (update = 1.0, contraction = 1.0, round = 0.0229) (update = 1.0, contraction = 1.0, round = 0.0231) ⋮ (update = 0.87, contraction = 0.88, round = 0.0247) (update = 0.477, contraction = 0.548, round = 0.0167) (update = 0.0915, contraction = 0.192, round = 0.00446) (update = 0.00266, contraction = 0.029, round = 0.000177) (update = 2.22e-6, contraction = 0.000837, round = 1.9e-7) (update = 1.6e-12, contraction = 7.21e-7, round = 1.69e-13)
With default solver settings, for this particular problem, Newton's method needs 93 iteration steps.
check(sol) = isapprox(sum(sol), 2554.7106586964906; rtol = 1.0e-12)
check (generic function with 1 method)
@test check(sol)
Test Passed
Damping
Try to use a damped version of Newton method. The damping scheme is rather simple: an initial damping value damp_initial
is increased by a growth factor damp_growth
in each iteration until it reaches 1.
begin
sol1 = solve(system; log = true, inival = 1, damp_initial = 0.15, damp_growth = 1.5)
hist1 = history(sol1)
end
28-element NewtonSolverHistory: 97.81584868964667 9.162760369739365 0.9999999998511512 0.9999999997401952 0.999999999441547 0.9999999983981809 0.9999999941837444 ⋮ 0.7769510871365355 0.50248357049323 0.17071846878363076 0.015586640226738904 0.00011698226343053463 6.5218450955795574e-9
VoronoiFVM.details(hist1)
28-element Vector{Any}: (update = 97.8, contraction = 1.0, round = 5.29) (update = 9.16, contraction = 0.0937, round = 0.0298) (update = 1.0, contraction = 0.109, round = 0.0446) (update = 1.0, contraction = 1.0, round = 0.0655) (update = 1.0, contraction = 1.0, round = 0.0959) (update = 1.0, contraction = 1.0, round = 0.123) (update = 1.0, contraction = 1.0, round = 0.119) ⋮ (update = 0.777, contraction = 0.85, round = 0.00032) (update = 0.502, contraction = 0.647, round = 0.000207) (update = 0.171, contraction = 0.34, round = 7.04e-5) (update = 0.0156, contraction = 0.0913, round = 6.43e-6) (update = 0.000117, contraction = 0.00751, round = 4.83e-8) (update = 6.52e-9, contraction = 5.58e-5, round = 2.69e-12)
summary(hist1)
(seconds = 0.475, tasm = 0.43, tlinsolve = 0.0443, iters = 28, absnorm = 6.52e-9, relnorm = 6.67e-11, roundoff = 2.69e-12, factorizations = 1, liniters = 0)
We see that the number of iterations decreased significantly.
@test check(sol1)
Test Passed
Embedding
Another possibility is the embedding (or homotopy) via a parameter: start with solving a simple problem and increase the level of complexity by increasing the parameter until the full problem is solved. This process is controlled by the parameters
Δp
: initial parameter step sizeΔp_min
: minimal parameter step sizeΔp_max
: maximum parameter step sizeΔp_grow
: maximum growth factorΔu_opt
: optimal difference of solutions between two embedding steps
After successful solution of a parameter, the new parameter step size is calculated as Δp_new=min(Δp_max, Δp_grow, Δp*Δu_opt/(|u-u_old|+1.0e-14))
and adjusted to the end of the parameter interval.
If the solution is unsuccessful, the parameter stepsize is halved and solution is retried, until the minimum step size is reached.
function pbc(y, u, node, data)
boundary_dirichlet!(y, u, node; region = 1, value = 100 * embedparam(node))
boundary_dirichlet!(y, u, node; region = 2, value = 0)
end;
system2 = VoronoiFVM.System(X;
flux = flux,
reaction = function (y, u, node, data)
reaction(y, u, node, data)
y[1] = y[1] * embedparam(node)
end,
bcondition = pbc,
species = 1,);
begin
sol2 = solve(system2;
inival = 0,
log = true,
embed = (0, 1),
Δp = 0.1,
Δp_grow = 1.2,
Δu_opt = 15)
history2 = history(sol2)
end
9-element TransientSolverHistory: [0.0] [9.989343055885163, 0.98148159332474, 0.8053267053718508, 0.34118506593594755, 0.03588566608807058, 0.00030501510006535184, 2.0656145341786982e-8, 3.5058194770979925e-15] [11.34146393762356, 0.9997315206431725, 0.9990808258671514, 0.9966384496135858, 0.9877041414138586, 0.9368309534589466, 0.7284028335581785, 0.3000679588018384, 0.034624159753369424, 0.00039224137194873325, 5.021311921932754e-8, 1.3972756883275171e-15] [1.5085724152572166, 0.4205506561491683, 0.10826945171253802, 0.005760884131907183, 1.5240635820497445e-5, 1.0633964457291663e-10] [0.32730498752976783, 0.04842614744099976, 0.001078162227433007, 4.82292409015689e-7, 8.958093218077431e-14] [0.18959174175719817, 0.01947525749984668, 0.00017208751333215453, 1.2283470243608807e-8, 2.6864715379034067e-15] [0.22151709368816527, 0.01350594149309612, 8.295254108143132e-5, 2.8747465562928216e-9, 1.5200200846066591e-15] [0.9999839855864412, 0.9999129411082036, 0.9996450908924625, 0.9987144967036693, 0.9956410853148298, 0.9858738195909444, 0.956086313069598, 0.8714057437368012, 0.6663892012877948, 0.32544591325884925, 0.05964842227578665, 0.001663844974785831, 1.2451922222765833e-6, 6.958364377298142e-13] [0.9999999999255985, 0.9999999995955116, 0.9999999983507291, 0.9999999940224222, 0.9999999796890737, 0.9999999337470156, 0.9999997898900223, 0.9999993472708903, 0.9999980039123316, 0.9999939712056518 … 0.9888382056759235, 0.9682851805189645, 0.912237194563469, 0.7725313381228643, 0.49505826421414995, 0.16481352791231269, 0.014468420882723345, 0.00010072395527063066, 4.834945474226937e-9, 1.0894463251846006e-15]
summary(history2)
(seconds = 0.982, tasm = 0.86, tlinsolve = 0.12, steps = 9, iters = 81, maxabsnorm = 1.06e-10, maxrelnorm = 7.05e-11, maxroundoff = 2.26e-13, iters_per_step = 10.1, facts_per_step = 0.0, liniters_per_step = 0.0)
plothistory(vcat(history2[2:end]...))
sol2.u[end]
1×1001 Matrix{Float64}: 79.6896 17.8835 14.5192 13.1759 12.3601 … 0.00695716 0.00347857 3.47857e-30
@test check(sol2.u[end])
Test Passed
For this particular problem, embedding uses less overall Newton steps than the default settings, but the damped method is faster.
SolverControl
Here we show the docsctring of SolverControl
(formerly NewtonControl
). This is a struct which can be passed to the solve
method. Alternatively, as shown in this notebook, keyword arguments named like its entries can be passed directly to the solve
method.
@doc VoronoiFVM.SolverControl
SolverControl
SolverControl(;kwargs...)
Solver control parameter for time stepping, embedding, Newton method and linear solver control. All field names can be used as keyword arguments for solve(system::VoronoiFVM.AbstractSystem; kwargs...)
Newton's method solves \(F(u)=0\) by the iterative procedure \(u_{i+1}=u_{i} - d_i F'(u_i)^{-1}F(u_i)\) starting with some initial value \(u_0\), where \(d_i\) is a damping parameter.
verbose::Union{Bool, String}
: Verbosity control. A collection of output categories is given in a string composed of the following letters:a: allocation warnings
d: deprecation warnings
e: time/parameter evolution log
n: newton solver log
l: linear solver log
Alternatively, a Bool value can be given, resulting in
true: "neda"
false: "da"
Switch off all output including deprecation warnings via
verbose=""
. In the output, corresponding messages are marked e.g. via '[n]',[a]
etc. (besides of '[l]')
abstol::Float64
: Tolerance (in terms of norm of Newton update): terminate if \(\Delta u_i=||u_{i+1}-u_i||_\infty <\)abstol
.
reltol::Float64
: Tolerance (relative to the size of the first update): terminate if \(\Delta u_i/\Delta u_1<\)reltol
.
maxiters::Int64
: Maximum number of newton iterations.
tol_round::Float64
: Tolerance for roundoff error detection: terminate if \(|\;||u_{i+1}||_1 - ||u_{i}||_1\;|/ ||u_{i}||_1<\)tol_round
occurredmax_round
times in a row.
tol_mono::Float64
: Tolerance for monotonicity test: terminate with error if \(\Delta u_i/\Delta u_{i-1}>\)1/tol_mono
.
damp_initial::Float64
: Initial damping parameter \(d_0\). To handle convergence problems, set this to a value less than 1.
damp_growth::Float64
: Damping parameter growth factor: \(d_{i+1}=\min(d_i\cdot\)max_growth
\(,1)\). It should be larger than 1.
max_round::Int64
: Maximum number of consecutive iterations within roundoff error tolerance The default effectively disables this criterion.
unorm::Function
: Calculation of Newton update norm
rnorm::Function
: Functional for roundoff error calculation
method_linear::Union{Nothing, LinearSolve.SciMLLinearSolveAlgorithm}
: Solver method for linear systems (see LinearSolve.jl). If givennothing
, as default are chosen:UMFPACKFactorization()
for sparse matrices with Float64SparspakFactorization()
for sparse matrices with general number types.Defaults from LinearSolve.jl for tridiagonal and banded matrices
Users should experiment with what works best for their problem.
For an overeview on possible alternatives, see
VoronoiFVM.LinearSolverStrategy
.
reltol_linear::Float64
: Relative tolerance of iterative linear solver.
abstol_linear::Float64
: Absolute tolerance of iterative linear solver.
maxiters_linear::Int64
: Maximum number of iterations of linear solver
precon_linear::Union{Nothing, Function, ExtendableSparse.AbstractFactorization, LinearSolve.SciMLLinearSolveAlgorithm, Type}
: Constructor for preconditioner for linear systems. This should work as a functionprecon_linear(A)
which takes an AbstractSparseMatrixCSC (e.g. an ExtendableSparseMatrix) and returns a preconditioner object in the sense ofLinearSolve.jl
, i.e. which has anldiv!(u,A,v)
method. Useful examples:ExtendableSparse.ILUZero
ExtendableSparse.Jacobi
For easy access to this functionality, see see also
VoronoiFVM.LinearSolverStrategy
.
keepcurrent_linear::Bool
: Update preconditioner in each Newton step ? Translates toreuse_precs=!keepcurrent_linear
for LinearSolve.
Δp::Float64
: Initial parameter step for embedding.
Δp_max::Float64
: Maximal parameter step size.
Δp_min::Float64
: Minimal parameter step size.
Δp_grow::Float64
: Maximal parameter step size growth.
Δp_decrease::Float64
: Parameter step decrease factor upon rejection
Δt::Float64
: Initial time step size.
Δt_max::Float64
: Maximal time step size.
Δt_min::Float64
: Minimal time step size.
Δt_grow::Float64
: Maximal time step size growth.
Δt_decrease::Float64
: Time step decrease factor upon rejection
Δu_opt::Float64
: Optimal size of update for time stepping and embedding. The algorithm tries to keep the difference in norm between "old" and "new" solutions approximately at this value.
Δu_max_factor::Float64
: Control maximum sice of updateΔu
for time stepping and embedding relative toΔu_opt
. Time steps withΔu > Δu_max_factor*Δu_opt
will be rejected.
force_first_step::Bool
: Force first timestep.
num_final_steps::Int64
: Number of final steps to adjust at end of time interval in order to prevent breakdown of step size.
handle_exceptions::Bool
: Handle exceptions during transient solver and parameter embedding. Iftrue
, exceptions in Newton solves are caught, the embedding resp. time step is lowered, and solution is retried. Moreover, if embedding or time stepping fails (e.g. due to reaching minimal step size), a warning is issued, and a solution is returned with all steps calculated so far.Otherwise (by default) errors are thrown.
store_all::Bool
: Store all steps of transient/embedding problem:
in_memory::Bool
: Store transient/embedding solution in memory
log::Any
: Record history
edge_cutoff::Float64
: Edge parameter cutoff for rectangular triangles.
pre::Function
: Functionpre(sol,t)
called before time/embedding step
post::Function
: Functionpost(sol,oldsol,t,Δt)
called after successful time/embedding step
sample::Function
: Functionsample(sol,t)
to be called for eacht in times[2:end]
delta::Function
: Time step error estimator. A functionΔu=delta(system,u,uold,t,Δt)
to calculateΔu
.
tol_absolute::Union{Nothing, Float64}
tol_relative::Union{Nothing, Float64}
damp::Union{Nothing, Float64}
damp_grow::Union{Nothing, Float64}
max_iterations::Union{Nothing, Int64}
tol_linear::Union{Nothing, Float64}
max_lureuse::Union{Nothing, Int64}
mynorm::Union{Nothing, Function}
myrnorm::Union{Nothing, Function}
Built with Julia 1.11.1 and
CairoMakie 0.12.16ExtendableGrids 1.9.2
GridVisualize 1.7.0
LinearAlgebra 1.11.0
Pkg 1.11.0
PlutoUI 0.7.60
Revise 3.5.18
Test 1.11.0
VoronoiFVM 1.25.1