RipQP
RipQP.ripqp
— Functionstats = ripqp(QM :: QuadraticModel{T0};
itol = InputTol(T0), scaling = true, ps = true,
normalize_rtol = true, kc = 0, mode = :mono, perturb = false,
early_multi_stop = true,
sp = (mode == :mono) ? K2LDLParams{T0}() : K2LDLParams{Float32}(),
sp2 = nothing, sp3 = nothing,
solve_method = PC(),
history = false, w = SystemWrite(), display = true) where {T0<:Real}
Minimize a convex quadratic problem. Algorithm stops when the criteria in pdd, rb, and rc are valid. Returns a GenericExecutionStats containing information about the solved problem.
QM :: QuadraticModel
: problem to solveitol :: InputTol{T, Int}
input Tolerances for the stopping criteria. SeeRipQP.InputTol
.scaling :: Bool
: activate/deactivate scaling of A and Q inQM
ps :: Bool
: activate/deactivate presolvenormalize_rtol :: Bool = true
: iftrue
, the primal and dual tolerance for the stopping criteria are normalized by the initial primal and dual residualskc :: Int
: number of centrality corrections (set to-1
for automatic computation)perturb :: Bool
: activate / deativate perturbation of the current point when μ is too smallmode :: Symbol
: should be:mono
to use the mono-precision mode,:multi
to use the multi-precision mode (start in single precision and gradually transitions toT0
),:zoom
to use the zoom procedure,:multizoom
to use the zoom procedure with multi-precision,ref
to use the QP refinement procedure, ormultiref
to use the QP refinement procedure with multi_precisionearly_multi_stop :: Bool
: stop the iterations in lower precision systems earlier in multi-precision mode, based on some quantities of the algorithmsp :: SolverParams
: choose a solver to solve linear systems that occurs at each iteration and during the initialization, seeRipQP.SolverParams
sp2 :: Union{Nothing, SolverParams}
andsp3 :: Union{Nothing, SolverParams}
: choose second and third solvers to solve linear systems that occurs at each iteration in the second and third solving phase. Whenmode != :mono
, leave tonothing
if you want to keep usingsp
. Ifsp2
is not nothing, thenmode
should be set to:multi
,:multiref
ormultizoom
.solve_method :: SolveMethod
: method used to solve the system at each iteration, usesolve_method = PC()
to use the Predictor-Corrector algorithm (default), and usesolve_method = IPF()
to use the Infeasible Path Following algorithm.solve_method2 :: Union{Nothing, SolveMethod}
andsolve_method3 :: Union{Nothing, SolveMethod}
should be used withsp2
andsp3
to choose their respective solve method. If they arenothing
, then the solve method used issolve_method
.history :: Bool
: set to true to return the primal and dual norm histories, the primal-dual relative difference history, and the number of products if using a Krylov method in thesolver_specific
field of the GenericExecutionStatsw :: SystemWrite
: configure writing of the systems to solve (no writing is done by default), seeRipQP.SystemWrite
display::Bool
: activate/deactivate iteration data display
You can also use ripqp
to solve a LLSModel:
stats = ripqp(LLS::LLSModel{T0}; mode = :mono,
sp = (mode == :mono) ? K2LDLParams{T0}() : K2LDLParams{Float32}(),
kwargs...) where {T0 <: Real}
Input Types
RipQP.InputTol
— TypeType to specify the tolerances used by RipQP.
max_iter :: Int
: maximum number of iterationsϵ_pdd
: relative primal-dual difference toleranceϵ_rb
: primal toleranceϵ_rc
: dual tolerancemax_iter1
,ϵ_pdd1
,ϵ_rb1
,ϵ_rc1
: same asmax_iter
,ϵ_pdd
,ϵ_rb
andϵ_rc
, but used for switching fromsp1
tosp2
(or from single to double precision ifsp2
isnothing
). They are only usefull whenmode=:multi
max_iter2
,ϵ_pdd2
,ϵ_rb2
,ϵ_rc2
: same asmax_iter
,ϵ_pdd
,ϵ_rb
andϵ_rc
, but used for switching fromsp2
tosp3
(or from double to quadruple precision ifsp3
isnothing
). They are only usefull whenmode=:multi
and/orT0=Float128
ϵ_rbz
: primal transition tolerance for the zoom procedure, (used only ifrefinement=:zoom
)ϵ_Δx
: step tolerance for the current point estimate (note: this criterion is currently disabled)ϵ_μ
: duality measure tolerance (note: this criterion is currently disabled)max_time
: maximum time to solve the QP
The constructor
itol = InputTol(::Type{T};
max_iter :: I = 200, max_iter1 :: I = 40, max_iter2 :: I = 180,
ϵ_pdd :: T = 1e-8, ϵ_pdd1 :: T = 1e-2, ϵ_pdd2 :: T = 1e-4,
ϵ_rb :: T = 1e-6, ϵ_rb1 :: T = 1e-4, ϵ_rb2 :: T = 1e-5, ϵ_rbz :: T = 1e-3,
ϵ_rc :: T = 1e-6, ϵ_rc1 :: T = 1e-4, ϵ_rc2 :: T = 1e-5,
ϵ_Δx :: T = 1e-16, ϵ_μ :: T = 1e-9) where {T<:Real, I<:Integer}
InputTol(; kwargs...) = InputTol(Float64; kwargs...)
returns a InputTol
struct that initializes the stopping criteria for RipQP. The 1 and 2 characters refer to the transitions between the chosen solvers in :multi
. If sp2
and sp3
are not precised when calling RipQP.ripqp
, they refer to transitions between floating-point systems.
RipQP.SystemWrite
— TypeType to write the matrix (.mtx format) and the right hand side (.rhs format) of the system to solve at each iteration.
write::Bool
: activate/deactivate writing of the systemname::String
: name of the sytem to solvekfirst::Int
: first iteration where a system should be writtenkgap::Int
: iteration gap between two problem writings
The constructor
SystemWrite(; write = false, name = "", kfirst = 0, kgap = 1)
returns a SystemWrite
structure that should be used to tell RipQP to save the system. See the tutorial for more information.
Solvers
RipQP.SolverParams
— TypeAbstract type for tuning the parameters of the different solvers. Each solver has its own SolverParams
type.
RipQP.K2LDLParams
— TypeType to use the K2 formulation with a LDLᵀ factorization. The package LDLFactorizations.jl
is used by default. The outer constructor
sp = K2LDLParams(; fact_alg = LDLFact(regul = :classic),
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = sqrt(eps()), δ_min = sqrt(eps()),
safety_dist_bnd = true)
creates a RipQP.SolverParams
. regul = :dynamic
uses a dynamic regularization (the regularization is only added if the LDLᵀ factorization encounters a pivot that has a small magnitude). regul = :none
uses no regularization (not recommended). When regul = :classic
, the parameters ρ0
and δ0
are used to choose the initial regularization values. fact_alg
should be a RipQP.AbstractFactorization
. safety_dist_bnd = true
: boolean used to determine if the regularization values should be updated (or if the algorithm should transition to another solver in multi mode) if the variables are too close from their bound.
RipQP.K2_5LDLParams
— TypeType to use the K2.5 formulation with a LDLᵀ factorization. The package LDLFactorizations.jl
is used by default. The outer constructor
sp = K2_5LDLParams(; fact_alg = LDLFact(regul = :classic), ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5)
creates a RipQP.SolverParams
. regul = :dynamic
uses a dynamic regularization (the regularization is only added if the LDLᵀ factorization encounters a pivot that has a small magnitude). regul = :none
uses no regularization (not recommended). When regul = :classic
, the parameters ρ0
and δ0
are used to choose the initial regularization values. fact_alg
should be a RipQP.AbstractFactorization
.
RipQP.K2KrylovParams
— TypeType to use the K2 formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K2KrylovParams(; uplo = :L, kmethod = :minres, preconditioner = Identity(),
rhs_scale = true, form_mat = false, equilibrate = false,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e2 * sqrt(eps()), δ_min = 1e2 * sqrt(eps()),
itmax = 0, memory = 20)
creates a RipQP.SolverParams
.
RipQP.K2_5KrylovParams
— TypeType to use the K2.5 formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K2_5KrylovParams(; uplo = :L, kmethod = :minres, preconditioner = Identity(),
rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e2 * sqrt(eps()), δ_min = 1e2 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:minres
:minres_qlp
:symmlq
RipQP.K2StructuredParams
— TypeType to use the K2 formulation with a structured Krylov method, using the package Krylov.jl
. This only works for solving Linear Problems. The outer constructor
K2StructuredParams(; uplo = :L, kmethod = :trimr, rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ_min = 1e2 * sqrt(eps()), δ_min = 1e2 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:tricg
:trimr
:gpmr
The mem
argument sould be used only with gpmr
.
RipQP.K2_5StructuredParams
— TypeType to use the K2.5 formulation with a structured Krylov method, using the package Krylov.jl
. This only works for solving Linear Problems. The outer constructor
K2_5StructuredParams(; uplo = :L, kmethod = :trimr, rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e2 * sqrt(eps()), δ_min = 1e2 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:tricg
:trimr
:gpmr
The mem
argument sould be used only with gpmr
.
RipQP.K3KrylovParams
— TypeType to use the K3 formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K3KrylovParams(; uplo = :L, kmethod = :qmr, preconditioner = Identity(),
rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e3 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:qmr
:bicgstab
:usymqr
RipQP.K3SKrylovParams
— TypeType to use the K3S formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K3SKrylovParams(; uplo = :L, kmethod = :minres, preconditioner = Identity(),
rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e3 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:minres
:minres_qlp
:symmlq
RipQP.K3SStructuredParams
— TypeType to use the K3S formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K3SStructuredParams(; uplo = :U, kmethod = :trimr, rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e3, δ0 = sqrt(eps()) * 1e4,
ρ_min = 1e4 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:tricg
:trimr
:gpmr
The mem
argument sould be used only with gpmr
.
RipQP.K3_5KrylovParams
— TypeType to use the K3.5 formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K3_5KrylovParams(; uplo = :L, kmethod = :minres, preconditioner = Identity(),
rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e3 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:minres
:minres_qlp
:symmlq
RipQP.K3_5StructuredParams
— TypeType to use the K3.5 formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K3_5StructuredParams(; uplo = :U, kmethod = :trimr, rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e3, δ0 = sqrt(eps()) * 1e4,
ρ_min = 1e4 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:tricg
:trimr
:gpmr
The mem
argument sould be used only with gpmr
.
RipQP.K1CholParams
— TypeType to use the K1 formulation with a Cholesky factorization. The outer constructor
sp = K1CholParams(; fact_alg = LDLFact(regul = :classic),
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = sqrt(eps()), δ_min = sqrt(eps()))
creates a RipQP.SolverParams
.
RipQP.K1KrylovParams
— TypeType to use the K1 formulation with a Krylov method, using the package Krylov.jl
. The outer constructor
K1KrylovParams(; uplo = :L, kmethod = :cg, preconditioner = Identity(),
rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5,
ρ_min = 1e2 * sqrt(eps()), δ_min = 1e2 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:cg
:cg_lanczos
:cr
:minres
:minres_qlp
:symmlq
RipQP.K1CholDenseParams
— TypeType to use the K1 formulation with a dense Cholesky factorization. The input QuadraticModel should have lcon .== ucon
. The outer constructor
sp = K1CholDenseParams(; ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5)
creates a RipQP.SolverParams
.
RipQP.K2LDLDenseParams
— TypeType to use the K2 formulation with a LDLᵀ factorization. The outer constructor
sp = K2LDLDenseParams(; ρ0 = sqrt(eps()) * 1e5, δ0 = sqrt(eps()) * 1e5)
creates a RipQP.SolverParams
.
RipQP.K1_1StructuredParams
— TypeType to use the K1.1 formulation with a structured Krylov method, using the package Krylov.jl
. This only works for solving Linear Problems. The outer constructor
K1_1StructuredParams(; uplo = :L, kmethod = :lsqr, rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ_min = 1e3 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:lslq
:lsqr
:lsmr
RipQP.K1_2StructuredParams
— TypeType to use the K1.2 formulation with a structured Krylov method, using the package Krylov.jl
. This only works for solving Linear Problems. The outer constructor
K1_2StructuredParams(; uplo = :L, kmethod = :craig, rhs_scale = true,
atol0 = 1.0e-4, rtol0 = 1.0e-4,
atol_min = 1.0e-10, rtol_min = 1.0e-10,
ρ_min = 1e3 * sqrt(eps()), δ_min = 1e4 * sqrt(eps()),
itmax = 0, mem = 20)
creates a RipQP.SolverParams
. The available methods are:
:lnlq
:craig
:craigmr
Preconditioners
RipQP.AbstractPreconditioner
— TypeAbstract type for the preconditioners used with a solver using a Krylov method.
RipQP.Identity
— Typepreconditioner = Identity()
Tells RipQP not to use a preconditioner.
RipQP.Jacobi
— Typepreconditioner = Jacobi()
Preconditioner using the inverse of the diagonal of the system to solve. Works with:
RipQP.Equilibration
— Typepreconditioner = Equilibration()
Preconditioner using the equilibration algorithm in infinity norm. Works with:
RipQP.LDL
— Typepreconditioner = LDL(; T = Float32, pos = :C, warm_start = true, fact_alg = LDLFact())
Preconditioner for K2KrylovParams
using a LDL factorization in precision T
. The pos
argument is used to choose the type of preconditioning with an unsymmetric Krylov method. It can be :C
(center), :L
(left) or :R
(right). The warm_start
argument tells RipQP to solve the system with the LDL factorization before using the Krylov method with the LDLFactorization as a preconditioner. fact_alg
should be a RipQP.AbstractFactorization
.
Factorizations
RipQP.AbstractFactorization
— TypeAbstract type to select a factorization algorithm.
RipQP.LDLFact
— Typefact_alg = LDLFact(; regul = :classic)
Choose LDLFactorizations.jl
to compute factorizations.
RipQP.CholmodFact
— Typefact_alg = CholmodFact(; regul = :classic)
Choose ldlt
from Cholmod to compute factorizations. using SuiteSparse
should be used before using RipQP
.
RipQP.QDLDLFact
— Typefact_alg = QDLDLFact(; regul = :classic)
Choose QDLDL.jl
to compute factorizations. using QDLDL
should be used before using RipQP
.
RipQP.HSLMA57Fact
— Typefact_alg = HSLMA57Fact(; regul = :classic)
Choose HSL.jl
MA57 to compute factorizations. using HSL
should be used before using RipQP
.
RipQP.HSLMA97Fact
— Typefact_alg = HSLMA97Fact(; regul = :classic)
Choose HSL.jl
MA57 to compute factorizations. using HSL
should be used before using RipQP
.
RipQP.LLDLFact
— Typefact_alg = LLDLFact(; regul = :classic, mem = 0, droptol = 0.0)
Choose LimitedLDLFactorizations.jl
to compute factorizations.