Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add a normal equations solver #93

Open
wants to merge 1 commit into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
59 changes: 40 additions & 19 deletions src/KKT/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -210,22 +210,25 @@ mutable struct CholmodSPD <: CholmodSolver
# Regularization
regP::Vector{Float64} # primal
regD::Vector{Float64} # dual
variant::Bool

# Factorization
F::CHOLMOD.Factor{Float64}

# Constructor and initial memory allocation
# TODO: symbolic only + allocation
function CholmodSPD(A::AbstractMatrix{Float64})
function CholmodSPD(A::AbstractMatrix{Float64}; variant::Bool=false)
m, n = size(A)
θ = ones(Float64, n)
θ = ones(n)
regP = zeros(n)
regD = ones(m)

S = sparse(A * A') + spdiagm(0 => ones(m))
S = !variant ? sparse(A * A') + spdiagm(0 => ones(m)) : sparse(A' * A) + spdiagm(0 => ones(n))

# TODO: PSD-ness checks
F = cholesky(Symmetric(S))

return new(m, n, A, θ, zeros(Float64, n), ones(Float64, m), F)
return new(m, n, A, θ, regP, regD, variant, F)
end

end
Expand Down Expand Up @@ -260,11 +263,15 @@ function update!(
kkt.regP .= regP # Primal regularization is disabled for normal equations
kkt.regD .= regD

# Re-compute factorization
# D = (Θ^{-1} + Rp)^{-1}
D = Diagonal(one(Float64) ./ (kkt.θ .+ kkt.regP))
Rd = spdiagm(0 => kkt.regD)
S = kkt.A * D * kkt.A' + Rd
if !kkt.variant
M⁻¹ = Diagonal(1.0 ./ (kkt.θ .+ kkt.regP))
N = spdiagm(0 => kkt.regD)
S = kkt.A * M⁻¹ * kkt.A' + N
else
N⁻¹ = Diagonal(1.0 ./ kkt.regD)
M = spdiagm(0 => kkt.θ .+ kkt.regP)
S = kkt.A' * N⁻¹ * kkt.A + M
end

# Update factorization
cholesky!(kkt.F, Symmetric(S), check=false)
Expand All @@ -284,18 +291,32 @@ function solve!(
ξp::Vector{Float64}, ξd::Vector{Float64}
)
m, n = kkt.m, kkt.n

d = one(Float64) ./ (kkt.θ .+ kkt.regP)
D = Diagonal(d)

# Set-up right-hand side
ξ_ = ξp .+ kkt.A * (D * ξd)
if !kkt.variant
# Compute M⁻¹
M⁻¹ = Diagonal(1.0 ./ (kkt.θ .+ kkt.regP))

# Solve augmented system
dy .= (kkt.F \ ξ_)
# Set-up right-hand side
ξ_ = ξp .+ kkt.A * (M⁻¹ * ξd)

# Solve augmented system
dy .= (kkt.F \ ξ_)

# Recover dx
dx .= D * (kkt.A' * dy - ξd)
# Recover dx
dx .= M⁻¹ * (kkt.A' * dy - ξd)
else
# Compute N⁻¹
N⁻¹ = Diagonal(1.0 ./ kkt.regD)

# Set-up right-hand side
ξ_ = kkt.A' * (N⁻¹ * ξp) .- ξd

# Solve augmented system
dx .= (kkt.F \ ξ_)

# Recover dx
dy .= N⁻¹ * (ξp - kkt.A * dx)
end

# TODO: Iterative refinement
# * Max number of refine steps
Expand All @@ -305,6 +326,6 @@ function solve!(
# resD = - D \ dx + kkt.A' * dy - ξd
# println("\n|resP| = $(norm(resP, Inf))\n|resD| = $(norm(resD, Inf))")


return nothing
end
7 changes: 6 additions & 1 deletion test/KKT/cholmod.jl
Original file line number Diff line number Diff line change
Expand Up @@ -11,7 +11,12 @@
end

@testset "Cholesky" begin
kkt = KKT.CholmodSolver(A, normal_equations=true)
kkt = KKT.CholmodSolver(A, normal_equations=true, variant=false)
KKT.run_ls_tests(A, kkt)
end

@testset "Cholesky" begin
kkt = KKT.CholmodSolver(A, normal_equations=true, variant=true)
KKT.run_ls_tests(A, kkt)
end

Expand Down