#include "linalg.h"
#include "linalg/shortclawrap.h"
Include dependency graph for linalg.cpp:
Functions | |
const complex | I (0, 1) |
Definition of the imaginary identity, I (declared in linalg.h). | |
void | ThrowError (char *errstring, bool fatal) |
A routine for easily printing an error string, and aborting if necessary. | |
void | Profile (char *methodname, int p1, int p2, double p3) |
Prints methodname, p1, p2, and p3 to "profile.log". | |
std::istream & | eatwhite (std::istream &input) |
Fast-forward a stream pointer past all whitespace. Used by FlexibleTable. | |
template<> | |
void | Copy (DenseVector< complex > &dest, const DenseVector< complex > &src) |
Copies one DenseVector<complex> into another. Uses BLAS if possible. | |
template<> | |
void | Copy (DenseVector< double > &dest, const DenseVector< double > &src) |
Copy one DenseVector<double> into another. Uses BLAS if possible. | |
template<> | |
void | Add (DenseVector< complex > &result, const DenseVector< complex > &v1, const DenseVector< complex > &v2, complex scale) |
Adds a scalar multiple of one DenseVector<complex> to another, and stores the result in another. Uses BLAS if possible. | |
template<> | |
void | Add (DenseVector< double > &result, const DenseVector< double > &v1, const DenseVector< double > &v2, double scale) |
Adds a scalar multiple of one DenseVector<double> to another, storing the result in a third. Uses BLAS if possible. | |
template<> | |
void | AddTo (DenseVector< complex > &result, const DenseVector< complex > &v, complex scale) |
Adds a scalar multiple of one DenseVector<complex> to another. Uses BLAS if possible. | |
template<> | |
void | AddTo (DenseVector< double > &result, const DenseVector< double > &v, double scale) |
Adds a scalar multiple of one DenseVector<double> to another in place. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseVector< complex > &result, const DenseVector< complex > &v, complex scalar) |
Scales a DenseVector<complex> by a scalar, storing the result in another. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseVector< double > &result, const DenseVector< double > &v, double scalar) |
Scales one DenseVector<double> by a scalar, storing the result in another. Uses BLAS if possible. | |
template<> | |
void | MultiplyBy (DenseVector< complex > &result, complex scalar) |
Scales a DenseVector<complex> by a scalar in place. Uses BLAS if possible. | |
template<> | |
void | MultiplyBy (DenseVector< double > &result, double scalar) |
Scales a DenseVector<double> by a scalar in place. Uses BLAS if possible. | |
template<> | |
double | dot (const DenseVector< double > &v1, const DenseVector< double > &v2) |
Computes the dot product between two DenseVector<double> objects. Uses BLAS if possible. | |
template<> | |
complex | dot (const DenseVector< complex > &v1, const DenseVector< complex > &v2) |
Computes the dot product between two DenseVector<complex> objects. Uses BLAS if possible. | |
template<> | |
double | c_dot (const DenseVector< double > &v1, const DenseVector< double > &v2) |
Computes the conjugate dot product between two DenseVector<double> objects. Uses BLAS if possible. | |
template<> | |
complex | c_dot (const DenseVector< complex > &v1, const DenseVector< complex > &v2) |
Computes the conjugate dot product between two DenseVector<complex> objects. Uses BLAS if possible. | |
DenseVector< complex > & | Randomize (DenseVector< complex > &v) |
Assigns a Gaussian random vector (all real and imaginary components are independent, univariate Gaussian random numbers). | |
DenseVector< double > | RealPart (const DenseVector< complex > &v) |
Returns a DenseVector<double> containing the real part of the argument. | |
DenseVector< double > | ImaginaryPart (const DenseVector< complex > &v) |
Returns a DenseVector<double> containing the imaginary part of the argument. | |
DenseVector< double > & | Randomize (DenseVector< double > &v) |
Assigns a Gaussian random vector (all components are independent, univariate Gaussian random numbers). | |
DenseMatrix< double > | RealPart (const DenseMatrix< complex > &m) |
Returns a DenseMatrix<double> containing the real part of the argument. | |
DenseMatrix< double > | ImaginaryPart (const DenseMatrix< complex > &m) |
Returns a DenseMatrix<double> containing the imaginary part of the argument. | |
bool | internalEigenFactor (DenseMatrix< complex > *matrix, DenseMatrix< complex > *L_ev, DenseMatrix< complex > *R_ev, DenseVector< complex > *evals) |
Private (internal) routine that interfaces to LAPACK complex eigenfactoring routines. | |
bool | internalEigenFactor (DenseMatrix< double > *matrix, DenseMatrix< complex > *L_ev, DenseMatrix< complex > *R_ev, DenseVector< complex > *evals) |
Private (internal) routine that interfaces to LAPACK double eigenfactoring routines. | |
bool | internalSymmetricEigenFactor (DenseMatrix< double > *matrix, DenseMatrix< double > *ev, DenseVector< double > *evals) |
Private (internal) routine that interfaces to LAPACK Hermitian-complex eigenfactoring routines. | |
bool | internalHermitianEigenFactor (DenseMatrix< complex > *matrix, DenseMatrix< complex > *ev, DenseVector< double > *evals) |
Private (internal) routine that interfaces to LAPACK symmetric-double eigenfactoring routines. | |
bool | EigenFactorSymmetric (DenseMatrix< double > &matrix, DenseMatrix< double > &ev, DenseVector< double > &evals) |
Computes the eigenvalues and eigenvectors of a symmetric real matrix. | |
bool | EigenValuesSymmetric (DenseMatrix< double > &matrix, DenseVector< double > &evals) |
bool | EigenFactorHermitian (DenseMatrix< complex > &matrix, DenseMatrix< complex > &ev, DenseVector< double > &evals) |
Computes the eigenvalues and eigenvectors of a Hermitian complex matrix. | |
bool | EigenValuesHermitian (DenseMatrix< complex > &matrix, DenseVector< double > &evals) |
Computes the eigenvalues (only) of a Hermitian complex matrix. | |
bool | CholeskyFactor (DenseMatrix< double > &M, bool UpperTriangle) |
Computes the Cholesky factorization of a symmetric positive definite matrix in place. | |
bool | CholeskyFactor (DenseMatrix< complex > &M, bool UpperTriangle) |
Computes the Cholesky factorization of a Hermitian positive definite matrix in place. | |
template<> | |
void | Copy (DenseMatrix< double > &dest, const DenseMatrix< double > &src) |
Copies one DenseMatrix<double> into another. Uses BLAS if possible. | |
template<> | |
void | Copy (DenseMatrix< complex > &dest, const DenseMatrix< complex > &src) |
Copies one DenseMatrix<complex> into another. Uses BLAS if possible. | |
template<> | |
void | Add (DenseMatrix< complex > &result, const DenseMatrix< complex > &m1, const DenseMatrix< complex > &m2, complex scale) |
Adds a scalar multiple of one DenseMatrix<complex> to another, storing the result in a third. Uses BLAS if possible. | |
template<> | |
void | Add (DenseMatrix< double > &result, const DenseMatrix< double > &m1, const DenseMatrix< double > &m2, double scale) |
Adds a scalar multiple of one DenseMatrix<double> to another, storing the result in a third. Uses BLAS if possible. | |
template<> | |
void | AddTo (DenseMatrix< complex > &result, const DenseMatrix< complex > &m, complex scale) |
Adds a scalar multiple of one DenseMatrix<complex> to another in place. Uses BLAS if possible. | |
template<> | |
void | AddTo (DenseMatrix< double > &result, const DenseMatrix< double > &m, double scale) |
Adds a scalar multiple of one DenseMatrix<double> to another in place. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseMatrix< double > &result, const DenseMatrix< double > &m1, const DenseMatrix< double > &m2) |
Multiplies one DenseMatrix<double> by another, storing the result in a third. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseMatrix< complex > &result, const DenseMatrix< complex > &m1, const DenseMatrix< complex > &m2) |
Multiplies one DenseMatrix<complex> by another, storing the result in a third. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseMatrix< complex > &result, const DenseMatrix< complex > &m, complex scalar) |
Scales a DenseMatrix<complex> by a scalar, storing the result in another. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseMatrix< double > &result, const DenseMatrix< double > &m, double scalar) |
Scales a DenseMatrix<double> by a scalar, storing the result in another. Uses BLAS if possible. | |
template<> | |
void | MultiplyBy (DenseMatrix< complex > &result, complex scalar) |
Scales a DenseMatrix<complex> by a scalar in place. Uses BLAS if possible. | |
template<> | |
void | MultiplyBy (DenseMatrix< double > &result, double scalar) |
Scales a DenseMatrix<double> by a scalar in place. Uses BLAS if possible. | |
template<> | |
complex | dot (const DenseMatrix< complex > &m1, const DenseMatrix< complex > &m2) |
Computes the dot product of two DenseMatrix<complex> objects. Uses BLAS if possible. | |
template<> | |
double | dot (const DenseMatrix< double > &m1, const DenseMatrix< double > &m2) |
Computes the Hilbert-Schmidt dot product of two DenseMatrix<double> objects. Uses BLAS if possible. | |
template<> | |
complex | c_dot (const DenseMatrix< complex > &m1, const DenseMatrix< complex > &m2) |
Computes the conjugate dot product of two DenseMatrix<complex> objects. Uses BLAS if possible. | |
template<> | |
double | c_dot (const DenseMatrix< double > &m1, const DenseMatrix< double > &m2) |
Computes the Hilbert-Schmidt conjugate dot product of two DenseMatrix<double> objects. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseVector< complex > &result, const DenseMatrix< complex > &m, const DenseVector< complex > &v) |
Multiplies a DenseVector<complex> by a DenseMatrix<complex> from the left and stores the result in another DenseVector<complex>. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseVector< double > &result, const DenseMatrix< double > &m, const DenseVector< double > &v) |
Multiplies a DenseVector<double> by a DenseMatrix<double> from the left and stores the result in another DenseVector<double>. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseVector< complex > &result, const DenseVector< complex > &v, const DenseMatrix< complex > &m) |
Multiplies a DenseVector<complex> by a DenseMatrix<complex> from the right and stores the result in another DenseVector<complex>. Uses BLAS if possible. | |
template<> | |
void | Multiply (DenseVector< double > &result, const DenseVector< double > &v, const DenseMatrix< double > &m) |
Multiplies a DenseVector<double> by a DenseMatrix<double> from the right and stores the result in another DenseVector<double>. Uses BLAS if possible. | |
Variables | |
const double | Pi = M_PI |
Defined in linalg.cpp as M_PI. | |
std::ofstream * | ProfileStream = NULL |
Output stream, used by Profile(), which opens it to "profile.log". |
|
Computes the Cholesky factorization of a symmetric positive definite matrix in place.
The input matrix The input matrix is assumed to be symmetric positive definite. Either the upper or lower triangle (depending on the value of UpperTriangle) will be ignored, and if the matrix turns out not to be positive, the routine will fail. The return value is "TRUE" unless the routine failed. IMPORTANT NOTE: The input matrix will be overwritten. ingroup LAPACK |
|
Computes the Cholesky factorization of a Hermitian positive definite matrix in place.
The input matrix The input matrix is assumed to be Hermitian positive definite. Either the upper or lower triangle (depending on the value of UpperTriangle) will be ignored, and if the matrix turns out not to be positive, the routine will fail. The return value is "TRUE" unless the routine failed. IMPORTANT NOTE: The input matrix will be overwritten. ingroup LAPACK |