Matrix.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: Matrix.hpp
4//
5// Created: Tue Jun 30 11:25:46 2009
6//
7// Author(s): Bård Skaflestad <bard.skaflestad@sintef.no>
8// Atgeirr F Rasmussen <atgeirr@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18 Copyright 2009, 2010 Statoil ASA.
19
20 This file is part of The Open Reservoir Simulator Project (OpenRS).
21
22 OpenRS is free software: you can redistribute it and/or modify
23 it under the terms of the GNU General Public License as published by
24 the Free Software Foundation, either version 3 of the License, or
25 (at your option) any later version.
26
27 OpenRS is distributed in the hope that it will be useful,
28 but WITHOUT ANY WARRANTY; without even the implied warranty of
29 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 GNU General Public License for more details.
31
32 You should have received a copy of the GNU General Public License
33 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34*/
35
36#ifndef OPENRS_MATRIX_HEADER
37#define OPENRS_MATRIX_HEADER
38
39#include <algorithm>
40#include <ostream>
41#include <vector>
42#include <boost/bind.hpp>
43
44#include <dune/common/fvector.hh>
45#include <opm/common/ErrorMacros.hpp>
46
49
50namespace Opm {
51
52 // ----------------------------------------------------------------------
53 // FullMatrix storage policies.
54 //
55
63 template<typename T>
64 class OwnData {
65 public:
73 T& operator[](int i) { return data_[i]; }
74 const T& operator[](int i) const { return data_[i]; }
75
76
80 int size() const { return data_.size(); }
81
85 T* data() { return &data_[0]; }
86 const T* data() const { return &data_[0]; }
87
88 protected:
100 OwnData(int sz, const T* data)
101 {
102 if (data) {
103 data_.assign(data, data + sz);
104 } else {
105 data_.resize(sz);
106 }
107 }
108
109 private:
110 std::vector<T> data_;
111 };
112
120 template<typename T>
122 public:
130 T& operator[](int i) { return data_[i]; }
131 const T& operator[](int i) const { return data_[i]; }
132
136 int size() const { return sz_; }
137
141 T* data() { return data_; }
142 const T* data() const { return data_; }
143
144 protected:
154 SharedData(int sz, T* data)
155 : sz_(sz), data_(data)
156 {
157 assert ((sz == 0) == (data == 0));
158 }
159
160 private:
161 int sz_;
162 T* data_;
163 };
164
172 template<typename T>
174 public:
182 const T& operator[](int i) const { return data_[i]; }
183
187 int size() const { return sz_; }
188
192 const T* data() const { return data_; }
193
194 protected:
203 ImmutableSharedData(int sz, const T* data)
204 : sz_(sz), data_(data)
205 {
206 assert (data_ != 0);
207 }
208
209 private:
210 int sz_;
211 const T* data_;
212 };
213
214
215
216
217
218 // ----------------------------------------------------------------------
219 // FullMatrix ordering policies.
220 //
221
224 class OrderingBase {
225 public:
231 int numRows() const { return rows_; }
232
238 int numCols() const { return cols_; }
239
240 protected:
243 OrderingBase()
244 : rows_(0), cols_(0)
245 {}
246
255 OrderingBase(int rows, int cols)
256 : rows_(rows), cols_(cols)
257 {}
258
259 private:
260 int rows_, cols_;
261 };
262
263
267 class COrdering : public OrderingBase {
268 public:
276 int leadingDimension() const { return numCols(); }
277
291 int idx(int row, int col) const
292 {
293 assert ((0 <= row) && (row < numRows()));
294 assert ((0 <= col) && (col < numCols()));
295
296 return row*numCols() + col;
297 }
298
299 protected:
301 COrdering()
302 : OrderingBase()
303 {}
304
313 COrdering(int rows, int cols)
314 : OrderingBase(rows, cols)
315 {}
316 };
317
318
322 class FortranOrdering : public OrderingBase {
323 public:
331 int leadingDimension() const { return numRows(); }
332
346 int idx(int row, int col) const
347 {
348 assert ((0 <= row) && (row < numRows()));
349 assert ((0 <= col) && (col < numCols()));
350
351 return row + col*numRows();
352 }
353
354 protected:
356 FortranOrdering()
357 : OrderingBase()
358 {}
359
368 FortranOrdering(int rows, int cols)
369 : OrderingBase(rows, cols)
370 {}
371 };
372
373
374
375
376 // ----------------------------------------------------------------------
377 // Class FullMatrix.
378 //
379
380
395 template<typename T,
396 template<typename> class StoragePolicy,
397 class OrderingPolicy>
398 class FullMatrix : private StoragePolicy<T>,
399 private OrderingPolicy
400 {
401 public:
404 FullMatrix()
405 : StoragePolicy<T>(0, 0),
406 OrderingPolicy()
407 {}
408
426 template <typename DataPointer>
427 FullMatrix(int rows, int cols, DataPointer data)
428 : StoragePolicy<T>(rows * cols, data),
429 OrderingPolicy(rows, cols)
430 {}
431
440 template <template<typename> class OtherSP>
441 explicit FullMatrix(const FullMatrix<T, OtherSP, OrderingPolicy>& m)
442 : StoragePolicy<T>(m.numRows()*m.numCols(), m.data()),
443 OrderingPolicy(m.numRows(), m.numCols())
444 {
445 }
446
458 template <template<typename> class OtherSP, class OtherOP>
459 FullMatrix& operator=(const FullMatrix<T, OtherSP, OtherOP>& m)
460 {
461 assert(numRows() == m.numRows());
462 assert(numCols() == m.numCols());
463 for (int r = 0; r < numRows(); ++r) {
464 for (int c = 0; c < numCols(); ++c) {
465 this->operator()(r, c) = m(r,c);
466 }
467 }
468 return *this;
469 }
470
482 template <template<typename> class OtherSP, class OtherOP>
483 bool operator==(const FullMatrix<T, OtherSP, OtherOP>& m)
484 {
485 if (numRows() != m.numRows() || numCols() != m.numCols()) {
486 return false;
487 }
488 for (int r = 0; r < numRows(); ++r) {
489 for (int c = 0; c < numCols(); ++c) {
490 if (this->operator()(r,c) != m(r,c)) {
491 return false;
492 }
493 }
494 }
495 return true;
496 }
497
506 template <template<typename> class OtherSP>
507 void operator+= (const FullMatrix<T, OtherSP, OrderingPolicy>& m)
508 {
509 assert(numRows() == m.numRows() && numCols() == m.numCols());
510 std::transform(data(), data() + this->size(),
511 m.data(), data(), std::plus<T>());
512 }
513
519 void operator*= (const T& scalar)
520 {
521 std::transform(data(), data() + this->size(),
522 data(), boost::bind(std::multiplies<T>(), _1, scalar));
523 }
524
527 typedef T value_type;
528
529 using StoragePolicy<T>::data;
530 using OrderingPolicy::numRows;
531 using OrderingPolicy::numCols;
532 using OrderingPolicy::leadingDimension;
533
545 value_type& operator()(int row, int col)
546 {
547 return this->operator[](this->idx(row, col));
548 }
549
561 const value_type& operator()(int row, int col) const
562 {
563 return this->operator[](this->idx(row, col));
564 }
565 };
566
567
568
569 // ----------------------------------------------------------------------
570 // FullMatrix operations.
571 //
572
573
574 // Convenience typedefs.
575
580 typedef FullMatrix<double, OwnData, COrdering> OwnCMatrix;
581 typedef FullMatrix<double, SharedData, COrdering> SharedCMatrix;
582 typedef const FullMatrix<double, ImmutableSharedData, COrdering> ImmutableCMatrix;
583
584
589 typedef FullMatrix<double, OwnData, FortranOrdering> OwnFortranMatrix;
590 typedef FullMatrix<double, SharedData, FortranOrdering> SharedFortranMatrix;
591 typedef const FullMatrix<double, ImmutableSharedData, FortranOrdering> ImmutableFortranMatrix;
592
593
602 template<class Matrix>
603 void zero(Matrix& A)
604 {
605 std::fill_n(A.data(), A.numRows() * A.numCols(),
606 typename Matrix::value_type(0.0));
607 }
608
617 template<class Matrix>
618 void eye(Matrix& A)
619 {
620 zero(A);
621 for (int i = 0; i < std::min(A.numRows(), A.numCols()); ++i) {
622 A(i, i) = 1.0;
623 }
624 }
625
637 template<class Matrix>
638 typename Matrix::value_type
639 trace(const Matrix& A)
640 {
641 typename Matrix::value_type ret(0);
642
643 for (int i = 0; i < std::min(A.numRows(), A.numCols()); ++i) {
644 ret += A(i,i);
645 }
646 return ret;
647 }
648
649
667 template<class Matrix, int rows>
668 Dune::FieldVector<typename Matrix::value_type, rows>
669 prod(const Matrix& A, const Dune::FieldVector<typename Matrix::value_type,rows>& x)
670 {
671 const int cols = rows;
672 assert (A.numRows() == rows);
673 assert (A.numCols() == cols);
674
675 Dune::FieldVector<typename Matrix::value_type, rows> res(0.0);
676 for (int c = 0; c < cols; ++c) {
677 for (int r = 0; r < rows; ++r) {
678 res[r] += A(r, c)*x[c];
679 }
680 }
681 return res;
682 }
683
691 template<class Matrix1, class Matrix2, class MutableMatrix>
692 void prod(const Matrix1& A, const Matrix2& B, MutableMatrix& C)
693 {
694 int result_rows = A.numRows();
695 int result_cols = B.numCols();
696 int inner_dim = A.numCols();
697 assert (inner_dim == B.numRows());
698 assert(C.numRows() == result_rows);
699 assert(C.numCols() == result_cols);
700
701 for (int c = 0; c < result_cols; ++c) {
702 for (int r = 0; r < result_rows; ++r) {
703 C(r,c) = 0.0;
704 for (int i = 0; i < inner_dim; ++i) {
705 C(r,c) += A(r,i)*B(i,c);
706 }
707 }
708 }
709 }
710
711
729 template<typename T, template<typename> class StoragePolicy>
730 int orthogonalizeColumns(FullMatrix<T,StoragePolicy,FortranOrdering>& A)
731 {
732 typedef typename FullMatrix<T,StoragePolicy,FortranOrdering>::value_type value_type;
733
734 static std::vector<value_type> tau;
735 static std::vector<value_type> work;
736
737 if (int(tau .size()) < A.numCols()) tau .resize( A.numCols());
738 if (int(work.size()) < 64 * A.numRows()) work.resize(64 * A.numRows()); // 64 from ILAENV
739
740 int info = 0;
741
742 // Generate factorization. Matrix Q is implicitly represented.
743 Opm::BLAS_LAPACK::GEQRF(A.numRows(), A.numCols() ,
744 A.data() , A.leadingDimension() ,
745 &tau[0] , &work[0], int(work.size()),
746 info);
747
748 if (info == 0)
749 {
750 // QR factorization successfully generated--extract
751 // explict representation of orthogonal matrix Q.
752 Opm::BLAS_LAPACK::ORGQR(A.numRows(), A.numCols(), A.numCols() ,
753 A.data() , A.leadingDimension() ,
754 &tau[0] , &work[0], int(work.size()),
755 info);
756 }
757
758 return info;
759 }
760
761
780 template<typename T, template<typename> class StoragePolicy, class OrderingPolicy>
781 int invert(FullMatrix<T,StoragePolicy,OrderingPolicy>& A)
782 {
783 typedef typename FullMatrix<T,StoragePolicy,OrderingPolicy>::value_type value_type;
784
785 assert (A.numRows() == A.numCols());
786
787 std::vector<int> ipiv(A.numRows());
788 int info = 0;
789
790 // Correct both for COrdering and FortranOrdering (inv(A)' == inv(A')).
791 Opm::BLAS_LAPACK::GETRF(A.numRows(), A.numCols(), A.data(),
792 A.leadingDimension(), &ipiv[0], info);
793
794 if (info == 0) {
795 std::vector<value_type> work(A.numRows());
796
797 Opm::BLAS_LAPACK::GETRI(A.numRows(), A.data(), A.leadingDimension(),
798 &ipiv[0], &work[0], int(work.size()), info);
799 }
800
801 return info;
802 }
803
804
829 template<typename T, template<typename> class StoragePolicy>
830 void symmetricUpdate(const T& a1,
831 const FullMatrix<T,StoragePolicy,FortranOrdering>& A ,
832 const T& a2,
833 FullMatrix<T,StoragePolicy,FortranOrdering>& C )
834 {
835 Opm::BLAS_LAPACK::SYRK("Upper" , "No transpose" ,
836 C.numRows() , A.numCols() ,
837 a1, A.data(), A.leadingDimension(),
838 a2, C.data(), C.leadingDimension());
839
840 // Account for SYRK (in this case) only updating the upper
841 // leading n-by-n submatrix.
842 //
843 for (int j = 0; j < C.numCols(); ++j) {
844 for (int i = j+1; i < C.numRows(); ++i) {
845 C(i,j) = C(j,i);
846 }
847 }
848 }
849
850
851 // B <- A*B*A'
852 // Assumes T is an arithmetic (floating point) type, and that A==A'.
857 template<typename T, template<typename> class StoragePolicy>
858 void symmetricUpdate(const FullMatrix<T,StoragePolicy,FortranOrdering>& A,
859 FullMatrix<T,StoragePolicy,FortranOrdering>& B)
860 {
861 typedef typename FullMatrix<T,StoragePolicy,FortranOrdering>::value_type value_type;
862
863 // B <- A*B
864 Opm::BLAS_LAPACK::TRMM("Left" , "Upper", "No transpose", "Non-unit",
865 B.numRows(), B.numCols(), value_type(1.0),
866 A.data(), A.leadingDimension(),
867 B.data(), B.leadingDimension());
868
869 // B <- B*A (== A * B_orig * A)
870 Opm::BLAS_LAPACK::TRMM("Right", "Upper", "No transpose", "Non-unit",
871 B.numRows(), B.numCols(), value_type(1.0),
872 A.data(), A.leadingDimension(),
873 B.data(), B.leadingDimension());
874
875 // Account for TRMM (in this case) only updating the upper
876 // leading n-by-n submatrix.
877 //
878 for (int j = 0; j < B.numCols(); ++j) {
879 for (int i = j+1; i < B.numRows(); ++i) {
880 B(i,j) = B(j,i);
881 }
882 }
883 }
884
885
912 template<typename T ,
913 template<typename> class SP>
914 void vecMulAdd_N(const T& a1,
915 const FullMatrix<T,SP,FortranOrdering>& A ,
916 const std::vector<T>& x ,
917 const T& a2,
918 std::vector<T>& y)
919 {
920 assert(A.numRows() == y.size());
921 assert(A.numCols() == x.size());
922
923 Opm::BLAS_LAPACK::GEMV("No Transpose",
924 A.numRows(), A.numCols(),
925 a1, A.data(), A.leadingDimension(),
926 &x[0], 1, a2, &y[0], 1);
927 }
928
929
956 template<typename T ,
957 template<typename> class SP>
958 void vecMulAdd_N(const T& a1,
959 const FullMatrix<T,SP,FortranOrdering>& A ,
960 const T* x ,
961 const T& a2,
962 T* y)
963 {
964 Opm::BLAS_LAPACK::GEMV("No Transpose",
965 A.numRows(), A.numCols(),
966 a1, A.data(), A.leadingDimension(),
967 x, 1, a2, y, 1);
968 }
969
970
998 template<typename T ,
999 template<typename> class SP>
1000 void vecMulAdd_T(const T& a1,
1001 const FullMatrix<T,SP,FortranOrdering>& A ,
1002 const std::vector<T>& x ,
1003 const T& a2,
1004 std::vector<T>& y)
1005 {
1006 assert (A.numCols() == y.size());
1007 assert (A.numRows() == x.size());
1008
1009 Opm::BLAS_LAPACK::GEMV("Transpose",
1010 A.numRows(), A.numCols(),
1011 a1, A.data(), A.leadingDimension(),
1012 &x[0], 1, a2, &y[0], 1);
1013 }
1014
1015
1043 template<typename T ,
1044 template<typename> class SP>
1045 void vecMulAdd_T(const T& a1,
1046 const FullMatrix<T,SP,FortranOrdering>& A ,
1047 const T* x ,
1048 const T& a2,
1049 T* y)
1050 {
1051 Opm::BLAS_LAPACK::GEMV("Transpose",
1052 A.numRows(), A.numCols(),
1053 a1, A.data(), A.leadingDimension(),
1054 x, 1, a2, y, 1);
1055 }
1056
1057
1085 template<typename T ,
1086 template<typename> class SP>
1087 void vecMulAdd_N(const T& a1,
1088 const FullMatrix<T,SP,COrdering>& A ,
1089 const T* x ,
1090 const T& a2,
1091 T* y)
1092 {
1093 typedef FullMatrix<T, ImmutableSharedData, FortranOrdering> FMAT;
1094
1095 const FMAT At(A.numCols(), A.numRows(), A.data());
1096
1097 vecMulAdd_T(a1, At, x, a2, y);
1098 }
1099
1100
1133 template<typename T ,
1134 template<typename> class SP1,
1135 template<typename> class SP2,
1136 template<typename> class SP3>
1137 void matMulAdd_NN(const T& a1,
1138 const FullMatrix<T,SP1,FortranOrdering>& A ,
1139 const FullMatrix<T,SP2,FortranOrdering>& B ,
1140 const T& a2,
1141 FullMatrix<T,SP3,FortranOrdering>& C)
1142 {
1143 assert(A.numRows() == C.numRows());
1144 assert(A.numCols() == B.numRows());
1145 assert(B.numCols() == C.numCols());
1146
1147 int m = A.numRows(); // Number of *rows* in A
1148 int n = B.numCols(); // Number of *cols* in B
1149 int k = A.numCols(); // Number of *cols* in A (== numer of *rows* in B)
1150
1151 Opm::BLAS_LAPACK::GEMM("No Transpose", "No Transpose", m, n, k,
1152 a1, A.data(), A.leadingDimension(),
1153 B.data(), B.leadingDimension(),
1154 a2, C.data(), C.leadingDimension());
1155 }
1156
1157
1191 template<typename T ,
1192 template<typename> class SP1,
1193 template<typename> class SP2,
1194 template<typename> class SP3>
1195 void matMulAdd_NT(const T& a1,
1196 const FullMatrix<T,SP1,FortranOrdering>& A ,
1197 const FullMatrix<T,SP2,FortranOrdering>& B ,
1198 const T& a2,
1199 FullMatrix<T,SP3,FortranOrdering>& C)
1200 {
1201 assert(A.numRows() == C.numRows());
1202 assert(B.numRows() == C.numCols());
1203 assert(A.numCols() == B.numCols());
1204
1205 int m = A.numRows(); // Number of *rows* in A
1206 int n = B.numRows(); // Number of *cols* in B'
1207 int k = A.numCols(); // Number of *cols* in A (== numer of *rows* in B')
1208
1209 Opm::BLAS_LAPACK::GEMM("No Transpose", "Transpose", m, n, k,
1210 a1, A.data(), A.leadingDimension(),
1211 B.data(), B.leadingDimension(),
1212 a2, C.data(), C.leadingDimension());
1213 }
1214
1215
1249 template<typename T ,
1250 template<typename> class SP1,
1251 template<typename> class SP2,
1252 template<typename> class SP3>
1253 void matMulAdd_TN(const T& a1,
1254 const FullMatrix<T,SP1,FortranOrdering>& A ,
1255 const FullMatrix<T,SP2,FortranOrdering>& B ,
1256 const T& a2,
1257 FullMatrix<T,SP3,FortranOrdering>& C)
1258 {
1259 assert (A.numCols() == C.numRows());
1260 assert (A.numRows() == B.numRows());
1261 assert (B.numCols() == C.numCols());
1262
1263 int m = A.numCols(); // Number of *rows* in A'
1264 int n = B.numCols(); // Number of *cols* in B
1265 int k = A.numRows(); // Number of *cols* in A' (== numer of *rows* in B)
1266
1267 Opm::BLAS_LAPACK::GEMM("Transpose", "No Transpose", m, n, k,
1268 a1, A.data(), A.leadingDimension(),
1269 B.data(), B.leadingDimension(),
1270 a2, C.data(), C.leadingDimension());
1271 }
1272
1273
1307 template<typename T ,
1308 template<typename> class SP1,
1309 template<typename> class SP2,
1310 template<typename> class SP3>
1311 void matMulAdd_NN(const T& a1,
1312 const FullMatrix<T,SP1,FortranOrdering>& A ,
1313 const FullMatrix<T,SP2,COrdering>& B ,
1314 const T& a2,
1315 FullMatrix<T,SP3,FortranOrdering>& C)
1316 {
1317 typedef typename FullMatrix<T,SP2,COrdering>::value_type value_type;
1318 typedef FullMatrix<value_type,ImmutableSharedData,FortranOrdering> FMat;
1319
1320 const FMat Bt(B.numCols(), B.numRows(), B.data());
1321
1322 matMulAdd_NT(a1, A, Bt, a2, C);
1323 }
1324
1325
1352 template<class charT, class traits,
1353 typename T, template<typename> class SP, class OP>
1354 std::basic_ostream<charT,traits>&
1355 operator<<(std::basic_ostream<charT,traits>& os,
1356 const FullMatrix<T,SP,OP>& A)
1357 {
1358 for (int i = 0; i < A.numRows(); ++i) {
1359 for (int j = 0; j < A.numCols(); ++j)
1360 os << A(i,j) << ' ';
1361 os << '\n';
1362 }
1363
1364 return os;
1365 }
1366} // namespace Opm
1367#endif // OPENRS_MATRIX_HEADER
FullMatrix StoragePolicy which provides immutable object sharing semantics.
Definition: Matrix.hpp:173
const T * data() const
Direct access to all data.
Definition: Matrix.hpp:192
ImmutableSharedData(int sz, const T *data)
Constructor.
Definition: Matrix.hpp:203
const T & operator[](int i) const
Storage element access.
Definition: Matrix.hpp:182
int size() const
Data size query.
Definition: Matrix.hpp:187
FullMatrix StoragePolicy which provides object owning semantics.
Definition: Matrix.hpp:64
OwnData(int sz, const T *data)
Constructor.
Definition: Matrix.hpp:100
T & operator[](int i)
Storage element access.
Definition: Matrix.hpp:73
T * data()
Direct access to all data.
Definition: Matrix.hpp:85
int size() const
Data size query.
Definition: Matrix.hpp:80
const T * data() const
Definition: Matrix.hpp:86
const T & operator[](int i) const
Definition: Matrix.hpp:74
FullMatrix StoragePolicy which provides object sharing semantics.
Definition: Matrix.hpp:121
const T & operator[](int i) const
Definition: Matrix.hpp:131
T & operator[](int i)
Storage element access.
Definition: Matrix.hpp:130
int size() const
Data size query.
Definition: Matrix.hpp:136
SharedData(int sz, T *data)
Constructor.
Definition: Matrix.hpp:154
T * data()
Direct access to all data.
Definition: Matrix.hpp:141
const T * data() const
Definition: Matrix.hpp:142
void GEMM(const char *transA, const char *transB, const int m, const int n, const int k, const T &a1, const T *A, const int ldA, const T *B, const int ldB, const T &a2, T *C, const int ldC)
GEneral Matrix Matrix product (Level 3 BLAS).
void GEMV(const char *transA, const int m, const int n, const T &a1, const T *A, const int ldA, const T *x, const int incX, const T &a2, T *y, const int incY)
GEneral Matrix Vector product (Level 2 BLAS).
void GEQRF(const int m, const int n, T *A, const int ld, T *tau, T *work, const int lwork, int &info)
GEneral matrix QR Factorization (LAPACK)
void SYRK(const char *uplo, const char *trans, const int n, const int k, const T &a1, const T *A, const int ldA, const T &a2, T *C, const int ldC)
SYmmetric Rank K update of symmetric matrix (Level 3 BLAS)
void GETRI(const int n, T *A, const int ld, const int *ipiv, T *work, int lwork, int &info)
GEneral matrix TRiangular Inversion (LAPACK).
void TRMM(const char *side, const char *uplo, const char *transA, const char *diag, const int m, const int n, const T &a, const T *A, const int ldA, T *B, const int ldB)
TRiangular Matrix Matrix product (Level 2 BLAS)
void GETRF(const int m, const int n, T *A, const int ld, int *ipiv, int &info)
GEneral matrix TRiangular Factorization (LAPACK).
void ORGQR(const int m, const int n, const int k, T *A, const int ld, const T *tau, T *work, const int lwork, int &info)
ORthogonal matrix Generator from QR factorization (LAPACK).
Definition: BlackoilFluid.hpp:32
void vecMulAdd_T(const T &a1, const FullMatrix< T, SP, FortranOrdering > &A, const std::vector< T > &x, const T &a2, std::vector< T > &y)
GEneral Matrix-Vector product (GAXPY operation). Specifically, .
Definition: Matrix.hpp:1000
void matMulAdd_NN(const T &a1, const FullMatrix< T, SP1, FortranOrdering > &A, const FullMatrix< T, SP2, FortranOrdering > &B, const T &a2, FullMatrix< T, SP3, FortranOrdering > &C)
GEneral Matrix-Matrix product update of other matrix. Specificlly, .
Definition: Matrix.hpp:1137
int orthogonalizeColumns(FullMatrix< T, StoragePolicy, FortranOrdering > &A)
Construct orthonormal basis for matrix range (i.e., column space). Based on a QR factorization of the...
Definition: Matrix.hpp:730
void symmetricUpdate(const T &a1, const FullMatrix< T, StoragePolicy, FortranOrdering > &A, const T &a2, FullMatrix< T, StoragePolicy, FortranOrdering > &C)
Symmetric, rank update of symmetric matrix. Specifically, .
Definition: Matrix.hpp:830
FullMatrix< double, SharedData, COrdering > SharedCMatrix
Definition: Matrix.hpp:581
const FullMatrix< double, ImmutableSharedData, FortranOrdering > ImmutableFortranMatrix
Definition: Matrix.hpp:591
Matrix::value_type trace(const Matrix &A)
Compute matrix trace (i.e., sum(diag(A))).
Definition: Matrix.hpp:639
void matMulAdd_NT(const T &a1, const FullMatrix< T, SP1, FortranOrdering > &A, const FullMatrix< T, SP2, FortranOrdering > &B, const T &a2, FullMatrix< T, SP3, FortranOrdering > &C)
GEneral Matrix-Matrix product update of other matrix. Specificlly, .
Definition: Matrix.hpp:1195
const FullMatrix< double, ImmutableSharedData, COrdering > ImmutableCMatrix
Definition: Matrix.hpp:582
void zero(Matrix &A)
Zero-fill a.
Definition: Matrix.hpp:603
FullMatrix< double, OwnData, FortranOrdering > OwnFortranMatrix
Convenience typedefs for Fortran-ordered.
Definition: Matrix.hpp:589
int invert(FullMatrix< T, StoragePolicy, OrderingPolicy > &A)
Matrix inversion, .
Definition: Matrix.hpp:781
std::basic_ostream< charT, traits > & operator<<(std::basic_ostream< charT, traits > &os, const BCBase< T > &bc)
Stream insertion for BCBase.
Definition: BoundaryConditions.hpp:110
void vecMulAdd_N(const T &a1, const FullMatrix< T, SP, FortranOrdering > &A, const std::vector< T > &x, const T &a2, std::vector< T > &y)
GEneral Matrix-Vector product (GAXPY operation). Specifically, .
Definition: Matrix.hpp:914
FullMatrix< double, SharedData, FortranOrdering > SharedFortranMatrix
Definition: Matrix.hpp:590
void matMulAdd_TN(const T &a1, const FullMatrix< T, SP1, FortranOrdering > &A, const FullMatrix< T, SP2, FortranOrdering > &B, const T &a2, FullMatrix< T, SP3, FortranOrdering > &C)
GEneral Matrix-Matrix product update of other matrix. Specificlly, .
Definition: Matrix.hpp:1253
void eye(Matrix &A)
Make an identity.
Definition: Matrix.hpp:618
FullMatrix< double, OwnData, COrdering > OwnCMatrix
Convenience typedefs for C-ordered.
Definition: Matrix.hpp:580
Dune::FieldVector< typename Matrix::value_type, rows > prod(const Matrix &A, const Dune::FieldVector< typename Matrix::value_type, rows > &x)
Matrix applied to a vector.
Definition: Matrix.hpp:669