Library documentation
Constants and types
MatrixAlgebraKit.AbstractAlgorithm — Type
abstract type AbstractAlgorithm endSupertype to dispatch on specific implementations of different the different functions. Concrete subtypes should represent both a way to dispatch to a given implementation, as well as the configuration of that implementation.
See also select_algorithm.
MatrixAlgebraKit.Algorithm — Type
MatrixAlgebraKit.CUSOLVER_DivideAndConquer — Type
CUSOLVER_DivideAndConquer(; fixgauge::Bool = true)Algorithm type to denote the CUSOLVER driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the Divide and Conquer algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.CUSOLVER_HouseholderQR — Type
CUSOLVER_HouseholderQR(; positive = false)Algorithm type to denote the standard CUSOLVER algorithm for computing the QR decomposition of a matrix using Householder reflectors. The keyword positive = true can be used to ensure that the diagonal elements of R are non-negative.
MatrixAlgebraKit.CUSOLVER_Jacobi — Type
CUSOLVER_Jacobi(; fixgauge::Bool = true)Algorithm type to denote the CUSOLVER driver for computing the singular value decomposition of a general matrix using the Jacobi algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the singular vectors, see also gaugefix!.
MatrixAlgebraKit.CUSOLVER_QRIteration — Type
CUSOLVER_QRIteration(; fixgauge::Bool = true)Algorithm type to denote the CUSOLVER driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the QR Iteration algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.CUSOLVER_Randomized — Type
CUSOLVER_Randomized(; k, p, niters)Algorithm type to denote the CUSOLVER driver for computing the singular value decomposition of a general matrix using the randomized SVD algorithm. Here, k denotes the number of singular values that should be computed, therefore requiring k <= min(size(A)). This method is accurate for small values of k compared to the size of the input matrix, where the accuracy can be improved by increasing p, the number of additional values used for oversampling, and niters, the number of iterations the solver uses, at the cost of increasing the runtime.
See also the CUSOLVER documentation for more information.
sourceMatrixAlgebraKit.CUSOLVER_SVDPolar — Type
CUSOLVER_SVDPolar(; fixgauge::Bool = true)Algorithm type to denote the CUSOLVER driver for computing the singular value decomposition of a general matrix by using Halley's iterative algorithm to compute the polar decompositon, followed by the hermitian eigenvalue decomposition of the positive definite factor. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the singular vectors, see also gaugefix!.
MatrixAlgebraKit.DiagonalAlgorithm — Type
DiagonalAlgorithm(; kwargs...)Algorithm type to denote a native Julia implementation of the decompositions making use of the diagonal structure of the input and outputs.
sourceMatrixAlgebraKit.GLA_HouseholderQR — Type
GLA_HouseholderQR(; positive = false)Algorithm type to denote the GenericLinearAlgebra.jl implementation for computing the QR decomposition of a matrix using Householder reflectors. Currently, only blocksize = 1 and pivoted == false are supported. The keyword positive = true can be used to ensure that the diagonal elements of R are non-negative.
MatrixAlgebraKit.GLA_QRIteration — Type
GLA_QRIteration(; fixgauge::Bool = true)Algorithm type to denote the GenericLinearAlgebra.jl implementation for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.GS_QRIteration — Type
GS_QRIteration()Algorithm type to denote the GenericSchur.jl implementation for computing the eigenvalue decomposition of a non-Hermitian matrix.
sourceMatrixAlgebraKit.LAPACK_Bisection — Type
LAPACK_Bisection(; fixgauge::Bool = true)Algorithm type to denote the LAPACK driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the Bisection algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.LAPACK_DivideAndConquer — Type
LAPACK_DivideAndConquer(; fixgauge::Bool = true)Algorithm type to denote the LAPACK driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the Divide and Conquer algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.LAPACK_HouseholderLQ — Type
LAPACK_HouseholderLQ(; blocksize, positive = false)Algorithm type to denote the standard LAPACK algorithm for computing the LQ decomposition of a matrix using Householder reflectors. The specific LAPACK function can be controlled using the keyword arugments, i.e. ?gelqt will be chosen if blocksize > 1 or ?gelqf will be chosen if blocksize == 1. The keyword positive = true can be used to ensure that the diagonal elements of L are non-negative.
MatrixAlgebraKit.LAPACK_HouseholderQR — Type
LAPACK_HouseholderQR(; blocksize, positive = false, pivoted = false)Algorithm type to denote the standard LAPACK algorithm for computing the QR decomposition of a matrix using Householder reflectors. The specific LAPACK function can be controlled using the keyword arugments, i.e. ?geqrt will be chosen if blocksize > 1. With blocksize == 1, ?geqrf will be chosen if pivoted == false and ?geqp3 will be chosen if pivoted == true. The keyword positive = true can be used to ensure that the diagonal elements of R are non-negative.
MatrixAlgebraKit.LAPACK_Jacobi — Type
LAPACK_Jacobi(; fixgauge::Bool = true)Algorithm type to denote the LAPACK driver for computing the singular value decomposition of a general matrix using the Jacobi algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the singular vectors, see also gaugefix!.
MatrixAlgebraKit.LAPACK_MultipleRelativelyRobustRepresentations — Type
LAPACK_MultipleRelativelyRobustRepresentations(; fixgauge::Bool = true)Algorithm type to denote the LAPACK driver for computing the eigenvalue decomposition of a Hermitian matrix using the Multiple Relatively Robust Representations algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigenvectors, see also gaugefix!.
MatrixAlgebraKit.LAPACK_QRIteration — Type
LAPACK_QRIteration(; fixgauge::Bool = true)Algorithm type to denote the LAPACK driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the QR Iteration algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.LQViaTransposedQR — Type
LQViaTransposedQR(qr_alg)Algorithm type to denote finding the LQ decomposition of A by computing the QR decomposition of Aᵀ. The qr_alg specifies which QR-decomposition implementation to use.
MatrixAlgebraKit.NativeBlocked — Type
NativeBlocked(; blocksize = 32)
Algorithm type to denote a native blocked algorithm with given blocksize for computing the hermitian or anti-hermitian part of a matrix.
MatrixAlgebraKit.PolarNewton — Type
PolarNewton(; maxiter = 10, tol = defaulttol(A))Algorithm for computing the polar decomposition of a matrix A via scaled Newton iteration, with a maximum of maxiter iterations and until convergence up to tolerance tol.
MatrixAlgebraKit.PolarViaSVD — Type
PolarViaSVD(svd_alg)Algorithm for computing the polar decomposition of a matrix A via the singular value decomposition (SVD) of A. The svd_alg argument specifies the SVD algorithm to use.
MatrixAlgebraKit.ROCSOLVER_Bisection — Type
ROCSOLVER_Bisection(; fixgauge::Bool = true)Algorithm type to denote the ROCSOLVER driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the Bisection algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.ROCSOLVER_DivideAndConquer — Type
ROCSOLVER_DivideAndConquer(; fixgauge::Bool = true)Algorithm type to denote the ROCSOLVER driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the Divide and Conquer algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.ROCSOLVER_HouseholderQR — Type
ROCSOLVER_HouseholderQR(; positive = false)Algorithm type to denote the standard ROCSOLVER algorithm for computing the QR decomposition of a matrix using Householder reflectors. The keyword positive=true can be used to ensure that the diagonal elements of R are non-negative.
MatrixAlgebraKit.ROCSOLVER_Jacobi — Type
ROCSOLVER_Jacobi(; fixgauge::Bool = true)Algorithm type to denote the ROCSOLVER driver for computing the singular value decomposition of a general matrix using the Jacobi algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the singular vectors, see also gaugefix!.
MatrixAlgebraKit.ROCSOLVER_QRIteration — Type
ROCSOLVER_QRIteration(; fixgauge::Bool = true)Algorithm type to denote the ROCSOLVER driver for computing the eigenvalue decomposition of a Hermitian matrix, or the singular value decomposition of a general matrix using the QR Iteration algorithm. The fixgauge keyword can be used to toggle whether or not to fix the gauge of the eigen or singular vectors, see also gaugefix!.
MatrixAlgebraKit.RightNullAlgorithm — Type
RightNullAlgorithm{Kind, Alg <: AbstractAlgorithm}(alg)Wrapper type to denote the Kind of factorization that is used as a backend for right_null. By default Kind is a symbol, which can be either :lq or :svd.
MatrixAlgebraKit.RightOrthAlgorithm — Type
RightOrthAlgorithm{Kind, Alg <: AbstractAlgorithm}(alg)Wrapper type to denote the Kind of factorization that is used as a backend for right_orth. By default Kind is a symbol, which can be either :lq, :polar or :svd.
MatrixAlgebraKit.TruncatedAlgorithm — Type
TruncatedAlgorithm(alg::AbstractAlgorithm, trunc::TruncationAlgorithm)Generic wrapper type for algorithms that consist of first using alg, followed by a truncation through trunc.
MatrixAlgebraKit.TruncationByError — Type
TruncationByError(; atol::Real, rtol::Real, p::Real)Truncation strategy to discard values until the error caused by the discarded values exceeds some tolerances. See also truncerror.
MatrixAlgebraKit.TruncationByFilter — Type
TruncationByFilter(filter::Function)Truncation strategy to keep the values for which filter returns true.
See also truncfilter.
MatrixAlgebraKit.TruncationIntersection — Type
TruncationIntersection(trunc::TruncationStrategy, truncs::TruncationStrategy...)Truncation strategy that composes multiple truncation strategies, keeping values that are common between them.
sourceMatrixAlgebraKit.TruncationStrategy — Method
TruncationStrategy(; kwargs...)Select a truncation strategy based on the provided keyword arguments.
Keyword arguments
The following keyword arguments are all optional, and their default value (nothing) will be ignored. It is also allowed to combine multiple of these, in which case the kept values will consist of the intersection of the different truncated strategies.
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
Functions
MatrixAlgebraKit.copy_input — Function
copy_input(f, A)Preprocess the input A for a given function, such that it may be handled correctly later. This may include a copy whenever the implementation would destroy the original matrix, or a change of element type to something that is supported.
MatrixAlgebraKit.default_algorithm — Function
MatrixAlgebraKit.default_algorithm(f, A; kwargs...)
MatrixAlgebraKit.default_algorithm(f, ::Type{TA}; kwargs...) where {TA}Select the default algorithm for a given factorization function f and input A. In general, this is called by select_algorithm if no algorithm is specified explicitly. New types should prefer to register their default algorithms in the type domain.
MatrixAlgebraKit.default_fixgauge — Function
default_fixgauge() -> current_value
default_fixgauge(new_value::Bool) -> previous_valueGlobal toggle for enabling or disabling the default behavior of gauge fixing the output of the eigen- and singular value decompositions.
sourceMatrixAlgebraKit.default_hermitian_tol — Method
default_hermitian_tol(A)Default tolerance for deciding to warn if the provided A is not hermitian.
MatrixAlgebraKit.default_pullback_degeneracy_atol — Method
default_pullback_degeneracy_atol(A)Default tolerance for deciding when values should be considered as degenerate.
sourceMatrixAlgebraKit.default_pullback_gauge_atol — Method
default_pullback_gauge_atol(ΔA...)Default tolerance for deciding to warn if incoming adjoints of a pullback rule has components that are not gauge-invariant.
sourceMatrixAlgebraKit.default_pullback_rank_atol — Method
default_pullback_rank_atol(A)Default tolerance for deciding what values should be considered equal to 0.
sourceMatrixAlgebraKit.defaulttol — Method
defaulttol(x)Default tolerance or precision for a given object, e.g. to decide when it can be considerd to be zero or ignored in some other way, or how accurate some quantity needs to be computed.
sourceMatrixAlgebraKit.does_truncate — Method
does_truncate(alg::AbstractAlgorithm) -> BoolIndicate whether or not an algorithm will compute a truncated decomposition (such that composing the factors only approximates the input up to some tolerance).
sourceMatrixAlgebraKit.eig_full — Function
eig_full(A; kwargs...) -> D, V
eig_full(A, alg::AbstractAlgorithm) -> D, V
eig_full!(A, [DV]; kwargs...) -> D, V
eig_full!(A, [DV], alg::AbstractAlgorithm) -> D, VCompute the full eigenvalue decomposition of the square matrix A, such that A * V = V * D, where the invertible matrix V contains the eigenvectors and the diagonal matrix D contains the associated eigenvalues.
The bang method eig_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
See also eig_vals(!) and eig_trunc(!).
MatrixAlgebraKit.eig_full! — Function
eig_full(A; kwargs...) -> D, V
eig_full(A, alg::AbstractAlgorithm) -> D, V
eig_full!(A, [DV]; kwargs...) -> D, V
eig_full!(A, [DV], alg::AbstractAlgorithm) -> D, VCompute the full eigenvalue decomposition of the square matrix A, such that A * V = V * D, where the invertible matrix V contains the eigenvectors and the diagonal matrix D contains the associated eigenvalues.
The bang method eig_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
See also eig_vals(!) and eig_trunc(!).
MatrixAlgebraKit.eig_pullback! — Function
eig_pullback!(
ΔA::AbstractMatrix, A, DV, ΔDV, [ind];
degeneracy_atol::Real = default_pullback_rank_atol(DV[1]),
gauge_atol::Real = default_pullback_gauge_atol(ΔDV[2])
)Adds the pullback from the full eigenvalue decomposition of A to ΔA, given the output DV of eig_full and the cotangent ΔDV of eig_full or eig_trunc.
In particular, it is assumed that A ≈ V * D * inv(V) with thus size(A) == size(V) == size(D) and D diagonal. For the cotangents, an arbitrary number of eigenvectors or eigenvalues can be missing, i.e. for a matrix A of size (n, n), ΔV can have size (n, pV) and diagview(ΔD) can have length pD. In those cases, additionally ind is required to specify which eigenvectors or eigenvalues are present in ΔV or ΔD. By default, it is assumed that all eigenvectors and eigenvalues are present.
A warning will be printed if the cotangents are not gauge-invariant, i.e. if the restriction of V' * ΔV to rows i and columns j for which abs(D[i] - D[j]) < degeneracy_atol, is not small compared to gauge_atol.
MatrixAlgebraKit.eig_trunc — Function
eig_trunc(A; [trunc], kwargs...) -> D, V, ϵ
eig_trunc(A, alg::AbstractAlgorithm) -> D, V, ϵ
eig_trunc!(A, [DV]; [trunc], kwargs...) -> D, V, ϵ
eig_trunc!(A, [DV], alg::AbstractAlgorithm) -> D, V, ϵCompute a partial or truncated eigenvalue decomposition of the matrix A, such that A * V ≈ V * D, where the (possibly rectangular) matrix V contains a subset of eigenvectors and the diagonal matrix D contains the associated eigenvalues, selected according to a truncation strategy.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded eigenvalues.
Truncation
The truncation strategy can be controlled via the trunc keyword argument. This can be either a NamedTuple or a TruncationStrategy. If trunc is not provided or nothing, all values will be kept.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword Arguments
Other keyword arguments are passed to the algorithm selection procedure. If no explicit alg is provided, these keywords are used to select and configure the algorithm through MatrixAlgebraKit.select_algorithm. The remaining keywords after algorithm selection are passed to the algorithm constructor. See MatrixAlgebraKit.default_algorithm for the default algorithm selection behavior.
When alg is a TruncatedAlgorithm, the trunc keyword cannot be specified as the truncation strategy is already embedded in the algorithm.
The bang method eig_trunc! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
See also eig_full(!), eig_vals(!), and Truncations for more information on truncation strategies.
MatrixAlgebraKit.eig_trunc! — Function
eig_trunc(A; [trunc], kwargs...) -> D, V, ϵ
eig_trunc(A, alg::AbstractAlgorithm) -> D, V, ϵ
eig_trunc!(A, [DV]; [trunc], kwargs...) -> D, V, ϵ
eig_trunc!(A, [DV], alg::AbstractAlgorithm) -> D, V, ϵCompute a partial or truncated eigenvalue decomposition of the matrix A, such that A * V ≈ V * D, where the (possibly rectangular) matrix V contains a subset of eigenvectors and the diagonal matrix D contains the associated eigenvalues, selected according to a truncation strategy.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded eigenvalues.
Truncation
The truncation strategy can be controlled via the trunc keyword argument. This can be either a NamedTuple or a TruncationStrategy. If trunc is not provided or nothing, all values will be kept.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword Arguments
Other keyword arguments are passed to the algorithm selection procedure. If no explicit alg is provided, these keywords are used to select and configure the algorithm through MatrixAlgebraKit.select_algorithm. The remaining keywords after algorithm selection are passed to the algorithm constructor. See MatrixAlgebraKit.default_algorithm for the default algorithm selection behavior.
When alg is a TruncatedAlgorithm, the trunc keyword cannot be specified as the truncation strategy is already embedded in the algorithm.
The bang method eig_trunc! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
See also eig_full(!), eig_vals(!), and Truncations for more information on truncation strategies.
MatrixAlgebraKit.eig_trunc_pullback! — Method
eig_trunc_pullback!(
ΔA::AbstractMatrix, ΔDV, A, DV;
degeneracy_atol::Real = default_pullback_rank_atol(DV[1]),
gauge_atol::Real = default_pullback_gauge_atol(ΔDV[2])
)Adds the pullback from the truncated eigenvalue decomposition of A to ΔA, given the output DV and the cotangent ΔDV of eig_trunc.
In particular, it is assumed that A * V ≈ V * D with V a rectangular matrix of eigenvectors and D diagonal. For the cotangents, it is assumed that if ΔV is not zero, then it has the same number of columns as V, and if ΔD is not zero, then it is a diagonal matrix of the same size as D.
For this method to work correctly, it is also assumed that the remaining eigenvalues (not included in D) are (sufficiently) separated from those in D.
A warning will be printed if the cotangents are not gauge-invariant, i.e. if the restriction of V' * ΔV to rows i and columns j for which abs(D[i] - D[j]) < degeneracy_atol, is not small compared to gauge_atol.
MatrixAlgebraKit.eig_vals — Function
eig_vals(A; kwargs...) -> D
eig_vals(A, alg::AbstractAlgorithm) -> D
eig_vals!(A, [D]; kwargs...) -> D
eig_vals!(A, [D], alg::AbstractAlgorithm) -> DCompute the list of eigenvalues of A.
The bang method eig_vals! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided D as output.
See also eig_full(!) and eig_trunc(!).
MatrixAlgebraKit.eig_vals! — Function
eig_vals(A; kwargs...) -> D
eig_vals(A, alg::AbstractAlgorithm) -> D
eig_vals!(A, [D]; kwargs...) -> D
eig_vals!(A, [D], alg::AbstractAlgorithm) -> DCompute the list of eigenvalues of A.
The bang method eig_vals! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided D as output.
See also eig_full(!) and eig_trunc(!).
MatrixAlgebraKit.eigh_full — Function
eigh_full(A; kwargs...) -> D, V, ϵ
eigh_full(A, alg::AbstractAlgorithm) -> D, V, ϵ
eigh_full!(A, [DV]; kwargs...) -> D, V, ϵ
eigh_full!(A, [DV], alg::AbstractAlgorithm) -> D, V, ϵCompute the full eigenvalue decomposition of the symmetric or hermitian matrix A, such that A * V = V * D, where the unitary matrix V contains the orthogonal eigenvectors and the real diagonal matrix D contains the associated eigenvalues.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded eigenvalues.
The bang method eigh_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
Note that eigh_full and its variants assume that the input matrix is hermitian, or thus symmetric if the input is real. The resulting algorithms exploit this structure, and return eigenvalues that are always real, and eigenvectors that are orthogonal and have the same eltype as the input matrix. If the input matrix does not have this structure, the generic eigenvalue decomposition provided by eig_full and its variants should be used instead.
See also eigh_vals(!) and eigh_trunc(!).
MatrixAlgebraKit.eigh_full! — Function
eigh_full(A; kwargs...) -> D, V, ϵ
eigh_full(A, alg::AbstractAlgorithm) -> D, V, ϵ
eigh_full!(A, [DV]; kwargs...) -> D, V, ϵ
eigh_full!(A, [DV], alg::AbstractAlgorithm) -> D, V, ϵCompute the full eigenvalue decomposition of the symmetric or hermitian matrix A, such that A * V = V * D, where the unitary matrix V contains the orthogonal eigenvectors and the real diagonal matrix D contains the associated eigenvalues.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded eigenvalues.
The bang method eigh_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
Note that eigh_full and its variants assume that the input matrix is hermitian, or thus symmetric if the input is real. The resulting algorithms exploit this structure, and return eigenvalues that are always real, and eigenvectors that are orthogonal and have the same eltype as the input matrix. If the input matrix does not have this structure, the generic eigenvalue decomposition provided by eig_full and its variants should be used instead.
See also eigh_vals(!) and eigh_trunc(!).
MatrixAlgebraKit.eigh_pullback! — Function
eigh_pullback!(
ΔA::AbstractMatrix, A, DV, ΔDV, [ind];
degeneracy_atol::Real = default_pullback_rank_atol(DV[1]),
gauge_atol::Real = default_pullback_gauge_atol(ΔDV[2])
)Adds the pullback from the Hermitian eigenvalue decomposition of A to ΔA, given the output DV of eigh_full and the cotangent ΔDV of eigh_full or eigh_trunc.
In particular, it is assumed that A ≈ V * D * V' with thus size(A) == size(V) == size(D) and D diagonal. For the cotangents, an arbitrary number of eigenvectors or eigenvalues can be missing, i.e. for a matrix A of size (n, n), ΔV can have size (n, pV) and diagview(ΔD) can have length pD. In those cases, additionally ind is required to specify which eigenvectors or eigenvalues are present in ΔV or ΔD. By default, it is assumed that all eigenvectors and eigenvalues are present.
A warning will be printed if the cotangents are not gauge-invariant, i.e. if the anti-hermitian part of V' * ΔV, restricted to rows i and columns j for which `abs(D[i]
- D[j]) < degeneracyatol
, is not small compared togaugeatol`.
MatrixAlgebraKit.eigh_trunc — Function
eigh_trunc(A; [trunc], kwargs...) -> D, V, ϵ
eigh_trunc(A, alg::AbstractAlgorithm) -> D, V, ϵ
eigh_trunc!(A, [DV]; [trunc], kwargs...) -> D, V, ϵ
eigh_trunc!(A, [DV], alg::AbstractAlgorithm) -> D, V, ϵCompute a partial or truncated eigenvalue decomposition of the symmetric or hermitian matrix A, such that A * V ≈ V * D, where the isometric matrix V contains a subset of the orthogonal eigenvectors and the real diagonal matrix D contains the associated eigenvalues, selected according to a truncation strategy.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded eigenvalues.
Truncation
The truncation strategy can be controlled via the trunc keyword argument. This can be either a NamedTuple or a TruncationStrategy. If trunc is not provided or nothing, all values will be kept.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
Other keyword arguments are passed to the algorithm selection procedure. If no explicit alg is provided, these keywords are used to select and configure the algorithm through MatrixAlgebraKit.select_algorithm. The remaining keywords after algorithm selection are passed to the algorithm constructor. See MatrixAlgebraKit.default_algorithm for the default algorithm selection behavior.
When alg is a TruncatedAlgorithm, the trunc keyword cannot be specified as the truncation strategy is already embedded in the algorithm.
The bang method eigh_trunc! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
Note that eigh_full and its variants assume that the input matrix is hermitian, or thus symmetric if the input is real. The resulting algorithms exploit this structure, and return eigenvalues that are always real, and eigenvectors that are orthogonal and have the same eltype as the input matrix. If the input matrix does not have this structure, the generic eigenvalue decomposition provided by eig_full and its variants should be used instead.
See also eigh_full(!), eigh_vals(!), and Truncations for more information on truncation strategies.
MatrixAlgebraKit.eigh_trunc! — Function
eigh_trunc(A; [trunc], kwargs...) -> D, V, ϵ
eigh_trunc(A, alg::AbstractAlgorithm) -> D, V, ϵ
eigh_trunc!(A, [DV]; [trunc], kwargs...) -> D, V, ϵ
eigh_trunc!(A, [DV], alg::AbstractAlgorithm) -> D, V, ϵCompute a partial or truncated eigenvalue decomposition of the symmetric or hermitian matrix A, such that A * V ≈ V * D, where the isometric matrix V contains a subset of the orthogonal eigenvectors and the real diagonal matrix D contains the associated eigenvalues, selected according to a truncation strategy.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded eigenvalues.
Truncation
The truncation strategy can be controlled via the trunc keyword argument. This can be either a NamedTuple or a TruncationStrategy. If trunc is not provided or nothing, all values will be kept.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
Other keyword arguments are passed to the algorithm selection procedure. If no explicit alg is provided, these keywords are used to select and configure the algorithm through MatrixAlgebraKit.select_algorithm. The remaining keywords after algorithm selection are passed to the algorithm constructor. See MatrixAlgebraKit.default_algorithm for the default algorithm selection behavior.
When alg is a TruncatedAlgorithm, the trunc keyword cannot be specified as the truncation strategy is already embedded in the algorithm.
The bang method eigh_trunc! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
Note that eigh_full and its variants assume that the input matrix is hermitian, or thus symmetric if the input is real. The resulting algorithms exploit this structure, and return eigenvalues that are always real, and eigenvectors that are orthogonal and have the same eltype as the input matrix. If the input matrix does not have this structure, the generic eigenvalue decomposition provided by eig_full and its variants should be used instead.
See also eigh_full(!), eigh_vals(!), and Truncations for more information on truncation strategies.
MatrixAlgebraKit.eigh_trunc_pullback! — Method
eigh_trunc_pullback!(
ΔA::AbstractMatrix, A, DV, ΔDV;
degeneracy_atol::Real = default_pullback_rank_atol(DV[1]),
gauge_atol::Real = default_pullback_gauge_atol(ΔDV[2])
)Adds the pullback from the truncated Hermitian eigenvalue decomposition of A to ΔA, given the output DV and the cotangent ΔDV of eig_trunc.
In particular, it is assumed that A * V ≈ V * D with V a rectangular matrix of eigenvectors and D diagonal. For the cotangents, it is assumed that if ΔV is not zero, then it has the same number of columns as V, and if ΔD is not zero, then it is a diagonal matrix of the same size as D.
For this method to work correctly, it is also assumed that the remaining eigenvalues (not included in D) are (sufficiently) separated from those in D.
A warning will be printed if the cotangents are not gauge-invariant, i.e. if the restriction of V' * ΔV to rows i and columns j for which abs(D[i] - D[j]) < degeneracy_atol, is not small compared to gauge_atol.
MatrixAlgebraKit.eigh_vals — Function
eigh_vals(A; kwargs...) -> D
eigh_vals(A, alg::AbstractAlgorithm) -> D
eigh_vals!(A, [D]; kwargs...) -> D
eigh_vals!(A, [D], alg::AbstractAlgorithm) -> DCompute the list of (real) eigenvalues of the symmetric or hermitian matrix A.
The bang method eigh_vals! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
Note that eigh_full and its variants assume that the input matrix is hermitian, or thus symmetric if the input is real. The resulting algorithms exploit this structure, and return eigenvalues that are always real, and eigenvectors that are orthogonal and have the same eltype as the input matrix. If the input matrix does not have this structure, the generic eigenvalue decomposition provided by eig_full and its variants should be used instead.
See also eigh_full(!) and eigh_trunc(!).
MatrixAlgebraKit.eigh_vals! — Function
eigh_vals(A; kwargs...) -> D
eigh_vals(A, alg::AbstractAlgorithm) -> D
eigh_vals!(A, [D]; kwargs...) -> D
eigh_vals!(A, [D], alg::AbstractAlgorithm) -> DCompute the list of (real) eigenvalues of the symmetric or hermitian matrix A.
The bang method eigh_vals! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided DV as output.
Note that eigh_full and its variants assume that the input matrix is hermitian, or thus symmetric if the input is real. The resulting algorithms exploit this structure, and return eigenvalues that are always real, and eigenvectors that are orthogonal and have the same eltype as the input matrix. If the input matrix does not have this structure, the generic eigenvalue decomposition provided by eig_full and its variants should be used instead.
See also eigh_full(!) and eigh_trunc(!).
MatrixAlgebraKit.findtruncated — Function
MatrixAlgebraKit.findtruncated(values::AbstractVector, strategy::TruncationStrategy)Generic interface for finding truncated values of the spectrum of a decomposition based on the strategy. The output should be a collection of indices specifying which values to keep. MatrixAlgebraKit.findtruncated is used inside of the default implementation of truncate to perform the truncation. It does not assume that the values are sorted. For a version that assumes the values are reverse sorted (which is the standard case for SVD) see MatrixAlgebraKit.findtruncated_svd.
MatrixAlgebraKit.findtruncated_svd — Function
MatrixAlgebraKit.findtruncated_svd(values::AbstractVector, strategy::TruncationStrategy)Like MatrixAlgebraKit.findtruncated but assumes that the values are real and sorted in descending order, as typically obtained by the SVD. This assumption is not checked, and this is used in the default implementation of svd_trunc!.
MatrixAlgebraKit.gaugefix! — Function
gaugefix!(f_eig, V)
gaugefix!(f_svd, U, Vᴴ)Fix the residual gauge degrees of freedom in the eigen or singular vectors, that are obtained from the decomposition performed by f_eig or f_svd. This is achieved by ensuring that the entry with the largest magnitude in V or U is real and positive.
MatrixAlgebraKit.gen_eig_full — Function
gen_eig_full(A, B; kwargs...) -> W, V
gen_eig_full(A, B, alg::AbstractAlgorithm) -> W, V
gen_eig_full!(A, B, [WV]; kwargs...) -> W, V
gen_eig_full!(A, B, [WV], alg::AbstractAlgorithm) -> W, VCompute the full generalized eigenvalue decomposition of the square matrices A and B, such that A * V = B * V * W, where the invertible matrix V contains the generalized eigenvectors and the diagonal matrix W contains the associated generalized eigenvalues.
The bang method gen_eig_full! optionally accepts the output structure and possibly destroys the input matrices A and B. Always use the return value of the function as it may not always be possible to use the provided WV as output.
Note that gen_eig_full and its variants do not assume additional structure on the inputs, and therefore will always return complex eigenvalues and eigenvectors. Real eigenvalues can be expected when both input matrices are hermitian and one of them is positive definite, but specialized methods that exploit this structure are not yet implemented or supported.
See also gen_eig_vals(!).
MatrixAlgebraKit.gen_eig_full! — Function
gen_eig_full(A, B; kwargs...) -> W, V
gen_eig_full(A, B, alg::AbstractAlgorithm) -> W, V
gen_eig_full!(A, B, [WV]; kwargs...) -> W, V
gen_eig_full!(A, B, [WV], alg::AbstractAlgorithm) -> W, VCompute the full generalized eigenvalue decomposition of the square matrices A and B, such that A * V = B * V * W, where the invertible matrix V contains the generalized eigenvectors and the diagonal matrix W contains the associated generalized eigenvalues.
The bang method gen_eig_full! optionally accepts the output structure and possibly destroys the input matrices A and B. Always use the return value of the function as it may not always be possible to use the provided WV as output.
Note that gen_eig_full and its variants do not assume additional structure on the inputs, and therefore will always return complex eigenvalues and eigenvectors. Real eigenvalues can be expected when both input matrices are hermitian and one of them is positive definite, but specialized methods that exploit this structure are not yet implemented or supported.
See also gen_eig_vals(!).
MatrixAlgebraKit.gen_eig_vals — Function
gen_eig_vals(A, B; kwargs...) -> W
gen_eig_vals(A, B, alg::AbstractAlgorithm) -> W
gen_eig_vals!(A, B, [W]; kwargs...) -> W
gen_eig_vals!(A, B, [W], alg::AbstractAlgorithm) -> WCompute the list of generalized eigenvalues of A and B.
The bang method gen_eig_vals! optionally accepts the output structure and possibly destroys the input matrices A and B. Always use the return value of the function as it may not always be possible to use the provided W as output.
Note that gen_eig_full and its variants do not assume additional structure on the inputs, and therefore will always return complex eigenvalues and eigenvectors. Real eigenvalues can be expected when both input matrices are hermitian and one of them is positive definite, but specialized methods that exploit this structure are not yet implemented or supported.
See also gen_eig_full(!).
MatrixAlgebraKit.gen_eig_vals! — Function
gen_eig_vals(A, B; kwargs...) -> W
gen_eig_vals(A, B, alg::AbstractAlgorithm) -> W
gen_eig_vals!(A, B, [W]; kwargs...) -> W
gen_eig_vals!(A, B, [W], alg::AbstractAlgorithm) -> WCompute the list of generalized eigenvalues of A and B.
The bang method gen_eig_vals! optionally accepts the output structure and possibly destroys the input matrices A and B. Always use the return value of the function as it may not always be possible to use the provided W as output.
Note that gen_eig_full and its variants do not assume additional structure on the inputs, and therefore will always return complex eigenvalues and eigenvectors. Real eigenvalues can be expected when both input matrices are hermitian and one of them is positive definite, but specialized methods that exploit this structure are not yet implemented or supported.
See also gen_eig_full(!).
MatrixAlgebraKit.initialize_output — Function
initialize_output(f, A, alg)Whenever possible, allocate the destination for applying a given algorithm in-place. If this is not possible, for example when the output size is not known a priori or immutable, this function may return nothing.
MatrixAlgebraKit.inv_regularized — Function
inv_regularized(a::Number, tol = defaulttol(a))
inv_regularized(A::Matrix, tol = defaulttol(A); isposdef = false, kwargs...)Compute a smooth regularised inverse (L2 Tikhonov regularisation) of a number or square matrix a.
For numbers, this is given by
inv(hypot(a, tol)).For matrices, this is computed using the singular value decomposition and aplying
inv_regularizedto the singular values. Ifisposdef = true, the singular value decomposition is equivalent to the (Hermitian) eigenvalue decomposition ofAand the latter is used instead.
MatrixAlgebraKit.inv_safe — Function
function inv_safe(a::Number, tol = defaulttol(a))Compute the inverse of a number a, but return zero if a is smaller than tol.
MatrixAlgebraKit.is_left_isometric — Function
is_left_isometric(A; isapprox_kwargs...) -> BoolTest whether a linear map is a (left) isometry, i.e. A' * A ≈ I. The isapprox_kwargs can be used to control the tolerances of the equality.
See also isisometric and MatrixAlgebraKit.is_right_isometric.
MatrixAlgebraKit.is_right_isometric — Function
is_right_isometric(A; isapprox_kwargs...) -> BoolTest whether a linear map is a (right) isometry, i.e. A * A' ≈ I. The isapprox_kwargs can be used to control the tolerances of the equality.
See also isisometric and MatrixAlgebraKit.is_left_isometric.
MatrixAlgebraKit.isantihermitian — Method
isantihermitian(A; isapprox_kwargs...)Test whether a linear map is anti-Hermitian, i.e. A = -A'. The isapprox_kwargs can be used to control the tolerances of the equality.
MatrixAlgebraKit.ishermitian — Method
ishermitian(A; isapprox_kwargs...)Test whether a linear map is Hermitian, i.e. A = A'. The isapprox_kwargs can be used to control the tolerances of the equality.
MatrixAlgebraKit.isisometric — Method
isisometric(A; side=:left, isapprox_kwargs...) -> BoolTest whether a linear map is an isometry, where the type of isometry is controlled by kind:
side = :left:A' * A ≈ I.side = :right:A * A' ≈ I.
The isapprox_kwargs are passed on to isapprox to control the tolerances.
New specializations should overload MatrixAlgebraKit.is_left_isometric and MatrixAlgebraKit.is_right_isometric.
See also isunitary.
MatrixAlgebraKit.isunitary — Method
isunitary(A; isapprox_kwargs...)Test whether a linear map is unitary, i.e. A * A' ≈ I ≈ A' * A. The isapprox_kwargs are passed on to isapprox to control the tolerances.
See also isisometric.
MatrixAlgebraKit.iszerotangent — Function
iszerotangent(x)Return true if x is of a type that the different AD engines use to communicate a (co)tangent that is identically zero. By overloading this method, and writing pullback definitions in term of it, we will be able to hook into different AD ecosystems
MatrixAlgebraKit.left_null — Function
left_null(A; [alg], [trunc], kwargs...) -> N
left_null!(A, [N], [alg]; [trunc], kwargs...) -> NCompute an orthonormal basis N for the cokernel of the matrix A, i.e. the nullspace of adjoint(A), such that adjoint(A) * N ≈ 0 and N' * N ≈ I.
This is a high-level wrapper where the keyword arguments can be used to specify and control the underlying orthogonal decomposition that should be used to find the null space of A', whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decomposition. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxnullity::Real: Maximal rank for the truncation
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Here notrunc has special meaning, and signifies keeping the values that correspond to the exact zeros determined from the additional columns of A.
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will be QR-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:qr: Factorize via QR nullspace, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
N = qr_null(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument. It is roughly equivalent to:
U, S, _ = svd_full(A; kwargs...)
N = truncate(left_null, (U, S), trunc)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a LeftNullAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with left_null_alg(alg).
The bang method left_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided N as output.
See also right_null(!), left_orth(!) and right_orth(!).
MatrixAlgebraKit.left_null! — Function
left_null(A; [alg], [trunc], kwargs...) -> N
left_null!(A, [N], [alg]; [trunc], kwargs...) -> NCompute an orthonormal basis N for the cokernel of the matrix A, i.e. the nullspace of adjoint(A), such that adjoint(A) * N ≈ 0 and N' * N ≈ I.
This is a high-level wrapper where the keyword arguments can be used to specify and control the underlying orthogonal decomposition that should be used to find the null space of A', whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decomposition. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxnullity::Real: Maximal rank for the truncation
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Here notrunc has special meaning, and signifies keeping the values that correspond to the exact zeros determined from the additional columns of A.
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will be QR-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:qr: Factorize via QR nullspace, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
N = qr_null(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument. It is roughly equivalent to:
U, S, _ = svd_full(A; kwargs...)
N = truncate(left_null, (U, S), trunc)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a LeftNullAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with left_null_alg(alg).
The bang method left_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided N as output.
See also right_null(!), left_orth(!) and right_orth(!).
MatrixAlgebraKit.left_null_alg — Method
left_null_alg(alg::AbstractAlgorithm) -> LeftNullAlgorithmConvert an algorithm to a LeftNullAlgorithm wrapper for use with left_null.
This function attempts to deduce the appropriate factorization kind (:qr or :svd) from the algorithm type and wraps it in a LeftNullAlgorithm. Custom algorithm types can be registered by defining:
MatrixAlgebraKit.left_null_alg(alg::CustomAlgorithm) = LeftNullAlgorithm{kind}(alg)where kind specifies the factorization backend to use.
See also LeftNullAlgorithm, left_null.
MatrixAlgebraKit.left_orth — Function
left_orth(A; [alg], [trunc], kwargs...) -> V, C
left_orth!(A, [VC], [alg]; [trunc], kwargs...) -> V, CCompute an orthonormal basis V for the image of the matrix A, as well as a matrix C (the corestriction) such that A factors as A = V * C.
This is a high-level wrapper where the keyword arguments can be used to specify and control the specific orthogonal decomposition that should be used to factor A, whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decompositions. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will typically be QR-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:qr: Factorize via QR decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
V, C = qr_compact(A; kwargs...):polar: Factorize via polar decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
V, C = left_polar(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument, and is roughly equivalent to:
V, S, C = svd_trunc(A; trunc, kwargs...)
C = lmul!(S, C)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a LeftOrthAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with left_orth_alg(alg).
The bang method left_orth! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided CV as output.
See also right_orth(!), left_null(!) and right_null(!).
MatrixAlgebraKit.left_orth! — Function
left_orth(A; [alg], [trunc], kwargs...) -> V, C
left_orth!(A, [VC], [alg]; [trunc], kwargs...) -> V, CCompute an orthonormal basis V for the image of the matrix A, as well as a matrix C (the corestriction) such that A factors as A = V * C.
This is a high-level wrapper where the keyword arguments can be used to specify and control the specific orthogonal decomposition that should be used to factor A, whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decompositions. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will typically be QR-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:qr: Factorize via QR decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
V, C = qr_compact(A; kwargs...):polar: Factorize via polar decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
V, C = left_polar(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument, and is roughly equivalent to:
V, S, C = svd_trunc(A; trunc, kwargs...)
C = lmul!(S, C)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a LeftOrthAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with left_orth_alg(alg).
The bang method left_orth! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided CV as output.
See also right_orth(!), left_null(!) and right_null(!).
MatrixAlgebraKit.left_orth_alg — Method
left_orth_alg(alg::AbstractAlgorithm) -> LeftOrthAlgorithmConvert an algorithm to a LeftOrthAlgorithm wrapper for use with left_orth.
This function attempts to deduce the appropriate factorization kind (:qr, :polar, or :svd) from the algorithm type and wraps it in a LeftOrthAlgorithm. Custom algorithm types can be registered by defining:
MatrixAlgebraKit.left_orth_alg(alg::CustomAlgorithm) = LeftOrthAlgorithm{kind}(alg)where kind specifies the factorization backend to use.
See also LeftOrthAlgorithm, left_orth.
MatrixAlgebraKit.left_polar — Function
left_polar(A; kwargs...) -> W, P
left_polar(A, alg::AbstractAlgorithm) -> W, P
left_polar!(A, [WP]; kwargs...) -> W, P
left_polar!(A, [WP], alg::AbstractAlgorithm) -> W, PCompute the full polar decomposition of the rectangular matrix A of size (m, n) with m >= n, such that A = W * P. Here, W is an isometric matrix (orthonormal columns) of size (m, n), whereas P is a positive (semi)definite matrix of size (n, n).
The bang method left_polar! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided WP as output.
See also right_polar(!).
MatrixAlgebraKit.left_polar! — Function
left_polar(A; kwargs...) -> W, P
left_polar(A, alg::AbstractAlgorithm) -> W, P
left_polar!(A, [WP]; kwargs...) -> W, P
left_polar!(A, [WP], alg::AbstractAlgorithm) -> W, PCompute the full polar decomposition of the rectangular matrix A of size (m, n) with m >= n, such that A = W * P. Here, W is an isometric matrix (orthonormal columns) of size (m, n), whereas P is a positive (semi)definite matrix of size (n, n).
The bang method left_polar! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided WP as output.
See also right_polar(!).
MatrixAlgebraKit.left_polar_pullback! — Method
left_polar_pullback!(ΔA, A, WP, ΔWP)Adds the pullback from the left polar decomposition of A to ΔA given the output WP and cotangent ΔWP of left_polar(A).
MatrixAlgebraKit.lq_compact — Function
lq_compact(A; kwargs...) -> L, Q
lq_compact(A, alg::AbstractAlgorithm) -> L, Q
lq_compact!(A, [LQ]; kwargs...) -> L, Q
lq_compact!(A, [LQ], alg::AbstractAlgorithm) -> L, QCompute the compact LQ decomposition of the rectangular matrix A of size (m,n), such that A = L * Q where the matrix Q of size (min(m,n), n) has orthogonal rows spanning the image of A', and the matrix L of size (m, min(m,n)) is lower triangular.
The bang method lq_compact! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided LQ as output.
The compact QR decomposition is equivalent to the full LQ decomposition when m >= n. Some algorithms may require m <= n.
See also lq_full(!).
MatrixAlgebraKit.lq_compact! — Function
lq_compact(A; kwargs...) -> L, Q
lq_compact(A, alg::AbstractAlgorithm) -> L, Q
lq_compact!(A, [LQ]; kwargs...) -> L, Q
lq_compact!(A, [LQ], alg::AbstractAlgorithm) -> L, QCompute the compact LQ decomposition of the rectangular matrix A of size (m,n), such that A = L * Q where the matrix Q of size (min(m,n), n) has orthogonal rows spanning the image of A', and the matrix L of size (m, min(m,n)) is lower triangular.
The bang method lq_compact! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided LQ as output.
The compact QR decomposition is equivalent to the full LQ decomposition when m >= n. Some algorithms may require m <= n.
See also lq_full(!).
MatrixAlgebraKit.lq_full — Function
lq_full(A; kwargs...) -> L, Q
lq_full(A, alg::AbstractAlgorithm) -> L, Q
lq_full!(A, [LQ]; kwargs...) -> L, Q
lq_full!(A, [LQ], alg::AbstractAlgorithm) -> L, QCompute the full LQ decomposition of the rectangular matrix A, such that A = L * Q where Q is a square unitary matrix with the same number of rows as A and L is a lower triangular matrix with the same size as A.
The bang method lq_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided LQ as output.
See also lq_compact(!).
MatrixAlgebraKit.lq_full! — Function
lq_full(A; kwargs...) -> L, Q
lq_full(A, alg::AbstractAlgorithm) -> L, Q
lq_full!(A, [LQ]; kwargs...) -> L, Q
lq_full!(A, [LQ], alg::AbstractAlgorithm) -> L, QCompute the full LQ decomposition of the rectangular matrix A, such that A = L * Q where Q is a square unitary matrix with the same number of rows as A and L is a lower triangular matrix with the same size as A.
The bang method lq_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided LQ as output.
See also lq_compact(!).
MatrixAlgebraKit.lq_null — Function
lq_null(A; kwargs...) -> Nᴴ
lq_null(A, alg::AbstractAlgorithm) -> Nᴴ
lq_null!(A, [Nᴴ]; kwargs...) -> Nᴴ
lq_null!(A, [Nᴴ], alg::AbstractAlgorithm) -> NᴴFor a (m, n) matrix A, compute the matrix Nᴴ corresponding the final n - min(m, n) rows oft the unitary Q factor in the full LQ decomposition of A, i.e. the rows that are not present in the Q factor of the compact LQ decomposition. The matrix Nᴴ is such that the isometric matrix N = adjoint(Nᴴ) contains an orthonormal basis for the kernel (null space) of A as its columns, i.e. A * N = 0 or thus A * adjoint(Nᴴ) = 0.
The bang method lq_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided Nᴴ argument as output.
See also qr_full(!) and qr_compact(!).
MatrixAlgebraKit.lq_null! — Function
lq_null(A; kwargs...) -> Nᴴ
lq_null(A, alg::AbstractAlgorithm) -> Nᴴ
lq_null!(A, [Nᴴ]; kwargs...) -> Nᴴ
lq_null!(A, [Nᴴ], alg::AbstractAlgorithm) -> NᴴFor a (m, n) matrix A, compute the matrix Nᴴ corresponding the final n - min(m, n) rows oft the unitary Q factor in the full LQ decomposition of A, i.e. the rows that are not present in the Q factor of the compact LQ decomposition. The matrix Nᴴ is such that the isometric matrix N = adjoint(Nᴴ) contains an orthonormal basis for the kernel (null space) of A as its columns, i.e. A * N = 0 or thus A * adjoint(Nᴴ) = 0.
The bang method lq_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided Nᴴ argument as output.
See also qr_full(!) and qr_compact(!).
MatrixAlgebraKit.lq_null_pullback! — Method
lq_null_pullback!(
ΔA::AbstractMatrix, A, Nᴴ, ΔNᴴ;
gauge_atol::Real = default_pullback_gauge_atol(ΔNᴴ)
)Adds the pullback from the left nullspace of A to ΔA, given the nullspace basis Nᴴ and its cotangent ΔNᴴ of lq_null(A).
See also lq_pullback!.
MatrixAlgebraKit.lq_pullback! — Method
lq_pullback!(
ΔA, A, LQ, ΔLQ;
rank_atol::Real = default_pullback_rank_atol(LQ[1]),
gauge_atol::Real = default_pullback_gauge_atol(ΔLQ[2])
)Adds the pullback from the LQ decomposition of A to ΔA given the output LQ and cotangent ΔLQ of lq_compact(A; positive = true) or lq_full(A; positive = true).
In the case where the rank r of the original matrix A ≈ L * Q (as determined by rank_atol) is less then the minimum of the number of rows and columns of the cotangents ΔL and ΔQ, only the first r columns of L and the first r rows of Q are well-defined, and also the adjoint variables ΔL and ΔQ should have nonzero values only in the first r columns and rows respectively. If nonzero values in the remaining columns or rows exceed gauge_atol, a warning will be printed.
MatrixAlgebraKit.notrunc — Method
MatrixAlgebraKit.null_truncation_strategy — Method
null_truncation_strategy(; kwargs...)Select a nullspace truncation strategy based on the provided keyword arguments.
Keyword arguments
The following keyword arguments are all optional, and their default value (nothing) will be ignored. It is also allowed to combine multiple of these, in which case the discarded values will consist of the intersection of the different truncated strategies.
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxnullity::Real: Maximal rank for the truncation
MatrixAlgebraKit.project_antihermitian — Function
project_antihermitian(A; kwargs...)
project_antihermitian(A, alg)
project_antihermitian!(A, [Aₐ]; kwargs...)
project_antihermitian!(A, [Aₐ], alg)Compute the anti-hermitian part of a (square) matrix A, defined as (A - A') / 2. For real matrices this corresponds to the antisymmetric part of A. In the bang method, the output storage can be provided via the optional argument Aₐ; by default it is equal toAand so the input matrixA` is replaced by its antihermitian projection.
See also project_hermitian.
MatrixAlgebraKit.project_antihermitian! — Function
project_antihermitian(A; kwargs...)
project_antihermitian(A, alg)
project_antihermitian!(A, [Aₐ]; kwargs...)
project_antihermitian!(A, [Aₐ], alg)Compute the anti-hermitian part of a (square) matrix A, defined as (A - A') / 2. For real matrices this corresponds to the antisymmetric part of A. In the bang method, the output storage can be provided via the optional argument Aₐ; by default it is equal toAand so the input matrixA` is replaced by its antihermitian projection.
See also project_hermitian.
MatrixAlgebraKit.project_hermitian — Function
project_hermitian(A; kwargs...)
project_hermitian(A, alg)
project_hermitian!(A, [Aₕ]; kwargs...)
project_hermitian!(A, [Aₕ], alg)Compute the hermitian part of a (square) matrix A, defined as (A + A') / 2. For real matrices this corresponds to the symmetric part of A. In the bang method, the output storage can be provided via the optional argument Aₕ; by default it is equal to A and so the input matrix A is replaced by its hermitian projection.
See also project_antihermitian.
MatrixAlgebraKit.project_hermitian! — Function
project_hermitian(A; kwargs...)
project_hermitian(A, alg)
project_hermitian!(A, [Aₕ]; kwargs...)
project_hermitian!(A, [Aₕ], alg)Compute the hermitian part of a (square) matrix A, defined as (A + A') / 2. For real matrices this corresponds to the symmetric part of A. In the bang method, the output storage can be provided via the optional argument Aₕ; by default it is equal to A and so the input matrix A is replaced by its hermitian projection.
See also project_antihermitian.
MatrixAlgebraKit.project_isometric — Function
project_isometric(A; kwargs...)
project_isometric(A, alg)
project_isometric!(A, [W]; kwargs...)
project_isometric!(A, [W], alg)Compute the projection of A onto the manifold of isometric matrices, i.e. matrices satisfying A' * A ≈ I. This projection is computed via the polar decomposition, i.e. W corresponds to the first return value of left_polar!, but avoids computing the positive definite factor explicitly.
The bang method project_isometric! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided W as output.
MatrixAlgebraKit.project_isometric! — Function
project_isometric(A; kwargs...)
project_isometric(A, alg)
project_isometric!(A, [W]; kwargs...)
project_isometric!(A, [W], alg)Compute the projection of A onto the manifold of isometric matrices, i.e. matrices satisfying A' * A ≈ I. This projection is computed via the polar decomposition, i.e. W corresponds to the first return value of left_polar!, but avoids computing the positive definite factor explicitly.
The bang method project_isometric! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided W as output.
MatrixAlgebraKit.qr_compact — Function
qr_compact(A; kwargs...) -> Q, R
qr_compact(A, alg::AbstractAlgorithm) -> Q, R
qr_compact!(A, [QR]; kwargs...) -> Q, R
qr_compact!(A, [QR], alg::AbstractAlgorithm) -> Q, RCompute the compact QR decomposition of the rectangular matrix A of size (m,n), such that A = Q * R where the isometric matrix Q of size (m, min(m,n)) has orthogonal columns spanning the image of A, and the matrix R of size (min(m,n), n) is upper triangular.
The bang method qr_compact! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided QR as output.
The compact QR decomposition is equivalent to the full QR decomposition when m >= n. Some algorithms may require m >= n.
See also qr_full(!).
MatrixAlgebraKit.qr_compact! — Function
qr_compact(A; kwargs...) -> Q, R
qr_compact(A, alg::AbstractAlgorithm) -> Q, R
qr_compact!(A, [QR]; kwargs...) -> Q, R
qr_compact!(A, [QR], alg::AbstractAlgorithm) -> Q, RCompute the compact QR decomposition of the rectangular matrix A of size (m,n), such that A = Q * R where the isometric matrix Q of size (m, min(m,n)) has orthogonal columns spanning the image of A, and the matrix R of size (min(m,n), n) is upper triangular.
The bang method qr_compact! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided QR as output.
The compact QR decomposition is equivalent to the full QR decomposition when m >= n. Some algorithms may require m >= n.
See also qr_full(!).
MatrixAlgebraKit.qr_full — Function
qr_full(A; kwargs...) -> Q, R
qr_full(A, alg::AbstractAlgorithm) -> Q, R
qr_full!(A, [QR]; kwargs...) -> Q, R
qr_full!(A, [QR], alg::AbstractAlgorithm) -> Q, RCompute the full QR decomposition of the rectangular matrix A, such that A = Q * R where Q is a square unitary matrix with the same number of rows as A and R is an upper triangular matrix with the same size as A.
The bang method qr_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided QR as output.
See also qr_compact(!).
MatrixAlgebraKit.qr_full! — Function
qr_full(A; kwargs...) -> Q, R
qr_full(A, alg::AbstractAlgorithm) -> Q, R
qr_full!(A, [QR]; kwargs...) -> Q, R
qr_full!(A, [QR], alg::AbstractAlgorithm) -> Q, RCompute the full QR decomposition of the rectangular matrix A, such that A = Q * R where Q is a square unitary matrix with the same number of rows as A and R is an upper triangular matrix with the same size as A.
The bang method qr_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided QR as output.
See also qr_compact(!).
MatrixAlgebraKit.qr_null — Function
qr_null(A; kwargs...) -> N
qr_null(A, alg::AbstractAlgorithm) -> N
qr_null!(A, [N]; kwargs...) -> N
qr_null!(A, [N], alg::AbstractAlgorithm) -> NFor a (m, n) matrix A, compute the matrix N corresponding the final m - min(m, n) columns of the unitary Q factor in the full QR decomposition of A, i.e. the columns that are not present in the Q factor of the compact QR decomposition. The isometric matrix N contains an orthonormal basis for the cokernel of A as its columns, i.e. adjoint(A) * N = 0.
The bang method qr_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided N argument as output.
See also lq_full(!) and lq_compact(!).
MatrixAlgebraKit.qr_null! — Function
qr_null(A; kwargs...) -> N
qr_null(A, alg::AbstractAlgorithm) -> N
qr_null!(A, [N]; kwargs...) -> N
qr_null!(A, [N], alg::AbstractAlgorithm) -> NFor a (m, n) matrix A, compute the matrix N corresponding the final m - min(m, n) columns of the unitary Q factor in the full QR decomposition of A, i.e. the columns that are not present in the Q factor of the compact QR decomposition. The isometric matrix N contains an orthonormal basis for the cokernel of A as its columns, i.e. adjoint(A) * N = 0.
The bang method qr_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided N argument as output.
See also lq_full(!) and lq_compact(!).
MatrixAlgebraKit.qr_null_pullback! — Method
qr_null_pullback!(
ΔA::AbstractMatrix, A, N, ΔN;
gauge_atol::Real = default_pullback_gauge_atol(ΔN)
)Adds the pullback from the right nullspace of A to ΔA, given the nullspace basis N and its cotangent ΔN of qr_null(A).
See also qr_pullback!.
MatrixAlgebraKit.qr_pullback! — Method
qr_pullback!(
ΔA, A, QR, ΔQR;
tol::Real = default_pullback_gaugetol(QR[2]),
rank_atol::Real = default_pullback_rank_atol(QR[2]),
gauge_atol::Real = default_pullback_gauge_atol(ΔQR[1])
)Adds the pullback from the QR decomposition of A to ΔA given the output QR and cotangent ΔQR of qr_compact(A; positive = true) or qr_full(A; positive = true).
In the case where the rank r of the original matrix A ≈ Q * R (as determined by rank_atol) is less then the minimum of the number of rows and columns, the cotangents ΔQ and ΔR, only the first r columns of Q and the first r rows of R are well-defined, and also the adjoint variables ΔQ and ΔR should have nonzero values only in the first r columns and rows respectively. If nonzero values in the remaining columns or rows exceed gauge_atol, a warning will be printed.
MatrixAlgebraKit.right_null — Function
right_null(A; [alg], [trunc], kwargs...) -> Nᴴ
right_null!(A, [Nᴴ], [alg]; [trunc], kwargs...) -> NᴴCompute an orthonormal basis N = adjoint(Nᴴ) for the kernel of the matrix A, i.e. the nullspace of A, such that A * Nᴴ' ≈ 0 and Nᴴ * Nᴴ' ≈ I.
This is a high-level wrapper where the keyword arguments can be used to specify and control the underlying orthogonal decomposition that should be used to find the null space of A, whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decomposition. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxnullity::Real: Maximal rank for the truncation
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Here notrunc has special meaning, and signifies keeping the values that correspond to the exact zeros determined from the additional rows of A.
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will be LQ-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:lq: Factorize via LQ nullspace, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
Nᴴ = lq_null(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument. It is roughly equivalent to:
_, S, Vᴴ = svd_full(A; kwargs...)
Nᴴ = truncate(right_null, (S, Vᴴ), trunc)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a RightNullAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with right_null_alg(alg).
The bang method right_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided Nᴴ as output.
See also left_null(!), left_orth(!) and right_orth(!).
MatrixAlgebraKit.right_null! — Function
right_null(A; [alg], [trunc], kwargs...) -> Nᴴ
right_null!(A, [Nᴴ], [alg]; [trunc], kwargs...) -> NᴴCompute an orthonormal basis N = adjoint(Nᴴ) for the kernel of the matrix A, i.e. the nullspace of A, such that A * Nᴴ' ≈ 0 and Nᴴ * Nᴴ' ≈ I.
This is a high-level wrapper where the keyword arguments can be used to specify and control the underlying orthogonal decomposition that should be used to find the null space of A, whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decomposition. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxnullity::Real: Maximal rank for the truncation
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Here notrunc has special meaning, and signifies keeping the values that correspond to the exact zeros determined from the additional rows of A.
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will be LQ-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:lq: Factorize via LQ nullspace, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
Nᴴ = lq_null(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument. It is roughly equivalent to:
_, S, Vᴴ = svd_full(A; kwargs...)
Nᴴ = truncate(right_null, (S, Vᴴ), trunc)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a RightNullAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with right_null_alg(alg).
The bang method right_null! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided Nᴴ as output.
See also left_null(!), left_orth(!) and right_orth(!).
MatrixAlgebraKit.right_null_alg — Method
right_null_alg(alg::AbstractAlgorithm) -> RightNullAlgorithmConvert an algorithm to a RightNullAlgorithm wrapper for use with right_null.
This function attempts to deduce the appropriate factorization kind (:lq or :svd) from the algorithm type and wraps it in a RightNullAlgorithm. Custom algorithm types can be registered by defining:
MatrixAlgebraKit.right_null_alg(alg::CustomAlgorithm) = RightNullAlgorithm{kind}(alg)where kind specifies the factorization backend to use.
See also RightNullAlgorithm, right_null.
MatrixAlgebraKit.right_orth — Function
right_orth(A; [alg], [trunc], kwargs...) -> C, Vᴴ
right_orth!(A, [CVᴴ], [alg]; [trunc], kwargs...) -> C, VᴴCompute an orthonormal basis V = adjoint(Vᴴ) for the coimage of the matrix A, i.e. for the image of adjoint(A), as well as a matrix C such that A factors as A = C * Vᴴ.
This is a high-level wrapper where the keyword arguments can be used to specify and control the specific orthogonal decomposition that should be used to factor A, whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decompositions. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will typicall be LQ-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:lq: Factorize via LQ decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
C, Vᴴ = lq_compact(A; kwargs...):polar: Factorize via polar decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
C, Vᴴ = right_polar(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument, and is roughly equivalent to:
C, S, Vᴴ = svd_trunc(A; trunc, kwargs...)
C = rmul!(C, S)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a RightOrthAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with right_orth_alg.
The bang method right_orth! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided CVᴴ as output.
See also left_orth(!), left_null(!) and right_null(!).
MatrixAlgebraKit.right_orth! — Function
right_orth(A; [alg], [trunc], kwargs...) -> C, Vᴴ
right_orth!(A, [CVᴴ], [alg]; [trunc], kwargs...) -> C, VᴴCompute an orthonormal basis V = adjoint(Vᴴ) for the coimage of the matrix A, i.e. for the image of adjoint(A), as well as a matrix C such that A factors as A = C * Vᴴ.
This is a high-level wrapper where the keyword arguments can be used to specify and control the specific orthogonal decomposition that should be used to factor A, whereas trunc can optionally be used to control the precision in determining the rank of A, typically via its singular values.
Truncation
The optional truncation strategy can be controlled via the trunc keyword argument, and any non-trivial strategy typically requires an SVD-based decompositions. This keyword can be either a NamedTuple or a TruncationStrategy.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
There are 3 major modes of operation, based on the alg keyword, with slightly different application purposes.
alg::Nothing
This default mode uses the presence of a truncation strategy trunc to determine an optimal decomposition type, which will typicall be LQ-based for no truncation, or SVD-based for truncation. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type.
alg::Symbol
Here, the driving selector is alg, which is used to select the kind of decomposition. The remaining keyword arguments are passed on directly to the algorithm selection procedure of the chosen decomposition type. By default, the supported kinds are:
:lq: Factorize via LQ decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
C, Vᴴ = lq_compact(A; kwargs...):polar: Factorize via polar decomposition, with further customizations through the other keywords. This mode requiresisnothing(trunc), and is equivalent to
C, Vᴴ = right_polar(A; kwargs...):svd: Factorize via SVD, with further customizations through the other keywords. This mode further allows truncation, which can be selected through thetruncargument, and is roughly equivalent to:
C, S, Vᴴ = svd_trunc(A; trunc, kwargs...)
C = rmul!(C, S)alg::AbstractAlgorithm
In this expert mode the algorithm is supplied directly, and the kind of decomposition is deduced from that. This is achieved either directly by providing a RightOrthAlgorithm{kind}, or automatically by attempting to deduce the decomposition kind with right_orth_alg.
The bang method right_orth! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided CVᴴ as output.
See also left_orth(!), left_null(!) and right_null(!).
MatrixAlgebraKit.right_orth_alg — Method
right_orth_alg(alg::AbstractAlgorithm) -> RightOrthAlgorithmConvert an algorithm to a RightOrthAlgorithm wrapper for use with right_orth.
This function attempts to deduce the appropriate factorization kind (:lq, :polar, or :svd) from the algorithm type and wraps it in a RightOrthAlgorithm. Custom algorithm types can be registered by defining:
MatrixAlgebraKit.right_orth_alg(alg::CustomAlgorithm) = RightOrthAlgorithm{kind}(alg)where kind specifies the factorization backend to use.
See also RightOrthAlgorithm, right_orth.
MatrixAlgebraKit.right_polar — Function
right_polar(A; kwargs...) -> P, Wᴴ
right_polar(A, alg::AbstractAlgorithm) -> P, Wᴴ
right_polar!(A, [PWᴴ]; kwargs...) -> P, Wᴴ
right_polar!(A, [PWᴴ], alg::AbstractAlgorithm) -> P, WᴴCompute the full polar decomposition of the rectangular matrix A of size (m, n) with n >= m, such that A = P * Wᴴ. Here, P is a positive (semi)definite matrix of size (m, m), whereas Wᴴ is a matrix with orthonormal rows (its adjoint is isometric) of size (n, m).
The bang method right_polar! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided WP as output.
See also left_polar(!).
MatrixAlgebraKit.right_polar! — Function
right_polar(A; kwargs...) -> P, Wᴴ
right_polar(A, alg::AbstractAlgorithm) -> P, Wᴴ
right_polar!(A, [PWᴴ]; kwargs...) -> P, Wᴴ
right_polar!(A, [PWᴴ], alg::AbstractAlgorithm) -> P, WᴴCompute the full polar decomposition of the rectangular matrix A of size (m, n) with n >= m, such that A = P * Wᴴ. Here, P is a positive (semi)definite matrix of size (m, m), whereas Wᴴ is a matrix with orthonormal rows (its adjoint is isometric) of size (n, m).
The bang method right_polar! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided WP as output.
See also left_polar(!).
MatrixAlgebraKit.right_polar_pullback! — Method
right_polar_pullback!(ΔA, A, PWᴴ, ΔPWᴴ)Adds the pullback from the left polar decomposition of A to ΔA given the output PWᴴ and cotangent ΔPWᴴ of right_polar(A).
MatrixAlgebraKit.schur_full — Function
schur_full(A; kwargs...) -> T, Z, vals
schur_full(A, alg::AbstractAlgorithm) -> T, Z, vals
schur_full!(A, [TZv]; kwargs...) -> T, Z, vals
schur_full!(A, [TZv], alg::AbstractAlgorithm) -> T, Z, valsCompute the full Schur decomposition of the square matrix A, such that A * Z = Z * T, where the orthogonal or unitary matrix Z contains the Schur vectors and the square matrix T is upper triangular (in the complex case) or quasi-upper triangular (in the real case). The list vals contains the (complex-valued) eigenvalues of A, as extracted from the (quasi-)diagonal of T.
The bang method schur_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided TZv as output.
MatrixAlgebraKit.schur_full! — Function
schur_full(A; kwargs...) -> T, Z, vals
schur_full(A, alg::AbstractAlgorithm) -> T, Z, vals
schur_full!(A, [TZv]; kwargs...) -> T, Z, vals
schur_full!(A, [TZv], alg::AbstractAlgorithm) -> T, Z, valsCompute the full Schur decomposition of the square matrix A, such that A * Z = Z * T, where the orthogonal or unitary matrix Z contains the Schur vectors and the square matrix T is upper triangular (in the complex case) or quasi-upper triangular (in the real case). The list vals contains the (complex-valued) eigenvalues of A, as extracted from the (quasi-)diagonal of T.
The bang method schur_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided TZv as output.
MatrixAlgebraKit.schur_vals — Function
schur_vals(A; kwargs...) -> vals
schur_vals(A, alg::AbstractAlgorithm) -> vals
schur_vals!(A, [vals]; kwargs...) -> vals
schur_vals!(A, [vals], alg::AbstractAlgorithm) -> valsCompute the list of eigenvalues of A by computing the Schur decomposition of A.
The bang method schur_vals! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided vals as output.
See also eig_full(!) and eig_trunc(!).
MatrixAlgebraKit.schur_vals! — Function
schur_vals(A; kwargs...) -> vals
schur_vals(A, alg::AbstractAlgorithm) -> vals
schur_vals!(A, [vals]; kwargs...) -> vals
schur_vals!(A, [vals], alg::AbstractAlgorithm) -> valsCompute the list of eigenvalues of A by computing the Schur decomposition of A.
The bang method schur_vals! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided vals as output.
See also eig_full(!) and eig_trunc(!).
MatrixAlgebraKit.select_algorithm — Function
MatrixAlgebraKit.select_algorithm(f, A, alg::AbstractAlgorithm)
MatrixAlgebraKit.select_algorithm(f, A, alg::Symbol; kwargs...)
MatrixAlgebraKit.select_algorithm(f, A, alg::Type; kwargs...)
MatrixAlgebraKit.select_algorithm(f, A; kwargs...)
MatrixAlgebraKit.select_algorithm(f, A, (; kwargs...))Decide on an algorithm to use for implementing the function f on inputs of type A. This can be obtained both for values A or types A.
If alg is an AbstractAlgorithm instance, it will be returned as-is.
If alg is a Symbol or a Type of algorithm, the return value is obtained by calling the corresponding algorithm constructor; keyword arguments in kwargs are passed along to this constructor.
If alg is not specified (or nothing), an algorithm will be selected automatically with MatrixAlgebraKit.default_algorithm and the keyword arguments in kwargs will be passed to the algorithm constructor. Finally, the same behavior is obtained when the keyword arguments are passed as the third positional argument in the form of a NamedTuple.
MatrixAlgebraKit.select_null_truncation — Function
MatrixAlgebraKit.select_null_truncation(trunc)Construct a TruncationStrategy from the given NamedTuple of keywords or input strategy, to implement a nullspace selection.
MatrixAlgebraKit.select_truncation — Function
MatrixAlgebraKit.select_truncation(trunc)Construct a TruncationStrategy from the given NamedTuple of keywords or input strategy.
MatrixAlgebraKit.sign_safe — Method
sign_safe(s::Number)Compute the sign of a number s, but return +1 if s is zero so that the result is always a number with modulus 1, i.e. an element of the unitary group U(1).
MatrixAlgebraKit.svd_compact — Function
svd_compact(A; kwargs...) -> U, S, Vᴴ
svd_compact(A, alg::AbstractAlgorithm) -> U, S, Vᴴ
svd_compact!(A, [USVᴴ]; kwargs...) -> U, S, Vᴴ
svd_compact!(A, [USVᴴ], alg::AbstractAlgorithm) -> U, S, VᴴCompute the compact singular value decomposition (SVD) of the rectangular matrix A of size (m, n), such that A = U * S * Vᴴ. Here, U is an isometric matrix (orthonormal columns) of size (m, k), whereas Vᴴ is a matrix of size (k, n) with orthonormal rows and S is a square diagonal matrix of size (k, k), with k = min(m, n).
The bang method svd_compact! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided USVᴴ as output.
See also svd_full(!), svd_vals(!) and svd_trunc(!).
MatrixAlgebraKit.svd_compact! — Function
svd_compact(A; kwargs...) -> U, S, Vᴴ
svd_compact(A, alg::AbstractAlgorithm) -> U, S, Vᴴ
svd_compact!(A, [USVᴴ]; kwargs...) -> U, S, Vᴴ
svd_compact!(A, [USVᴴ], alg::AbstractAlgorithm) -> U, S, VᴴCompute the compact singular value decomposition (SVD) of the rectangular matrix A of size (m, n), such that A = U * S * Vᴴ. Here, U is an isometric matrix (orthonormal columns) of size (m, k), whereas Vᴴ is a matrix of size (k, n) with orthonormal rows and S is a square diagonal matrix of size (k, k), with k = min(m, n).
The bang method svd_compact! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided USVᴴ as output.
See also svd_full(!), svd_vals(!) and svd_trunc(!).
MatrixAlgebraKit.svd_full — Function
svd_full(A; kwargs...) -> U, S, Vᴴ
svd_full(A, alg::AbstractAlgorithm) -> U, S, Vᴴ
svd_full!(A, [USVᴴ]; kwargs...) -> U, S, Vᴴ
svd_full!(A, [USVᴴ], alg::AbstractAlgorithm) -> U, S, VᴴCompute the full singular value decomposition (SVD) of the rectangular matrix A of size (m, n), such that A = U * S * Vᴴ. Here, U and Vᴴ are unitary matrices of size (m, m) and (n, n) respectively, and S is a diagonal matrix of size (m, n).
The bang method svd_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided USVᴴ as output.
See also svd_compact(!), svd_vals(!) and svd_trunc(!).
MatrixAlgebraKit.svd_full! — Function
svd_full(A; kwargs...) -> U, S, Vᴴ
svd_full(A, alg::AbstractAlgorithm) -> U, S, Vᴴ
svd_full!(A, [USVᴴ]; kwargs...) -> U, S, Vᴴ
svd_full!(A, [USVᴴ], alg::AbstractAlgorithm) -> U, S, VᴴCompute the full singular value decomposition (SVD) of the rectangular matrix A of size (m, n), such that A = U * S * Vᴴ. Here, U and Vᴴ are unitary matrices of size (m, m) and (n, n) respectively, and S is a diagonal matrix of size (m, n).
The bang method svd_full! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided USVᴴ as output.
See also svd_compact(!), svd_vals(!) and svd_trunc(!).
MatrixAlgebraKit.svd_pullback! — Function
svd_pullback!(
ΔA, A, USVᴴ, ΔUSVᴴ, [ind];
rank_atol::Real = default_pullback_rank_atol(USVᴴ[2]),
degeneracy_atol::Real = default_pullback_rank_atol(USVᴴ[2]),
gauge_atol::Real = default_pullback_gauge_atol(ΔUSVᴴ[1], ΔUSVᴴ[3])
)Adds the pullback from the SVD of A to ΔA given the output USVᴴ of svd_compact or svd_full and the cotangent ΔUSVᴴ of svd_compact, svd_full or svd_trunc.
In particular, it is assumed that A ≈ U * S * Vᴴ, or thus, that no singular values with magnitude less than rank_atol are missing from S. For the cotangents, an arbitrary number of singular vectors or singular values can be missing, i.e. for a matrix A with size (m, n), ΔU and ΔVᴴ can have sizes (m, pU) and (pV, n) respectively, whereas diagview(ΔS) can have length pS. In those cases, additionally ind is required to specify which singular vectors and values are present in ΔU, ΔS and ΔVᴴ.
A warning will be printed if the cotangents are not gauge-invariant, i.e. if the anti-hermitian part of U' * ΔU + Vᴴ * ΔVᴴ', restricted to rows i and columns j for which abs(S[i] - S[j]) < degeneracy_atol, is not small compared to gauge_atol.
MatrixAlgebraKit.svd_trunc — Function
svd_trunc(A; [trunc], kwargs...) -> U, S, Vᴴ, ϵ
svd_trunc(A, alg::AbstractAlgorithm) -> U, S, Vᴴ, ϵ
svd_trunc!(A, [USVᴴ]; [trunc], kwargs...) -> U, S, Vᴴ, ϵ
svd_trunc!(A, [USVᴴ], alg::AbstractAlgorithm) -> U, S, Vᴴ, ϵCompute a partial or truncated singular value decomposition (SVD) of A, such that A * (Vᴴ)' ≈ U * S. Here, U is an isometric matrix (orthonormal columns) of size (m, k), whereas Vᴴ is a matrix of size (k, n) with orthonormal rows and S is a square diagonal matrix of size (k, k), with k is set by the truncation strategy.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded singular values.
Truncation
The truncation strategy can be controlled via the trunc keyword argument. This can be either a NamedTuple or a TruncationStrategy. If trunc is not provided or nothing, all values will be kept.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
Other keyword arguments are passed to the algorithm selection procedure. If no explicit alg is provided, these keywords are used to select and configure the algorithm through MatrixAlgebraKit.select_algorithm. The remaining keywords after algorithm selection are passed to the algorithm constructor. See MatrixAlgebraKit.default_algorithm for the default algorithm selection behavior.
When alg is a TruncatedAlgorithm, the trunc keyword cannot be specified as the truncation strategy is already embedded in the algorithm.
The bang method svd_trunc! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided USVᴴ as output.
See also svd_full(!), svd_compact(!), svd_vals(!), and Truncations for more information on truncation strategies.
MatrixAlgebraKit.svd_trunc! — Function
svd_trunc(A; [trunc], kwargs...) -> U, S, Vᴴ, ϵ
svd_trunc(A, alg::AbstractAlgorithm) -> U, S, Vᴴ, ϵ
svd_trunc!(A, [USVᴴ]; [trunc], kwargs...) -> U, S, Vᴴ, ϵ
svd_trunc!(A, [USVᴴ], alg::AbstractAlgorithm) -> U, S, Vᴴ, ϵCompute a partial or truncated singular value decomposition (SVD) of A, such that A * (Vᴴ)' ≈ U * S. Here, U is an isometric matrix (orthonormal columns) of size (m, k), whereas Vᴴ is a matrix of size (k, n) with orthonormal rows and S is a square diagonal matrix of size (k, k), with k is set by the truncation strategy.
The function also returns ϵ, the truncation error defined as the 2-norm of the discarded singular values.
Truncation
The truncation strategy can be controlled via the trunc keyword argument. This can be either a NamedTuple or a TruncationStrategy. If trunc is not provided or nothing, all values will be kept.
trunc::NamedTuple
The supported truncation keyword arguments are:
atol::Real: Absolute tolerance for the truncationrtol::Real: Relative tolerance for the truncationmaxrank::Real: Maximal rank for the truncationmaxerror::Real: Maximal truncation error.filter: Custom filter to select truncated values.
trunc::TruncationStrategy
For more control, a truncation strategy can be supplied directly. By default, MatrixAlgebraKit supplies the following:
Keyword arguments
Other keyword arguments are passed to the algorithm selection procedure. If no explicit alg is provided, these keywords are used to select and configure the algorithm through MatrixAlgebraKit.select_algorithm. The remaining keywords after algorithm selection are passed to the algorithm constructor. See MatrixAlgebraKit.default_algorithm for the default algorithm selection behavior.
When alg is a TruncatedAlgorithm, the trunc keyword cannot be specified as the truncation strategy is already embedded in the algorithm.
The bang method svd_trunc! optionally accepts the output structure and possibly destroys the input matrix A. Always use the return value of the function as it may not always be possible to use the provided USVᴴ as output.
See also svd_full(!), svd_compact(!), svd_vals(!), and Truncations for more information on truncation strategies.
MatrixAlgebraKit.svd_trunc_pullback! — Method
svd_trunc_pullback!(
ΔA, A, USVᴴ, ΔUSVᴴ;
rank_atol::Real = default_pullback_rank_atol(USVᴴ[2]),
degeneracy_atol::Real = default_pullback_rank_atol(USVᴴ[2]),
gauge_atol::Real = default_pullback_gauge_atol(ΔUSVᴴ[1], ΔUSVᴴ[3])
)Adds the pullback from the truncated SVD of A to ΔA, given the output USVᴴ and the cotangent ΔUSVᴴ of svd_trunc.
In particular, it is assumed that A * Vᴴ' ≈ U * S and U' * A = S * Vᴴ, with U and Vᴴ rectangular matrices of left and right singular vectors, and S diagonal. For the cotangents, it is assumed that if ΔU and ΔVᴴ are not zero, then they have the same size as U and Vᴴ (respectively), and if ΔS is not zero, then it is a diagonal matrix of the same size as S. For this method to work correctly, it is also assumed that the remaining singular values (not included in S) are (sufficiently) separated from those in S.
A warning will be printed if the cotangents are not gauge-invariant, i.e. if the anti-hermitian part of U' * ΔU + Vᴴ * ΔVᴴ', restricted to rows i and columns j for which abs(S[i] - S[j]) < degeneracy_atol, is not small compared to gauge_atol.
MatrixAlgebraKit.svd_vals — Function
svd_vals(A; kwargs...) -> S
svd_vals(A, alg::AbstractAlgorithm) -> S
svd_vals!(A, [S]; kwargs...) -> S
svd_vals!(A, [S], alg::AbstractAlgorithm) -> SCompute the vector of singular values of A, such that for an M×N matrix A, S is a vector of size K = min(M, N), the number of kept singular values.
See also svd_full(!), svd_compact(!) and svd_trunc(!).
MatrixAlgebraKit.svd_vals! — Function
svd_vals(A; kwargs...) -> S
svd_vals(A, alg::AbstractAlgorithm) -> S
svd_vals!(A, [S]; kwargs...) -> S
svd_vals!(A, [S], alg::AbstractAlgorithm) -> SCompute the vector of singular values of A, such that for an M×N matrix A, S is a vector of size K = min(M, N), the number of kept singular values.
See also svd_full(!), svd_compact(!) and svd_trunc(!).
MatrixAlgebraKit.truncate — Function
truncate(::typeof(f), F, strategy::TruncationStrategy) -> F′, indGiven a factorization function f and truncation strategy, truncate the factors F such that the rows or columns at the indices ind are kept.
See also findtruncated and findtruncated_svd for determining the indices.
MatrixAlgebraKit.truncation_error — Function
truncation_error(values, ind)Compute the truncation error as the 2-norm of the values that are not kept by ind.
MatrixAlgebraKit.truncation_error! — Function
truncation_error(values, ind)Compute the truncation error as the 2-norm of the values that are not kept by ind.
MatrixAlgebraKit.truncerror — Method
truncerror(; atol::Real=0, rtol::Real=0, p::Real=2)Truncation strategy for truncating values such that the error in the factorization is smaller than max(atol, rtol * norm), where the error is determined using the p-norm.
MatrixAlgebraKit.truncfilter — Method
MatrixAlgebraKit.truncrank — Method
truncrank(howmany::Integer; by=abs, rev::Bool=true)Truncation strategy to keep the first howmany values when sorted according to by or the last howmany if rev is true.
MatrixAlgebraKit.trunctol — Method
trunctol(; atol::Real=0, rtol::Real=0, p::Real=2, by=abs, keep_below::Bool=false)Truncation strategy to keep the values that satisfy by(val) > max(atol, rtol * norm(values, p). If keep_below = true, discard these values instead.
Other
MatrixAlgebraKit.@algdef — Macro
@algdef AlgorithmNameConvenience macro to define an algorithm AlgorithmName that accepts generic keywords. This defines an exported alias for Algorithm{:AlgorithmName} along with some utility methods.
MatrixAlgebraKit.@check_scalar — Macro
@check_scalar(x, y, [op], [eltype])Check if eltype(x) == op(eltype(y)) and throw an error if not. By default op = identity and `eltype = eltype'.
MatrixAlgebraKit.@check_size — Macro
@check_size(x, sz, [size])Check if size(x) == sz and throw an error if not. By default, size = size.
MatrixAlgebraKit.@functiondef — Macro
@functiondef [n_args=1] fConvenience macro to define the boilerplate code that dispatches between several versions of f and f!. By default, f accepts a single argument A. This enables the following signatures to be defined in terms of the final f!(A, out, alg::Algorithm).
f(A; kwargs...)
f(A, alg::Algorithm)
f!(A, [out]; kwargs...)
f!(A, alg::Algorithm)The number of inputs can be set with the n_args keyword argument, so that
@functiondef n_args=2 fwould create
f(A, B; kwargs...)
f(A, B, alg::Algorithm)
f!(A, B, [out]; kwargs...)
f!(A, B, alg::Algorithm)See also copy_input, select_algorithm and initialize_output.