queso-0.53.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 
55  const Map& map,
56  unsigned int numCols);
57 
59  GslMatrix(const BaseEnvironment& env,
60  const Map& map,
61  double diagValue); // MATLAB eye
62 
64  GslMatrix(const GslVector& v,
65  double diagValue); // MATLAB eye
66 
68 
69  GslMatrix(const GslVector& v); // MATLAB diag
70 
72 
73  GslMatrix(const GslMatrix& B);
74 
76  ~GslMatrix();
78 
79 
81 
82  GslMatrix& operator= (const GslMatrix& rhs);
84 
86  GslMatrix& operator*=(double a);
87 
89  GslMatrix& operator/=(double a);
90 
92  GslMatrix& operator+=(const GslMatrix& rhs);
93 
95  GslMatrix& operator-=(const GslMatrix& rhs);
97 
98 
100 
101  double& operator()(unsigned int i, unsigned int j);
103 
105  const double& operator()(unsigned int i, unsigned int j) const;
106 
107 
109 
111 
112  unsigned int numRowsLocal () const;
114 
116  unsigned int numRowsGlobal () const;
117 
119  unsigned int numCols () const;
120 
122  double max () const;
123 
125 
126  unsigned int rank (double absoluteZeroThreshold, double relativeZeroThreshold) const;
127 
129  GslMatrix transpose () const;
130 
132  GslMatrix inverse () const;
133 
134 
136  double determinant () const;
137 
139  double lnDeterminant () const;
140 
142 
144 
145  double normFrob () const;
147 
149  double normMax () const;
151 
152 
153 
154 
155 
157 
158 
160 
161  int chol ();
162 
164  int svd (GslMatrix& matU, GslVector& vecS, GslMatrix& matVt) const;
165 
167 
168  const GslMatrix& svdMatU () const;
169 
171 
172  const GslMatrix& svdMatV () const;
173 
175  int svdSolve (const GslVector& rhsVec, GslVector& solVec) const;
176 
178  int svdSolve (const GslMatrix& rhsMat, GslMatrix& solMat) const;
179 
180 
181 
183  GslVector multiply (const GslVector& x) const;
184 
186 
187  GslVector invertMultiply (const GslVector& b) const;
188 
190 
192  void invertMultiply (const GslVector& b, GslVector& x) const;
193 
195 
196  GslMatrix invertMultiply (const GslMatrix& B) const;
197 
199 
201  void invertMultiply (const GslMatrix& B, GslMatrix& X) const;
202 
204 
206  GslVector invertMultiplyForceLU (const GslVector& b) const;
207 
209 
210  void invertMultiplyForceLU (const GslVector& b, GslVector& x) const;
211 
213  void getColumn (const unsigned int column_num, GslVector& column) const;
214 
216  GslVector getColumn (const unsigned int column_num) const;
217 
219  void setColumn (const unsigned int column_num, const GslVector& column);
220 
222  void getRow (const unsigned int row_num, GslVector& row) const;
223 
225  GslVector getRow (const unsigned int row_num) const;
226 
228  void setRow (const unsigned int row_num, const GslVector& row);
229 
231  void eigen (GslVector& eigenValues, GslMatrix* eigenVectors) const;
232 
234  void largestEigen (double& eigenValue, GslVector& eigenVector) const;
235 
237  void smallestEigen (double& eigenValue, GslVector& eigenVector) const;
238 
239 
241 
243 
244 
246  void cwSet (double value);
247 
249  void cwSet (unsigned int rowId, unsigned int colId, const GslMatrix& mat);
250 
251  void cwExtract (unsigned int rowId, unsigned int colId, GslMatrix& mat) const;
252 
254 
256  void zeroLower (bool includeDiagonal = false);
257 
259 
261  void zeroUpper (bool includeDiagonal = false);
262 
264 
265  void filterSmallValues (double thresholdValue);
266 
268 
269  void filterLargeValues (double thresholdValue);
270 
272  void fillWithTranspose (unsigned int rowId,
273  unsigned int colId,
274  const GslMatrix& mat,
275  bool checkForExactNumRowsMatching,
276  bool checkForExactNumColsMatching);
277 
279  void fillWithBlocksDiagonally (unsigned int rowId,
280  unsigned int colId,
281  const std::vector<const GslMatrix* >& matrices,
282  bool checkForExactNumRowsMatching,
283  bool checkForExactNumColsMatching);
284 
286  void fillWithBlocksDiagonally (unsigned int rowId,
287  unsigned int colId,
288  const std::vector< GslMatrix* >& matrices,
289  bool checkForExactNumRowsMatching,
290  bool checkForExactNumColsMatching);
291 
293  void fillWithBlocksHorizontally(unsigned int rowId,
294  unsigned int colId,
295  const std::vector<const GslMatrix* >& matrices,
296  bool checkForExactNumRowsMatching,
297  bool checkForExactNumColsMatching);
298 
300  void fillWithBlocksHorizontally(unsigned int rowId,
301  unsigned int colId,
302  const std::vector< GslMatrix* >& matrices,
303  bool checkForExactNumRowsMatching,
304  bool checkForExactNumColsMatching);
305 
307  void fillWithBlocksVertically (unsigned int rowId,
308  unsigned int colId,
309  const std::vector<const GslMatrix* >& matrices,
310  bool checkForExactNumRowsMatching,
311  bool checkForExactNumColsMatching);
312 
314  void fillWithBlocksVertically (unsigned int rowId,
315  unsigned int colId,
316  const std::vector< GslMatrix* >& matrices,
317  bool checkForExactNumRowsMatching,
318  bool checkForExactNumColsMatching);
319 
321  void fillWithTensorProduct (unsigned int rowId,
322  unsigned int colId,
323  const GslMatrix& mat1,
324  const GslMatrix& mat2,
325  bool checkForExactNumRowsMatching,
326  bool checkForExactNumColsMatching);
327 
329  void fillWithTensorProduct (unsigned int rowId,
330  unsigned int colId,
331  const GslMatrix& mat1,
332  const GslVector& vec2,
333  bool checkForExactNumRowsMatching,
334  bool checkForExactNumColsMatching);
336 
337 
339 
340  gsl_matrix* data ();
342 
343  void mpiSum (const MpiComm& comm, GslMatrix& M_global) const;
344 
345  void matlabLinearInterpExtrap (const GslVector& x1Vec, const GslMatrix& y1Mat, const GslVector& x2Vec);
347 
348 
349 
351 
352 
354  void print (std::ostream& os) const;
355 
357  void subWriteContents (const std::string& varNamePrefix,
358  const std::string& fileName,
359  const std::string& fileType,
360  const std::set<unsigned int>& allowedSubEnvIds) const;
361 
363  void subReadContents (const std::string& fileName,
364  const std::string& fileType,
365  const std::set<unsigned int>& allowedSubEnvIds);
367 
368 private:
370 
371  GslMatrix();
372 
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
const GslMatrix & svdMatV() const
This function calls private member GslMatrix::internalSvd() to set a N-by-N orthogonal square matrix ...
Definition: GslMatrix.C:547
void zeroUpper(bool includeDiagonal=false)
This function sets all the entries above the main diagonal of this matrix to zero.
Definition: GslMatrix.C:644
GslMatrix operator*(double a, const GslMatrix &mat)
Definition: GslMatrix.C:1972
double m_determinant
The determinant of this matrix.
Definition: GslMatrix.h:423
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:1791
GslMatrix & operator=(const GslMatrix &rhs)
Copies values from matrix rhs to this.
Definition: GslMatrix.C:169
bool m_isSingular
Indicates whether or not this matrix is singular.
Definition: GslMatrix.h:441
unsigned int numCols() const
Returns the column dimension of this matrix.
Definition: GslMatrix.C:312
int m_signum
m_signum stores the sign of the permutation of the LU decomposition PA = LU.
Definition: GslMatrix.h:438
GslVector invertMultiply(const GslVector &b) const
This function calculates the inverse of this matrix and multiplies it with vector b...
Definition: GslMatrix.C:1212
Class for matrix operations using GSL library.
Definition: GslMatrix.h:47
GslMatrix & operator/=(double a)
Stores in this the coordinate-wise division of this by a.
Definition: GslMatrix.C:187
void eigen(GslVector &eigenValues, GslMatrix *eigenVectors) const
This function computes the eigenvalues of a real symmetric matrix.
Definition: GslMatrix.C:1410
double determinant() const
Calculates the determinant of this matrix.
Definition: GslMatrix.C:1061
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:481
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:460
GslMatrix leftDiagScaling(const GslVector &vec, const GslMatrix &mat)
Definition: GslMatrix.C:2042
double & operator()(unsigned int i, unsigned int j)
Element access method (non-const).
Definition: GslMatrix.C:219
Class for matrix operations (virtual).
Definition: Matrix.h:46
void mpiSum(const MpiComm &comm, GslMatrix &M_global) const
Definition: GslMatrix.C:1688
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
GslMatrix & operator*=(double a)
Stores in this the coordinate-wise multiplication of this and a.
Definition: GslMatrix.C:177
GslMatrix()
Default Constructor.
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:1671
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:1037
gsl_matrix * m_mat
GSL matrix, also referred to as this matrix.
Definition: GslMatrix.h:386
std::ostream & operator<<(std::ostream &os, const BaseEnvironment &obj)
GslMatrix operator+(const GslMatrix &m1, const GslMatrix &m2)
Definition: GslMatrix.C:2012
double max() const
Returns the maximum element value of the matrix.
Definition: GslMatrix.C:354
The QUESO MPI Communicator Class.
Definition: MpiComm.h:75
gsl_permutation * m_permutation
The permutation matrix of a LU decomposition.
Definition: GslMatrix.h:432
void cwExtract(unsigned int rowId, unsigned int colId, GslMatrix &mat) const
Definition: GslMatrix.C:410
void smallestEigen(double &eigenValue, GslVector &eigenVector) const
This function finds smallest eigenvalue, namely eigenValue, of this matrix and its corresponding eige...
Definition: GslMatrix.C:1505
GslVector invertMultiplyForceLU(const GslVector &b) const
This function calculates the inverse of this matrix and multiplies it with vector b...
Definition: GslMatrix.C:1363
unsigned int numRowsGlobal() const
Returns the global row dimension of this matrix.
Definition: GslMatrix.C:306
void zeroLower(bool includeDiagonal=false)
This function sets all the entries bellow the main diagonal of this matrix to zero.
Definition: GslMatrix.C:616
unsigned int numRowsLocal() const
Returns the local row dimension of this matrix.
Definition: GslMatrix.C:300
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:912
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:1762
gsl_matrix * data()
Returns this matrix.
Definition: GslMatrix.C:1756
double m_lnDeterminant
The natural logarithm of the determinant of this matrix.
Definition: GslMatrix.h:426
Class for vector operations using GSL library.
Definition: GslVector.h:48
GslMatrix inverse() const
This function calculated the inverse of this matrix (square).
Definition: GslMatrix.C:736
GslMatrix rightDiagScaling(const GslMatrix &mat, const GslVector &vec)
Definition: GslMatrix.C:2063
GslMatrix & operator+=(const GslMatrix &rhs)
Stores in this the coordinate-wise addition of this and rhs.
Definition: GslMatrix.C:196
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:1602
GslVector multiply(const GslVector &x) const
This function multiplies this matrix by vector x and returns the resulting vector.
Definition: GslMatrix.C:1178
GslMatrix matrixProduct(const GslVector &v1, const GslVector &v2)
Definition: GslMatrix.C:2026
double lnDeterminant() const
Calculates the ln(determinant) of this matrix.
Definition: GslMatrix.C:1102
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
Definition: Environment.h:193
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:557
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:976
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
void matlabLinearInterpExtrap(const GslVector &x1Vec, const GslMatrix &y1Mat, const GslVector &x2Vec)
Definition: GslMatrix.C:1731
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
const BaseEnvironment & env() const
Definition: Matrix.C:47
void copy(const GslMatrix &src)
In this function this matrix receives a copy of matrix src.
Definition: GslMatrix.C:245
const Map & map() const
Definition: Matrix.C:54
void resetLU()
In this function resets the LU decomposition of this matrix, as well as deletes the private member po...
Definition: GslMatrix.C:257
unsigned int rank(double absoluteZeroThreshold, double relativeZeroThreshold) const
This function returns the number of singular values of this matrix (rank).
Definition: GslMatrix.C:1143
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:1576
gsl_matrix * m_LU
GSL matrix for the LU decomposition of m_mat.
Definition: GslMatrix.h:389
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 filterLargeValues(double thresholdValue)
This function sets to zero (filters) all entries of this matrix which are greater than thresholdValue...
Definition: GslMatrix.C:695
void largestEigen(double &eigenValue, GslVector &eigenVector) const
This function finds largest eigenvalue, namely eigenValue, of this matrix and its corresponding eigen...
Definition: GslMatrix.C:1436
GslMatrix operator-(const GslMatrix &m1, const GslMatrix &m2)
Definition: GslMatrix.C:2019
GslMatrix * m_inverse
Inverse matrix of this.
Definition: GslMatrix.h:392
A class for partitioning vectors and matrices.
Definition: Map.h:49
~GslMatrix()
Destructor.
Definition: GslMatrix.C:161
GslMatrix transpose() const
This function calculated the transpose of this matrix (square).
Definition: GslMatrix.C:718
void cwSet(double value)
Component-wise set all values to this with value.
Definition: GslMatrix.C:372
int chol()
Computes Cholesky factorization of a real symmetric positive definite matrix this.
Definition: GslMatrix.C:433
double normMax() const
Returns the Frobenius norm of this matrix.
Definition: GslMatrix.C:336
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:776
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:537
double normFrob() const
Returns the Frobenius norm of this matrix.
Definition: GslMatrix.C:318
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:1654
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:1831
GslMatrix & operator-=(const GslMatrix &rhs)
Stores in this the coordinate-wise subtraction of this by rhs.
Definition: GslMatrix.C:207
void filterSmallValues(double thresholdValue)
This function sets to zero (filters) all entries of this matrix which are smaller than thresholdValue...
Definition: GslMatrix.C:672
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:848
Map * m_svdColMap
Mapping for matrices involved in the singular value decomposition (svd) routine.
Definition: GslMatrix.h:395

Generated on Thu Jun 11 2015 13:52:31 for queso-0.53.0 by  doxygen 1.8.5