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>
42class SymmSCMatrixdouble;
43class DiagSCMatrixdouble;
45class SCMatrixBlockIter;
46class SCMatrixRectBlock;
47class SCMatrixLTriBlock;
48class SCMatrixDiagBlock;
49class SCVectorSimpleBlock;
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;
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;
211 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;
339 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>
The SymmSCMatrix class is the abstract base class for diagonal double valued matrices.
Definition abstract.h:555
static std::ostream & out0()
Return an ostream that writes from node 0.
The RefDiagSCMatrix class is a smart pointer to an DiagSCMatrix specialization.
Definition matrix.h:389
RefDiagSCMatrix()
Initializes the matrix pointer to 0.
RefDiagSCMatrix(const RefDiagSCMatrix &m)
Make this and m refer to the same SCMatrix.
RefDiagSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this.
RefDiagSCMatrix operator+(const RefDiagSCMatrix &) const
Matrix addition and subtraction.
RefSCMatrix operator*(const RefSCMatrix &) const
Multiply this by a matrix and return a matrix.
RefDiagSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
RefDiagSCMatrix i() const
Return the inverse of this.
RefDiagSCMatrix(const RefSCDimension &, const Ref< SCMatrixKit > &)
Create a diagonal matrix with dimension d by d.
RefDiagSCMatrix & operator=(DiagSCMatrix *m)
Make this refer to m.
int nblock() const
If this matrix is blocked return the number of blocks.
RefDiagSCMatrix & operator=(const RefDiagSCMatrix &m)
Make this and m refer to the same matrix.
double determ() const
Returns the determinant of the referenced matrix.
RefDiagSCMatrix 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.
RefDiagSCMatrix(DiagSCMatrix *m)
Make this refer to m.
DiagSCMatrixdouble operator()(int i) const
Assign and examine matrix elements.
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition dim.h:152
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition matrix.h:135
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().
RefSCMatrix(const RefSCMatrix &m)
Make this and m refer to the same SCMatrix.
double determ() const
Returns the determinant of the referenced matrix.
RefSCMatrix operator*(const RefSCMatrix &) const
Multiply this by a matrix and return a matrix.
RefSCMatrix operator-(const RefSCMatrix &) const
Matrix subtraction.
double solve_lin(const RefSCVector &v) const
Solves this x = v.
RefSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this using SVD decomposition.
RefSCMatrix operator+(const RefSCMatrix &) const
Matrix addition.
RefSCMatrix i() const
Return the inverse of this.
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
SCMatrixdouble operator()(int i, int j) const
Assign and examine matrix elements.
RefSCMatrix t() const
Return the transpose of this.
RefSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
RefSCVector operator*(const RefSCVector &) const
Multiply this by a vector and return a vector.
RefSCMatrix operator*(double) const
Multiply this by a scalar and return the result.
RefSCMatrix block(int i) const
If this matrix is blocked return block i.
RefSCMatrix(SCMatrix *m)
Make this refer to m.
RefSCMatrix & operator=(const RefSCMatrix &m)
Make this and m refer to the same matrix.
int nblock() const
If this matrix is blocked return the number of blocks.
RefSCMatrix(const RefSCDimension &d1, const RefSCDimension &d2, const Ref< SCMatrixKit > &)
Create a vector with dimension d1 by d2.
RefSCMatrix & operator=(SCMatrix *m)
Make this refer to m.
The RefSCVector class is a smart pointer to an SCVector specialization.
Definition matrix.h:55
RefSCVector(const RefSCDimension &dim, const Ref< SCMatrixKit > &)
Create a vector with dimension dim.
SCVectordouble operator[](int) const
Return an l-value that can be used to assign or retrieve an element.
RefSCVector(const RefSCVector &v)
Make this and v refer to the same SCVector.
RefSCVector operator-(const RefSCVector &a) const
Subtract two vectors.
RefSCVector & operator=(const RefSCVector &v)
Make this and v refer to the same SCVector.
SCVectordouble operator()(int) const
Return an l-value that can be used to assign or retrieve an element.
RefSymmSCMatrix symmetric_outer_product() const
The outer product of this with itself is a symmetric matrix.
RefSCMatrix outer_product(const RefSCVector &v) const
Return the outer product between this and v.
RefSCVector(SCVector *v)
Make this refer to v.
void restore(StateIn &)
Restores the matrix from StateIn object. The vector must have been initialized already.
RefSCVector & operator=(SCVector *v)
Make this refer to v.
RefSCVector operator+(const RefSCVector &a) const
Add two vectors.
RefSCVector operator*(double) const
Scale a vector.
RefSCVector()
Initializes the vector pointer to 0.
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization.
Definition matrix.h:265
RefSymmSCMatrix i() const
Return the inverse of this.
RefSCMatrix convert2RefSCMat()
convert RefSymmSCMatrix to RefSCMatrix
RefSymmSCMatrix & operator=(const RefSymmSCMatrix &m)
Make this and m refer to the same matrix.
RefSymmSCMatrix operator+(const RefSymmSCMatrix &) const
Matrix addition and subtraction.
double solve_lin(const RefSCVector &) const
Solves this x = v.
RefSCMatrix eigvecs() const
Returns the eigenvectors of the reference matrix.
RefSymmSCMatrix(SymmSCMatrix *m)
Make this refer to m.
RefSymmSCMatrix(const RefSCDimension &d, const Ref< SCMatrixKit > &)
Create a vector with dimension d by d.
RefSCVector operator*(const RefSCVector &a) const
Multiply this by a vector and return a vector.
int nblock() const
If this matrix is blocked return the number of blocks.
SymmSCMatrixdouble operator()(int i, int j) const
Assign and examine matrix elements.
RefSymmSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this.
RefSymmSCMatrix & operator=(SymmSCMatrix *m)
Make this refer to m.
RefSymmSCMatrix()
Initializes the matrix pointer to 0.
RefSymmSCMatrix(const RefSymmSCMatrix &m)
Make this and m refer to the same SCMatrix.
RefSCMatrix operator*(const RefSCMatrix &) const
Multiply this by a matrix and return a matrix.
RefDiagSCMatrix eigvals() const
Returns the eigenvalues of the reference matrix.
void accumulate_transform(const RefSCMatrix &a, const RefSymmSCMatrix &b, SCMatrix::Transform=SCMatrix::NormalTransform) const
Add a * b * a.t() to this.
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...
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.
RefSymmSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
double determ() const
Returns the determinant of the referenced matrix.
A template class that maintains references counts.
Definition ref.h:361
SCMatrix & operator*() const
Definition ref.h:420
T * pointer() const
Returns a pointer the reference counted object.
Definition ref.h:413
The SCMatrix class is the abstract base class for general double valued n by m matrices.
Definition abstract.h:199
Transform
types of matrix transforms. Only real-valued matrices are assumed.
Definition abstract.h:205
The SCVector class is the abstract base class for double valued vectors.
Definition abstract.h:97
Restores fundamental and user-defined types from images created with StateOut.
Definition statein.h:79
Serializes fundamental and user-defined types.
Definition stateout.h:71
The SymmSCMatrix class is the abstract base class for symmetric double valued matrices.
Definition abstract.h:389
Contains all MPQC code up to version 3.
Definition mpqcin.h:14
void FromStateIn< sc::RefSymmSCMatrix >(sc::RefSymmSCMatrix &t, StateIn &so, int &count)
specialization for RefSymmSCMatrix
void FromStateIn< sc::RefDiagSCMatrix >(sc::RefDiagSCMatrix &t, StateIn &so, int &count)
specialization for RefDiagSCMatrix
void FromStateIn< sc::RefSCMatrix >(sc::RefSCMatrix &t, StateIn &so, int &count)
specialization for RefSCMatrix
void ToStateOut< sc::RefSCVector >(const sc::RefSCVector &t, StateOut &so, int &count)
specialization for RefSCVector
void FromStateIn< sc::RefSCVector >(sc::RefSCVector &t, StateIn &so, int &count)
specialization for RefSCVector
void ToStateOut< sc::RefSymmSCMatrix >(const sc::RefSymmSCMatrix &t, StateOut &so, int &count)
specialization for RefSymmSCMatrix
void ToStateOut< sc::RefDiagSCMatrix >(const sc::RefDiagSCMatrix &t, StateOut &so, int &count)
specialization for RefDiagSCMatrix
void ToStateOut< sc::RefSCMatrix >(const sc::RefSCMatrix &t, StateOut &so, int &count)
specialization for RefSCMatrix
this functor compares RefSymmSCMatrix objects.
Definition matrix.h:537