Nonlinear solvers

FESolvers.NewtonSolverType
NewtonSolver(;
    linsolver=BackslashSolver(), linesearch=NoLineSearch(), 
    maxiter=10, tolerance=1.e-6,
    update_jac_first=true, update_jac_each=true)

Use the standard Newton-Raphson solver to solve the nonlinear problem r(x) = 0 with tolerance within the maximum number of iterations maxiter.

Quasi-Newton methods

Linesearch: The linsolver argument determines the used linear solver whereas the linesearch can be set currently between NoLineSearch or ArmijoGoldstein. The latter globalizes the Newton strategy.

**Jacobian updates: ** The keyword update_jac_first decides if the jacobian from the previously converged time step should be updated after calling update_to_next_step!, or to use the old. Setting update_jac_each implies that the jacobian will not be updated during the iterations. If both update_jac_each and update_jac_first are false, the initial jacobian will be used throughout. Note that these keywords require that the problem respects the update_jacobian keyword given to update_problem!. For time-independent problems or time-depdent problems with constant time steps, update_jac_first=false is often a good choice. However, for time-dependent problems with changing time step length, the standard solver (default), may work better.

source
FESolvers.AdaptiveNewtonSolverType
AdaptiveNewtonSolver(;
    update_types, switch_criterion,
    linsolver=BackslashSolver(), linesearch=NoLineSearch(),
    maxiter=10, tolerance=1e-6, update_jac_first=true)

Define an adaptive newton solver, where the update type (given to UpdateSpec) changes during the iterations. The different options are given to update_types, and the criterion to switch between them is given to switch_criterion. Remaining options are similar to the standard newton solver.

source
FESolvers.SteepestDescentType
SteepestDescent(;maxiter=10, tolerance=1.e-6)

Use a steepest descent solver to solve the nonlinear problem r(x) = 0, which minimizes a potential $\Pi$ with tolerance and the maximum number of iterations maxiter.

This method is second derivative free and is not as locally limited as a Newton-Raphson scheme. Thus, it is especially suited for strongly nonlinear behavior with potentially vanishing tangent stiffnesses. For this method, it is required to implement getdescentpreconditioner or alternatively getsystemmatrix with SteepestDescent.

source
FESolvers.LinearProblemSolverType
LinearProblemSolver(;linsolver=BackslashSolver())

This is a special type of "Nonlinear solver", which actually only solves linear problems, but allows all other features (i.e. time stepping and postprocessing) of the FESolvers package to be used. In particular, it allows you to maintain all other parts of your problem exactly the same as for a nonlinear problem, but it is possible to get better performance as it is, in principle, not necessary to assemble twice in each time step.

This solver is specialized for linear problems of the form

\[\boldsymbol{r}(\boldsymbol{x}(t),t)=\boldsymbol{K}(t) \boldsymbol{x}(t) - \boldsymbol{f}(t)\]

where $\boldsymbol{K}=\partial \boldsymbol{r}/\partial \boldsymbol{x}$. It expects that $\boldsymbol{x}(t)$ and $\boldsymbol{r}(t)$ have been updated to $\boldsymbol{x}_\mathrm{bc}$ and $\boldsymbol{r}_\mathrm{bc}=\boldsymbol{K}(t)\boldsymbol{x}_\mathrm{bc}-\boldsymbol{f}(t)$, such that

\[\boldsymbol{x}(t) = \boldsymbol{x}_\mathrm{bc} - \boldsymbol{K}^{-1}(t)\boldsymbol{r}_\mathrm{bc} = \boldsymbol{x}_\mathrm{bc} - \boldsymbol{K}^{-1}(t)\left[\boldsymbol{K}(t)\boldsymbol{x}_\mathrm{bc}-\boldsymbol{f}(t)\right] = \boldsymbol{K}^{-1}(t)\boldsymbol{f}(t)\]

is the solution to the current time step. This normally implies when using Ferrite the same procedure as for nonlinear problems, i.e. that the boundary conditions are applied in update_to_next_step! and update_problem!, as well as the calculation of the residual according to above and that apply_zero!(K,r,ch) is called (on both the residual and the stiffness matrix).

If you have strange results when running the LinearProblemSolver, please ensure that the problem converges in one iteration for the NewtonSolver

source
FESolvers.DynamicSolverType
DynamicSolver(nlsolver, updater)

DynamicSolver contains a base nlsolver that must support set_update_type! such that the value of type in UpdateSpec can be changed dynamically. This is typically used to change regularization factors for the jacobian calculation. updater should be a function with the signature type, reset, finished = updater(nlsolver, num_attempts), where type is given to UpdateSpec in nlsolver, reset says if the problem should be reset, and finished tells if this is the last update that can be done (and if not converged then, it will fail). The inputs are the current nlsolver as well as the number of attempts.

source

Further details on selected solvers

AdaptiveNewtonSolver

For the adaptive newton solver, a few basic switchers have been implemented.

FESolvers.ToleranceSwitchType
ToleranceSwitch(;switch_at)

Use update_types[1] when the convergence measure is larger than switch_at. When the convergence measure is below switch_at, use update_types[2]. If the convergence measure increases above switch_at again, reset the problem and change to update_types[1] iterations with update_types[1].

source
FESolvers.IncreaseSwitchType
IncreaseSwitch(;num_slow)

Use the (typically) fast but less stable update_types[1] as long as the convergence measure is decreasing. If it starts increasing, switch to the (typically) slower but more stable update_types[2]. Use this for num_slow iterations after the convergence measure starts to decrease again.

