queso-0.53.0
GslMatrix.C
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 #include <queso/GslMatrix.h>
26 #include <queso/GslVector.h>
27 #include <queso/Defines.h>
28 #include <gsl/gsl_linalg.h>
29 #include <gsl/gsl_eigen.h>
30 #include <sys/time.h>
31 #include <cmath>
32 
33 namespace QUESO {
34 
35 GslMatrix::GslMatrix( // can be a rectangular matrix
36  const BaseEnvironment& env,
37  const Map& map,
38  unsigned int nCols)
39  :
40  Matrix (env,map),
41  m_mat (gsl_matrix_calloc(map.NumGlobalElements(),nCols)),
42  m_LU (NULL),
43  m_inverse (NULL),
44  m_svdColMap (NULL),
45  m_svdUmat (NULL),
46  m_svdSvec (NULL),
47  m_svdVmat (NULL),
48  m_svdVTmat (NULL),
49  m_determinant (-INFINITY),
50  m_lnDeterminant(-INFINITY),
51  m_permutation (NULL),
52  m_signum (0),
53  m_isSingular (false)
54 {
55  queso_require_msg(m_mat, "null matrix generated");
56 }
57 
58 // Shaped constructor --------------------------------------------------
59 GslMatrix::GslMatrix( // square matrix
60  const BaseEnvironment& env,
61  const Map& map,
62  double diagValue)
63  :
64  Matrix (env,map),
65  m_mat (gsl_matrix_calloc(map.NumGlobalElements(),map.NumGlobalElements())),
66  m_LU (NULL),
67  m_inverse (NULL),
68  m_svdColMap (NULL),
69  m_svdUmat (NULL),
70  m_svdSvec (NULL),
71  m_svdVmat (NULL),
72  m_svdVTmat (NULL),
73  m_determinant (-INFINITY),
74  m_lnDeterminant(-INFINITY),
75  m_permutation (NULL),
76  m_signum (0),
77  m_isSingular (false)
78 {
79  queso_require_msg(m_mat, "null matrix generated");
80 
81  for (unsigned int i = 0; i < m_mat->size1; ++i) {
82  (*this)(i,i) = diagValue;
83  }
84 }
85 // Shaped constructor --------------------------------------------------
86 GslMatrix::GslMatrix( // square matrix
87  const GslVector& v,
88  double diagValue)
89  :
90  Matrix (v.env(),v.map()),
91  m_mat (gsl_matrix_calloc(v.sizeLocal(),v.sizeLocal())),
92  m_LU (NULL),
93  m_inverse (NULL),
94  m_svdColMap (NULL),
95  m_svdUmat (NULL),
96  m_svdSvec (NULL),
97  m_svdVmat (NULL),
98  m_svdVTmat (NULL),
99  m_determinant (-INFINITY),
100  m_lnDeterminant(-INFINITY),
101  m_permutation (NULL),
102  m_signum (0),
103  m_isSingular (false)
104 {
105  queso_require_msg(m_mat, "null matrix generated");
106 
107  for (unsigned int i = 0; i < m_mat->size1; ++i) {
108  (*this)(i,i) = diagValue;
109  }
110 }
111 // Shaped constructor --------------------------------------------------
112 GslMatrix::GslMatrix(const GslVector& v) // square matrix
113  :
114  Matrix (v.env(),v.map()),
115  m_mat (gsl_matrix_calloc(v.sizeLocal(),v.sizeLocal())),
116  m_LU (NULL),
117  m_inverse (NULL),
118  m_svdColMap (NULL),
119  m_svdUmat (NULL),
120  m_svdSvec (NULL),
121  m_svdVmat (NULL),
122  m_svdVTmat (NULL),
123  m_determinant (-INFINITY),
124  m_lnDeterminant(-INFINITY),
125  m_permutation (NULL),
126  m_signum (0),
127  m_isSingular (false)
128 {
129  queso_require_msg(m_mat, "null matrix generated");
130 
131  unsigned int dim = v.sizeLocal();
132  for (unsigned int i = 0; i < dim; ++i) {
133  (*this)(i,i) = v[i];
134  }
135 }
136 
137 // Shaped constructor --------------------------------------------------
138 GslMatrix::GslMatrix(const GslMatrix& B) // can be a rectangular matrix
139  :
140  Matrix (B.env(),B.map()),
141  m_mat (gsl_matrix_calloc(B.numRowsLocal(),B.numCols())),
142  m_LU (NULL),
143  m_inverse (NULL),
144  m_svdColMap (NULL),
145  m_svdUmat (NULL),
146  m_svdSvec (NULL),
147  m_svdVmat (NULL),
148  m_svdVTmat (NULL),
149  m_determinant (-INFINITY),
150  m_lnDeterminant(-INFINITY),
151  m_permutation (NULL),
152  m_signum (0),
153  m_isSingular (false)
154 {
155  queso_require_msg(m_mat, "null vector generated");
156  this->Matrix::base_copy(B);
157  this->copy(B);
158 }
159 
160 // Destructor ----------------------------------------------------------
162 {
163  this->resetLU();
164  if (m_mat) gsl_matrix_free(m_mat);
165 }
166 
167 // Set methods (operators) ---------------------------------------------
168 GslMatrix&
170 {
171  this->resetLU();
172  this->copy(obj);
173  return *this;
174 }
175 
176 GslMatrix&
178 {
179  this->resetLU();
180  int iRC;
181  iRC = gsl_matrix_scale(m_mat,a);
182  queso_require_msg(!(iRC), "scaling failed");
183  return *this;
184 }
185 
186 GslMatrix&
188 {
189  this->resetLU();
190  *this *= (1./a);
191 
192  return *this;
193 }
194 
195 GslMatrix&
197 {
198  this->resetLU();
199  int iRC;
200  iRC = gsl_matrix_add(m_mat,rhs.m_mat);
201  queso_require_msg(!(iRC), "failed");
202 
203  return *this;
204 }
205 
206 GslMatrix&
208 {
209  this->resetLU();
210  int iRC;
211  iRC = gsl_matrix_sub(m_mat,rhs.m_mat);
212  queso_require_msg(!(iRC), "failed");
213 
214  return *this;
215 }
216 
217 // Acessor methods (operators) ----------------------------------------
218 double&
219 GslMatrix::operator()(unsigned int i, unsigned int j)
220 {
221  this->resetLU();
222  if ((i >= m_mat->size1) ||
223  (j >= m_mat->size2)) {
224  std::cerr << "In GslMatrix::operator(i,j)"
225  << ": i = " << i
226  << ", j = " << j
227  << ", m_mat->size1 = " << m_mat->size1
228  << ", m_mat->size2 = " << m_mat->size2
229  << std::endl;
230  queso_require_less_msg(i, m_mat->size1, "i is too large");
231  queso_require_less_msg(j, m_mat->size2, "j is too large");
232  }
233  return *gsl_matrix_ptr(m_mat,i,j);
234 }
235 
236 const double&
237 GslMatrix::operator()(unsigned int i, unsigned int j) const
238 {
239  queso_require_less_msg(i, m_mat->size1, "i is too large");
240  queso_require_less_msg(j, m_mat->size2, "j is too large");
241  return *gsl_matrix_const_ptr(m_mat,i,j);
242 }
243 
244 void
246 {
247  // FIXME - should we be calling Matrix::base_copy here? - RHS
248  this->resetLU();
249  int iRC;
250  iRC = gsl_matrix_memcpy(this->m_mat, src.m_mat);
251  queso_require_msg(!(iRC), "failed");
252 
253  return;
254 }
255 
256 void
258 {
259  if (m_LU) {
260  gsl_matrix_free(m_LU);
261  m_LU = NULL;
262  }
263  if (m_inverse) {
264  delete m_inverse;
265  m_inverse = NULL;
266  }
267  if (m_svdColMap) {
268  delete m_svdColMap;
269  m_svdColMap = NULL;
270  }
271  if (m_svdUmat) {
272  delete m_svdUmat;
273  m_svdUmat = NULL;
274  }
275  if (m_svdSvec) {
276  delete m_svdSvec;
277  m_svdSvec = NULL;
278  }
279  if (m_svdVmat) {
280  delete m_svdVmat;
281  m_svdVmat = NULL;
282  }
283  if (m_svdVTmat) {
284  delete m_svdVTmat;
285  m_svdVTmat = NULL;
286  }
287  m_determinant = -INFINITY;
288  m_lnDeterminant = -INFINITY;
289  if (m_permutation) {
290  gsl_permutation_free(m_permutation);
291  m_permutation = NULL;
292  }
293  m_signum = 0;
294  m_isSingular = false;
295 
296  return;
297 }
298 
299 unsigned int
301 {
302  return m_mat->size1;
303 }
304 
305 unsigned int
307 {
308  return m_mat->size1;
309 }
310 
311 unsigned int
313 {
314  return m_mat->size2;
315 }
316 
317 double
319 {
320  double value = 0.;
321 
322  unsigned int nRows = this->numRowsLocal();
323  unsigned int nCols = this->numCols();
324  double aux = 0.;
325  for (unsigned int i = 0; i < nRows; i++) {
326  for (unsigned int j = 0; j < nCols; j++) {
327  aux = (*this)(i,j);
328  value += aux*aux;
329  }
330  }
331 
332  return sqrt(value);
333 }
334 
335 double
337 {
338  double value = 0.;
339 
340  unsigned int nRows = this->numRowsLocal();
341  unsigned int nCols = this->numCols();
342  double aux = 0.;
343  for (unsigned int i = 0; i < nRows; i++) {
344  for (unsigned int j = 0; j < nCols; j++) {
345  aux = fabs((*this)(i,j));
346  if (aux > value) value = aux;
347  }
348  }
349 
350  return value;
351 }
352 
353 double
355 {
356  double value = -INFINITY;
357 
358  unsigned int nRows = this->numRowsLocal();
359  unsigned int nCols = this->numCols();
360  double aux = 0.;
361  for (unsigned int i = 0; i < nRows; i++) {
362  for (unsigned int j = 0; j < nCols; j++) {
363  aux = (*this)(i,j);
364  if (aux > value) value = aux;
365  }
366  }
367 
368  return value;
369 }
370 
371 void
372 GslMatrix::cwSet(double value)
373 {
374  unsigned int nRows = this->numRowsLocal();
375  unsigned int nCols = this->numCols();
376 
377  for (unsigned int row = 0; row < nRows; ++row) {
378  for (unsigned int col = 0; col < nCols; ++col) {
379  *gsl_matrix_ptr(m_mat,row,col) = value;
380  }
381  }
382 
383  return;
384 }
385 
386 void
388  unsigned int initialTargetRowId,
389  unsigned int initialTargetColId,
390  const GslMatrix& mat)
391 {
392  queso_require_less_msg(initialTargetRowId, this->numRowsLocal(), "invalid initialTargetRowId");
393 
394  queso_require_less_equal_msg((initialTargetRowId + mat.numRowsLocal()), this->numRowsLocal(), "invalid vec.numRowsLocal()");
395 
396  queso_require_less_msg(initialTargetColId, this->numCols(), "invalid initialTargetColId");
397 
398  queso_require_less_equal_msg((initialTargetColId + mat.numCols()), this->numCols(), "invalid vec.numCols()");
399 
400  for (unsigned int i = 0; i < mat.numRowsLocal(); ++i) {
401  for (unsigned int j = 0; j < mat.numCols(); ++j) {
402  (*this)(initialTargetRowId+i,initialTargetColId+j) = mat(i,j);
403  }
404  }
405 
406  return;
407 }
408 
409 void
411  unsigned int initialTargetRowId,
412  unsigned int initialTargetColId,
413  GslMatrix& mat) const
414 {
415  queso_require_less_msg(initialTargetRowId, this->numRowsLocal(), "invalid initialTargetRowId");
416 
417  queso_require_less_equal_msg((initialTargetRowId + mat.numRowsLocal()), this->numRowsLocal(), "invalid vec.numRowsLocal()");
418 
419  queso_require_less_msg(initialTargetColId, this->numCols(), "invalid initialTargetColId");
420 
421  queso_require_less_equal_msg((initialTargetColId + mat.numCols()), this->numCols(), "invalid vec.numCols()");
422 
423  for (unsigned int i = 0; i < mat.numRowsLocal(); ++i) {
424  for (unsigned int j = 0; j < mat.numCols(); ++j) {
425  mat(i,j) = (*this)(initialTargetRowId+i,initialTargetColId+j) ;
426  }
427  }
428 
429  return;
430 }
431 
432 int
434 {
435  int iRC;
436  //std::cout << "Calling gsl_linalg_cholesky_decomp()..." << std::endl;
437  gsl_error_handler_t* oldHandler;
438  oldHandler = gsl_set_error_handler_off();
439  iRC = gsl_linalg_cholesky_decomp(m_mat);
440  if (iRC != 0) {
441  std::cerr << "In GslMatrix::chol()"
442  << ": iRC = " << iRC
443  << ", gsl error message = " << gsl_strerror(iRC)
444  << std::endl;
445  std::cerr << "Here is the offending matrix: " << std::endl;
446  std::cerr << *this << std::endl;
447  }
448  gsl_set_error_handler(oldHandler);
449  //std::cout << "Returned from gsl_linalg_cholesky_decomp() with iRC = " << iRC << std::endl;
450  UQ_RC_MACRO(iRC, // Yes, *not* a fatal check on RC
451  m_env.worldRank(),
452  "GslMatrix::chol()",
453  "matrix is not positive definite",
455 
456  return iRC;
457 }
458 
459 int
460 GslMatrix::svd(GslMatrix& matU, GslVector& vecS, GslMatrix& matVt) const
461 {
462  unsigned int nRows = this->numRowsLocal();
463  unsigned int nCols = this->numCols();
464 
465  queso_require_msg(!((matU.numRowsLocal() != nRows) || (matU.numCols() != nCols)), "invalid matU");
466 
467  queso_require_equal_to_msg(vecS.sizeLocal(), nCols, "invalid vecS");
468 
469  queso_require_msg(!((matVt.numRowsLocal() != nCols) || (matVt.numCols() != nCols)), "invalid matVt");
470 
471  int iRC = internalSvd();
472 
473  matU = *m_svdUmat;
474  vecS = *m_svdSvec;
475  matVt = *m_svdVTmat;
476 
477  return iRC;
478 }
479 
480 int
481 GslMatrix::svdSolve(const GslVector& rhsVec, GslVector& solVec) const
482 {
483  unsigned int nRows = this->numRowsLocal();
484  unsigned int nCols = this->numCols();
485 
486  queso_require_equal_to_msg(rhsVec.sizeLocal(), nRows, "invalid rhsVec");
487 
488  queso_require_equal_to_msg(solVec.sizeLocal(), nCols, "invalid solVec");
489 
490  int iRC = internalSvd();
491 
492  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 5)) {
493  *m_env.subDisplayFile() << "In GslMatrix::svdSolve():"
494  << "\n this->numRowsLocal() = " << this->numRowsLocal()
495  << ", this->numCols() = " << this->numCols()
496  << "\n m_svdUmat->numRowsLocal() = " << m_svdUmat->numRowsLocal()
497  << ", m_svdUmat->numCols() = " << m_svdUmat->numCols()
498  << "\n m_svdVmat->numRowsLocal() = " << m_svdVmat->numRowsLocal()
499  << ", m_svdVmat->numCols() = " << m_svdVmat->numCols()
500  << "\n m_svdSvec->sizeLocal() = " << m_svdSvec->sizeLocal()
501  << "\n rhsVec.sizeLocal() = " << rhsVec.sizeLocal()
502  << "\n solVec.sizeLocal() = " << solVec.sizeLocal()
503  << std::endl;
504  }
505 
506  if (iRC == 0) iRC = gsl_linalg_SV_solve(m_svdUmat->data(), m_svdVmat->data(), m_svdSvec->data(), rhsVec.data(), solVec.data());
507 
508  return iRC;
509 }
510 
511 int
512 GslMatrix::svdSolve(const GslMatrix& rhsMat, GslMatrix& solMat) const
513 {
514  unsigned int nRows = this->numRowsLocal();
515  unsigned int nCols = this->numCols();
516 
517  queso_require_equal_to_msg(rhsMat.numRowsLocal(), nRows, "invalid rhsMat");
518 
519  queso_require_equal_to_msg(solMat.numRowsLocal(), nCols, "invalid solMat");
520 
521  queso_require_equal_to_msg(rhsMat.numCols(), solMat.numCols(), "rhsMat and solMat are not compatible");
522 
523  GslVector rhsVec(m_env,rhsMat.map());
524  GslVector solVec(m_env,solMat.map());
525  int iRC = 0;
526  for (unsigned int j = 0; j < rhsMat.numCols(); ++j) {
527  rhsVec = rhsMat.getColumn(j);
528  iRC = this->svdSolve(rhsVec, solVec);
529  if (iRC) break;
530  solMat.setColumn(j,solVec);
531  }
532 
533  return iRC;
534 }
535 
536 const GslMatrix&
538 {
539  int iRC = 0;
540  iRC = internalSvd();
541  if (iRC) {}; // just to remove compiler warning
542 
543  return *m_svdUmat;
544 }
545 
546 const GslMatrix&
548 {
549  int iRC = 0;
550  iRC = internalSvd();
551  if (iRC) {}; // just to remove compiler warning
552 
553  return *m_svdVmat;
554 }
555 
556 int
558 {
559  int iRC = 0;
560 
561  if (m_svdColMap == NULL) {
562  unsigned int nRows = this->numRowsLocal();
563  unsigned int nCols = this->numCols();
564  queso_require_greater_equal_msg(nRows, nCols, "GSL only supports cases where nRows >= nCols");
565 
566  m_svdColMap = new Map(this->numCols(),0,this->map().Comm()); // see 'VectorSpace<.,.>::newMap()' in src/basic/src/GslVectorSpace.C
567  m_svdUmat = new GslMatrix(*this); // Yes, 'this'
569  m_svdVmat = new GslMatrix(*m_svdSvec);
571 
572  //std::cout << "In GslMatrix::internalSvd()"
573  // << ", calling gsl_linalg_SV_decomp_jacobi()..."
574  // << ": nRows = " << nRows
575  // << ", nCols = " << nCols
576  // << std::endl;
577  struct timeval timevalBegin;
578  gettimeofday(&timevalBegin, NULL);
579  gsl_error_handler_t* oldHandler;
580  oldHandler = gsl_set_error_handler_off();
581 #if 1
582  iRC = gsl_linalg_SV_decomp_jacobi(m_svdUmat->data(), m_svdVmat->data(), m_svdSvec->data());
583 #else
584  GslVector vecWork(*m_svdSvec );
585  iRC = gsl_linalg_SV_decomp(m_svdUmat->data(), m_svdVmat->data(), m_svdSvec->data(), vecWork.data());
586 #endif
587  if (iRC != 0) {
588  std::cerr << "In GslMatrix::internalSvd()"
589  << ": iRC = " << iRC
590  << ", gsl error message = " << gsl_strerror(iRC)
591  << std::endl;
592  }
593  gsl_set_error_handler(oldHandler);
594 
595  struct timeval timevalNow;
596  gettimeofday(&timevalNow, NULL);
597  //std::cout << "In GslMatrix::internalSvd()"
598  // << ": returned from gsl_linalg_SV_decomp_jacobi() with iRC = " << iRC
599  // << " after " << timevalNow.tv_sec - timevalBegin.tv_sec
600  // << " seconds"
601  // << std::endl;
602  UQ_RC_MACRO(iRC, // Yes, *not* a fatal check on RC
603  m_env.worldRank(),
604  "GslMatrix::internalSvd()",
605  "matrix svd failed",
608  }
609 
610  return iRC;
611 }
612 
613 
614 
615 void
616 GslMatrix::zeroLower(bool includeDiagonal)
617 {
618  unsigned int nRows = this->numRowsLocal();
619  unsigned int nCols = this->numCols();
620 
621  queso_require_equal_to_msg(nRows, nCols, "routine works only for square matrices");
622 
623  this->resetLU();
624 
625  if (includeDiagonal) {
626  for (unsigned int i = 0; i < nRows; i++) {
627  for (unsigned int j = 0; j <= i; j++) {
628  (*this)(i,j) = 0.;
629  }
630  }
631  }
632  else {
633  for (unsigned int i = 0; i < nRows; i++) {
634  for (unsigned int j = 0; j < i; j++) {
635  (*this)(i,j) = 0.;
636  }
637  }
638  }
639 
640  return;
641 }
642 
643 void
644 GslMatrix::zeroUpper(bool includeDiagonal)
645 {
646  unsigned int nRows = this->numRowsLocal();
647  unsigned int nCols = this->numCols();
648 
649  queso_require_equal_to_msg(nRows, nCols, "routine works only for square matrices");
650 
651  this->resetLU();
652 
653  if (includeDiagonal) {
654  for (unsigned int i = 0; i < nRows; i++) {
655  for (unsigned int j = i; j < nCols; j++) {
656  (*this)(i,j) = 0.;
657  }
658  }
659  }
660  else {
661  for (unsigned int i = 0; i < nRows; i++) {
662  for (unsigned int j = (i+1); j < nCols; j++) {
663  (*this)(i,j) = 0.;
664  }
665  }
666  }
667 
668  return;
669 }
670 
671 void
672 GslMatrix::filterSmallValues(double thresholdValue)
673 {
674  unsigned int nRows = this->numRowsLocal();
675  unsigned int nCols = this->numCols();
676  for (unsigned int i = 0; i < nRows; ++i) {
677  for (unsigned int j = 0; j < nCols; ++j) {
678  double aux = (*this)(i,j);
679  // If 'thresholdValue' is negative, no values will be filtered
680  if ((aux < 0. ) &&
681  (-thresholdValue < aux)) {
682  (*this)(i,j) = 0.;
683  }
684  if ((aux > 0. ) &&
685  (thresholdValue > aux)) {
686  (*this)(i,j) = 0.;
687  }
688  }
689  }
690 
691  return;
692 }
693 
694 void
695 GslMatrix::filterLargeValues(double thresholdValue)
696 {
697  unsigned int nRows = this->numRowsLocal();
698  unsigned int nCols = this->numCols();
699  for (unsigned int i = 0; i < nRows; ++i) {
700  for (unsigned int j = 0; j < nCols; ++j) {
701  double aux = (*this)(i,j);
702  // If 'thresholdValue' is negative, no values will be filtered
703  if ((aux < 0. ) &&
704  (-thresholdValue > aux)) {
705  (*this)(i,j) = 0.;
706  }
707  if ((aux > 0. ) &&
708  (thresholdValue < aux)) {
709  (*this)(i,j) = 0.;
710  }
711  }
712  }
713 
714  return;
715 }
716 
717 GslMatrix
719 {
720  unsigned int nRows = this->numRowsLocal();
721  unsigned int nCols = this->numCols();
722 
723  queso_require_equal_to_msg(nRows, nCols, "routine works only for square matrices");
724 
725  GslMatrix mat(m_env,m_map,nCols);
726  for (unsigned int row = 0; row < nRows; ++row) {
727  for (unsigned int col = 0; col < nCols; ++col) {
728  mat(row,col) = (*this)(col,row);
729  }
730  }
731 
732  return mat;
733 }
734 
735 GslMatrix
737 {
738  unsigned int nRows = this->numRowsLocal();
739  unsigned int nCols = this->numCols();
740 
741  queso_require_equal_to_msg(nRows, nCols, "matrix is not square");
742 
743  if (m_inverse == NULL) {
744  m_inverse = new GslMatrix(m_env,m_map,nCols);
745  GslVector unitVector(m_env,m_map);
746  unitVector.cwSet(0.);
747  GslVector multVector(m_env,m_map);
748  for (unsigned int j = 0; j < nCols; ++j) {
749  if (j > 0) unitVector[j-1] = 0.;
750  unitVector[j] = 1.;
751  this->invertMultiply(unitVector, multVector);
752  for (unsigned int i = 0; i < nRows; ++i) {
753  (*m_inverse)(i,j) = multVector[i];
754  }
755  }
756  }
757  if (m_env.checkingLevel() >= 1) {
758  *m_env.subDisplayFile() << "CHECKING In GslMatrix::inverse()"
759  << ": M.lnDet = " << this->lnDeterminant()
760  << ", M^{-1}.lnDet = " << m_inverse->lnDeterminant()
761  << std::endl;
762  }
763  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 5)) {
764  *m_env.subDisplayFile() << "In GslMatrix::inverse():"
765  << "\n M = " << *this
766  << "\n M^{-1} = " << *m_inverse
767  << "\n M*M^{-1} = " << (*this)*(*m_inverse)
768  << "\n M^{-1}*M = " << (*m_inverse)*(*this)
769  << std::endl;
770  }
771 
772  return *m_inverse;
773 }
774 
775 void
777  unsigned int initialTargetRowId,
778  unsigned int initialTargetColId,
779  const std::vector<const GslMatrix* >& matrices,
780  bool checkForExactNumRowsMatching,
781  bool checkForExactNumColsMatching)
782 {
783  unsigned int sumNumRowsLocals = 0;
784  unsigned int sumNumCols = 0;
785  for (unsigned int i = 0; i < matrices.size(); ++i) {
786  sumNumRowsLocals += matrices[i]->numRowsLocal();
787  sumNumCols += matrices[i]->numCols();
788  }
789  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRowsLocals), "too big number of rows");
790  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRowsLocals), "inconsistent number of rows");
791  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + sumNumCols), "too big number of cols");
792  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + sumNumCols), "inconsistent number of cols");
793 
794  unsigned int cumulativeRowId = 0;
795  unsigned int cumulativeColId = 0;
796  for (unsigned int i = 0; i < matrices.size(); ++i) {
797  unsigned int nRows = matrices[i]->numRowsLocal();
798  unsigned int nCols = matrices[i]->numCols();
799  for (unsigned int rowId = 0; rowId < nRows; ++rowId) {
800  for (unsigned int colId = 0; colId < nCols; ++colId) {
801  (*this)(initialTargetRowId + cumulativeRowId + rowId, initialTargetColId + cumulativeColId + colId) = (*(matrices[i]))(rowId,colId);
802  }
803  }
804  cumulativeRowId += nRows;
805  cumulativeColId += nCols;
806  }
807 
808  return;
809 }
810 
811 void
813  unsigned int initialTargetRowId,
814  unsigned int initialTargetColId,
815  const std::vector<GslMatrix* >& matrices,
816  bool checkForExactNumRowsMatching,
817  bool checkForExactNumColsMatching)
818 {
819  unsigned int sumNumRowsLocals = 0;
820  unsigned int sumNumCols = 0;
821  for (unsigned int i = 0; i < matrices.size(); ++i) {
822  sumNumRowsLocals += matrices[i]->numRowsLocal();
823  sumNumCols += matrices[i]->numCols();
824  }
825  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRowsLocals), "too big number of rows");
826  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRowsLocals), "inconsistent number of rows");
827  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + sumNumCols), "too big number of cols");
828  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + sumNumCols), "inconsistent number of cols");
829 
830  unsigned int cumulativeRowId = 0;
831  unsigned int cumulativeColId = 0;
832  for (unsigned int i = 0; i < matrices.size(); ++i) {
833  unsigned int nRows = matrices[i]->numRowsLocal();
834  unsigned int nCols = matrices[i]->numCols();
835  for (unsigned int rowId = 0; rowId < nRows; ++rowId) {
836  for (unsigned int colId = 0; colId < nCols; ++colId) {
837  (*this)(initialTargetRowId + cumulativeRowId + rowId, initialTargetColId + cumulativeColId + colId) = (*(matrices[i]))(rowId,colId);
838  }
839  }
840  cumulativeRowId += nRows;
841  cumulativeColId += nCols;
842  }
843 
844  return;
845 }
846 
847 void
849  unsigned int initialTargetRowId,
850  unsigned int initialTargetColId,
851  const std::vector<const GslMatrix* >& matrices,
852  bool checkForExactNumRowsMatching,
853  bool checkForExactNumColsMatching)
854 {
855  unsigned int sumNumCols = 0;
856  for (unsigned int i = 0; i < matrices.size(); ++i) {
857  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + matrices[i]->numRowsLocal()), "too big number of rows");
858  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + matrices[i]->numRowsLocal()), "inconsistent number of rows");
859  sumNumCols += matrices[i]->numCols();
860  }
861  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + sumNumCols), "too big number of cols");
862  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + sumNumCols), "inconsistent number of cols");
863 
864  unsigned int cumulativeColId = 0;
865  for (unsigned int i = 0; i < matrices.size(); ++i) {
866  unsigned int nRows = matrices[i]->numRowsLocal();
867  unsigned int nCols = matrices[i]->numCols();
868  for (unsigned int rowId = 0; rowId < nRows; ++rowId) {
869  for (unsigned int colId = 0; colId < nCols; ++colId) {
870  (*this)(initialTargetRowId + rowId, initialTargetColId + cumulativeColId + colId) = (*(matrices[i]))(rowId,colId);
871  }
872  }
873  cumulativeColId += nCols;
874  }
875 
876  return;
877 }
878 
879 void
881  unsigned int initialTargetRowId,
882  unsigned int initialTargetColId,
883  const std::vector<GslMatrix* >& matrices,
884  bool checkForExactNumRowsMatching,
885  bool checkForExactNumColsMatching)
886 {
887  unsigned int sumNumCols = 0;
888  for (unsigned int i = 0; i < matrices.size(); ++i) {
889  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + matrices[i]->numRowsLocal()), "too big number of rows");
890  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + matrices[i]->numRowsLocal()), "inconsistent number of rows");
891  sumNumCols += matrices[i]->numCols();
892  }
893  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + sumNumCols), "too big number of cols");
894  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + sumNumCols), "inconsistent number of cols");
895 
896  unsigned int cumulativeColId = 0;
897  for (unsigned int i = 0; i < matrices.size(); ++i) {
898  unsigned int nRows = matrices[i]->numRowsLocal();
899  unsigned int nCols = matrices[i]->numCols();
900  for (unsigned int rowId = 0; rowId < nRows; ++rowId) {
901  for (unsigned int colId = 0; colId < nCols; ++colId) {
902  (*this)(initialTargetRowId + rowId, initialTargetColId + cumulativeColId + colId) = (*(matrices[i]))(rowId,colId);
903  }
904  }
905  cumulativeColId += nCols;
906  }
907 
908  return;
909 }
910 
911 void
913  unsigned int initialTargetRowId,
914  unsigned int initialTargetColId,
915  const std::vector<const GslMatrix* >& matrices,
916  bool checkForExactNumRowsMatching,
917  bool checkForExactNumColsMatching)
918 {
919  unsigned int sumNumRows = 0;
920  for (unsigned int i = 0; i < matrices.size(); ++i) {
921  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + matrices[i]->numCols()), "too big number of cols");
922  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + matrices[i]->numCols()), "inconsistent number of cols");
923  sumNumRows += matrices[i]->numRowsLocal();
924  }
925  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRows), "too big number of rows");
926  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRows), "inconsistent number of rows");
927 
928  unsigned int cumulativeRowId = 0;
929  for (unsigned int i = 0; i < matrices.size(); ++i) {
930  unsigned int nRows = matrices[i]->numRowsLocal();
931  unsigned int nCols = matrices[i]->numCols();
932  for (unsigned int rowId = 0; rowId < nRows; ++rowId) {
933  for (unsigned int colId = 0; colId < nCols; ++colId) {
934  (*this)(initialTargetRowId + cumulativeRowId + rowId, initialTargetColId + colId) = (*(matrices[i]))(rowId,colId);
935  }
936  }
937  cumulativeRowId += nRows;
938  }
939 
940  return;
941 }
942 
943 void
945  unsigned int initialTargetRowId,
946  unsigned int initialTargetColId,
947  const std::vector<GslMatrix* >& matrices,
948  bool checkForExactNumRowsMatching,
949  bool checkForExactNumColsMatching)
950 {
951  unsigned int sumNumRows = 0;
952  for (unsigned int i = 0; i < matrices.size(); ++i) {
953  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + matrices[i]->numCols()), "too big number of cols");
954  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + matrices[i]->numCols()), "inconsistent number of cols");
955  sumNumRows += matrices[i]->numRowsLocal();
956  }
957  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRows), "too big number of rows");
958  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + sumNumRows), "inconsistent number of rows");
959 
960  unsigned int cumulativeRowId = 0;
961  for (unsigned int i = 0; i < matrices.size(); ++i) {
962  unsigned int nRows = matrices[i]->numRowsLocal();
963  unsigned int nCols = matrices[i]->numCols();
964  for (unsigned int rowId = 0; rowId < nRows; ++rowId) {
965  for (unsigned int colId = 0; colId < nCols; ++colId) {
966  (*this)(initialTargetRowId + cumulativeRowId + rowId, initialTargetColId + colId) = (*(matrices[i]))(rowId,colId);
967  }
968  }
969  cumulativeRowId += nRows;
970  }
971 
972  return;
973 }
974 
975 void
977  unsigned int initialTargetRowId,
978  unsigned int initialTargetColId,
979  const GslMatrix& mat1,
980  const GslMatrix& mat2,
981  bool checkForExactNumRowsMatching,
982  bool checkForExactNumColsMatching)
983 {
984  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + (mat1.numRowsLocal() * mat2.numRowsLocal())), "too big number of rows");
985  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + (mat1.numRowsLocal() * mat2.numRowsLocal())), "inconsistent number of rows");
986  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + (mat1.numCols() * mat2.numCols())), "too big number of columns");
987  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + (mat1.numCols() * mat2.numCols())), "inconsistent number of columns");
988 
989  for (unsigned int rowId1 = 0; rowId1 < mat1.numRowsLocal(); ++rowId1) {
990  for (unsigned int colId1 = 0; colId1 < mat1.numCols(); ++colId1) {
991  double multiplicativeFactor = mat1(rowId1,colId1);
992  unsigned int targetRowId = rowId1 * mat2.numRowsLocal();
993  unsigned int targetColId = colId1 * mat2.numCols();
994  for (unsigned int rowId2 = 0; rowId2 < mat2.numRowsLocal(); ++rowId2) {
995  for (unsigned int colId2 = 0; colId2 < mat2.numCols(); ++colId2) {
996  (*this)(initialTargetRowId + targetRowId + rowId2, initialTargetColId + targetColId + colId2) = multiplicativeFactor * mat2(rowId2,colId2);
997  }
998  }
999  }
1000  }
1001 
1002  return;
1003 }
1004 
1005 void
1007  unsigned int initialTargetRowId,
1008  unsigned int initialTargetColId,
1009  const GslMatrix& mat1,
1010  const GslVector& vec2,
1011  bool checkForExactNumRowsMatching,
1012  bool checkForExactNumColsMatching)
1013 {
1014  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + (mat1.numRowsLocal() * vec2.sizeLocal())), "too big number of rows");
1015  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + (mat1.numRowsLocal() * vec2.sizeLocal())), "inconsistent number of rows");
1016  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + (mat1.numCols() * 1)), "too big number of columns");
1017  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + (mat1.numCols() * 1)), "inconsistent number of columns");
1018 
1019  for (unsigned int rowId1 = 0; rowId1 < mat1.numRowsLocal(); ++rowId1) {
1020  for (unsigned int colId1 = 0; colId1 < mat1.numCols(); ++colId1) {
1021  double multiplicativeFactor = mat1(rowId1,colId1);
1022  unsigned int targetRowId = rowId1 * vec2.sizeLocal();
1023  unsigned int targetColId = colId1 * 1;
1024  for (unsigned int rowId2 = 0; rowId2 < vec2.sizeLocal(); ++rowId2) {
1025  for (unsigned int colId2 = 0; colId2 < 1; ++colId2) {
1026  (*this)(initialTargetRowId + targetRowId + rowId2, initialTargetColId + targetColId + colId2) = multiplicativeFactor * vec2[rowId2];
1027  }
1028  }
1029  }
1030  }
1031 
1032 
1033  return;
1034 }
1035 
1036 void
1038  unsigned int initialTargetRowId,
1039  unsigned int initialTargetColId,
1040  const GslMatrix& mat,
1041  bool checkForExactNumRowsMatching,
1042  bool checkForExactNumColsMatching)
1043 {
1044  unsigned int nRows = mat.numRowsLocal();
1045  unsigned int nCols = mat.numCols();
1046  queso_require_greater_equal_msg(this->numRowsLocal(), (initialTargetRowId + nCols), "too big number of rows");
1047  if (checkForExactNumRowsMatching) queso_require_equal_to_msg(this->numRowsLocal(), (initialTargetRowId + nCols), "inconsistent number of rows");
1048  queso_require_greater_equal_msg(this->numCols(), (initialTargetColId + nRows), "too big number of cols");
1049  if (checkForExactNumColsMatching) queso_require_equal_to_msg(this->numCols(), (initialTargetColId + nRows), "inconsistent number of cols");
1050 
1051  for (unsigned int row = 0; row < nRows; ++row) {
1052  for (unsigned int col = 0; col < nCols; ++col) {
1053  (*this)(initialTargetRowId + col, initialTargetColId + row) = mat(row,col);
1054  }
1055  }
1056 
1057  return;
1058 }
1059 
1060 double
1062 {
1063  if (m_determinant == -INFINITY) {
1064  if (m_LU == NULL) {
1065  GslVector tmpB(m_env,m_map);
1066  GslVector tmpX(m_env,m_map);
1067  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1068  *m_env.subDisplayFile() << "In GslMatrix::determinant()"
1069  << ": before 'this->invertMultiply()'"
1070  << std::endl;
1071  }
1072  this->invertMultiply(tmpB,tmpX);
1073  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1074  *m_env.subDisplayFile() << "In GslMatrix::determinant()"
1075  << ": after 'this->invertMultiply()'"
1076  << std::endl;
1077  }
1078  }
1079  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1080  *m_env.subDisplayFile() << "In GslMatrix::determinant()"
1081  << ": before 'gsl_linalg_LU_det()'"
1082  << std::endl;
1083  }
1084  m_determinant = gsl_linalg_LU_det(m_LU,m_signum);
1085  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1086  *m_env.subDisplayFile() << "In GslMatrix::determinant()"
1087  << ": after 'gsl_linalg_LU_det()'"
1088  << std::endl;
1089  }
1090  m_lnDeterminant = gsl_linalg_LU_lndet(m_LU);
1091  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1092  *m_env.subDisplayFile() << "In GslMatrix::determinant()"
1093  << ": after 'gsl_linalg_LU_lndet()'"
1094  << std::endl;
1095  }
1096  }
1097 
1098  return m_determinant;
1099 }
1100 
1101 double
1103 {
1104  if (m_lnDeterminant == -INFINITY) {
1105  if (m_LU == NULL) {
1106  GslVector tmpB(m_env,m_map);
1107  GslVector tmpX(m_env,m_map);
1108  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1109  *m_env.subDisplayFile() << "In GslMatrix::lnDeterminant()"
1110  << ": before 'this->invertMultiply()'"
1111  << std::endl;
1112  }
1113  this->invertMultiply(tmpB,tmpX);
1114  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1115  *m_env.subDisplayFile() << "In GslMatrix::lnDeterminant()"
1116  << ": after 'this->invertMultiply()'"
1117  << std::endl;
1118  }
1119  }
1120  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1121  *m_env.subDisplayFile() << "In GslMatrix::lnDeterminant()"
1122  << ": before 'gsl_linalg_LU_det()'"
1123  << std::endl;
1124  }
1125  m_determinant = gsl_linalg_LU_det(m_LU,m_signum);
1126  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1127  *m_env.subDisplayFile() << "In GslMatrix::lnDeterminant()"
1128  << ": after 'gsl_linalg_LU_det()'"
1129  << std::endl;
1130  }
1131  m_lnDeterminant = gsl_linalg_LU_lndet(m_LU);
1132  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1133  *m_env.subDisplayFile() << "In GslMatrix::lnDeterminant()"
1134  << ": before 'gsl_linalg_LU_lndet()'"
1135  << std::endl;
1136  }
1137  }
1138 
1139  return m_lnDeterminant;
1140 }
1141 
1142 unsigned int
1143 GslMatrix::rank(double absoluteZeroThreshold, double relativeZeroThreshold) const
1144 {
1145  int iRC = 0;
1146  iRC = internalSvd();
1147  if (iRC) {}; // just to remove compiler warning
1148 
1149  GslVector relativeVec(*m_svdSvec);
1150  if (relativeVec[0] > 0.) {
1151  relativeVec = (1./relativeVec[0])*relativeVec;
1152  }
1153 
1154  unsigned int rankValue = 0;
1155  for (unsigned int i = 0; i < relativeVec.sizeLocal(); ++i) {
1156  if (( (*m_svdSvec)[i] >= absoluteZeroThreshold ) &&
1157  ( relativeVec [i] >= relativeZeroThreshold )) {
1158  rankValue += 1;
1159  }
1160  }
1161 
1162  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 3)) {
1163  *m_env.subDisplayFile() << "In GslMatrix::rank()"
1164  << ": this->numRowsLocal() = " << this->numRowsLocal()
1165  << ", this->numCols() = " << this->numCols()
1166  << ", absoluteZeroThreshold = " << absoluteZeroThreshold
1167  << ", relativeZeroThreshold = " << relativeZeroThreshold
1168  << ", rankValue = " << rankValue
1169  << ", m_svdSvec = " << *m_svdSvec
1170  << ", relativeVec = " << relativeVec
1171  << std::endl;
1172  }
1173 
1174  return rankValue;
1175 }
1176 
1177 GslVector
1179  const GslVector& x) const
1180 {
1181  queso_require_equal_to_msg(this->numCols(), x.sizeLocal(), "matrix and vector have incompatible sizes");
1182 
1183  GslVector y(m_env,m_map);
1184  this->multiply(x,y);
1185 
1186  return y;
1187 }
1188 
1189 void
1191  const GslVector& x,
1192  GslVector& y) const
1193 {
1194  queso_require_equal_to_msg(this->numCols(), x.sizeLocal(), "matrix and x have incompatible sizes");
1195 
1196  queso_require_equal_to_msg(this->numRowsLocal(), y.sizeLocal(), "matrix and y have incompatible sizes");
1197 
1198  unsigned int sizeX = this->numCols();
1199  unsigned int sizeY = this->numRowsLocal();
1200  for (unsigned int i = 0; i < sizeY; ++i) {
1201  double value = 0.;
1202  for (unsigned int j = 0; j < sizeX; ++j) {
1203  value += (*this)(i,j)*x[j];
1204  }
1205  y[i] = value;
1206  }
1207 
1208  return;
1209 }
1210 
1211 GslVector
1213  const GslVector& b) const
1214 {
1215  queso_require_equal_to_msg(this->numCols(), b.sizeLocal(), "matrix and rhs have incompatible sizes");
1216 
1217  GslVector x(m_env,m_map);
1218  this->invertMultiply(b,x);
1219 
1220  return x;
1221 }
1222 
1223 void
1225  const GslVector& b,
1226  GslVector& x) const
1227 {
1228  queso_require_equal_to_msg(this->numCols(), b.sizeLocal(), "matrix and rhs have incompatible sizes");
1229 
1230  queso_require_equal_to_msg(x.sizeLocal(), b.sizeLocal(), "solution and rhs have incompatible sizes");
1231 
1232  int iRC;
1233  if (m_LU == NULL) {
1234  queso_require_msg(!(m_permutation), "m_permutation should be NULL");
1235 
1236  m_LU = gsl_matrix_calloc(this->numRowsLocal(),this->numCols());
1237  queso_require_msg(m_LU, "gsl_matrix_calloc() failed");
1238 
1239  iRC = gsl_matrix_memcpy(m_LU, m_mat);
1240  queso_require_msg(!(iRC), "gsl_matrix_memcpy() failed");
1241 
1242  m_permutation = gsl_permutation_calloc(numCols());
1243  queso_require_msg(m_permutation, "gsl_permutation_calloc() failed");
1244 
1245  if (m_inDebugMode) {
1246  std::cout << "In GslMatrix::invertMultiply()"
1247  << ": before LU decomposition, m_LU = ";
1248  gsl_matrix_fprintf(stdout, m_LU, "%f");
1249  std::cout << std::endl;
1250  }
1251 
1252  gsl_error_handler_t* oldHandler;
1253  oldHandler = gsl_set_error_handler_off();
1254  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1255  *m_env.subDisplayFile() << "In GslMatrix::invertMultiply()"
1256  << ": before 'gsl_linalg_LU_decomp()'"
1257  << std::endl;
1258  }
1259  iRC = gsl_linalg_LU_decomp(m_LU,m_permutation,&m_signum);
1260  if (iRC != 0) {
1261  std::cerr << "In GslMatrix::invertMultiply()"
1262  << ", after gsl_linalg_LU_decomp()"
1263  << ": iRC = " << iRC
1264  << ", gsl error message = " << gsl_strerror(iRC)
1265  << std::endl;
1266  }
1267  gsl_set_error_handler(oldHandler);
1268  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1269  *m_env.subDisplayFile() << "In GslMatrix::invertMultiply()"
1270  << ": after 'gsl_linalg_LU_decomp()'"
1271  << ", IRC = " << iRC
1272  << std::endl;
1273  }
1274  queso_require_msg(!(iRC), "gsl_linalg_LU_decomp() failed");
1275 
1276  if (m_inDebugMode) {
1277  std::cout << "In GslMatrix::invertMultiply()"
1278  << ": after LU decomposition, m_LU = ";
1279  gsl_matrix_fprintf(stdout, m_LU, "%f");
1280  std::cout << std::endl;
1281  }
1282  }
1283 
1284  gsl_error_handler_t* oldHandler;
1285  oldHandler = gsl_set_error_handler_off();
1286  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1287  *m_env.subDisplayFile() << "In GslMatrix::invertMultiply()"
1288  << ": before 'gsl_linalg_LU_solve()'"
1289  << std::endl;
1290  }
1291  iRC = gsl_linalg_LU_solve(m_LU,m_permutation,b.data(),x.data());
1292  if (iRC != 0) {
1293  m_isSingular = true;
1294  std::cerr << "In GslMatrix::invertMultiply()"
1295  << ", after gsl_linalg_LU_solve()"
1296  << ": iRC = " << iRC
1297  << ", gsl error message = " << gsl_strerror(iRC)
1298  << std::endl;
1299  }
1300  gsl_set_error_handler(oldHandler);
1301  if ((m_env.subDisplayFile()) && (m_env.displayVerbosity() >= 99)) {
1302  *m_env.subDisplayFile() << "In GslMatrix::invertMultiply()"
1303  << ": after 'gsl_linalg_LU_solve()'"
1304  << ", IRC = " << iRC
1305  << std::endl;
1306  }
1307 
1308 
1309  // m_env.worldRank(),
1310  // "GslMatrix::invertMultiply()",
1311  // "gsl_linalg_LU_solve() failed");
1312 
1313  if (m_inDebugMode) {
1314  GslVector tmpVec(b - (*this)*x);
1315  std::cout << "In GslMatrix::invertMultiply()"
1316  << ": ||b - Ax||_2 = " << tmpVec.norm2()
1317  << ": ||b - Ax||_2/||b||_2 = " << tmpVec.norm2()/b.norm2()
1318  << std::endl;
1319  }
1320 
1321  return;
1322 }
1323 
1324 GslMatrix
1326 {
1327  GslMatrix X(m_env,m_map,B.numCols());
1328  this->invertMultiply(B,X);
1329 
1330  return X;
1331 }
1332 
1333 void
1335 {
1336  // Sanity Checks
1338  "Matrices B and X are incompatible");
1340  "Matrices B and X are incompatible");
1342  "This and X matrices are incompatible");
1343 
1344  // Some local variables used within the loop.
1345  GslVector b(m_env, m_map);
1346  GslVector x(m_env, m_map);
1347 
1348  for( unsigned int j = 0; j < B.numCols(); ++j )
1349  {
1350  b = B.getColumn( j );
1351 
1352  //invertMultiply will only do the LU once and store it. So we don't
1353  //need to worry about it doing LU multiple times.
1354  x = this->invertMultiply( b );
1355 
1356  X.setColumn( j, x );
1357  }
1358 
1359  return;
1360 }
1361 
1362 GslVector
1364  const GslVector& b) const
1365 {
1366  queso_require_equal_to_msg(this->numCols(), b.sizeLocal(), "matrix and rhs have incompatible sizes");
1367 
1368  GslVector x(m_env,m_map);
1369  this->invertMultiplyForceLU(b,x);
1370 
1371  return x;
1372 }
1373 
1374 void
1376  const GslVector& b,
1377  GslVector& x) const
1378 {
1379  queso_require_equal_to_msg(this->numCols(), b.sizeLocal(), "matrix and rhs have incompatible sizes");
1380 
1381  queso_require_equal_to_msg(x.sizeLocal(), b.sizeLocal(), "solution and rhs have incompatible sizes");
1382 
1383  int iRC;
1384 
1385  if ( m_LU == NULL ) {
1386  queso_require_msg(!(m_permutation), "m_permutation should be NULL");
1387  m_LU = gsl_matrix_calloc(this->numRowsLocal(),this->numCols());
1388  }
1389  queso_require_msg(m_LU, "gsl_matrix_calloc() failed");
1390 
1391  iRC = gsl_matrix_memcpy(m_LU, m_mat);
1392  queso_require_msg(!(iRC), "gsl_matrix_memcpy() failed");
1393 
1394  if( m_permutation == NULL ) m_permutation = gsl_permutation_calloc(numCols());
1395  queso_require_msg(m_permutation, "gsl_permutation_calloc() failed");
1396 
1397  iRC = gsl_linalg_LU_decomp(m_LU,m_permutation,&m_signum);
1398  queso_require_msg(!(iRC), "gsl_linalg_LU_decomp() failed");
1399 
1400  iRC = gsl_linalg_LU_solve(m_LU,m_permutation,b.data(),x.data());
1401  if (iRC != 0) {
1402  m_isSingular = true;
1403  }
1404  queso_require_msg(!(iRC), "gsl_linalg_LU_solve() failed");
1405 
1406  return;
1407 }
1408 
1409 void
1410 GslMatrix::eigen(GslVector& eigenValues, GslMatrix* eigenVectors) const
1411 {
1412  unsigned int n = eigenValues.sizeLocal();
1413 
1414  queso_require_not_equal_to_msg(n, 0, "invalid input vector size");
1415 
1416  if (eigenVectors) {
1417  queso_require_equal_to_msg(eigenValues.sizeLocal(), eigenVectors->numRowsLocal(), "different input vector sizes");
1418  }
1419 
1420  if (eigenVectors == NULL) {
1421  gsl_eigen_symm_workspace* w = gsl_eigen_symm_alloc((size_t) n);
1422  gsl_eigen_symm(m_mat,eigenValues.data(),w);
1423  gsl_eigen_symm_free(w);
1424  }
1425  else {
1426  gsl_eigen_symmv_workspace* w = gsl_eigen_symmv_alloc((size_t) n);
1427  gsl_eigen_symmv(m_mat,eigenValues.data(),eigenVectors->m_mat,w);
1428  gsl_eigen_symmv_sort(eigenValues.data(),eigenVectors->m_mat,GSL_EIGEN_SORT_VAL_ASC);
1429  gsl_eigen_symmv_free(w);
1430  }
1431 
1432  return;
1433 }
1434 
1435 void
1436 GslMatrix::largestEigen(double& eigenValue, GslVector& eigenVector) const
1437 {
1438 
1439  // Sanity Checks
1440  unsigned int n = eigenVector.sizeLocal();
1441 
1442  queso_require_not_equal_to_msg(n, 0, "invalid input vector size");
1443 
1444  /* The following notation is used:
1445  z = vector used in iteration that ends up being the eigenvector corresponding to the
1446  largest eigenvalue
1447  w = vector used in iteration that we extract the largest eigenvalue from. */
1448 
1449  // Some parameters associated with the algorithm
1450  // TODO: Do we want to add the ability to have these set by the user?
1451  const unsigned int max_num_iterations = 10000;
1452  const double tolerance = 1.0e-13;
1453 
1454  // Create temporary working vectors.
1455  // TODO: We shouldn't have to use these - we ought to be able to work directly
1456  // TODO: with eigenValue and eigenVector. Optimize this once we have regression
1457  // TODO: tests.
1458  GslVector z(m_env, m_map, 1.0 ); // Needs to be initialized to 1.0
1459  GslVector w(m_env, m_map);
1460 
1461  // Some variables we use inside the loop.
1462  int index;
1463  double residual;
1464  double lambda;
1465 
1466  for( unsigned int k = 0; k < max_num_iterations; ++k )
1467  {
1468  w = (*this) * z;
1469 
1470  // For this algorithm, it's crucial to get the maximum in
1471  // absolute value, but then to normalize by the actual value
1472  // and *not* the absolute value.
1473  index = (w.abs()).getMaxValueIndex();
1474 
1475  lambda = w[index];
1476 
1477  z = (1.0/lambda) * w;
1478 
1479  // Here we use the norm of the residual as our convergence check:
1480  // norm( A*x - \lambda*x )
1481  residual = ( (*this)*z - lambda*z ).norm2();
1482 
1483  if( residual < tolerance )
1484  {
1485  eigenValue = lambda;
1486 
1487  // TODO: Do we want to normalize this so eigenVector.norm2() = 1?
1488  eigenVector = z;
1489  return;
1490  }
1491 
1492  }
1493 
1494  // If we reach this point, then we didn't converge. Print error message
1495  // to this effect.
1496  // TODO: We know we failed at this point - need a way to not need a test
1497  // TODO: Just specify the error.
1498  queso_require_less_msg(residual, tolerance, "Maximum num iterations exceeded");
1499 
1500 
1501  return;
1502 }
1503 
1504 void
1505 GslMatrix::smallestEigen(double& eigenValue, GslVector& eigenVector) const
1506 {
1507  // Sanity Checks
1508  unsigned int n = eigenVector.sizeLocal();
1509 
1510  queso_require_not_equal_to_msg(n, 0, "invalid input vector size");
1511 
1512  /* The following notation is used:
1513  z = vector used in iteration that ends up being the eigenvector corresponding to the
1514  largest eigenvalue
1515  w = vector used in iteration that we extract the largest eigenvalue from. */
1516 
1517  // Some parameters associated with the algorithm
1518  // TODO: Do we want to add the ability to have these set by the user?
1519  const unsigned int max_num_iterations = 1000;
1520  const double tolerance = 1.0e-13;
1521 
1522  // Create temporary working vectors.
1523  // TODO: We shouldn't have to use these - we ought to be able to work directly
1524  // TODO: with eigenValue and eigenVector. Optimize this once we have regression
1525  // TODO: tests.
1526  GslVector z(m_env, m_map, 1.0 ); // Needs to be initialized to 1.0
1527  GslVector w(m_env, m_map);
1528 
1529  // Some variables we use inside the loop.
1530  int index;
1531  double residual;
1532  double one_over_lambda;
1533  double lambda;
1534 
1535  for( unsigned int k = 0; k < max_num_iterations; ++k )
1536  {
1537  w = (*this).invertMultiplyForceLU( z );
1538 
1539  // For this algorithm, it's crucial to get the maximum in
1540  // absolute value, but then to normalize by the actual value
1541  // and *not* the absolute value.
1542  index = (w.abs()).getMaxValueIndex();
1543 
1544  // Remember: Inverse power method solves for max 1/lambda ==> lambda smallest
1545  one_over_lambda = w[index];
1546 
1547  lambda = 1.0/one_over_lambda;
1548 
1549  z = lambda * w;
1550 
1551  // Here we use the norm of the residual as our convergence check:
1552  // norm( A*x - \lambda*x )
1553  residual = ( (*this)*z - lambda*z ).norm2();
1554 
1555  if( residual < tolerance )
1556  {
1557  eigenValue = lambda;
1558 
1559  // TODO: Do we want to normalize this so eigenVector.norm2() = 1?
1560  eigenVector = z;
1561  return;
1562  }
1563 
1564  }
1565 
1566  // If we reach this point, then we didn't converge. Print error message
1567  // to this effect.
1568  // TODO: We know we failed at this point - need a way to not need a test
1569  // TODO: Just specify the error.
1570  queso_require_less_msg(residual, tolerance, "Maximum num iterations exceeded");
1571 
1572  return;
1573 }
1574 
1575 void
1576 GslMatrix::getColumn(unsigned int column_num, GslVector& column) const
1577 {
1578  // Sanity checks
1579  queso_require_less_msg(column_num, this->numCols(), "Specified row number not within range");
1580 
1581  queso_require_equal_to_msg(column.sizeLocal(), this->numRowsLocal(), "column vector not same size as this matrix");
1582 
1583  // Temporary working vector
1584  gsl_vector* gsl_column = gsl_vector_alloc( column.sizeLocal() );
1585 
1586  int error_code = gsl_matrix_get_col( gsl_column, m_mat, column_num );
1587  queso_require_equal_to_msg(error_code, 0, "gsl_matrix_get_col failed");
1588 
1589  // Copy column from gsl matrix into our GslVector object
1590  for( unsigned int i = 0; i < column.sizeLocal(); ++i )
1591  {
1592  column[i] = gsl_vector_get( gsl_column, i );
1593  }
1594 
1595  // Clean up gsl temporaries
1596  gsl_vector_free( gsl_column );
1597 
1598  return;
1599 }
1600 
1601 void
1602 GslMatrix::getRow(unsigned int row_num, GslVector& row) const
1603 {
1604  // Sanity checks
1605  queso_require_less_msg(row_num, this->numRowsLocal(), "Specified row number not within range");
1606 
1607  queso_require_equal_to_msg(row.sizeLocal(), this->numCols(), "row vector not same size as this matrix");
1608 
1609  // Temporary working vector
1610  gsl_vector* gsl_row = gsl_vector_alloc( row.sizeLocal() );
1611 
1612  int error_code = gsl_matrix_get_row( gsl_row, m_mat, row_num );
1613  queso_require_equal_to_msg(error_code, 0, "gsl_matrix_get_row failed");
1614 
1615  // Copy row from gsl matrix into our GslVector object
1616  for( unsigned int i = 0; i < row.sizeLocal(); ++i )
1617  {
1618  row[i] = gsl_vector_get( gsl_row, i );
1619  }
1620 
1621  // Clean up gsl temporaries
1622  gsl_vector_free( gsl_row );
1623 
1624  return;
1625 }
1626 
1627 GslVector
1628 GslMatrix::getRow(unsigned int row_num ) const
1629 {
1630  // We rely on the void getRow's to do sanity checks
1631  // So we don't do them here.
1632 
1633  GslVector row(m_env, m_map);
1634 
1635  this->getRow( row_num, row );
1636 
1637  return row;
1638 }
1639 
1640 GslVector
1641 GslMatrix::getColumn(unsigned int column_num ) const
1642 {
1643  // We rely on the void getRow's to do sanity checks
1644  // So we don't do them here.
1645 
1646  GslVector column(m_env, m_map);
1647 
1648  this->getColumn( column_num, column );
1649 
1650  return column;
1651 }
1652 
1653 void
1654 GslMatrix::setRow(unsigned int row_num, const GslVector& row)
1655 {
1656  this->resetLU();
1657  // Sanity checks
1658  queso_require_less_msg(row_num, this->numRowsLocal(), "Specified row number not within range");
1659 
1660  queso_require_equal_to_msg(row.sizeLocal(), this->numCols(), "row vector not same size as this matrix");
1661 
1662  gsl_vector* gsl_row = row.data();
1663 
1664  int error_code = gsl_matrix_set_row( m_mat, row_num, gsl_row );
1665  queso_require_equal_to_msg(error_code, 0, "gsl_matrix_set_row failed");
1666 
1667  return;
1668 }
1669 
1670 void
1671 GslMatrix::setColumn(unsigned int column_num, const GslVector& column)
1672 {
1673  this->resetLU();
1674  // Sanity checks
1675  queso_require_less_msg(column_num, this->numCols(), "Specified column number not within range");
1676 
1677  queso_require_equal_to_msg(column.sizeLocal(), this->numRowsLocal(), "column vector not same size as this matrix");
1678 
1679  gsl_vector* gsl_column = column.data();
1680 
1681  int error_code = gsl_matrix_set_col( m_mat, column_num, gsl_column );
1682  queso_require_equal_to_msg(error_code, 0, "gsl_matrix_set_col failed");
1683 
1684  return;
1685 }
1686 
1687 void
1688 GslMatrix::mpiSum( const MpiComm& comm, GslMatrix& M_global ) const
1689 {
1690  // Sanity Checks
1691  UQ_FATAL_RC_MACRO(((this->numRowsLocal() != M_global.numRowsLocal()) ||
1692  (this->numCols() != M_global.numCols() )),
1693  env().fullRank(),
1694  "GslMatrix::mpiSum()",
1695  "local and global matrices incompatible");
1696 
1697  /* TODO: Probably a better way to handle this unpacking/packing of data */
1698  int size = M_global.numRowsLocal()*M_global.numCols();
1699  std::vector<double> local( size, 0.0 );
1700  std::vector<double> global( size, 0.0 );
1701 
1702  int k;
1703  for( unsigned int i = 0; i < this->numRowsLocal(); i++ )
1704  {
1705  for( unsigned int j = 0; j < this->numCols(); j++ )
1706  {
1707  k = i + j*M_global.numCols();
1708 
1709  local[k] = (*this)(i,j);
1710  }
1711  }
1712 
1713  comm.Allreduce((void*) &local[0], (void*) &global[0], size, RawValue_MPI_DOUBLE, RawValue_MPI_SUM,
1714  "GslMatrix::mpiSum()",
1715  "failed MPI.Allreduce()");
1716 
1717  for( unsigned int i = 0; i < this->numRowsLocal(); i++ )
1718  {
1719  for( unsigned int j = 0; j < this->numCols(); j++ )
1720  {
1721  k = i + j*M_global.numCols();
1722 
1723  M_global(i,j) = global[k];
1724  }
1725  }
1726 
1727  return;
1728 }
1729 
1730 void
1732  const GslVector& x1Vec,
1733  const GslMatrix& y1Mat,
1734  const GslVector& x2Vec)
1735 {
1736  queso_require_greater_msg(x1Vec.sizeLocal(), 1, "invalid 'x1' size");
1737 
1738  queso_require_equal_to_msg(x1Vec.sizeLocal(), y1Mat.numRowsLocal(), "invalid 'x1' and 'y1' sizes");
1739 
1740  queso_require_equal_to_msg(x2Vec.sizeLocal(), this->numRowsLocal(), "invalid 'x2' and 'this' sizes");
1741 
1742  queso_require_equal_to_msg(y1Mat.numCols(), this->numCols(), "invalid 'y1' and 'this' sizes");
1743 
1744  GslVector y1Vec(x1Vec);
1745  GslVector y2Vec(x2Vec);
1746  for (unsigned int colId = 0; colId < y1Mat.numCols(); ++colId) {
1747  y1Mat.getColumn(colId,y1Vec);
1748  y2Vec.matlabLinearInterpExtrap(x1Vec,y1Vec,x2Vec);
1749  this->setColumn(colId,y2Vec);
1750  }
1751 
1752  return;
1753 }
1754 
1755 gsl_matrix*
1757 {
1758  return m_mat;
1759 }
1760 
1761 void
1762 GslMatrix::print(std::ostream& os) const
1763 {
1764  unsigned int nRows = this->numRowsLocal();
1765  unsigned int nCols = this->numCols();
1766 
1767  if (m_printHorizontally) {
1768  for (unsigned int i = 0; i < nRows; ++i) {
1769  for (unsigned int j = 0; j < nCols; ++j) {
1770  os << (*this)(i,j)
1771  << " ";
1772  }
1773  if (i != (nRows-1)) os << "; ";
1774  }
1775  //os << std::endl;
1776  }
1777  else {
1778  for (unsigned int i = 0; i < nRows; ++i) {
1779  for (unsigned int j = 0; j < nCols; ++j) {
1780  os << (*this)(i,j)
1781  << " ";
1782  }
1783  os << std::endl;
1784  }
1785  }
1786 
1787  return;
1788 }
1789 
1790 void
1792  const std::string& varNamePrefix,
1793  const std::string& fileName,
1794  const std::string& fileType,
1795  const std::set<unsigned int>& allowedSubEnvIds) const
1796 {
1797  queso_require_greater_equal_msg(m_env.subRank(), 0, "unexpected subRank");
1798 
1799  queso_require_less_equal_msg(this->numOfProcsForStorage(), 1, "implemented just for sequential vectors for now");
1800 
1801  FilePtrSetStruct filePtrSet;
1802  if (m_env.openOutputFile(fileName,
1803  fileType, // "m or hdf"
1804  allowedSubEnvIds,
1805  false,
1806  filePtrSet)) {
1807  unsigned int nRows = this->numRowsLocal();
1808  unsigned int nCols = this->numCols();
1809  *filePtrSet.ofsVar << varNamePrefix << "_sub" << m_env.subIdString() << " = zeros(" << nRows
1810  << "," << nCols
1811  << ");"
1812  << std::endl;
1813  *filePtrSet.ofsVar << varNamePrefix << "_sub" << m_env.subIdString() << " = [";
1814 
1815  for (unsigned int i = 0; i < nRows; ++i) {
1816  for (unsigned int j = 0; j < nCols; ++j) {
1817  *filePtrSet.ofsVar << (*this)(i,j)
1818  << " ";
1819  }
1820  *filePtrSet.ofsVar << "\n";
1821  }
1822  *filePtrSet.ofsVar << "];\n";
1823 
1824  m_env.closeFile(filePtrSet,fileType);
1825  }
1826 
1827  return;
1828 }
1829 
1830 void
1832  const std::string& fileName,
1833  const std::string& fileType,
1834  const std::set<unsigned int>& allowedSubEnvIds)
1835 {
1836  queso_require_greater_equal_msg(m_env.subRank(), 0, "unexpected subRank");
1837 
1838  queso_require_less_equal_msg(this->numOfProcsForStorage(), 1, "implemented just for sequential vectors for now");
1839 
1840  FilePtrSetStruct filePtrSet;
1841  if (m_env.openInputFile(fileName,
1842  fileType, // "m or hdf"
1843  allowedSubEnvIds,
1844  filePtrSet)) {
1845 
1846  // palms
1847  unsigned int nRowsLocal = this->numRowsLocal();
1848 
1849  // In the logic below, the id of a line' begins with value 0 (zero)
1850  unsigned int idOfMyFirstLine = 1;
1851  unsigned int idOfMyLastLine = nRowsLocal;
1852  unsigned int nCols = this->numCols();
1853 
1854  // Read number of matrix rows in the file by taking care of the first line,
1855  // which resembles something like 'variable_name = zeros(n_rows,n_cols);'
1856  std::string tmpString;
1857 
1858  // Read 'variable name' string
1859  *filePtrSet.ifsVar >> tmpString;
1860  //std::cout << "Just read '" << tmpString << "'" << std::endl;
1861 
1862  // Read '=' sign
1863  *filePtrSet.ifsVar >> tmpString;
1864  //std::cout << "Just read '" << tmpString << "'" << std::endl;
1865  queso_require_equal_to_msg(tmpString, "=", "string should be the '=' sign");
1866 
1867  // Read 'zeros(n_rows,n_cols)' string
1868  *filePtrSet.ifsVar >> tmpString;
1869  //std::cout << "Just read '" << tmpString << "'" << std::endl;
1870  unsigned int posInTmpString = 6;
1871 
1872  // Isolate 'n_rows' in a string
1873  char nRowsString[tmpString.size()-posInTmpString+1];
1874  unsigned int posInRowsString = 0;
1875  do {
1876  queso_require_less_msg(posInTmpString, tmpString.size(), "symbol ',' not found in first line of file");
1877  nRowsString[posInRowsString++] = tmpString[posInTmpString++];
1878  } while (tmpString[posInTmpString] != ',');
1879  nRowsString[posInRowsString] = '\0';
1880 
1881  // Isolate 'n_cols' in a string
1882  posInTmpString++; // Avoid reading ',' char
1883  char nColsString[tmpString.size()-posInTmpString+1];
1884  unsigned int posInColsString = 0;
1885  do {
1886  queso_require_less_msg(posInTmpString, tmpString.size(), "symbol ')' not found in first line of file");
1887  nColsString[posInColsString++] = tmpString[posInTmpString++];
1888  } while (tmpString[posInTmpString] != ')');
1889  nColsString[posInColsString] = '\0';
1890 
1891  // Convert 'n_rows' and 'n_cols' strings to numbers
1892  unsigned int numRowsInFile = (unsigned int) strtod(nRowsString,NULL);
1893  unsigned int numColsInFile = (unsigned int) strtod(nColsString,NULL);
1894  if (m_env.subDisplayFile()) {
1895  *m_env.subDisplayFile() << "In GslMatrix::subReadContents()"
1896  << ": fullRank " << m_env.fullRank()
1897  << ", numRowsInFile = " << numRowsInFile
1898  << ", numColsInFile = " << numColsInFile
1899  << ", nRowsLocal = " << nRowsLocal
1900  << ", nCols = " << nCols
1901  << std::endl;
1902  }
1903 
1904  // Check if [num of rows in file] == [requested matrix row size]
1905  queso_require_equal_to_msg(numRowsInFile, nRowsLocal, "size of vec in file is not big enough");
1906 
1907  // Check if [num of cols in file] == [num cols in current matrix]
1908  queso_require_equal_to_msg(numColsInFile, nCols, "number of parameters of vec in file is different than number of parameters in this vec object");
1909 
1910  // Code common to any core in a communicator
1911  unsigned int maxCharsPerLine = 64*nCols; // Up to about 60 characters to represent each parameter value
1912 
1913  unsigned int lineId = 0;
1914  while (lineId < idOfMyFirstLine) {
1915  filePtrSet.ifsVar->ignore(maxCharsPerLine,'\n');
1916  lineId++;
1917  };
1918 
1919  if (m_env.subDisplayFile()) {
1920  *m_env.subDisplayFile() << "In GslMatrix::subReadContents()"
1921  << ": beginning to read input actual data"
1922  << std::endl;
1923  }
1924 
1925  // Take care of initial part of the first data line,
1926  // which resembles something like 'variable_name = [value1 value2 ...'
1927  // Read 'variable name' string
1928  *filePtrSet.ifsVar >> tmpString;
1929  //std::cout << "Core 0 just read '" << tmpString << "'" << std::endl;
1930 
1931  // Read '=' sign
1932  *filePtrSet.ifsVar >> tmpString;
1933  //std::cout << "Core 0 just read '" << tmpString << "'" << std::endl;
1934  queso_require_equal_to_msg(tmpString, "=", "in core 0, string should be the '=' sign");
1935 
1936  // Take into account the ' [' portion
1937  std::streampos tmpPos = filePtrSet.ifsVar->tellg();
1938  filePtrSet.ifsVar->seekg(tmpPos+(std::streampos)2);
1939 
1940  if (m_env.subDisplayFile()) {
1941  *m_env.subDisplayFile() << "In GslMatrix::subReadContents()"
1942  << ": beginning to read lines with numbers only"
1943  << ", lineId = " << lineId
1944  << ", idOfMyFirstLine = " << idOfMyFirstLine
1945  << ", idOfMyLastLine = " << idOfMyLastLine
1946  << std::endl;
1947  }
1948 
1949  double tmpRead;
1950  while (lineId <= idOfMyLastLine) {
1951  for (unsigned int j = 0; j < nCols; ++j) {
1952  *filePtrSet.ifsVar >> tmpRead;
1953  (*this)(lineId-idOfMyFirstLine,j) = tmpRead;
1954  }
1955  lineId++;
1956  };
1957 
1958  m_env.closeFile(filePtrSet,fileType);
1959  }
1960 
1961  return;
1962 }
1963 
1964 std::ostream&
1965 operator<<(std::ostream& os, const GslMatrix& obj)
1966 {
1967  obj.print(os);
1968 
1969  return os;
1970 }
1971 
1972 GslMatrix operator*(double a, const GslMatrix& mat)
1973 {
1974  GslMatrix answer(mat);
1975  answer *= a;
1976  return answer;
1977 }
1978 
1979 GslVector operator*(const GslMatrix& mat, const GslVector& vec)
1980 {
1981  return mat.multiply(vec);
1982 }
1983 
1985 {
1986  unsigned int m1Rows = m1.numRowsLocal();
1987  unsigned int m1Cols = m1.numCols();
1988  unsigned int m2Rows = m2.numRowsLocal();
1989  unsigned int m2Cols = m2.numCols();
1990 
1991  queso_require_equal_to_msg(m1Cols, m2Rows, "different sizes m1Cols and m2Rows");
1992 
1993  GslMatrix mat(m1.env(),m1.map(),m2Cols);
1994 
1995  //std::cout << "In GslMatrix(mat * mat): m1Cols = " << m1Cols << std::endl;
1996 
1997  unsigned int commonSize = m1Cols;
1998  for (unsigned int row1 = 0; row1 < m1Rows; ++row1) {
1999  for (unsigned int col2 = 0; col2 < m2Cols; ++col2) {
2000  double result = 0.;
2001  for (unsigned int k = 0; k < commonSize; ++k) {
2002  result += m1(row1,k)*m2(k,col2);
2003  }
2004  mat(row1,col2) = result;
2005  }
2006  //std::cout << "In GslMatrix(mat * mat): ended row " << row1 << std::endl;
2007  }
2008 
2009  return mat;
2010 }
2011 
2013 {
2014  GslMatrix answer(m1);
2015  answer += m2;
2016  return answer;
2017 }
2018 
2020 {
2021  GslMatrix answer(m1);
2022  answer -= m2;
2023  return answer;
2024 }
2025 
2027 {
2028  unsigned int nRows = v1.sizeLocal();
2029  unsigned int nCols = v2.sizeLocal();
2030  GslMatrix answer(v1.env(),v1.map(),nCols);
2031 
2032  for (unsigned int i = 0; i < nRows; ++i) {
2033  double value1 = v1[i];
2034  for (unsigned int j = 0; j < nCols; ++j) {
2035  answer(i,j) = value1*v2[j];
2036  }
2037  }
2038 
2039  return answer;
2040 }
2041 
2043 {
2044  unsigned int vSize = vec.sizeLocal();
2045  unsigned int mRows = mat.numRowsLocal();
2046  unsigned int mCols = mat.numCols();
2047 
2048  queso_require_equal_to_msg(vSize, mRows, "size of vector is different from the number of rows in matrix");
2049 
2050  queso_require_equal_to_msg(mCols, mRows, "routine currently works for square matrices only");
2051 
2052  GslMatrix answer(mat);
2053  for (unsigned int i = 0; i < mRows; ++i) {
2054  double vecValue = vec[i];
2055  for (unsigned int j = 0; j < mCols; ++j) {
2056  answer(i,j) *= vecValue;
2057  }
2058  }
2059 
2060  return answer;
2061 }
2062 
2064 {
2065  unsigned int vSize = vec.sizeLocal();
2066  unsigned int mRows = mat.numRowsLocal();
2067  unsigned int mCols = mat.numCols();
2068 
2069  queso_require_equal_to_msg(vSize, mCols, "size of vector is different from the number of cols in matrix");
2070 
2071  queso_require_equal_to_msg(mCols, mRows, "routine currently works for square matrices only");
2072 
2073  GslMatrix answer(mat);
2074  for (unsigned int j = 0; j < mCols; ++j) {
2075  double vecValue = vec[j];
2076  for (unsigned int i = 0; i < mRows; ++i) {
2077  answer(i,j) *= vecValue;
2078  }
2079  }
2080 
2081  return answer;
2082 }
2083 
2084 } // End namespace QUESO
const GslMatrix & svdMatV() const
This function calls private member GslMatrix::internalSvd() to set a N-by-N orthogonal square matrix ...
Definition: GslMatrix.C:547
unsigned int displayVerbosity() const
Definition: Environment.C:396
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
virtual void base_copy(const Matrix &src)
Copies base data from matrix src to this matrix.
Definition: Matrix.C:99
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
#define RawValue_MPI_SUM
Definition: MpiComm.h:52
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
int subRank() const
Access function for sub-rank.
Definition: Environment.C:241
GslMatrix & operator/=(double a)
Stores in this the coordinate-wise division of this by a.
Definition: GslMatrix.C:187
GslVector abs() const
This function returns absolute value of elements in this.
Definition: GslVector.C:1094
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
const int UQ_MATRIX_SVD_FAILED_RC
Definition: Defines.h:87
GslMatrix leftDiagScaling(const GslVector &vec, const GslMatrix &mat)
Definition: GslMatrix.C:2042
const BaseEnvironment & m_env
QUESO environment variable.
Definition: Matrix.h:116
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
std::ofstream * ofsVar
Provides a stream interface to write data to files.
Definition: Environment.h:81
int fullRank() const
Returns the process full rank.
Definition: Environment.C:222
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
#define queso_require_not_equal_to_msg(expr1, expr2, msg)
Definition: asserts.h:86
Struct for handling data input and output from files.
Definition: Environment.h:72
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
#define queso_require_less_equal_msg(expr1, expr2, msg)
Definition: asserts.h:89
unsigned int sizeLocal() const
Returns the length of this vector.
Definition: GslVector.C:253
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
const std::string & subIdString() const
Access to the attribute m_subIdString; which stores the string for the sub-environment, and it will be used, for instance, to create the output files for each sub-environment.
Definition: Environment.C:301
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 closeFile(FilePtrSetStruct &filePtrSet, const std::string &fileType) const
Closes the file.
Definition: Environment.C:1020
#define queso_require_equal_to_msg(expr1, expr2, msg)
Definition: asserts.h:85
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
double norm2() const
Returns the 2-norm (Euclidean norm) of the vector.
Definition: GslVector.C:293
const int UQ_MATRIX_IS_NOT_POS_DEFINITE_RC
Definition: Defines.h:84
#define queso_require_msg(asserted, msg)
Definition: asserts.h:69
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
const BaseEnvironment & env() const
Definition: Vector.C:54
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
std::ofstream * subDisplayFile() const
Access function for m_subDisplayFile (displays file on stream).
Definition: Environment.C:274
GslMatrix inverse() const
This function calculated the inverse of this matrix (square).
Definition: GslMatrix.C:736
#define queso_require_less_msg(expr1, expr2, msg)
Definition: asserts.h:87
GslMatrix rightDiagScaling(const GslMatrix &mat, const GslVector &vec)
Definition: GslMatrix.C:2063
void cwSet(double value)
Component-wise sets all values to this with value.
Definition: GslVector.C:339
#define UQ_RC_MACRO(macroIRc, givenRank, where, what, retValue)
Definition: Defines.h:123
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
#define UQ_FATAL_RC_MACRO(macroIRc, givenRank, where, what)
Definition: Defines.h:154
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
bool m_printHorizontally
Flag for either or not print this matrix.
Definition: Matrix.h:131
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
const Map & m_map
Mapping variable.
Definition: Matrix.h:126
void copy(const GslMatrix &src)
In this function this matrix receives a copy of matrix src.
Definition: GslMatrix.C:245
unsigned int numOfProcsForStorage() const
Definition: Matrix.C:62
bool m_inDebugMode
Flag for either or not QUESO is in debug mode.
Definition: Matrix.h:134
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
const Map & map() const
Definition: Vector.C:61
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
bool openOutputFile(const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds, bool writeOver, FilePtrSetStruct &filePtrSet) const
Opens an output file for each sub-environment that was chosen to send data to the file...
Definition: Environment.C:467
gsl_vector * data() const
Definition: GslVector.C:982
A class for partitioning vectors and matrices.
Definition: Map.h:49
#define queso_require_greater_equal_msg(expr1, expr2, msg)
Definition: asserts.h:90
std::ifstream * ifsVar
Provides a stream interface to read data from files.
Definition: Environment.h:84
~GslMatrix()
Destructor.
Definition: GslMatrix.C:161
void matlabLinearInterpExtrap(const GslVector &x1Vec, const GslVector &y1Vec, const GslVector &x2Vec)
Definition: GslVector.C:569
unsigned int checkingLevel() const
Access function to private attribute m_checkingLevel.
Definition: Environment.C:410
int dim
Definition: ann2fig.cpp:81
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
int worldRank() const
Returns the process world rank.
Definition: Environment.C:216
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
#define queso_require_greater_msg(expr1, expr2, msg)
Definition: asserts.h:88
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
#define RawValue_MPI_DOUBLE
Definition: MpiComm.h:48
void Allreduce(void *sendbuf, void *recvbuf, int count, RawType_MPI_Datatype datatype, RawType_MPI_Op op, const char *whereMsg, const char *whatMsg) const
Combines values from all processes and distributes the result back to all processes.
Definition: MpiComm.C:113
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
int k
Definition: ann_sample.cpp:53
Map * m_svdColMap
Mapping for matrices involved in the singular value decomposition (svd) routine.
Definition: GslMatrix.h:395
bool openInputFile(const std::string &fileName, const std::string &fileType, const std::set< unsigned int > &allowedSubEnvIds, FilePtrSetStruct &filePtrSet) const
Opens an input file.
Definition: Environment.C:834

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