queso-0.52.0
GslMatrix.h
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008-2015 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #ifndef UQ_GSL_MATRIX_H
26 #define UQ_GSL_MATRIX_H
27 
32 #include <queso/Matrix.h>
33 #include <gsl/gsl_matrix.h>
34 #include <gsl/gsl_permutation.h>
35 #include <queso/GslVector.h>
36 
37 namespace QUESO {
38 
47 class GslMatrix : public Matrix
48 {
49 public:
51 
52 
54 
55  GslMatrix();
56 
59  const Map& map,
60  unsigned int numCols);
61 
63  GslMatrix(const BaseEnvironment& env,
64  const Map& map,
65  double diagValue); // MATLAB eye
66 
68  GslMatrix(const GslVector& v,
69  double diagValue); // MATLAB eye
70 
72 
73  GslMatrix(const GslVector& v); // MATLAB diag
74 
76 
77  GslMatrix(const GslMatrix& B);
78 
80  ~GslMatrix();
82 
83 
85 
86  GslMatrix& operator= (const GslMatrix& rhs);
88 
90  GslMatrix& operator*=(double a);
91 
93  GslMatrix& operator/=(double a);
94 
96  GslMatrix& operator+=(const GslMatrix& rhs);
97 
99  GslMatrix& operator-=(const GslMatrix& rhs);
101 
102 
104 
105  double& operator()(unsigned int i, unsigned int j);
107 
109  const double& operator()(unsigned int i, unsigned int j) const;
110 
111 
113 
115 
116  unsigned int numRowsLocal () const;
118 
120  unsigned int numRowsGlobal () const;
121 
123  unsigned int numCols () const;
124 
126  double max () const;
127 
129 
130  unsigned int rank (double absoluteZeroThreshold, double relativeZeroThreshold) const;
131 
133  GslMatrix transpose () const;
134 
136  GslMatrix inverse () const;
137 
138 
140  double determinant () const;
141 
143  double lnDeterminant () const;
144 
146 
148 
149  double normFrob () const;
151 
153  double normMax () const;
155 
156 
157 
158 
159 
161 
162 
164 
165  int chol ();
166 
168  int svd (GslMatrix& matU, GslVector& vecS, GslMatrix& matVt) const;
169 
171 
172  const GslMatrix& svdMatU () const;
173 
175 
176  const GslMatrix& svdMatV () const;
177 
179  int svdSolve (const GslVector& rhsVec, GslVector& solVec) const;
180 
182  int svdSolve (const GslMatrix& rhsMat, GslMatrix& solMat) const;
183 
184 
185 
187  GslVector multiply (const GslVector& x) const;
188 
190 
191  GslVector invertMultiply (const GslVector& b) const;
192 
194 
196  void invertMultiply (const GslVector& b, GslVector& x) const;
197 
199 
200  GslMatrix invertMultiply (const GslMatrix& B) const;
201 
203 
205  void invertMultiply (const GslMatrix& B, GslMatrix& X) const;
206 
208 
210  GslVector invertMultiplyForceLU (const GslVector& b) const;
211 
213 
214  void invertMultiplyForceLU (const GslVector& b, GslVector& x) const;
215 
217  void getColumn (const unsigned int column_num, GslVector& column) const;
218 
220  GslVector getColumn (const unsigned int column_num) const;
221 
223  void setColumn (const unsigned int column_num, const GslVector& column);
224 
226  void getRow (const unsigned int row_num, GslVector& row) const;
227 
229  GslVector getRow (const unsigned int row_num) const;
230 
232  void setRow (const unsigned int row_num, const GslVector& row);
233 
235  void eigen (GslVector& eigenValues, GslMatrix* eigenVectors) const;
236 
238  void largestEigen (double& eigenValue, GslVector& eigenVector) const;
239 
241  void smallestEigen (double& eigenValue, GslVector& eigenVector) const;
242 
243 
245 
247 
248 
250  void cwSet (double value);
251 
253  void cwSet (unsigned int rowId, unsigned int colId, const GslMatrix& mat);
254 
255  void cwExtract (unsigned int rowId, unsigned int colId, GslMatrix& mat) const;
256 
258 
260  void zeroLower (bool includeDiagonal = false);
261 
263 
265  void zeroUpper (bool includeDiagonal = false);
266 
268 
269  void filterSmallValues (double thresholdValue);
270 
272 
273  void filterLargeValues (double thresholdValue);
274 
276  void fillWithTranspose (unsigned int rowId,
277  unsigned int colId,
278  const GslMatrix& mat,
279  bool checkForExactNumRowsMatching,
280  bool checkForExactNumColsMatching);
281 
283  void fillWithBlocksDiagonally (unsigned int rowId,
284  unsigned int colId,
285  const std::vector<const GslMatrix* >& matrices,
286  bool checkForExactNumRowsMatching,
287  bool checkForExactNumColsMatching);
288 
290  void fillWithBlocksDiagonally (unsigned int rowId,
291  unsigned int colId,
292  const std::vector< GslMatrix* >& matrices,
293  bool checkForExactNumRowsMatching,
294  bool checkForExactNumColsMatching);
295 
297  void fillWithBlocksHorizontally(unsigned int rowId,
298  unsigned int colId,
299  const std::vector<const GslMatrix* >& matrices,
300  bool checkForExactNumRowsMatching,
301  bool checkForExactNumColsMatching);
302 
304  void fillWithBlocksHorizontally(unsigned int rowId,
305  unsigned int colId,
306  const std::vector< GslMatrix* >& matrices,
307  bool checkForExactNumRowsMatching,
308  bool checkForExactNumColsMatching);
309 
311  void fillWithBlocksVertically (unsigned int rowId,
312  unsigned int colId,
313  const std::vector<const GslMatrix* >& matrices,
314  bool checkForExactNumRowsMatching,
315  bool checkForExactNumColsMatching);
316 
318  void fillWithBlocksVertically (unsigned int rowId,
319  unsigned int colId,
320  const std::vector< GslMatrix* >& matrices,
321  bool checkForExactNumRowsMatching,
322  bool checkForExactNumColsMatching);
323 
325  void fillWithTensorProduct (unsigned int rowId,
326  unsigned int colId,
327  const GslMatrix& mat1,
328  const GslMatrix& mat2,
329  bool checkForExactNumRowsMatching,
330  bool checkForExactNumColsMatching);
331 
333  void fillWithTensorProduct (unsigned int rowId,
334  unsigned int colId,
335  const GslMatrix& mat1,
336  const GslVector& vec2,
337  bool checkForExactNumRowsMatching,
338  bool checkForExactNumColsMatching);
340 
341 
343 
344  gsl_matrix* data ();
346 
347  void mpiSum (const MpiComm& comm, GslMatrix& M_global) const;
348 
349  void matlabLinearInterpExtrap (const GslVector& x1Vec, const GslMatrix& y1Mat, const GslVector& x2Vec);
351 
352 
353 
355 
356 
358  void print (std::ostream& os) const;
359 
361  void subWriteContents (const std::string& varNamePrefix,
362  const std::string& fileName,
363  const std::string& fileType,
364  const std::set<unsigned int>& allowedSubEnvIds) const;
365 
367  void subReadContents (const std::string& fileName,
368  const std::string& fileType,
369  const std::set<unsigned int>& allowedSubEnvIds);
371 
372 private:
374  void copy (const GslMatrix& src);
375 
377  void resetLU ();
378 
380  void multiply (const GslVector& x, GslVector& y) const;
381 
383  int internalSvd () const;
384 
386  gsl_matrix* m_mat;
387 
389  mutable gsl_matrix* m_LU;
390 
393 
395  mutable Map* m_svdColMap;
396 
398 
402 
404 
409 
411 
415 
417 
421 
423  mutable double m_determinant;
424 
426  mutable double m_lnDeterminant;
427 
429 
432  mutable gsl_permutation* m_permutation;
433 
435 
438  mutable int m_signum;
439 
441  mutable bool m_isSingular;
442 };
443 
444 GslMatrix operator* (double a, const GslMatrix& mat);
445 GslVector operator* (const GslMatrix& mat, const GslVector& vec);
446 GslMatrix operator* (const GslMatrix& m1, const GslMatrix& m2 );
447 GslMatrix operator+ (const GslMatrix& m1, const GslMatrix& m2 );
448 GslMatrix operator- (const GslMatrix& m1, const GslMatrix& m2 );
449 GslMatrix matrixProduct (const GslVector& v1, const GslVector& v2 );
450 
451 // Row \c i of the returned matrix is equal to row \c i of \c mat multiplied by element \c i of \c v.
452 GslMatrix leftDiagScaling (const GslVector& vec, const GslMatrix& mat);
453 
454 // Column \c j of the returned matrix is equal to column \c j of \c mat multiplied by element \c j of \c v.
455 GslMatrix rightDiagScaling(const GslMatrix& mat, const GslVector& vec);
456 std::ostream& operator<< (std::ostream& os, const GslMatrix& obj);
457 
458 } // End namespace QUESO
459 
460 #endif // UQ_GSL_MATRIX_H
void fillWithBlocksHorizontally(unsigned int rowId, unsigned int colId, const std::vector< const GslMatrix * > &matrices, bool checkForExactNumRowsMatching, bool checkForExactNumColsMatching)
This function fills this matrix horizontally with const block matrices.
Definition: GslMatrix.C:984
GslMatrix & operator/=(double a)
Stores in this the coordinate-wise division of this by a.
Definition: GslMatrix.C:216
void setColumn(const unsigned int column_num, const GslVector &column)
This function copies vector column into the column_num-th column of this matrix.
Definition: GslMatrix.C:2002
GslMatrix & operator+=(const GslMatrix &rhs)
Stores in this the coordinate-wise addition of this and rhs.
Definition: GslMatrix.C:225
double max() const
Returns the maximum element value of the matrix.
Definition: GslMatrix.C:403
Class for matrix operations using GSL library.
Definition: GslMatrix.h:47
GslMatrix operator*(double a, const GslMatrix &mat)
Definition: GslMatrix.C:2354
gsl_matrix * m_LU
GSL matrix for the LU decomposition of m_mat.
Definition: GslMatrix.h:389
GslMatrix operator-(const GslMatrix &m1, const GslMatrix &m2)
Definition: GslMatrix.C:2404
GslMatrix operator+(const GslMatrix &m1, const GslMatrix &m2)
Definition: GslMatrix.C:2397
GslMatrix rightDiagScaling(const GslMatrix &mat, const GslVector &vec)
Definition: GslMatrix.C:2454
void zeroUpper(bool includeDiagonal=false)
This function sets all the entries above the main diagonal of this matrix to zero.
Definition: GslMatrix.C:747
GslMatrix matrixProduct(const GslVector &v1, const GslVector &v2)
Definition: GslMatrix.C:2411
void fillWithTensorProduct(unsigned int rowId, unsigned int colId, const GslMatrix &mat1, const GslMatrix &mat2, bool checkForExactNumRowsMatching, bool checkForExactNumColsMatching)
This function calculates the tensor product of matrices mat1 and mat2 and stores it in this matrix...
Definition: GslMatrix.C:1160
unsigned int numRowsLocal() const
Returns the local row dimension of this matrix.
Definition: GslMatrix.C:349
~GslMatrix()
Destructor.
Definition: GslMatrix.C:187
std::ostream & operator<<(std::ostream &os, const BaseEnvironment &obj)
void fillWithBlocksVertically(unsigned int rowId, unsigned int colId, const std::vector< const GslMatrix * > &matrices, bool checkForExactNumRowsMatching, bool checkForExactNumColsMatching)
This function fills this matrix vertically with const block matrices.
Definition: GslMatrix.C:1072
void eigen(GslVector &eigenValues, GslMatrix *eigenVectors) const
This function computes the eigenvalues of a real symmetric matrix.
Definition: GslMatrix.C:1696
GslMatrix leftDiagScaling(const GslVector &vec, const GslMatrix &mat)
Definition: GslMatrix.C:2427
GslVector invertMultiply(const GslVector &b) const
This function calculates the inverse of this matrix and multiplies it with vector b...
Definition: GslMatrix.C:1441
double & operator()(unsigned int i, unsigned int j)
Element access method (non-const).
Definition: GslMatrix.C:254
double normFrob() const
Returns the Frobenius norm of this matrix.
Definition: GslMatrix.C:367
GslVector * m_svdSvec
m_svdSvec stores the diagonal of the N-by-N diagonal matrix of singular values S after the singular v...
Definition: GslMatrix.h:408
void getRow(const unsigned int row_num, GslVector &row) const
This function gets the row_num-th column of this matrix and stores it into vector row...
Definition: GslMatrix.C:1915
GslMatrix * m_svdUmat
m_svdUmat stores the M-by-N orthogonal matrix U after the singular value decomposition of a matrix...
Definition: GslMatrix.h:401
void filterLargeValues(double thresholdValue)
This function sets to zero (filters) all entries of this matrix which are greater than thresholdValue...
Definition: GslMatrix.C:801
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
Definition: Environment.h:187
void print(std::ostream &os) const
Print method. Defines the behavior of the ostream &lt;&lt; operator inherited from the Object class...
Definition: GslMatrix.C:2114
void mpiSum(const MpiComm &comm, GslMatrix &M_global) const
Definition: GslMatrix.C:2028
A class for partitioning vectors and matrices.
Definition: Map.h:49
GslVector invertMultiplyForceLU(const GslVector &b) const
This function calculates the inverse of this matrix and multiplies it with vector b...
Definition: GslMatrix.C:1622
double m_lnDeterminant
The natural logarithm of the determinant of this matrix.
Definition: GslMatrix.h:426
GslMatrix & operator-=(const GslMatrix &rhs)
Stores in this the coordinate-wise subtraction of this by rhs.
Definition: GslMatrix.C:239
void getColumn(const unsigned int column_num, GslVector &column) const
This function gets the column_num-th column of this matrix and stores it into vector column...
Definition: GslMatrix.C:1880
void zeroLower(bool includeDiagonal=false)
This function sets all the entries bellow the main diagonal of this matrix to zero.
Definition: GslMatrix.C:716
void subWriteContents(const std::string &varNamePrefix, const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds) const
Write contents of subenvironment in file fileName.
Definition: GslMatrix.C:2143
GslMatrix transpose() const
This function calculated the transpose of this matrix (square).
Definition: GslMatrix.C:824
int svdSolve(const GslVector &rhsVec, GslVector &solVec) const
This function solves the system A x = b using the singular value decomposition (U, S, V) of A which must have been computed previously with GslMatrix::svd (x=solVec, b=rhsVec).
Definition: GslMatrix.C:563
Class for vector operations using GSL library.
Definition: GslVector.h:48
GslMatrix * m_svdVTmat
m_svdVmatT stores the transpose of N-by-N orthogonal square matrix V, namely V^T, after the singular ...
Definition: GslMatrix.h:420
void cwSet(double value)
Component-wise set all values to this with value.
Definition: GslMatrix.C:421
GslMatrix & operator=(const GslMatrix &rhs)
Copies values from matrix rhs to this.
Definition: GslMatrix.C:195
void largestEigen(double &eigenValue, GslVector &eigenVector) const
This function finds largest eigenvalue, namely eigenValue, of this matrix and its corresponding eigen...
Definition: GslMatrix.C:1728
void fillWithBlocksDiagonally(unsigned int rowId, unsigned int colId, const std::vector< const GslMatrix * > &matrices, bool checkForExactNumRowsMatching, bool checkForExactNumColsMatching)
This function fills this matrix diagonally with const block matrices.
Definition: GslMatrix.C:888
GslMatrix()
Default Constructor.
Definition: GslMatrix.C:36
Map * m_svdColMap
Mapping for matrices involved in the singular value decomposition (svd) routine.
Definition: GslMatrix.h:395
const GslMatrix & svdMatV() const
This function calls private member GslMatrix::internalSvd() to set a N-by-N orthogonal square matrix ...
Definition: GslMatrix.C:644
void setRow(const unsigned int row_num, const GslVector &row)
This function copies vector row into the row_num-th column of this matrix.
Definition: GslMatrix.C:1976
The QUESO MPI Communicator Class.
Definition: MpiComm.h:75
gsl_matrix * data()
Returns this matrix.
Definition: GslMatrix.C:2108
GslMatrix inverse() const
This function calculated the inverse of this matrix (square).
Definition: GslMatrix.C:845
unsigned int rank(double absoluteZeroThreshold, double relativeZeroThreshold) const
This function returns the number of singular values of this matrix (rank).
Definition: GslMatrix.C:1363
Class for matrix operations (virtual).
Definition: Matrix.h:46
void copy(const GslMatrix &src)
In this function this matrix receives a copy of matrix src.
Definition: GslMatrix.C:292
void smallestEigen(double &eigenValue, GslVector &eigenVector) const
This function finds smallest eigenvalue, namely eigenValue, of this matrix and its corresponding eige...
Definition: GslMatrix.C:1803
GslVector multiply(const GslVector &x) const
This function multiplies this matrix by vector x and returns the resulting vector.
Definition: GslMatrix.C:1398
void subReadContents(const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds)
Read contents of subenvironment from file fileName.
Definition: GslMatrix.C:2189
void fillWithTranspose(unsigned int rowId, unsigned int colId, const GslMatrix &mat, bool checkForExactNumRowsMatching, bool checkForExactNumColsMatching)
This function stores the transpose of this matrix into this matrix.
Definition: GslMatrix.C:1245
GslMatrix & operator*=(double a)
Stores in this the coordinate-wise multiplication of this and a.
Definition: GslMatrix.C:203
void filterSmallValues(double thresholdValue)
This function sets to zero (filters) all entries of this matrix which are smaller than thresholdValue...
Definition: GslMatrix.C:778
const Map & map() const
Definition: Matrix.C:123
int svd(GslMatrix &matU, GslVector &vecS, GslMatrix &matVt) const
Checks for the dimension of this matrix, matU, VecS and matVt, and calls the protected routine intern...
Definition: GslMatrix.C:533
gsl_permutation * m_permutation
The permutation matrix of a LU decomposition.
Definition: GslMatrix.h:432
double normMax() const
Returns the Frobenius norm of this matrix.
Definition: GslMatrix.C:385
void matlabLinearInterpExtrap(const GslVector &x1Vec, const GslMatrix &y1Mat, const GslVector &x2Vec)
Definition: GslMatrix.C:2071
gsl_matrix * m_mat
GSL matrix, also referred to as this matrix.
Definition: GslMatrix.h:386
void resetLU()
In this function resets the LU decomposition of this matrix, as well as deletes the private member po...
Definition: GslMatrix.C:306
int m_signum
m_signum stores the sign of the permutation of the LU decomposition PA = LU.
Definition: GslMatrix.h:438
unsigned int numRowsGlobal() const
Returns the global row dimension of this matrix.
Definition: GslMatrix.C:355
int internalSvd() const
This function factorizes the M-by-N matrix A into the singular value decomposition A = U S V^T for M ...
Definition: GslMatrix.C:654
double lnDeterminant() const
Calculates the ln(determinant) of this matrix.
Definition: GslMatrix.C:1322
double determinant() const
Calculates the determinant of this matrix.
Definition: GslMatrix.C:1281
bool m_isSingular
Indicates whether or not this matrix is singular.
Definition: GslMatrix.h:441
double m_determinant
The determinant of this matrix.
Definition: GslMatrix.h:423
GslMatrix * m_svdVmat
m_svdVmat stores the N-by-N orthogonal square matrix V after the singular value decomposition of a ma...
Definition: GslMatrix.h:414
unsigned int numCols() const
Returns the column dimension of this matrix.
Definition: GslMatrix.C:361
GslMatrix * m_inverse
Inverse matrix of this.
Definition: GslMatrix.h:392
void cwExtract(unsigned int rowId, unsigned int colId, GslMatrix &mat) const
Definition: GslMatrix.C:471
const GslMatrix & svdMatU() const
This function calls private member GslMatrix::internalSvd() to set a M-by-N orthogonal matrix U of th...
Definition: GslMatrix.C:634
const BaseEnvironment & env() const
Definition: Matrix.C:116
int chol()
Computes Cholesky factorization of a real symmetric positive definite matrix this.
Definition: GslMatrix.C:506

Generated on Thu Apr 23 2015 19:30:54 for queso-0.52.0 by  doxygen 1.8.5