Nonlinear Solvers

One of our goals is to make adding new and hacking existing nonlinear solvers easy!

Input File Syntax Example

Below is an example input file block for a trust region solver that leverages a direct solver with an LDLT factorization. The nonlinear solver is also utilizing warm start to improve the first few iterations of each load step.

linear solvers:
  direct:
    type: DirectLinearSolver
    factorization method: ldl

nonlinear solvers:
  trs:
    type: TrustRegionSolver
    linear solver: direct
    warm start: on

General methods and abstract types

Cthonios.AbstractNonlinearSolverType
abstract type AbstractNonlinearSolver{L, O, U, T}

a nonlinear solver needs to define the following methods

  1. check_convergence - returns a bool
  2. logger - @info information
  3. step! - step for the solver

it needs the following types

  1. a linear solver
  2. an objective
  3. an unknown vector
  4. an int called max_iter
source
Cthonios.create_unknownsMethod
create_unknowns(
    solver::Cthonios.AbstractNonlinearSolver
) -> Any

Creates a set of unknowns for the nonlinear solver

source
Cthonios.solve!Method
solve!(solver::Cthonios.AbstractNonlinearSolver, Uu, p)

Generic method to fall back on if step! is defined

source

Newton Solver

Cthonios.NewtonSolverType
struct NewtonSolver{L, O, U, T} <: Cthonios.AbstractNonlinearSolver{L, O, U, T}
  • linear_solver::Any

  • objective::Any

  • ΔUu::Any

  • timer::Any

  • max_iter::Int64

  • abs_tol::Float64

  • rel_tol::Float64

  • use_warm_start::Bool

source
Cthonios.NewtonSolverMethod
NewtonSolver(
    objective::Objective,
    p,
    linear_solver_type,
    timer
) -> NewtonSolver{_A, O} where {_A, O<:Objective}
source

Trust Region Solver

Cthonios.TrustRegionSolverType
struct TrustRegionSolver{L, O, U<:(AbstractVector), T<:TimerOutputs.TimerOutput, W, F<:FiniteElementContainers.NodalField, S<:Cthonios.TrustRegionSolverSettings} <: Cthonios.AbstractNonlinearSolver{L, O, U<:(AbstractVector), T<:TimerOutputs.TimerOutput}
  • preconditioner::Any

  • objective::Any

  • ΔUu::AbstractVector

  • timer::TimerOutputs.TimerOutput

  • warm_start::Any

  • o::AbstractVector

  • g::FiniteElementContainers.NodalField

  • Hv::FiniteElementContainers.NodalField

  • cauchy_point::AbstractVector

  • q_newton_point::AbstractVector

  • d::AbstractVector

  • y_scratch_1::AbstractVector

  • y_scratch_2::AbstractVector

  • y_scratch_3::AbstractVector

  • y_scratch_4::AbstractVector

  • use_warm_start::Bool

  • settings::Cthonios.TrustRegionSolverSettings

source
Cthonios.TrustRegionSolverMethod
TrustRegionSolver(
    objective::Objective,
    p,
    timer;
    preconditioner,
    use_warm_start
) -> TrustRegionSolver{L, O, Vector{Float64}, TimerOutputs.TimerOutput, W, F, Cthonios.TrustRegionSolverSettings} where {L<:(CholeskyPreconditioner{A, SparseArrays.CHOLMOD.Factor{Float64, Int64}} where A<:(FiniteElementContainers.StaticAssembler{Float64, Int64, _A, _B, _C, _D, _E, _F, _G, _H, _I, _J, _K, Vector{Int64}, Vector{Int64}, Vector{Float64}} where {_A<:AbstractVector{Int64}, _B<:AbstractVector{Int64}, _C<:AbstractVector{Int64}, _D<:AbstractVector{Int64}, _E<:AbstractVector{Int64}, _F<:FiniteElementContainers.NodalField, _G<:AbstractVector{Float64}, _H, _I, _J, _K})), O<:Objective, W<:Cthonios.WarmStart, F<:FiniteElementContainers.NodalField}

TODO figure out which scratch arrays can be nixed

source
Cthonios.TrustRegionSolverSettingsType
  • t_1::Float64

  • t_2::Float64

  • η_1::Float64

  • η_2::Float64

  • η_3::Float64

  • max_trust_iters::Int64

  • tol::Float64

  • max_cg_iters::Int64

  • max_cumulative_cg_iters::Int64

  • cg_tol::Float64

  • cg_inexact_solve_ratio::Float64

  • tr_size::Float64

  • min_tr_size::Float64

  • check_stability::Bool

  • use_preconditioned_inner_product_for_cg::Bool

  • use_incremental_objective::Bool

  • debug_info::Bool

  • over_iter::Int64

source
Cthonios.solve!Method
solve!(solver::TrustRegionSolver, Uu, p)

TODO Eventually map this to the interface

source
Cthonios.update_tr_sizeMethod
update_tr_size(
    solver,
    model_objective,
    real_objective,
    step_type,
    tr_size,
    real_res_norm,
    g_norm
) -> Tuple{Any, Any}
source