source

To implement a custom switcher for AdaptiveNewtonSolver, define a new struct and the function switch_information for that struct:

FESolvers.switch_informationFunction
switch_information(switch_criterion, nlsolver)

Create a custom switch_criterion by overloading this function, which given the defined switch_criterion and the nlsolver, should return reset_problem::Bool and new_nr::Int, which determines if the problem should be reset to the state before its last update and which update_type should be used next, respectively.

source

Problem update specification

An UpdateSpec, which can be querried for information, is passed to update_problem! to give instructions on how to update the problem.

FESolvers.UpdateSpecType
UpdateSpec(;jacobian, residual, type=nothing)

An UpdateSpec is sent to update_problem! to pass the problem information about how it should be updated. The following methods should be used by the problem to request information

  • should_update_jacobian(::UpdateSpec)::Bool: Self explanatory
  • should_update_residual(::UpdateSpec)::Bool: Self explanatory
  • get_update_type(::UpdateSpec): How should the problem be updated? This is used by special solvers where, for example, different approximations of the stiffness is available, see e.g. AdaptiveNewtonSolver
source

Custom solvers

A custom nonlinear solver can be written by tapping into the existing functions at different levels. For example, LinearProblemSolver defines a custom solve_nonlinear! that uses the default calculate_update!.

All nonlinear solvers are expected to implement the following methods,

FESolvers.get_solver_stateFunction
get_solver_state(nlsolver)

All nonlinear solvers that contain a SolverState should normally overload this function to make many other functions work automatically.

source
FESolvers.reset_solver_state!Function
reset_solver_state!(nlsolver)

Called at the beginning of each new time step, and resets the solver's status and potentially convergence history etc.

source

get_solver_state

If get_solver_state is not implemented, the following methods must be implemented instead

FESolvers.get_convergence_measureFunction
get_convergence_measure(nlsolver)

Get the last convergence measure for nlsolver.

get_convergence_measure(nlsolver, k::Integer)

Get the kth convergence measure in the current iteration

get_convergence_measure(nlsolver, inds)

Get a view to the vector of the convergence measures for iterations inds, i.e. view(residuals, inds)

get_convergence_measure(nlsolver, ::Colon)

Get a view to the vector of all convergence measures, i.e. get_convergence_measure(nlsolver, 1:get_num_iter(nlsolver))

source

Methods required by solve_nonlinear!

After implementing the methods above, one can either implement solve_nonlinear!, or a set of method described below to use the default implementation.

FESolvers.solve_nonlinear!Function
solve_nonlinear!(problem, nlsolver)

Solve the current time step in the nonlinear problem, (r(x) = 0), by using the nonlinear solver nlsolver. In many cases, it suffices to overload calculate_update! for a custom nonlinear solver.

source

To use the default implementation of solve_nonlinear!, calculate_update! and the other methods in the list below must be implemented for the specific nonlinear solver. The default implementation of calculate_update! may be used as well, see the description below.

FESolvers.calculate_update!Function
function calculate_update!(Δx, problem, nlsolver)

According to the nonlinear solver, nlsolver, at iteration iter, calculate the update, Δx to the unknowns x.

source
FESolvers.get_update_specFunction
get_update_spec(nlsolver)::UpdateSpec

Get the update specification during regular iterations. It is the nlsolver's job to keep track of any state required for deciding changes to the update specification during iterations.

source
FESolvers.should_reset_problemFunction
should_reset_problem(nlsolver)

Note: Custom nonlinear solvers may rely on the default false return value.

If this function returns true, the problem will be reset as update_problem!(problem, -Δa) to reset the problem to the state at the last iteration (e.g. when switching the jacobian calculation)

source

Methods required by calculate_update!

To support the default calculate_update! implementation, the following methods must be implemented for the given solver.

Additional methods that usually don't require specialization

FESolvers.should_do_initial_updateFunction
should_do_initial_update(nlsolver)

Should the problem be updated initially, before starting time stepping? Normally, not required to overload as the update specification from get_initial_update_spec is used to decide how if an update is required.

source
FESolvers.reset_problem!Function
reset_problem!(problem, nlsolver; x, Δx_old)

Reset the problem, either by giving the new vector of unknown values, x, or the last increment to be reset, Δx_old. x will not be modified if given, but Δx_old will be modified to -Δx_old if given.

source

Linesearch

Some nonlinear solvers can use linesearch as a complement, and the following linesearches are included.

FESolvers.ArmijoGoldsteinType
Armijo-Goldstein{T}(;β=0.9,μ=0.01,τ0=1.0,τmin=1e-4)

Backtracking line search based on the Armijo-Goldstein condition

\[\Pi(\boldsymbol{u} + \tau \Delta\boldsymbol{u}) \leq \Pi(\boldsymbol{u}) - \mu\tau\delta\Pi(\boldsymbol{u})[\Delta \boldsymbol{u}]\]

where $\Pi$ is the potential, $\tau$ the stepsize, and $\delta\Pi$ the residuum.

#Fields

  • β::T = 0.9 constant factor that changes the steplength τ in each iteration
  • μ::T = 0.01 second constant factor that determines how much the potential needs to decrease additionally
  • τ0::T = 1.0 start stepsize
  • τmin::T = 1e-4 minimal stepsize
source

Custom linesearch

A custom linesearch should implement the following function

FESolvers.linesearch!Function
linesearch!(Δx, problem, ls::AbstractLineSearch)

Search along Δx to find the minimum of the potential. Return the modified Δx.

source