MPQC 3.0.0-alpha
Loading...
Searching...
No Matches
matrix.h
1//
2// matrix.h
3//
4// Copyright (C) 1996 Limit Point Systems, Inc.
5//
6// Author: Curtis Janssen <cljanss@limitpt.com>
7// Maintainer: LPS
8//
9// This file is part of the SC Toolkit.
10//
11// The SC Toolkit is free software; you can redistribute it and/or modify
12// it under the terms of the GNU Library General Public License as published by
13// the Free Software Foundation; either version 2, or (at your option)
14// any later version.
15//
16// The SC Toolkit is distributed in the hope that it will be useful,
17// but WITHOUT ANY WARRANTY; without even the implied warranty of
18// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19// GNU Library General Public License for more details.
20//
21// You should have received a copy of the GNU Library General Public License
22// along with the SC Toolkit; see the file COPYING.LIB. If not, write to
23// the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
24//
25// The U.S. Government is granted a limited license as per AL 91-7.
26//
27
28#ifndef _math_scmat_matrix_h
29#define _math_scmat_matrix_h
30
31#include <iostream>
32#include <cfloat>
33
34#include <math/scmat/abstract.h>
35#include <util/state/statein.h>
36#include <util/state/stateout.h>
37
38namespace sc {
39
40class SCVectordouble;
41class SCMatrixdouble;
42class SymmSCMatrixdouble;
43class DiagSCMatrixdouble;
44
45class SCMatrixBlockIter;
46class SCMatrixRectBlock;
47class SCMatrixLTriBlock;
48class SCMatrixDiagBlock;
49class SCVectorSimpleBlock;
50
51class RefSCMatrix;
52class RefSymmSCMatrix;
55class RefSCVector: public Ref<SCVector> {
56 // standard overrides
57 public:
65 // don't allow automatic conversion from any reference to a
66 // described class
72
73 // vector specific members
74 public:
78
88 RefSCVector operator*(double) const;
93
94 void set_element(int i,double val) const;
95 void accumulate_element(int i,double val) const;
96 double get_element(int) const;
97 int n() const;
98 RefSCDimension dim() const;
99 Ref<SCMatrixKit> kit() const;
100 RefSCVector clone() const;
101 RefSCVector copy() const;
102 double maxabs() const;
103 double scalar_product(const RefSCVector&) const;
104 double dot(const RefSCVector&) const;
105 void normalize() const;
106 void randomize() const;
107 void assign(const RefSCVector& v) 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;
112 void accumulate(const RefSCVector& v) const;
113 void accumulate_product(const RefSymmSCMatrix&, const RefSCVector&);
114 void accumulate_product(const RefSCMatrix&, const RefSCVector&);
115 void element_op(const Ref<SCElementOp>& op) const;
116 void element_op(const Ref<SCElementOp2>&,
117 const RefSCVector&) const;
118 void element_op(const Ref<SCElementOp3>&,
119 const RefSCVector&,
120 const RefSCVector&) 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;
124 void save(StateOut&);
127};
128RefSCVector operator*(double,const RefSCVector&);
129
130class RefSymmSCMatrix;
131class RefDiagSCMatrix;
135class RefSCMatrix: public Ref<SCMatrix> {
136 // standard overrides
137 public:
138 typedef double value_type;
139 typedef double element_type;
140
148 ~RefSCMatrix();
153
154 // matrix specific members
155 public:
159 const Ref<SCMatrixKit>&);
160
163
168
170 RefSCMatrix operator*(double) const;
171
176
178 RefSCMatrix t() const;
180 RefSCMatrix i() const;
183 RefSCMatrix gi(double condition_number_threshold = 1e8) const;
184
188 RefSCMatrix copy() const;
189
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);
195 RefSCVector get_row(int) const;
196 RefSCVector get_column(int) const;
197 void assign_row(const RefSCVector&, int) const;
198 void assign_column(const RefSCVector&, int) const;
199 void accumulate_row(const RefSCVector&, int) const;
200 void accumulate_column(const RefSCVector&, int) const;
201
202 void accumulate_outer_product(const RefSCVector&,const RefSCVector&) const;
203 void accumulate_product(const RefSCMatrix&,const RefSCMatrix&) const;
204 void assign(const RefSCMatrix&) 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;
212 void accumulate(const RefSCMatrix&) const;
213 void accumulate(const RefSymmSCMatrix&) const;
214 void accumulate(const RefDiagSCMatrix&) const;
215 void element_op(const Ref<SCElementOp>&) const;
216 void element_op(const Ref<SCElementOp2>&,
217 const RefSCMatrix&) const;
218 void element_op(const Ref<SCElementOp3>&,
219 const RefSCMatrix&,
220 const RefSCMatrix&) const;
221 int nrow() const;
222 int ncol() const;
223 RefSCDimension rowdim() const;
224 RefSCDimension coldim() const;
225 Ref<SCMatrixKit> kit() 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,
231 std::ostream&out=ExEnv::out0(), int =10) const;
232 double trace() const;
233 void save(StateOut&);
236
241 void svd(const RefSCMatrix &U,
242 const RefDiagSCMatrix &sigma,
243 const RefSCMatrix &V);
245 double solve_lin(const RefSCVector& v) const;
247 double determ() const;
249 SCMatrixdouble operator()(int i,int j) const;
250
254 int nblock() const;
258 RefSCMatrix block(int i) const;
259};
261RefSCMatrix operator*(double,const RefSCMatrix&);
262
265class RefSymmSCMatrix: public Ref<SymmSCMatrix> {
266 // standard overrides
267 public:
280 void copyRefSCMatrix(const RefSCMatrix& m);// the same dimension, do elementwise copying
281
282 // matrix specific members
283 public:
293 RefSymmSCMatrix operator*(double) const;
296 RefSymmSCMatrix operator-(const RefSymmSCMatrix&) const;
300 RefSymmSCMatrix gi(double condition_number_threshold = 1e8) const;
304 RefSymmSCMatrix copy() const;
305 void set_element(int,int,double) const;
306 void accumulate_element(int,int,double) const;
307 double get_element(int,int) const;
308
309 RefSCMatrix get_subblock(int br, int er, int bc, int ec);
310 RefSymmSCMatrix get_subblock(int br, int er);
311 void assign_subblock(const RefSCMatrix&, int br, int er, int bc, int ec);
312 void assign_subblock(const RefSymmSCMatrix&, int br, int er);
313 void accumulate_subblock(const RefSCMatrix&, int, int, int, int);
314 void accumulate_subblock(const RefSymmSCMatrix&, int, int);
315 RefSCVector get_row(int);
316 void assign_row(const RefSCVector&, int);
317 void accumulate_row(const RefSCVector&, int);
318
319 void accumulate_symmetric_outer_product(const RefSCVector&) const;
320 double scalar_product(const RefSCVector&) const;
321 void accumulate_symmetric_product(const RefSCMatrix&) const;
322 void accumulate_symmetric_sum(const RefSCMatrix&) const;
325 SCMatrix::Transform = SCMatrix::NormalTransform) const;
327 SCMatrix::Transform = SCMatrix::NormalTransform) const;
329 const RefSymmSCMatrix&b) const;
330
331 void randomize() const;
332 void assign(const RefSymmSCMatrix&) 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;
340 void accumulate(const RefSymmSCMatrix&) const;
341 void element_op(const Ref<SCElementOp>&) const;
342 void element_op(const Ref<SCElementOp2>&,
343 const RefSymmSCMatrix&) const;
344 void element_op(const Ref<SCElementOp3>&,
345 const RefSymmSCMatrix&,
346 const RefSymmSCMatrix&) const;
347 double trace() const;
348 int n() const;
349 RefSCDimension dim() const;
350 Ref<SCMatrixKit> kit() const;
351 void print(std::ostream&) const;
352 void print(const char*title=0,
353 std::ostream&out=ExEnv::out0(), int =10) const;
354 void save(StateOut&);
357
359 double solve_lin(const RefSCVector&) const;
361 double determ() const;
370 const RefSCMatrix& eigvecs) const;
376 int nblock() const;
383};
385RefSymmSCMatrix operator*(double,const RefSymmSCMatrix&);
386
389class RefDiagSCMatrix: public Ref<DiagSCMatrix> {
390 // standard overrides
391 public:
404
405 // matrix specific members
406 public:
414 RefDiagSCMatrix operator*(double) const;
417 RefDiagSCMatrix operator-(const RefDiagSCMatrix&) const;
421 RefDiagSCMatrix gi(double condition_number_threshold = 1e8) const;
425 RefDiagSCMatrix copy() const;
426 void set_element(int,double) const;
427 void accumulate_element(int,double) const;
428 double get_element(int) const;
429 void randomize() const;
430 void assign(const RefDiagSCMatrix&) const;
431 void scale(double) const;
432 void assign(double) const;
433 void assign(const double*) const;
434 void convert(double*) const;
435 void accumulate(const RefDiagSCMatrix&) const;
436 void element_op(const Ref<SCElementOp>&) const;
437 void element_op(const Ref<SCElementOp2>&,
438 const RefDiagSCMatrix&) const;
439 void element_op(const Ref<SCElementOp3>&,
440 const RefDiagSCMatrix&,
441 const RefDiagSCMatrix&) const;
442 int n() const;
443 RefSCDimension dim() const;
444 Ref<SCMatrixKit> kit() const;
445 double trace() const;
446 void print(std::ostream&) const;
447 void print(const char*title=0,
448 std::ostream&out=ExEnv::out0(), int =10) const;
449 void save(StateOut&);
453 double determ() const;
459 int nblock() const;
464};
466RefDiagSCMatrix operator*(double,const RefDiagSCMatrix&);
467
469 friend class RefSCVector;
470 private:
471 RefSCVector vector;
472 int i;
473
475 public:
478 double operator=(double a);
479 double operator=(const SCVectordouble&);
480 operator double();
481 double val() const;
482};
483
485 friend class RefSCMatrix;
486 private:
487 RefSCMatrix matrix;
488 int i;
489 int j;
490
491 SCMatrixdouble(SCMatrix*,int,int);
492 public:
495 double operator=(double a);
496 double operator=(const SCMatrixdouble&);
497 operator double() const;
498 double val() const;
499};
500
502 friend class RefSymmSCMatrix;
503 private:
504 RefSymmSCMatrix matrix;
505 int i;
506 int j;
507
509 public:
512 double operator=(double a);
513 double operator=(const SymmSCMatrixdouble&);
514 operator double();
515 double val() const;
516};
517
519 friend class RefDiagSCMatrix;
520 private:
521 RefDiagSCMatrix matrix;
522 int i;
523 int j;
524
526 public:
529 double operator=(double a);
530 double operator=(const DiagSCMatrixdouble&);
531 operator double();
532 double val() const;
533};
534
538 bool operator()(const RefSymmSCMatrix& mat1, const RefSymmSCMatrix& mat2) {
539 if (mat1.pointer() == mat2.pointer()) return true;
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)
545 return false;
546 return true;
547 }
548};
549
551template <> void FromStateIn<sc::RefSCVector>(sc::RefSCVector& t, StateIn& so, int& count);
553template <> void FromStateIn<sc::RefSCMatrix>(sc::RefSCMatrix& t, StateIn& so, int& count);
559template <> void ToStateOut<sc::RefSCVector>(const sc::RefSCVector& t, StateOut& so, int& count);
561template <> void ToStateOut<sc::RefSCMatrix>(const sc::RefSCMatrix& t, StateOut& so, int& count);
563template <> void ToStateOut<sc::RefSymmSCMatrix>(const sc::RefSymmSCMatrix& t, StateOut& so, int& count);
565template <> void ToStateOut<sc::RefDiagSCMatrix>(const sc::RefDiagSCMatrix& t, StateOut& so, int& count);
566
567} // namespace sc
568
569#ifdef INLINE_FUNCTIONS
570#include <math/scmat/matrix_i.h>
571#endif
572
573#endif
574
575// Local Variables:
576// mode: c++
577// c-file-style: "CLJ"
578// End:
The SymmSCMatrix class is the abstract base class for diagonal double valued matrices.
Definition abstract.h:555
Definition matrix.h:518
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
Definition matrix.h:484
The SCVector class is the abstract base class for double valued vectors.
Definition abstract.h:97
Definition matrix.h:468
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
Definition matrix.h:501
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

Generated at Wed Sep 25 2024 02:45:30 for MPQC 3.0.0-alpha using the documentation package Doxygen 1.12.0.