Level 1¶
The prototypes for the following routines can be found at include/elemental/blas-like/decl.hpp, while the implementations are in include/elemental/blas-like/level1/.
Adjoint¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
\(B := A^H\).
-
void
Adjoint
(const DistMatrix<T, U, V> &A, DistMatrix<T, W, Z> &B)¶
Axpy¶
Performs \(Y := \alpha X + Y\) (hence the name axpy).
-
void
Axpy
(T alpha, const DistMatrix<T, U1, V1> &X, DistMatrix<T, U2, V2> &Y)¶
Conjugate¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
\(A := \bar A\). For real datatypes, this is a no-op.
-
void
Conjugate
(DistMatrix<T, U, V> &A)¶
\(B := \bar A\).
-
void
Conjugate
(const DistMatrix<T, U, V> &A, DistMatrix<T, W, Z> &B)¶
DiagonalScale¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
Performs either \(X := \mbox{op}(D) X\) or \(X := X \mbox{op}(D)\), where \(op(D)\) equals \(D=D^T\), or \(D^H=\bar D\), where \(D = \mbox{diag}(d)\) and \(d\) is a column vector.
-
void
DiagonalScale
(LeftOrRight side, Orientation orientation, const Matrix<T> &d, Matrix<T> &X)¶
-
void
DiagonalScale
(LeftOrRight side, Orientation orientation, const DistMatrix<T, U, V> &d, DistMatrix<T, W, Z> &X)¶
DiagonalScaleTrapezoid¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
Performs either \(A := \mbox{op}(D) A\) or \(A := A \mbox{op}(D)\), where \(A\) is trapezoidal (upper or lower with the boundary diagonal of given offset), \(op(D)\) equals \(D=D^T\), or \(D^H=\bar D\), where \(D = \mbox{diag}(d)\) and \(d\) is a column vector.
-
void
DiagonalScaleTrapezoid
(LeftOrRight side, UpperOrLower uplo, Orientation orientation, const Matrix<T> &d, Matrix<T> &A, int offset = 0)¶
-
void
DiagonalScaleTrapezoid
(LeftOrRight side, UpperOrLower uplo, Orientation orientation, const DistMatrix<T, U, V> &d, DistMatrix<T, W, Z> &A, int offset = 0)¶
DiagonalSolve¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
Performs either \(X := \mbox{op}(D)^{-1} X\) or \(X := X \mbox{op}(D)^{-1}\), where \(D = \mbox{diag}(d)\) and \(d\) is a column vector.
-
void
DiagonalSolve
(LeftOrRight side, Orientation orientation, const Matrix<F> &d, Matrix<F> &X, bool checkIfSingular = false)¶
-
void
DiagonalSolve
(LeftOrRight side, Orientation orientation, const DistMatrix<F, U, V> &d, DistMatrix<F, W, Z> &X, bool checkIfSingular = false)¶
Dot¶
Returns \((x,y) = x^H y\). \(x\) and \(y\) are both allowed to be stored as column or row vectors, but will be interpreted as column vectors.
-
T
Dot
(const DistMatrix<T, U, V> &x, const DistMatrix<T, U, V> &y)¶
Dotc¶
Same as Dot
. This routine name is provided since it is the usual
BLAS naming convention.
-
T
Dotc
(const DistMatrix<T, U, V> &x, const DistMatrix<T, U, V> &y)¶
Dotu¶
Returns \(x^T y\), which is not an inner product.
-
T
Dotu
(const DistMatrix<T, U, V> &x, const DistMatrix<T, U, V> &y)¶
EntrywiseMap¶
-
void
EntrywiseMap
(DistMatrix<T, U, V> &A, Function func)¶
-
void
EntrywiseMap
(BlockDistMatrix<T, U, V> &A, Function func)¶ Replaces each entry of the passed in matrix with a specified function of the existing entry.
func
will typically be a lambda function which accepts a single argument of type T and returns a value of type T.
Hadamard¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
The Hadamard product of two \(m \times n\) matrices \(A\) and \(B\) is given entrywise by \(\alpha_{i,j} \beta_{i,j}\) and denoted by \(C = A \circ B\).
-
void
Hadamard
(const DistMatrix<F, U, V> &A, const DistMatrix<F, U, V> &B, DistMatrix<F, U, V> &C)¶
HilbertSchmidt¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
The Hilbert-Schmidt inner-product of two \(m \times n\) matrices \(A\) and \(B\) is \(\mbox{tr}(A^H B)\).
-
F
HilbertSchmidt
(const DistMatrix<F, U, V> &A, const DistMatrix<F, U, V> &B)¶
MakeTrapezoidal¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
Sets all entries outside of the specified trapezoidal submatrix to zero.
Whether or not the trapezoid is upper or lower
(analogous to an upper or lower-triangular matrix) is determined by the
uplo
parameter, and the last diagonal is defined with the offset
integer.
-
void
MakeTrapezoidal
(UpperOrLower uplo, Matrix<T> &A, int offset = 0)¶
-
void
MakeTrapezoidal
(UpperOrLower uplo, DistMatrix<T, U, V> &A, int offset = 0)¶
Nrm2¶
Returns \(||x||_2 = \sqrt{(x,x)} = \sqrt{x^H x}\). As with most other routines, even if \(x\) is stored as a row vector, it will be interpreted as a column vector.
-
Base<F>
Nrm2
(const DistMatrix<F> &x)¶
ScaleTrapezoid¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
Scales the entries within the specified trapezoid of a general matrix.
The parameter conventions follow those of MakeTrapezoidal
described above.
-
void
ScaleTrapezoid
(T alpha, UpperOrLower uplo, Matrix<T> &A, int offset = 0)¶
-
void
ScaleTrapezoid
(T alpha, UpperOrLower uplo, DistMatrix<T, U, V> &A, int offset = 0)¶
Transpose¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
\(B := A^T\) or \(B := A^H\).
-
void
Transpose
(const DistMatrix<T, U, V> &A, DistMatrix<T, W, Z> &B)¶
Zero¶
Note
This is not a standard BLAS routine, but it is BLAS-like.
Sets all of the entries of the input matrix to zero.
-
void
Zero
(DistMatrix<T, U, V> &A)¶
SetDiagonal¶
Note
This is not a standard BLAS routine.
Sets all of the diagonal entries of a matrix to a given value.
-
void
SetDiagonal
(DistMatrix<T, U, V> &A, T alpha)¶
-
void
SetDiagonal
(Matrix<T> &A, T alpha, int offset = 0, LeftOrRight side = LEFT)¶
-
void
SetDiagonal
(DistMatrix<T, U, V> &A, T alpha, int offset = 0, LeftOrRight side = LEFT)¶
Swap¶
-
void
Swap
(Orientation orientation, Matrix<T> &A, Matrix<T> &B)¶
-
void
Swap
(Orientation orientation, DistMatrix<T, U1, V1> &A, DistMatrix<T, U2, V2> &B)¶ Replace \(A\) and \(B\) with each other, their transpose, or their adjoint.
-
void
RowSwap
(DistMatrix<T, U, V> &A, int to, int from)¶ Swap rows to and from in the matrix.
-
void
ColumnSwap
(DistMatrix<T, U, V> &A, int to, int from)¶ Swap columns to and from in the matrix.
-
void
SymmetricSwap
(UpperOrLower uplo, Matrix<T> &A, int to, int from, bool conjugate = false)¶
-
void
SymmetricSwap
(UpperOrLower uplo, DistMatrix<T> &A, int to, int from, bool conjugate = false)¶ Symmetrically permute the to and from degrees of freedom within the implicitly symmetric (Hermitian) matrix \(A\) which stores its data in the specified triangle.
QuasiDiagonalScale¶
Note
This is not a standard BLAS routine.
-
void
QuasiDiagonalScale
(LeftOrRight side, UpperOrLower uplo, const Matrix<FMain> &d, const Matrix<F> &dSub, Matrix<F> &X, bool conjugate = false)¶
-
void
QuasiDiagonalScale
(LeftOrRight side, UpperOrLower uplo, const DistMatrix<FMain, U, V> &d, const DistMatrix<F, U, V> &dSub, DistMatrix<F> &X, bool conjugate = false)¶ Apply a symmetric (Hermitian) quasi-diagonal matrix to the matrix X.
QuasiDiagonalSolve¶
Note
This is not a standard BLAS routine.
-
void
QuasiDiagonalSolve
(LeftOrRight side, UpperOrLower uplo, const Matrix<FMain> &d, const Matrix<F> &dSub, Matrix<F> &X, bool conjugate = false)¶
-
void
QuasiDiagonalSolve
(LeftOrRight side, UpperOrLower uplo, const DistMatrix<FMain, U, V> &d, const DistMatrix<F, U, V> &dSub, DistMatrix<F> &X, bool conjugate = false)¶ Apply the inverse of a symmetric (Hermitian) quasi-diagonal matrix to the matrix X.
Symmetric2x2Scale¶
Note
This is not a standard BLAS routine.
-
void
Symmetric2x2Scale
(LeftOrRight side, UpperOrLower uplo, const Matrix<F> &D, Matrix<F> &A, bool conjugate = false)¶
-
void
Symmetric2x2Scale
(LeftOrRight side, UpperOrLower uplo, const DistMatrix<F, STAR, STAR> &D, DistMatrix<F> &A, bool conjugate = false)¶ Apply a 2x2 symmetric (Hermitian) matrix to the matrix A.
Symmetric2x2Solve¶
Note
This is not a standard BLAS routine.
-
void
Symmetric2x2Solve
(LeftOrRight side, UpperOrLower uplo, const Matrix<F> &D, Matrix<F> &A, bool conjugate = false)¶
-
void
Symmetric2x2Solve
(LeftOrRight side, UpperOrLower uplo, const DistMatrix<F, STAR, STAR> &D, DistMatrix<F> &A, bool conjugate = false)¶ Apply the inverse of a 2x2 symmetric (Hermitian) matrix to the matrix A.
UpdateDiagonal¶
Note
This is not a standard BLAS routine.
Adds a given value to the diagonal of a matrix.
-
void
UpdateDiagonal
(DistMatrix<T, U, V> &A, T alpha)¶
-
void
UpdateDiagonal
(Matrix<T> &A, T alpha, int offset = 0, LeftOrRight side = LEFT)¶
-
void
UpdateDiagonal
(DistMatrix<T, U, V> &A, T alpha, int offset = 0, LeftOrRight side = LEFT)¶