Dory.jl Documentation
Dory is a package to extend the functionality of Hecke.jl.
The main features are:
- Convenience functions, such as various matrix get/set index methods for matrices supporting various types.
- Generic eigenvector methods.
- p-adic Linear algebra.
Broadcasting
Broadcasting is enabled for the AbstractAlgebra.Generic.MatElem{T}
type. the syntax is your_function.(A)
.
The return type is either of type AbstractAlgebra.Generic.MatElem{T}
or an array if the output type of your_function
is not a subtype of NCRingElem
.
pAdic linear algebra.
Dory.padic_qr
— Function.padic_qr(A :: Hecke.Generic.MatElem{padic} ; col_pivot :: Union{Val{true},Val{false}}) --> F :: QRPadicPivoted
The return type of F
is a QRPadicPivoted, with fields F.Q, F.R, F.p, F.q
described below.
Compute the p-adic QR factorization of A
. More precisely, compute matrices Q
,R
, and an arrays p
, q
such that
A[F.p,F.q] = Q*R
If col_pivot=Val(false)
, then F.q = [1,2,...,size(A,2)]
.
#–––––––––-
INPUTS: A
– a matrix over Qp col_pivot
– a type, either Val(true)
or Val(false)
, indicating whether column permutations should be used to move p-adically large entries to the pivot position.
padic_qr(A; col_pivot) --> F
The return type of F
is a QRPadicPivoted, with fields F.Q, F.R, F.p, F.q
described below.
Compute the p-adic QR factorization of A. More precisely, compute matrices Q
,R
, and an arrays p
, q
such that
A[F.p,F.q] = Q*R
If col_pivot=Val(false), then F.q = [1,2,...,size(A,2)].
#–––––––––-
INPUTS: A – a matrix over Qp, A::Hecke.Generic.MatElem{padic} col_pivot – a type, either Val(true) or Val(false), indicating whether column permutations should be used to move p-adically large entries to the pivot position.
Dory.rectangular_solve
— Function.rectangular_solve(A::Hecke.MatElem{T}, b::Hecke.MatElem{T}; suppress_error=false) where T
--> x ::Hecke.MatElem{T}
Solve the possibly overdetermined linear equation Ax = b. If no solution exists returns an error. Generally not intended for use.
WARNINGS:
This function is very unsafe. It does not do basic sanity checks and will fail if the top nxn block is singular.
rectangular_solve(A::Hecke.MatElem{padic}, b_input::Hecke.MatElem{padic}; stable::Bool=false)
--> (nu :: Int64,N::Hecke.MatElem{padic})
Solves the linear system A*N = b. The output nu
is the dimension of the nullspace. Parameter stable
determines whether padic_qr
or svd
method is used. Default is qr (for speed).
WARNINGS: If A,b_input
have different precisions, maximal precision output is not guarenteed. Underdetermined solve not implemented.
Dory.eigen
— Function.eigen(A::nmod_mat)
Computes the Eigenvalue decomposition of A
. Requires factorization of polynomials implemented over the base ring.
(Depreciated. eigspaces
is better to use.)
Dory.eigvecs
— Function.eigvecs( A :: Hecke.Generic.MatElem{T}) where T -> A :: Hecke.Generic.MatElem{T}
Return a matrix M
whose columns are the eigenvectors of A
. (The kth eigenvector can be obtained from the slice M[:, k]
.)
Dory.eigspaces
— Function.eigspaces(A::Hecke.Generic.MatElem{T}) where T --> EigenSpaceDec{T}
Computes the eigen spaces of a generic matrix, and returns a list of matrices whose columns are generators for the eigen spaces.
Dory.singular_values
— Function.singular_values(A::Hecke.MatElem{padic}) -> Array{padic, 1}
Returns the list of diagonal elements in the singular value decomposition of the matrix A
.