Reference

Contents

Index

LimitedLDLFactorizations.LimitedLDLFactorizationsModule

A Pure Julia Version of limited-memory LDLᵀ factorization. A left-looking implementation of the sparse LDLᵀ factorization of a symmetric matrix with the possibility to compute a limited-memory incomplete factorization.

Dominique Orban <dominique.orban@gmail.com> Montreal, April 2015, December 2015, July 2017.

This code is strongly inspired by Lin and Moré's ICFS [1,2]. The modified version is described in [3,4].

References

[1] C.-J. Lin and J. J. Moré. Incomplete Cholesky factorizations with limited memory. SIAM Journal on Scientific Computing, 21(1):24–45, 1999. [2] http://www.mcs.anl.gov/~more/icfs [3] D. Orban. Limited-Memory LDLᵀ Factorization of Symmetric Quasi-Definite Matrices with Application to Constrained Optimization. Numerical Algorithms 70(1):9–41, 2015. DOI 10.1007/s11075-014-9933-x [4] https://github.com/optimizers/lldl

source
LimitedLDLFactorizations.LimitedLDLFactorizationMethod
LLDL = LimitedLDLFactorization(T; P = amd(T), memory = 0, α = 0, α_increase_factor = 10)
LLDL = LimitedLDLFactorization(T, ::Type{Tf}; P = amd(T), memory = 0, α = 0, α_increase_factor = 10)

Perform the allocations for the LLDL factorization of symmetric matrix whose lower triangle is T with the permutation vector P.

Arguments

  • T::SparseMatrixCSC{Tv,Ti}: lower triangle of the matrix to factorize;
  • ::Type{Tf}: type used for the factorization, by default the type of the elements of A.

Keyword arguments

  • P::AbstractVector{<:Integer} = amd(T): permutation vector;
  • memory::Int=0: extra amount of memory to allocate for the incomplete factor L. The total memory allocated is nnz(T) + n * memory, where T is the strict lower triangle of A and n is the size of A;
  • α::Number=0: initial value of the shift in case the incomplete LDLᵀ factorization of A is found to not exist. The shift will be gradually increased from this initial value until success;
  • α_increase_factor::Number=10: value by which the shift will be increased after the incomplete LDLᵀ factorization of T is found to not exist.

Example

A = sprand(Float64, 10, 10, 0.2)
T = tril(A * A' + I)
LLDL = LimitedLDLFactorization(T) # Float64 factorization
LLDL = LimitedLDLFactorization(T, Float32) # Float32 factorization
source
LimitedLDLFactorizations.abspermute!Method

Permute the elements of keys in place so that abs(x[keys[i]]) ≤ abs(x[keys[k]]) for i = 1, ..., k abs(x[keys[k]]) ≤ abs(x[keys[i]]) for i = k, ..., n, where n is the length of keys. The length of x should be at least n. Only keys is modified. From the MINPACK2 function dsel2 by Kastak, Lin and Moré.

source
LimitedLDLFactorizations.lldlMethod
lldl(A; P = amd(A), memory = 0, α = 0, droptol = 0, check_tril = true)
lldl(A, ::Type{Tf}; P = amd(A), memory = 0, α = 0, droptol = 0, check_tril = true)

Compute the limited-memory LDLᵀ factorization of A. A should be a lower triangular matrix.

Arguments

  • A::SparseMatrixCSC{Tv,Ti}: matrix to factorize (its strict lower triangle and diagonal will be extracted);
  • ::Type{Tf}: type used for the factorization, by default the type of the elements of A.

Keyword arguments

  • P::AbstractVector{<:Integer} = amd(A): permutation vector.
  • memory::Int=0: extra amount of memory to allocate for the incomplete factor L. The total memory allocated is nnz(T) + n * memory, where T is the strict lower triangle of A and n is the size of A;
  • α::Number=0: initial value of the shift in case the incomplete LDLᵀ factorization of A is found to not exist. The shift will be gradually increased from this initial value until success;
  • α_increase_factor::Number = 10: value by which the shift will be increased after the incomplete LDLᵀ factorization of T is found to not exist.
  • droptol::Tv=Tv(0): to further sparsify L, all elements with magnitude smaller than droptol are dropped;
  • check_tril::Bool = true: check if A is a lower triangular matrix.

Example

A = sprand(Float64, 10, 10, 0.2)
As = A * A' + I
LLDL = lldl(As) # lower triangle is extracted
T = tril(As)
LLDL = lldl(T) # Float64 factorization
LLDL = lldl(T, Float32) # Float32 factorization
source
LimitedLDLFactorizations.lldl_factorize!Method
lldl_factorize!(S, T; droptol = 0.0)

Perform the in-place factorization of a symmetric matrix whose lower triangle is T with the permutation vector.

Arguments

  • S::LimitedLDLFactorization{Tf, Ti};
  • T::SparseMatrixCSC{Tv,Ti}: lower triangle of the matrix to factorize.

T should keep the same nonzero pattern and the sign of its diagonal elements.

Keyword arguments

  • droptol::Tf=Tf(0): to further sparsify L, all elements with magnitude smaller than droptol are dropped.
source