28 #ifndef _math_scmat_matrix_h
29 #define _math_scmat_matrix_h
34 #include <math/scmat/abstract.h>
35 #include <util/state/statein.h>
36 #include <util/state/stateout.h>
42 class SymmSCMatrixdouble;
43 class DiagSCMatrixdouble;
45 class SCMatrixBlockIter;
46 class SCMatrixRectBlock;
47 class SCMatrixLTriBlock;
48 class SCMatrixDiagBlock;
49 class SCVectorSimpleBlock;
52 class RefSymmSCMatrix;
94 void set_element(
int i,
double val)
const;
95 void accumulate_element(
int i,
double val)
const;
96 double get_element(
int)
const;
102 double maxabs()
const;
105 void normalize()
const;
106 void randomize()
const;
108 void assign(
double val)
const;
109 void assign(
const double* v)
const;
110 void convert(
double*)
const;
111 void scale(
double val)
const;
121 void print(std::ostream&out)
const;
122 void print(
const char*title=0,
123 std::ostream&out=
ExEnv::out0(),
int precision=10)
const;
138 typedef double value_type;
139 typedef double element_type;
183 RefSCMatrix gi(
double condition_number_threshold = 1e8)
const;
190 RefSCMatrix get_subblock(
int br,
int er,
int bc,
int ec);
191 void assign_subblock(
const RefSCMatrix&,
int br,
int er,
int bc,
int ec,
192 int source_br = 0,
int source_bc = 0);
193 void accumulate_subblock(
const RefSCMatrix&,
int,
int,
int,
int,
194 int source_br = 0,
int source_bc = 0);
199 void accumulate_row(
const RefSCVector&,
int)
const;
200 void accumulate_column(
const RefSCVector&,
int)
const;
205 void scale(
double)
const;
206 void randomize()
const;
207 void assign(
double)
const;
208 void assign(
const double*)
const;
209 void assign(
const double**)
const;
210 void convert(
double*)
const;
226 void set_element(
int,
int,
double)
const;
227 void accumulate_element(
int,
int,
double)
const;
228 double get_element(
int,
int)
const;
229 void print(std::ostream&)
const;
230 void print(
const char*title=0,
232 double trace()
const;
305 void set_element(
int,
int,
double)
const;
306 void accumulate_element(
int,
int,
double)
const;
307 double get_element(
int,
int)
const;
309 RefSCMatrix get_subblock(
int br,
int er,
int bc,
int ec);
311 void assign_subblock(
const RefSCMatrix&,
int br,
int er,
int bc,
int ec);
313 void accumulate_subblock(
const RefSCMatrix&,
int,
int,
int,
int);
319 void accumulate_symmetric_outer_product(
const RefSCVector&)
const;
321 void accumulate_symmetric_product(
const RefSCMatrix&)
const;
322 void accumulate_symmetric_sum(
const RefSCMatrix&)
const;
331 void randomize()
const;
333 void scale(
double)
const;
334 void scale_diagonal(
double)
const;
335 void assign(
double)
const;
336 void assign(
const double*)
const;
337 void assign(
const double**)
const;
338 void convert(
double*)
const;
347 double trace()
const;
351 void print(std::ostream&)
const;
352 void print(
const char*title=0,
426 void set_element(
int,
double)
const;
427 void accumulate_element(
int,
double)
const;
428 double get_element(
int)
const;
429 void randomize()
const;
431 void scale(
double)
const;
432 void assign(
double)
const;
433 void assign(
const double*)
const;
434 void convert(
double*)
const;
445 double trace()
const;
446 void print(std::ostream&)
const;
447 void print(
const char*title=0,
478 double operator=(
double a);
495 double operator=(
double a);
497 operator double()
const;
512 double operator=(
double a);
529 double operator=(
double a);
540 const int n = mat1.n();
541 if (n != mat2.n())
return false;
542 for(
int r=0; r<n; ++r)
543 for(
int c=0; c<=r; ++c)
544 if ( fabs(mat1(r,c) - mat2(r,c)) > DBL_EPSILON)
569 #ifdef INLINE_FUNCTIONS
570 #include <math/scmat/matrix_i.h>
RefSCMatrix i() const
Return the inverse of this.
double determ() const
Returns the determinant of the referenced matrix.
RefSCMatrix operator+(const RefSCMatrix &) const
Matrix addition.
void diagonalize(const RefDiagSCMatrix &eigvals, const RefSCMatrix &eigvecs) const
Sets eigvals to the eigenvalues and eigvecs to the eigenvalues and eigenvectors of the referenced mat...
RefDiagSCMatrix i() const
Return the inverse of this.
void sigma(const CI< Type, Index > &ci, const mpqc::Vector &h, const Matrix &V, ci::Vector &C, ci::Vector &S)
Computes sigma 1,2,3 contributions.
Definition: sigma.hpp:30
The SCVector class is the abstract base class for double valued vectors.
Definition: abstract.h:97
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization.
Definition: matrix.h:265
RefSymmSCMatrix i() const
Return the inverse of this.
int nblock() const
If this matrix is blocked return the number of blocks.
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition: matrix.h:135
RefSymmSCMatrix operator+(const RefSymmSCMatrix &) const
Matrix addition and subtraction.
A template class that maintains references counts.
Definition: ref.h:361
Matrix class derived from Eigen::Matrix with additional MPQC integration.
Definition: matrix.hpp:21
RefDiagSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
RefSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
int nblock() const
If this matrix is blocked return the number of blocks.
RefDiagSCMatrix block(int i) const
If this matrix is blocked return block i.
T * pointer() const
Returns a pointer the reference counted object.
Definition: ref.h:413
SCVector & operator*() const
Returns a C++ reference to the reference counted object.
Definition: ref.h:420
RefSCMatrix eigvecs() const
Returns the eigenvectors of the reference matrix.
The RefDiagSCMatrix class is a smart pointer to an DiagSCMatrix specialization.
Definition: matrix.h:389
RefSymmSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
The SymmSCMatrix class is the abstract base class for diagonal double valued matrices.
Definition: abstract.h:551
double solve_lin(const RefSCVector &v) const
Solves this x = v.
RefSCVector()
Initializes the vector pointer to 0.
std::vector< double > convert(const RefDiagSCMatrix &A)
Converts RefDiagSCMatrix to std::vector<double>
double solve_lin(const RefSCVector &) const
Solves this x = v.
RefSCVector & operator=(SCVector *v)
Make this refer to v.
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition: dim.h:152
SCMatrixdouble operator()(int i, int j) const
Assign and examine matrix elements.
RefDiagSCMatrix eigvals() const
Returns the eigenvalues of the reference matrix.
RefSymmSCMatrix block(int i) const
If this matrix is blocked return block i.
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
RefSCMatrix operator-(const RefSCMatrix &) const
Matrix subtraction.
RefSymmSCMatrix & operator=(SymmSCMatrix *m)
Make this refer to m.
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
this functor compares RefSymmSCMatrix objects.
Definition: matrix.h:537
RefDiagSCMatrix()
Initializes the matrix pointer to 0.
RefDiagSCMatrix operator+(const RefDiagSCMatrix &) const
Matrix addition and subtraction.
RefSCMatrix()
Initializes the matrix pointer to 0.
void svd(const RefSCMatrix &U, const RefDiagSCMatrix &sigma, const RefSCMatrix &V)
Compute the singular value decomposition, this = U sigma V.t().
int nblock() const
If this matrix is blocked return the number of blocks.
SCVectordouble operator()(int) const
Return an l-value that can be used to assign or retrieve an element.
void accumulate_transform(const RefSCMatrix &a, const RefSymmSCMatrix &b, SCMatrix::Transform=SCMatrix::NormalTransform) const
Add a * b * a.t() to this.
RefSCMatrix convert2RefSCMat()
convert RefSymmSCMatrix to RefSCMatrix
The RefSCVector class is a smart pointer to an SCVector specialization.
Definition: matrix.h:55
RefSymmSCMatrix symmetric_outer_product() const
The outer product of this with itself is a symmetric matrix.
Definition: stateout.h:71
RefSCMatrix & operator=(SCMatrix *m)
Make this refer to m.
RefSCMatrix t() const
Return the transpose of this.
The SymmSCMatrix class is the abstract base class for symmetric double valued matrices.
Definition: abstract.h:385
Transform
types of matrix transforms. Only real-valued matrices are assumed.
Definition: abstract.h:205
DiagSCMatrixdouble operator()(int i) const
Assign and examine matrix elements.
SCVectordouble operator[](int) const
Return an l-value that can be used to assign or retrieve an element.
static std::ostream & out0()
Return an ostream that writes from node 0.
RefSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this using SVD decomposition.
double determ() const
Returns the determinant of the referenced matrix.
The SCMatrix class is the abstract base class for general double valued n by m matrices.
Definition: abstract.h:195
SymmSCMatrixdouble operator()(int i, int j) const
Assign and examine matrix elements.
RefSCMatrix outer_product(const RefSCVector &v) const
Return the outer product between this and v.
RefDiagSCMatrix & operator=(DiagSCMatrix *m)
Make this refer to m.
Vector class derived from Eigen::Matrix with additional MPQC integration.
Definition: matrix.hpp:133
RefSCVector operator+(const RefSCVector &a) const
Add two vectors.
double determ() const
Returns the determinant of the referenced matrix.
RefSCMatrix block(int i) const
If this matrix is blocked return block i.
Contains all MPQC code up to version 3.
Definition: mpqcin.h:14
RefSymmSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this.
void restore(StateIn &)
Restores the matrix from StateIn object. The vector must have been initialized already.
RefSymmSCMatrix()
Initializes the matrix pointer to 0.
RefSCVector operator-(const RefSCVector &a) const
Subtract two vectors.
RefDiagSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this.
Generated at Sun Jan 26 2020 23:24:00 for MPQC
3.0.0-alpha using the documentation package Doxygen
1.8.16.