Dune Core Modules (2.6.0)

fmatrix.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_FMATRIX_HH
4 #define DUNE_FMATRIX_HH
5 
6 #include <cmath>
7 #include <cstddef>
8 #include <iostream>
9 #include <algorithm>
10 #include <initializer_list>
11 
14 #include <dune/common/fvector.hh>
16 #include <dune/common/precision.hh>
18 
19 namespace Dune
20 {
21 
33  template< class K, int ROWS, int COLS > class FieldMatrix;
34 
35  template< class K, int ROWS, int COLS >
36  struct DenseMatVecTraits< FieldMatrix<K,ROWS,COLS> >
37  {
38  typedef FieldMatrix<K,ROWS,COLS> derived_type;
39 
40  // each row is implemented by a field vector
41  typedef FieldVector<K,COLS> row_type;
42 
43  typedef row_type &row_reference;
44  typedef const row_type &const_row_reference;
45 
46  typedef std::array<row_type,ROWS> container_type;
47  typedef K value_type;
48  typedef typename container_type::size_type size_type;
49  };
50 
51  template< class K, int ROWS, int COLS >
52  struct FieldTraits< FieldMatrix<K,ROWS,COLS> >
53  {
54  typedef typename FieldTraits<K>::field_type field_type;
55  typedef typename FieldTraits<K>::real_type real_type;
56  };
57 
66  template<class K, int ROWS, int COLS>
67  class FieldMatrix : public DenseMatrix< FieldMatrix<K,ROWS,COLS> >
68  {
69  std::array< FieldVector<K,COLS>, ROWS > _data;
71  public:
72 
74  enum {
76  rows = ROWS,
78  cols = COLS
79  };
80 
81  typedef typename Base::size_type size_type;
82  typedef typename Base::row_type row_type;
83 
84  typedef typename Base::row_reference row_reference;
85  typedef typename Base::const_row_reference const_row_reference;
86 
87  //===== constructors
91 
94  FieldMatrix(std::initializer_list<Dune::FieldVector<K, cols> > const &l) {
95  assert(l.size() == rows); // Actually, this is not needed any more!
96  std::copy_n(l.begin(), std::min(static_cast<std::size_t>(ROWS),
97  l.size()),
98  _data.begin());
99  }
100 
101  template <class T,
102  typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>>
103  FieldMatrix(T const& rhs)
104  {
105  *this = rhs;
106  }
107 
108  using Base::operator=;
109 
110  // Specialisation: FieldMatrix assignment (compile-time bounds checking)
111  template <typename T, int rows, int cols>
112  FieldMatrix& operator=(FieldMatrix<T,rows,cols> const &rhs)
113  {
114  static_assert(rows == ROWS, "Size mismatch in matrix assignment (rows)");
115  static_assert(cols == COLS, "Size mismatch in matrix assignment (columns)");
116  _data = rhs._data;
117  return *this;
118  }
119 
121  template<int l>
123  {
125 
126  for (size_type i=0; i<l; i++) {
127  for (size_type j=0; j<cols; j++) {
128  C[i][j] = 0;
129  for (size_type k=0; k<rows; k++)
130  C[i][j] += M[i][k]*(*this)[k][j];
131  }
132  }
133  return C;
134  }
135 
136  using Base::rightmultiply;
137 
139  template <int r, int c>
141  {
142  static_assert(r == c, "Cannot rightmultiply with non-square matrix");
143  static_assert(r == cols, "Size mismatch");
144  FieldMatrix<K,rows,cols> C(*this);
145 
146  for (size_type i=0; i<rows; i++)
147  for (size_type j=0; j<cols; j++) {
148  (*this)[i][j] = 0;
149  for (size_type k=0; k<cols; k++)
150  (*this)[i][j] += C[i][k]*M[k][j];
151  }
152  return *this;
153  }
154 
156  template<int l>
158  {
160 
161  for (size_type i=0; i<rows; i++) {
162  for (size_type j=0; j<l; j++) {
163  C[i][j] = 0;
164  for (size_type k=0; k<cols; k++)
165  C[i][j] += (*this)[i][k]*M[k][j];
166  }
167  }
168  return C;
169  }
170 
171  // make this thing a matrix
172  constexpr size_type mat_rows() const { return ROWS; }
173  constexpr size_type mat_cols() const { return COLS; }
174 
175  row_reference mat_access ( size_type i )
176  {
177  DUNE_ASSERT_BOUNDS(i < ROWS);
178  return _data[i];
179  }
180 
181  const_row_reference mat_access ( size_type i ) const
182  {
183  DUNE_ASSERT_BOUNDS(i < ROWS);
184  return _data[i];
185  }
186  };
187 
188 #ifndef DOXYGEN // hide specialization
191  template<class K>
192  class FieldMatrix<K,1,1> : public DenseMatrix< FieldMatrix<K,1,1> >
193  {
194  FieldVector<K,1> _data;
195  typedef DenseMatrix< FieldMatrix<K,1,1> > Base;
196  public:
197  // standard constructor and everything is sufficient ...
198 
199  //===== type definitions and constants
200 
202  typedef typename Base::size_type size_type;
203 
205  enum {
208  blocklevel = 1
209  };
210 
211  typedef typename Base::row_type row_type;
212 
213  typedef typename Base::row_reference row_reference;
214  typedef typename Base::const_row_reference const_row_reference;
215 
217  enum {
220  rows = 1,
223  cols = 1
224  };
225 
226  //===== constructors
229  FieldMatrix () {}
230 
233  FieldMatrix(std::initializer_list<Dune::FieldVector<K, 1>> const &l)
234  {
235  std::copy_n(l.begin(), std::min(static_cast< std::size_t >( 1 ), l.size()), &_data);
236  }
237 
238  template <class T,
239  typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>>
240  FieldMatrix(T const& rhs)
241  {
242  *this = rhs;
243  }
244 
245  using Base::operator=;
246 
247  //===== solve
248 
250  template<int l>
251  FieldMatrix<K,l,1> leftmultiplyany (const FieldMatrix<K,l,1>& M) const
252  {
253  FieldMatrix<K,l,1> C;
254  for (size_type j=0; j<l; j++)
255  C[j][0] = M[j][0]*(*this)[0][0];
256  return C;
257  }
258 
261  {
262  _data[0] *= M[0][0];
263  return *this;
264  }
265 
267  template<int l>
268  FieldMatrix<K,1,l> rightmultiplyany (const FieldMatrix<K,1,l>& M) const
269  {
270  FieldMatrix<K,1,l> C;
271 
272  for (size_type j=0; j<l; j++)
273  C[0][j] = M[0][j]*_data[0];
274  return C;
275  }
276 
277  // make this thing a matrix
278  constexpr size_type mat_rows() const { return 1; }
279  constexpr size_type mat_cols() const { return 1; }
280 
281  row_reference mat_access ( size_type i )
282  {
284  DUNE_ASSERT_BOUNDS(i == 0);
285  return _data;
286  }
287 
288  const_row_reference mat_access ( size_type i ) const
289  {
291  DUNE_ASSERT_BOUNDS(i == 0);
292  return _data;
293  }
294 
296  FieldMatrix& operator+= (const K& k)
297  {
298  _data[0] += k;
299  return (*this);
300  }
301 
303  FieldMatrix& operator-= (const K& k)
304  {
305  _data[0] -= k;
306  return (*this);
307  }
308 
310  FieldMatrix& operator*= (const K& k)
311  {
312  _data[0] *= k;
313  return (*this);
314  }
315 
317  FieldMatrix& operator/= (const K& k)
318  {
319  _data[0] /= k;
320  return (*this);
321  }
322 
323  //===== conversion operator
324 
325  operator const K& () const { return _data[0]; }
326 
327  };
328 
330  template<typename K>
331  std::ostream& operator<< (std::ostream& s, const FieldMatrix<K,1,1>& a)
332  {
333  s << a[0][0];
334  return s;
335  }
336 
337 #endif // DOXYGEN
338 
339  namespace FMatrixHelp {
340 
342  template <typename K>
343  static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
344  {
345  inverse[0][0] = 1.0/matrix[0][0];
346  return matrix[0][0];
347  }
348 
350  template <typename K>
351  static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
352  {
353  return invertMatrix(matrix,inverse);
354  }
355 
356 
358  template <typename K>
359  static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
360  {
361  // code generated by maple
362  K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
363  K det_1 = 1.0/det;
364  inverse[0][0] = matrix[1][1] * det_1;
365  inverse[0][1] = - matrix[0][1] * det_1;
366  inverse[1][0] = - matrix[1][0] * det_1;
367  inverse[1][1] = matrix[0][0] * det_1;
368  return det;
369  }
370 
373  template <typename K>
374  static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
375  {
376  // code generated by maple
377  K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
378  K det_1 = 1.0/det;
379  inverse[0][0] = matrix[1][1] * det_1;
380  inverse[1][0] = - matrix[0][1] * det_1;
381  inverse[0][1] = - matrix[1][0] * det_1;
382  inverse[1][1] = matrix[0][0] * det_1;
383  return det;
384  }
385 
387  template <typename K>
388  static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
389  {
390  // code generated by maple
391  K t4 = matrix[0][0] * matrix[1][1];
392  K t6 = matrix[0][0] * matrix[1][2];
393  K t8 = matrix[0][1] * matrix[1][0];
394  K t10 = matrix[0][2] * matrix[1][0];
395  K t12 = matrix[0][1] * matrix[2][0];
396  K t14 = matrix[0][2] * matrix[2][0];
397 
398  K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
399  t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
400  K t17 = 1.0/det;
401 
402  inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
403  inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
404  inverse[0][2] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
405  inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
406  inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
407  inverse[1][2] = -(t6-t10) * t17;
408  inverse[2][0] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
409  inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
410  inverse[2][2] = (t4-t8) * t17;
411 
412  return det;
413  }
414 
416  template <typename K>
417  static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
418  {
419  // code generated by maple
420  K t4 = matrix[0][0] * matrix[1][1];
421  K t6 = matrix[0][0] * matrix[1][2];
422  K t8 = matrix[0][1] * matrix[1][0];
423  K t10 = matrix[0][2] * matrix[1][0];
424  K t12 = matrix[0][1] * matrix[2][0];
425  K t14 = matrix[0][2] * matrix[2][0];
426 
427  K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
428  t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
429  K t17 = 1.0/det;
430 
431  inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
432  inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
433  inverse[2][0] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
434  inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
435  inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
436  inverse[2][1] = -(t6-t10) * t17;
437  inverse[0][2] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
438  inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
439  inverse[2][2] = (t4-t8) * t17;
440 
441  return det;
442  }
443 
445  template< class K, int m, int n, int p >
446  static inline void multMatrix ( const FieldMatrix< K, m, n > &A,
447  const FieldMatrix< K, n, p > &B,
449  {
450  typedef typename FieldMatrix< K, m, p > :: size_type size_type;
451 
452  for( size_type i = 0; i < m; ++i )
453  {
454  for( size_type j = 0; j < p; ++j )
455  {
456  ret[ i ][ j ] = K( 0 );
457  for( size_type k = 0; k < n; ++k )
458  ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
459  }
460  }
461  }
462 
464  template <typename K, int rows, int cols>
466  {
467  typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
468 
469  for(size_type i=0; i<cols; i++)
470  for(size_type j=0; j<cols; j++)
471  {
472  ret[i][j]=0.0;
473  for(size_type k=0; k<rows; k++)
474  ret[i][j]+=matrix[k][i]*matrix[k][j];
475  }
476  }
477 
478  using Dune::DenseMatrixHelp::multAssign;
479 
481  template <typename K, int rows, int cols>
482  static inline void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
483  {
484  typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
485 
486  for(size_type i=0; i<cols; ++i)
487  {
488  ret[i] = 0.0;
489  for(size_type j=0; j<rows; ++j)
490  ret[i] += matrix[j][i]*x[j];
491  }
492  }
493 
495  template <typename K, int rows, int cols>
497  {
499  multAssign(matrix,x,ret);
500  return ret;
501  }
502 
504  template <typename K, int rows, int cols>
506  {
508  multAssignTransposed( matrix, x, ret );
509  return ret;
510  }
511 
512  } // end namespace FMatrixHelp
513 
516 } // end namespace
517 
518 #include "fmatrixev.hh"
519 #endif
Macro for wrapping boundary checks.
A dense n x m matrix.
Definition: densematrix.hh:160
derived_type & operator-=(const DenseMatrix< Other > &y)
vector space subtraction
Definition: densematrix.hh:322
derived_type & operator+=(const DenseMatrix< Other > &y)
vector space addition
Definition: densematrix.hh:312
FieldMatrix< K, ROWS, COLS > & rightmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the right to this matrix.
Definition: densematrix.hh:634
size_type M() const
number of columns
Definition: densematrix.hh:692
derived_type & operator/=(const field_type &k)
vector space division by scalar
Definition: densematrix.hh:339
derived_type & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: densematrix.hh:331
Traits::row_type row_type
The type used to represent a row (must fulfill the Dune::DenseVector interface)
Definition: densematrix.hh:186
Traits::size_type size_type
The type used for the index access and size operation.
Definition: densematrix.hh:183
Traits::const_row_reference const_row_reference
The type used to represent a reference to a constant row (usually const row_type &)
Definition: densematrix.hh:192
Traits::row_reference row_reference
The type used to represent a reference to a row (usually row_type &)
Definition: densematrix.hh:189
@ blocklevel
The number of block levels we contain. This is 1.
Definition: densematrix.hh:197
A dense n x m matrix.
Definition: fmatrix.hh:68
@ rows
The number of rows.
Definition: fmatrix.hh:76
@ cols
The number of columns.
Definition: fmatrix.hh:78
FieldMatrix< K, l, cols > leftmultiplyany(const FieldMatrix< K, l, rows > &M) const
Multiplies M from the left to this matrix, this matrix is not modified.
Definition: fmatrix.hh:122
FieldMatrix< K, rows, l > rightmultiplyany(const FieldMatrix< K, cols, l > &M) const
Multiplies M from the right to this matrix, this matrix is not modified.
Definition: fmatrix.hh:157
FieldMatrix()
Default constructor.
Definition: fmatrix.hh:90
FieldMatrix(std::initializer_list< Dune::FieldVector< K, cols > > const &l)
Constructor initializing the matrix from a list of vector.
Definition: fmatrix.hh:94
FieldMatrix & rightmultiply(const FieldMatrix< K, r, c > &M)
Multiplies M from the right to this matrix.
Definition: fmatrix.hh:140
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Implements a matrix constructed from a given type representing a field and a compile-time given numbe...
A few common exception classes.
static K invertMatrix_retTransposed(const FieldMatrix< K, 1, 1 > &matrix, FieldMatrix< K, 1, 1 > &inverse)
invert scalar without changing the original matrix
Definition: fmatrix.hh:351
static FieldVector< K, rows > mult(const FieldMatrix< K, rows, cols > &matrix, const FieldVector< K, cols > &x)
calculates ret = matrix * x
Definition: fmatrix.hh:496
static void multMatrix(const FieldMatrix< K, m, n > &A, const FieldMatrix< K, n, p > &B, FieldMatrix< K, m, p > &ret)
calculates ret = A * B
Definition: fmatrix.hh:446
static K invertMatrix(const FieldMatrix< K, 1, 1 > &matrix, FieldMatrix< K, 1, 1 > &inverse)
invert scalar without changing the original matrix
Definition: fmatrix.hh:343
static FieldVector< K, cols > multTransposed(const FieldMatrix< K, rows, cols > &matrix, const FieldVector< K, rows > &x)
calculates ret = matrix^T * x
Definition: fmatrix.hh:505
static void multTransposedMatrix(const FieldMatrix< K, rows, cols > &matrix, FieldMatrix< K, cols, cols > &ret)
calculates ret= A_t*A
Definition: fmatrix.hh:465
static void multAssignTransposed(const FieldMatrix< K, rows, cols > &matrix, const FieldVector< K, rows > &x, FieldVector< K, cols > &ret)
calculates ret = matrix^T * x
Definition: fmatrix.hh:482
Eigenvalue computations for the FieldMatrix class.
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_ASSERT_BOUNDS(cond)
If DUNE_CHECK_BOUNDS is defined: check if condition cond holds; otherwise, do nothing.
Definition: boundschecking.hh:28
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
Dune namespace.
Definition: alignedallocator.hh:10
Various precision settings for calculations with FieldMatrix and FieldVector.
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 24, 22:30, 2024)