Detailed Docs

AJD.addrandomnoiseMethod
addrandomnoise(
    A::Vector{M};
    σ::AbstractFloat = 0.5,
    same_noise::Bool = true,
) where {T<:Number,M<:AbstractMatrix{T}}

Add ranom noise to a vector of matrices A. If same noise selected same random matrix R used for all (!) matrices. Completley random noise added to each matrix in A with same_noise=false.

source
AJD.check_inputMethod
check_input(
    A::Vector{<:AbstractMatrix{<:Number}},
    max_iter::Int,
    threshold::AbstractFloat,
    )

Check input of diagonalize(). Validate matrices, threshold and max iteration as selected.

source
AJD.diagonalizeMethod
(1)
    diagonalize(
        M::Vector{<:AbstractMatrix{<:Number}};
        algorithm::AbstractDiagonalization = JDiagEdourdPineau(),
        max_iter::Int = 1000,
        threshold::AbstractFloat = eps())


    (2) diagonalize(
        M::Vector{<:AbstractMatrix{<:Number}},
        only_plot::Symbol;
        algorithm::AbstractDiagonalization = JDiagGabrielDernbach(),
        max_iter::Int = 1000,
        threshold::AbstractFloat = eps())


    (3) diagonalize(
        benchmark::Symbol,
        n_dims::Int,
        n_matrices::Int)

Calculate joint diagonalization of multiple input matrices $M_k$.

Main function of the AJD package. Implemented algorithms are JDiag and FFDiag. Input of matrices $M_k$ need to be a vector of matrices. The matrices can be of types Float64 or Complex.

See the Getting Started Guide for information on the algorithms.

Dispatch (1)

Inputs:

  • M: Vector of matrices (requiered).
  • algorithm: Selected algorithm from JDiagGabrielDernbach(), JDiagEdourdPineau(), JDiagCardoso() or FFDiag().
  • max_iter: Maximum iteration step as integer.
  • threshold: Desired threshold minimizing the off-diagonal elements.

Output:

  • Return LinearFilter object. Filter fil contains filter matrix fil.F and the inverse fil.iF.

Dispatch (2) - Benchmark Extension to AJD.jl

Inputs:

  • Additional symbol used to generate overview plot of the result. Use diagonalize(M, :plot).

Output:

  • Overview plot. Shows heatmap of the filter matrix, heatmap of the mean of the diagonalized matrices and the vonvergence behaviour of the algorithm.

Dispatch (3) - Benchmark Extension to AJD.jl

Inputs:

  • Symbol :benchmark used as diagonalize(:benchmark, 10, 10).
  • Automatic benchmark is run comparing JDiag and FFDiag algorithms. Using n_dims $\times$ n_dims matrices of count n_matrices.

Output:

  • BenchmarkGroup of package BenchmarkTools comparing the algorithms using diffrent types of test data.
source
AJD.ffdMethod
function ffd(A::Vector{M}; <keyword_arguments>) where {T <: Number, M<:AbstractMatrix{T}}
  • A: Vector of matrices of dimension $n × n × k$

  • threshold: absolute threshold until calculation should stop. default is eps().

  • rel_threshold`: relative threshold for stopping calculation. default is 1e-3.

  • max_iter: max number of iterations. default is 100.

  • norm_: norm by which the update matrix is divided by. can either be :frobenius or :inf. default is :frobenius

  • θ: normation criterion. should be smaller than 1. default is 0.99. See paper on ffdiag for more information.

*plot_convergence: whether the convergence plot should be shown or not. default is false.

  • initial_guess: is the matrix by which the matrices can be diagonalized. if close to the solution the calculation gets faster.

default is identity matrix with size $n × n$.

Calculates the diagonalization of a set of matrices proposed in 2.

Might be faster than the jade algorithms for certain matrices if initial guess is close to the diagonalized solution.

Won't work if input A is of size $n × n × 1$ since this will lead to NaN values due to how the update matrix is calculated or if all of the entries of A are the same or have the same symmetric offdiagonal elements.

Will lead to a warning and stops the calculation.

source
AJD.frobenius_offdiag_normMethod
frobenius_offdiag_norm(A::AbstractArray{T,3})::Real where {T<:Number}
  • A: Vector of matrices with sizen × n × k`

Takes the offdiagonal elements of an Array of matrices $A^k$ and applies the frobenius norm ($\sum |a_{i,j}|^{2}$).

source
AJD.get_diag_elementsMethod
get_diag_elements(A::Array)
  • A: Vector of matrices

Takes an array of matrices and returns the diagonal elements as a diagonal matrix D.

