Operators
LinearOperators.LinearOperator
— TypeBase type to represent a linear operator. The usual arithmetic operations may be applied to operators to combine or otherwise alter them. They can be combined with other operators, with matrices and with scalars. Operators may be transposed and conjugate-transposed using the usual Julia syntax.
LinearOperators.PreallocatedLinearOperator
— TypeType to represent a linear operator with preallocation. Implicit modifications may happen if used without care:
op = PreallocatedLinearOperator(rand(5, 5))
v = rand(5)
x = op * v # Uses internal storage and passes pointer to x
y = op * ones(5) # Overwrites the same memory as x.
y === x # true. op * v is lost
x = op * v # Uses internal storage and passes pointer to x
y = op * x # Silently overwrite x to zeros! Equivalent to mul!(x, A, x).
y == zeros(5) # true. op * v and op * x are lost
LinearOperators.opEye
— TypeopEye()
Identity operator.
opI = opEye()
v = rand(5)
@assert opI * v === v
LinearOperators.opOnes
— FunctionopOnes(T, nrow, ncol)
opOnes(nrow, ncol)
Operator of all ones of size nrow
-by-ncol
and of data type T
(defaults to Float64
).
LinearOperators.opZeros
— FunctionopZeros(T, nrow, ncol)
opZeros(nrow, ncol)
Zero operator of size nrow
-by-ncol
and of data type T
(defaults to Float64
).
LinearOperators.opDiagonal
— FunctionopDiagonal(d)
Diagonal operator with the vector d
on its main diagonal.
opDiagonal(nrow, ncol, d)
Rectangular diagonal operator of size nrow
-by-ncol
with the vector d
on its main diagonal.
LinearOperators.opInverse
— FunctionopInverse(M; symmetric=false, hermitian=false)
Inverse of a matrix as a linear operator using \
. Useful for triangular matrices. Note that each application of this operator applies \
.
LinearOperators.opCholesky
— FunctionopCholesky(M, [check=false])
Inverse of a Hermitian and positive definite matrix as a linear operator using its Cholesky factorization. The factorization is computed only once. The optional check
argument will perform cheap hermicity and definiteness checks.
LinearOperators.opLDL
— FunctionopLDL(M, [check=false])
Inverse of a symmetric matrix as a linear operator using its LDL' factorization if it exists. The factorization is computed only once. The optional check
argument will perform a cheap hermicity check.
LinearOperators.opHouseholder
— FunctionopHouseholder(h)
Apply a Householder transformation defined by the vector h
. The result is x -> (I - 2 h h') x
.
LinearOperators.opHermitian
— FunctionopHermitian(d, A)
A symmetric/hermitian operator based on the diagonal d
and lower triangle of A
.
opHermitian(A)
A symmetric/hermitian operator based on a matrix.
LinearOperators.opRestriction
— FunctionZ = opRestriction(I, ncol)
Z = opRestriction(:, ncol)
Creates a LinearOperator restricting a ncol
-sized vector to indices I
. The operation Z * v
is equivalent to v[I]
. I
can be :
.
Z = opRestriction(k, ncol)
Alias for opRestriction([k], ncol)
.
LinearOperators.opExtension
— FunctionZ = opExtension(I, ncol)
Z = opExtension(:, ncol)
Creates a LinearOperator extending a vector of size length(I)
to size ncol
, where the position of the elements on the new vector are given by the indices I
. The operation w = Z * v
is equivalent to w = zeros(ncol); w[I] = v
.
Z = opExtension(k, ncol)
Alias for opExtension([k], ncol)
.
LinearOperators.LBFGSOperator
— TypeA type for limited-memory BFGS approximations.
LinearOperators.InverseLBFGSOperator
— FunctionInverseLBFGSOperator(T, n, [mem=5; scaling=true])
InverseLBFGSOperator(n, [mem=5; scaling=true])
Construct a limited-memory BFGS approximation in inverse form. If the type T
is omitted, then Float64
is used.
LinearOperators.LSR1Operator
— TypeA type for limited-memory SR1 approximations.
Base.kron
— Functionkron(A, B)
Kronecker tensor product of A and B in linear operator form, if either or both are linear operators. If both A and B are matrices, then Base.kron
is used.
Utility functions
LinearOperators.check_ctranspose
— Functioncheck_ctranspose(op)
Cheap check that the operator and its conjugate transposed are related.
LinearOperators.check_hermitian
— Functioncheck_hermitian(op)
Cheap check that the operator is Hermitian.
LinearOperators.check_positive_definite
— Functioncheck_positive_definite(op; semi=false)
Cheap check that the operator is positive (semi-)definite.
LinearAlgebra.diag
— Functiondiag(op)
Extract the diagonal of a L-BFGS operator in forward mode.
diag(op)
Extract the diagonal of a L-SR1 operator in forward mode.
Base.Matrix
— TypeA = Matrix(op)
Materialize an operator as a dense array using op.ncol
products.
LinearOperators.hermitian
— Functionhermitian(op)
ishermitian(op)
Determine whether the operator is Hermitian.
Base.push!
— Functionpush!(op, s, y)
push!(op, s, y, α, g)
Push a new {s,y} pair into a L-BFGS operator. The second calling sequence is used in inverse LBFGS updating in conjunction with damping, where α is the most recent steplength and g the gradient used when solving d=-Hg
. In forward updating with damping, it is not necessary to supply α and g.
push!(op, s, y)
Push a new {s,y} pair into a L-SR1 operator.
LinearOperators.reset!
— Functionreset!(op)
Reset the product counters of a linear operator.
reset!(data)
Resets the given LBFGS data.
reset!(op)
Resets the LBFGS data of the given operator.
reset!(data)
Reset the given LSR1 data.
reset!(op)
Resets the LSR1 data of the given operator.
LinearOperators.shape
— Functionm, n = shape(op)
An alias for size.
Base.show
— Functionshow(io, op)
Display basic information about a linear operator.
show(io, op)
Display basic information about a linear operator.
Base.size
— Functionm, n = size(op)
Return the size of a linear operator as a tuple.
m = size(op, d)
Return the size of a linear operator along dimension d
.
LinearOperators.symmetric
— Functionsymmetric(op)
issymmetric(op)
Determine whether the operator is symmetric.
Internal
LinearOperators.LBFGSData
— TypeA data type to hold information relative to LBFGS operators.
LinearOperators.LSR1Data
— TypeA data type to hold information relative to LSR1 operators.