MPQC  3.0.0-alpha
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 
38 namespace sc {
39 
40 class SCVectordouble;
41 class SCMatrixdouble;
42 class SymmSCMatrixdouble;
43 class DiagSCMatrixdouble;
44 
45 class SCMatrixBlockIter;
46 class SCMatrixRectBlock;
47 class SCMatrixLTriBlock;
48 class SCMatrixDiagBlock;
49 class SCVectorSimpleBlock;
50 
51 class RefSCMatrix;
52 class RefSymmSCMatrix;
55 class RefSCVector: public Ref<SCVector> {
56  // standard overrides
57  public:
60  RefSCVector();
62  RefSCVector(const RefSCVector& v);
65  // don't allow automatic conversion from any reference to a
66  // described class
67  ~RefSCVector();
72 
73  // vector specific members
74  public:
77  RefSCVector(const RefSCDimension& dim,const Ref<SCMatrixKit>&);
78 
80  SCVectordouble operator()(int) const;
82  SCVectordouble operator[](int) const;
84  RefSCVector operator+(const RefSCVector&a) const;
86  RefSCVector operator-(const RefSCVector&a) const;
88  RefSCVector operator*(double) const;
90  RefSCMatrix outer_product(const RefSCVector& v) 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&);
126  void restore(StateIn&);
127 };
128 RefSCVector operator*(double,const RefSCVector&);
129 
130 class RefSymmSCMatrix;
131 class RefDiagSCMatrix;
135 class RefSCMatrix: public Ref<SCMatrix> {
136  // standard overrides
137  public:
138  typedef double value_type;
139  typedef double element_type;
140 
143  RefSCMatrix();
145  RefSCMatrix(const RefSCMatrix& m);
147  RefSCMatrix(SCMatrix* m);
148  ~RefSCMatrix();
152  RefSCMatrix& operator=(const RefSCMatrix& m);
153 
154  // matrix specific members
155  public:
158  RefSCMatrix(const RefSCDimension& d1,const RefSCDimension& d2,
159  const Ref<SCMatrixKit>&);
160 
162  RefSCVector operator*(const RefSCVector&) const;
163 
165  RefSCMatrix operator*(const RefSCMatrix&) const;
166  RefSCMatrix operator*(const RefSymmSCMatrix&) const;
167  RefSCMatrix operator*(const RefDiagSCMatrix&) const;
168 
170  RefSCMatrix operator*(double) const;
171 
173  RefSCMatrix operator+(const RefSCMatrix&) const;
175  RefSCMatrix operator-(const RefSCMatrix&) const;
176 
178  RefSCMatrix t() const;
180  RefSCMatrix i() const;
183  RefSCMatrix gi(double condition_number_threshold = 1e8) const;
184 
187  RefSCMatrix clone() const;
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&);
235  void restore(StateIn&);
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 };
261 RefSCMatrix operator*(double,const RefSCMatrix&);
262 
265 class RefSymmSCMatrix: public Ref<SymmSCMatrix> {
266  // standard overrides
267  public:
270  RefSymmSCMatrix();
275  ~RefSymmSCMatrix();
280  void copyRefSCMatrix(const RefSCMatrix& m);// the same dimension, do elementwise copying
281 
282  // matrix specific members
283  public:
288  RefSCMatrix operator*(const RefSCMatrix&) const;
289  RefSCMatrix operator*(const RefSymmSCMatrix&) const;
290  RefSCMatrix operator*(const RefDiagSCMatrix&) const;
292  RefSCVector operator*(const RefSCVector&a) const;
293  RefSymmSCMatrix operator*(double) const;
296  RefSymmSCMatrix operator-(const RefSymmSCMatrix&) const;
298  RefSymmSCMatrix i() const;
300  RefSymmSCMatrix gi(double condition_number_threshold = 1e8) const;
303  RefSymmSCMatrix clone() 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;
324  void accumulate_transform(const RefSCMatrix&a,const RefSymmSCMatrix&b,
325  SCMatrix::Transform = SCMatrix::NormalTransform) const;
326  void accumulate_transform(const RefSCMatrix&a,const RefDiagSCMatrix&b,
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&);
356  void restore(StateIn&);
357 
359  double solve_lin(const RefSCVector&) const;
361  double determ() const;
363  RefDiagSCMatrix eigvals() const;
365  RefSCMatrix eigvecs() const;
369  void diagonalize(const RefDiagSCMatrix& eigvals,
370  const RefSCMatrix& eigvecs) const;
372  SymmSCMatrixdouble operator()(int i,int j) const;
376  int nblock() const;
380  RefSymmSCMatrix block(int i) const;
383 };
385 RefSymmSCMatrix operator*(double,const RefSymmSCMatrix&);
386 
389 class RefDiagSCMatrix: public Ref<DiagSCMatrix> {
390  // standard overrides
391  public:
394  RefDiagSCMatrix();
399  ~RefDiagSCMatrix();
404 
405  // matrix specific members
406  public:
411  RefSCMatrix operator*(const RefSCMatrix&) const;
412  RefSCMatrix operator*(const RefSymmSCMatrix&) const;
414  RefDiagSCMatrix operator*(double) const;
417  RefDiagSCMatrix operator-(const RefDiagSCMatrix&) const;
419  RefDiagSCMatrix i() const;
421  RefDiagSCMatrix gi(double condition_number_threshold = 1e8) const;
424  RefDiagSCMatrix clone() 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&);
451  void restore(StateIn&);
453  double determ() const;
455  DiagSCMatrixdouble operator()(int i) const;
459  int nblock() const;
463  RefDiagSCMatrix block(int i) const;
464 };
466 RefDiagSCMatrix operator*(double,const RefDiagSCMatrix&);
467 
469  friend class RefSCVector;
470  private:
472  int i;
473 
474  SCVectordouble(SCVector*,int);
475  public:
477  ~SCVectordouble();
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:
488  int i;
489  int j;
490 
491  SCMatrixdouble(SCMatrix*,int,int);
492  public:
494  ~SCMatrixdouble();
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:
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:
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 
551 template <> void FromStateIn<sc::RefSCVector>(sc::RefSCVector& t, StateIn& so, int& count);
553 template <> void FromStateIn<sc::RefSCMatrix>(sc::RefSCMatrix& t, StateIn& so, int& count);
555 template <> void FromStateIn<sc::RefSymmSCMatrix>(sc::RefSymmSCMatrix& t, StateIn& so, int& count);
557 template <> void FromStateIn<sc::RefDiagSCMatrix>(sc::RefDiagSCMatrix& t, StateIn& so, int& count);
559 template <> void ToStateOut<sc::RefSCVector>(const sc::RefSCVector& t, StateOut& so, int& count);
561 template <> void ToStateOut<sc::RefSCMatrix>(const sc::RefSCMatrix& t, StateOut& so, int& count);
563 template <> void ToStateOut<sc::RefSymmSCMatrix>(const sc::RefSymmSCMatrix& t, StateOut& so, int& count);
565 template <> 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:
sc::RefSCMatrix::i
RefSCMatrix i() const
Return the inverse of this.
sc::RefSymmSCMatrix::determ
double determ() const
Returns the determinant of the referenced matrix.
sc::RefSCMatrix::operator+
RefSCMatrix operator+(const RefSCMatrix &) const
Matrix addition.
sc::RefSymmSCMatrix::diagonalize
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...
sc::RefDiagSCMatrix::i
RefDiagSCMatrix i() const
Return the inverse of this.
mpqc::ci::sigma
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
sc::SCVector
The SCVector class is the abstract base class for double valued vectors.
Definition: abstract.h:97
sc::RefSymmSCMatrix::restore
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
sc::SymmSCMatrixdouble
Definition: matrix.h:501
sc::RefSymmSCMatrix
The RefSymmSCMatrix class is a smart pointer to an SCSymmSCMatrix specialization.
Definition: matrix.h:265
sc::RefSymmSCMatrix::i
RefSymmSCMatrix i() const
Return the inverse of this.
sc::RefSCMatrix::nblock
int nblock() const
If this matrix is blocked return the number of blocks.
sc::RefSCMatrix
The RefSCMatrix class is a smart pointer to an SCMatrix specialization.
Definition: matrix.h:135
sc::RefSymmSCMatrix::operator+
RefSymmSCMatrix operator+(const RefSymmSCMatrix &) const
Matrix addition and subtraction.
sc::Ref
A template class that maintains references counts.
Definition: ref.h:361
mpqc::matrix
Matrix class derived from Eigen::Matrix with additional MPQC integration.
Definition: matrix.hpp:21
sc::RefDiagSCMatrix::clone
RefDiagSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
sc::RefSCMatrix::clone
RefSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
sc::DiagSCMatrixdouble
Definition: matrix.h:518
sc::RefSymmSCMatrix::nblock
int nblock() const
If this matrix is blocked return the number of blocks.
sc::RefDiagSCMatrix::block
RefDiagSCMatrix block(int i) const
If this matrix is blocked return block i.
sc::Ref::pointer
T * pointer() const
Returns a pointer the reference counted object.
Definition: ref.h:413
sc::Ref< SCVector >::operator*
SCVector & operator*() const
Returns a C++ reference to the reference counted object.
Definition: ref.h:420
sc::RefSymmSCMatrix::eigvecs
RefSCMatrix eigvecs() const
Returns the eigenvectors of the reference matrix.
sc::RefDiagSCMatrix
The RefDiagSCMatrix class is a smart pointer to an DiagSCMatrix specialization.
Definition: matrix.h:389
sc::RefSymmSCMatrix::clone
RefSymmSCMatrix clone() const
These call the SCMatrix members of the same name after checking for references to 0.
sc::DiagSCMatrix
The SymmSCMatrix class is the abstract base class for diagonal double valued matrices.
Definition: abstract.h:551
sc::RefSCMatrix::solve_lin
double solve_lin(const RefSCVector &v) const
Solves this x = v.
sc::RefSCVector::RefSCVector
RefSCVector()
Initializes the vector pointer to 0.
sc::convert
std::vector< double > convert(const RefDiagSCMatrix &A)
Converts RefDiagSCMatrix to std::vector<double>
sc::RefSymmSCMatrix::solve_lin
double solve_lin(const RefSCVector &) const
Solves this x = v.
sc::RefSCVector::operator=
RefSCVector & operator=(SCVector *v)
Make this refer to v.
sc::StateIn
Definition: statein.h:79
sc::RefSCDimension
The RefSCDimension class is a smart pointer to an SCDimension specialization.
Definition: dim.h:152
sc::RefSCMatrix::operator()
SCMatrixdouble operator()(int i, int j) const
Assign and examine matrix elements.
sc::RefSymmSCMatrix::eigvals
RefDiagSCMatrix eigvals() const
Returns the eigenvalues of the reference matrix.
sc::RefSymmSCMatrix::block
RefSymmSCMatrix block(int i) const
If this matrix is blocked return block i.
sc::SCVectordouble
Definition: matrix.h:468
sc::RefDiagSCMatrix::restore
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
sc::RefSCMatrix::operator-
RefSCMatrix operator-(const RefSCMatrix &) const
Matrix subtraction.
sc::RefSymmSCMatrix::operator=
RefSymmSCMatrix & operator=(SymmSCMatrix *m)
Make this refer to m.
sc::RefSCMatrix::restore
void restore(StateIn &)
Restores the matrix from StateIn object. The matrix must have been initialized already.
sc::RefSymmSCMatrixEqual
this functor compares RefSymmSCMatrix objects.
Definition: matrix.h:537
sc::RefDiagSCMatrix::RefDiagSCMatrix
RefDiagSCMatrix()
Initializes the matrix pointer to 0.
sc::RefDiagSCMatrix::operator+
RefDiagSCMatrix operator+(const RefDiagSCMatrix &) const
Matrix addition and subtraction.
sc::RefSCMatrix::RefSCMatrix
RefSCMatrix()
Initializes the matrix pointer to 0.
sc::RefSCMatrix::svd
void svd(const RefSCMatrix &U, const RefDiagSCMatrix &sigma, const RefSCMatrix &V)
Compute the singular value decomposition, this = U sigma V.t().
sc::RefDiagSCMatrix::nblock
int nblock() const
If this matrix is blocked return the number of blocks.
sc::RefSCVector::operator()
SCVectordouble operator()(int) const
Return an l-value that can be used to assign or retrieve an element.
sc::RefSymmSCMatrix::accumulate_transform
void accumulate_transform(const RefSCMatrix &a, const RefSymmSCMatrix &b, SCMatrix::Transform=SCMatrix::NormalTransform) const
Add a * b * a.t() to this.
sc::RefSymmSCMatrix::convert2RefSCMat
RefSCMatrix convert2RefSCMat()
convert RefSymmSCMatrix to RefSCMatrix
sc::RefSCVector
The RefSCVector class is a smart pointer to an SCVector specialization.
Definition: matrix.h:55
sc::RefSCVector::symmetric_outer_product
RefSymmSCMatrix symmetric_outer_product() const
The outer product of this with itself is a symmetric matrix.
sc::StateOut
Definition: stateout.h:71
sc::RefSCMatrix::operator=
RefSCMatrix & operator=(SCMatrix *m)
Make this refer to m.
sc::RefSCMatrix::t
RefSCMatrix t() const
Return the transpose of this.
sc::SCMatrixdouble
Definition: matrix.h:484
sc::SymmSCMatrix
The SymmSCMatrix class is the abstract base class for symmetric double valued matrices.
Definition: abstract.h:385
sc::SCMatrix::Transform
Transform
types of matrix transforms. Only real-valued matrices are assumed.
Definition: abstract.h:205
sc::RefDiagSCMatrix::operator()
DiagSCMatrixdouble operator()(int i) const
Assign and examine matrix elements.
sc::RefSCVector::operator[]
SCVectordouble operator[](int) const
Return an l-value that can be used to assign or retrieve an element.
sc::ExEnv::out0
static std::ostream & out0()
Return an ostream that writes from node 0.
sc::RefSCMatrix::gi
RefSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this using SVD decomposition.
sc::RefSCMatrix::determ
double determ() const
Returns the determinant of the referenced matrix.
sc::SCMatrix
The SCMatrix class is the abstract base class for general double valued n by m matrices.
Definition: abstract.h:195
sc::RefSymmSCMatrix::operator()
SymmSCMatrixdouble operator()(int i, int j) const
Assign and examine matrix elements.
sc::RefSCVector::outer_product
RefSCMatrix outer_product(const RefSCVector &v) const
Return the outer product between this and v.
sc::RefDiagSCMatrix::operator=
RefDiagSCMatrix & operator=(DiagSCMatrix *m)
Make this refer to m.
mpqc::vector
Vector class derived from Eigen::Matrix with additional MPQC integration.
Definition: matrix.hpp:133
sc::RefSCVector::operator+
RefSCVector operator+(const RefSCVector &a) const
Add two vectors.
sc::RefDiagSCMatrix::determ
double determ() const
Returns the determinant of the referenced matrix.
sc::RefSCMatrix::block
RefSCMatrix block(int i) const
If this matrix is blocked return block i.
sc
Contains all MPQC code up to version 3.
Definition: mpqcin.h:14
sc::RefSymmSCMatrix::gi
RefSymmSCMatrix gi(double condition_number_threshold=1e8) const
Return the generalized inverse of this.
sc::RefSCVector::restore
void restore(StateIn &)
Restores the matrix from StateIn object. The vector must have been initialized already.
sc::RefSymmSCMatrix::RefSymmSCMatrix
RefSymmSCMatrix()
Initializes the matrix pointer to 0.
sc::RefSCVector::operator-
RefSCVector operator-(const RefSCVector &a) const
Subtract two vectors.
sc::RefDiagSCMatrix::gi
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.