source
AJD.get_diagonalizationMethod
get_diagonalization(
    A::Vector{<:AbstractMatrix{<:Number}};
    algorithm::AbstractDiagonalization = JDiagGabrielDernbach(),
    max_iter::Int = 1000,
    threshold::AbstractFloat = eps(),
    only_plot::Symbol = :no_plot
    )

Get the actual diagonalization. Function is seperated from diagonalize() to facilitate plotting functionality in the REPL and Pluto. All implemented algorithms are called from this function. To generate the error histories of the algorithm runs, as used for the plots, select only_plot=:plot. Input is checked here as well.

source
AJD.get_offdiag_elementsMethod
get_offdiag_elements(A::Array{<:Number,3})
  • A: Vector of matrices

Takes an array of matrices and returns the offdiagonal elements of A.

source
AJD.get_y_fdiagMethod
get_y_fdiag(D::AbstractArray{<:Number}, E::AbstractArray{<:Number}, i::Int, j::Int)
  • D: Diagonal Matrix with offdiagonal elements set to zero
  • E: Diagonal Matrix with diagonal elements set to zero
  • i,j: Denotes the indexes of the matrices D and E

Calculates the factor $y_{ij}$ which is defined by: $∑_{k} D_{j,j}^{k}E_{j,i}^{k}$

source
AJD.get_z_fdiagMethod
get_z_fdiag(D::AbstractArray{<:Number}, i::Int, j::Int)
  • D: Diagonal Matrix with offdiagonal elements set to zero
  • i,j: Denotes the indexes of matrix D

Calculates the factor $z_{ij}$ which is defined by: $∑_{k} D_{i,i}^{k}D_{j,j}^{k}$

source
AJD.is_commutingMethod
is_commuting(A::AbstractMatrix, B::AbstractMatrix)
  • A: AbstractMatrix of dimension $n × n$
  • B: AbstractMatrix of dimension $n × n$

Check if two matrices A, B are commuting. $AB = BA$ must hold.

source
AJD.jdiag_cardosoMethod
jdiag_cardoso(M,threshhold,max_iter = 800)

Only works for matrix with real valued entries. Based on Matlab Code by Cardoso.

  • A: a $m × m × n$ matrix,($A_1,...,A_n$), each with dimension $m × m$
  • thresh: absolute threshold for approximation stops.default is 10e-8.
  • max_iter: number of iterations before algorithm stops. default is 800.

Output:

  • V : is a $m × m$ matrix, which accumulates givens rotations G in each iteration.
  • A : is a $m × m × n$ matrix, which contains [$VA_1V'$,...,$VA_nV'$]
  • iter: accumulates the iteration numbers
source
AJD.jdiag_edourdpineauMethod
jdiag_edourdpineau(X::Vector{M}; iter=100,  rtol = 1e-3, atol = eps())
    where {T<:Union{Real,Complex},M<:AbstractMatrix{T}}
  • X: set of matrices with dimension $n × n × k$
  • rtol: relative error before iteration stops
  • atol: absolute error before iteration stops

Diagonalize a set of matrices using the Jacobi method ("Jacobi Angles for Simultaneous Diagonalization"). Code adapted from Edouardpineaus Python implementation

source
AJD.jdiag_gabrieldernbach!Method
(1) jdiag_gabrieldernbach!(
    A::Vector{M};
    threshold::AbstractFloat = eps(),
    rel_threshold = 1e-3, #not used in function but easier for multiple dispatch
    max_iter = 1000,
    plot_convergence::Bool = false) where {T<:Real, M<:AbstractMatrix{T}}

JDiag algorithm based on the implementation by Gabrieldernbach in Python. Source: Algorithm

(2) jdiag_gabrieldernbach!(
A::Vector{M};
threshold = eps(),
rel_threshold = 1e-3,
max_iter = 1000,
plot_convergence::Bool = false) where {T<:Complex, M<:AbstractMatrix{T}}
  • A: set of matrices with dimension $n × n × k$
  • threshhold: absolute error used for stopping the iteration
  • rel_threshold: relative error used for stopping the iteration

JDiag algorithm for complex matrices based on the implementation by Gabrieldernbach in Python, the Cardoso Paper and the code of Algorithm

source
AJD.random_normal_commuting_matricesMethod
random_normal_commuting_matrices(n::Int, m::Int; complex::Bool=false)

Generate m random normal commuting matrices of size $n × n$ These can be exactly diagonalized

$M_i M_j = M_j M_i$ for all i,j $M_i M_i' = M_i' M_i$ for all i

source