Dune Core Modules (2.9.1)

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// SPDX-FileCopyrightInfo: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_FMATRIX_HH
6#define DUNE_FMATRIX_HH
7
8#include <cmath>
9#include <cstddef>
10#include <iostream>
11#include <algorithm>
12#include <initializer_list>
13
21#include <dune/common/matrixconcepts.hh>
22
23namespace Dune
24{
25
26 namespace Impl
27 {
28
29 template<class M>
30 class ColumnVectorView
31 {
32 public:
33
34 using value_type = typename M::value_type;
35 using size_type = typename M::size_type;
36
37 constexpr ColumnVectorView(M& matrix, size_type col) :
38 matrix_(matrix),
39 col_(col)
40 {}
41
42 constexpr size_type N () const {
43 return matrix_.N();
44 }
45
46 template<class M_ = M,
47 std::enable_if_t<std::is_same_v<M_,M> and not std::is_const_v<M_>, int> = 0>
48 constexpr value_type& operator[] (size_type row) {
49 return matrix_[row][col_];
50 }
51
52 constexpr const value_type& operator[] (size_type row) const {
53 return matrix_[row][col_];
54 }
55
56 protected:
57 M& matrix_;
58 const size_type col_;
59 };
60
61 }
62
63 template<typename M>
64 struct FieldTraits< Impl::ColumnVectorView<M> >
65 {
66 using field_type = typename FieldTraits<M>::field_type;
67 using real_type = typename FieldTraits<M>::real_type;
68 };
69
81 template< class K, int ROWS, int COLS = ROWS > class FieldMatrix;
82
83
84 template< class K, int ROWS, int COLS >
85 struct DenseMatVecTraits< FieldMatrix<K,ROWS,COLS> >
86 {
87 typedef FieldMatrix<K,ROWS,COLS> derived_type;
88
89 // each row is implemented by a field vector
90 typedef FieldVector<K,COLS> row_type;
91
92 typedef row_type &row_reference;
93 typedef const row_type &const_row_reference;
94
95 typedef std::array<row_type,ROWS> container_type;
96 typedef K value_type;
97 typedef typename container_type::size_type size_type;
98 };
99
100 template< class K, int ROWS, int COLS >
101 struct FieldTraits< FieldMatrix<K,ROWS,COLS> >
102 {
103 typedef typename FieldTraits<K>::field_type field_type;
104 typedef typename FieldTraits<K>::real_type real_type;
105 };
106
115 template<class K, int ROWS, int COLS>
116 class FieldMatrix : public DenseMatrix< FieldMatrix<K,ROWS,COLS> >
117 {
118 std::array< FieldVector<K,COLS>, ROWS > _data;
120 public:
121
123 constexpr static int rows = ROWS;
125 constexpr static int cols = COLS;
126
127 typedef typename Base::size_type size_type;
128 typedef typename Base::row_type row_type;
129
130 typedef typename Base::row_reference row_reference;
132
133 //===== constructors
136 constexpr FieldMatrix() = default;
137
140 FieldMatrix(std::initializer_list<Dune::FieldVector<K, cols> > const &l) {
141 assert(l.size() == rows); // Actually, this is not needed any more!
142 std::copy_n(l.begin(), std::min(static_cast<std::size_t>(ROWS),
143 l.size()),
144 _data.begin());
145 }
146
147 template <class T,
148 typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>>
149 FieldMatrix(T const& rhs)
150 {
151 *this = rhs;
152 }
153
154 using Base::operator=;
155
158
160 template<typename T>
162 {
163 _data = x._data;
164 return *this;
165 }
166
168 template <typename T, int rows, int cols>
170
173 {
175 for( int i = 0; i < ROWS; ++i )
176 for( int j = 0; j < COLS; ++j )
177 AT[j][i] = (*this)[i][j];
178 return AT;
179 }
180
182 template <class OtherScalar>
183 friend auto operator+ ( const FieldMatrix& matrixA,
185 {
187
188 for (size_type i = 0; i < ROWS; ++i)
189 for (size_type j = 0; j < COLS; ++j)
190 result[i][j] = matrixA[i][j] + matrixB[i][j];
191
192 return result;
193 }
194
196 template <class OtherScalar>
197 friend auto operator- ( const FieldMatrix& matrixA,
199 {
201
202 for (size_type i = 0; i < ROWS; ++i)
203 for (size_type j = 0; j < COLS; ++j)
204 result[i][j] = matrixA[i][j] - matrixB[i][j];
205
206 return result;
207 }
208
210 template <class Scalar,
211 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
212 friend auto operator* ( const FieldMatrix& matrix, Scalar scalar)
213 {
215
216 for (size_type i = 0; i < ROWS; ++i)
217 for (size_type j = 0; j < COLS; ++j)
218 result[i][j] = matrix[i][j] * scalar;
219
220 return result;
221 }
222
224 template <class Scalar,
225 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
226 friend auto operator* ( Scalar scalar, const FieldMatrix& matrix)
227 {
229
230 for (size_type i = 0; i < ROWS; ++i)
231 for (size_type j = 0; j < COLS; ++j)
232 result[i][j] = scalar * matrix[i][j];
233
234 return result;
235 }
236
238 template <class Scalar,
239 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
240 friend auto operator/ ( const FieldMatrix& matrix, Scalar scalar)
241 {
243
244 for (size_type i = 0; i < ROWS; ++i)
245 for (size_type j = 0; j < COLS; ++j)
246 result[i][j] = matrix[i][j] / scalar;
247
248 return result;
249 }
250
253 template <class OtherScalar, int otherCols>
254 friend auto operator* ( const FieldMatrix& matrixA,
256 {
258
259 for (size_type i = 0; i < matrixA.mat_rows(); ++i)
260 for (size_type j = 0; j < matrixB.mat_cols(); ++j)
261 {
262 result[i][j] = 0;
263 for (size_type k = 0; k < matrixA.mat_cols(); ++k)
264 result[i][j] += matrixA[i][k] * matrixB[k][j];
265 }
266
267 return result;
268 }
269
276 template <class OtherMatrix, std::enable_if_t<
277 Impl::IsStaticSizeMatrix_v<OtherMatrix>
278 and not Impl::IsFieldMatrix_v<OtherMatrix>
279 , int> = 0>
280 friend auto operator* ( const FieldMatrix& matrixA,
281 const OtherMatrix& matrixB)
282 {
283 using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
285 for (std::size_t j=0; j<rows; ++j)
286 matrixB.mtv(matrixA[j], result[j]);
287 return result;
288 }
289
296 template <class OtherMatrix, std::enable_if_t<
297 Impl::IsStaticSizeMatrix_v<OtherMatrix>
298 and not Impl::IsFieldMatrix_v<OtherMatrix>
299 , int> = 0>
300 friend auto operator* ( const OtherMatrix& matrixA,
301 const FieldMatrix& matrixB)
302 {
303 using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
305 for (std::size_t j=0; j<cols; ++j)
306 {
307 auto B_j = Impl::ColumnVectorView(matrixB, j);
308 auto result_j = Impl::ColumnVectorView(result, j);
309 matrixA.mv(B_j, result_j);
310 }
311 return result;
312 }
313
315 template<int l>
317 {
319
320 for (size_type i=0; i<l; i++) {
321 for (size_type j=0; j<cols; j++) {
322 C[i][j] = 0;
323 for (size_type k=0; k<rows; k++)
324 C[i][j] += M[i][k]*(*this)[k][j];
325 }
326 }
327 return C;
328 }
329
331
333 template <int r, int c>
335 {
336 static_assert(r == c, "Cannot rightmultiply with non-square matrix");
337 static_assert(r == cols, "Size mismatch");
339
340 for (size_type i=0; i<rows; i++)
341 for (size_type j=0; j<cols; j++) {
342 (*this)[i][j] = 0;
343 for (size_type k=0; k<cols; k++)
344 (*this)[i][j] += C[i][k]*M[k][j];
345 }
346 return *this;
347 }
348
350 template<int l>
352 {
354
355 for (size_type i=0; i<rows; i++) {
356 for (size_type j=0; j<l; j++) {
357 C[i][j] = 0;
358 for (size_type k=0; k<cols; k++)
359 C[i][j] += (*this)[i][k]*M[k][j];
360 }
361 }
362 return C;
363 }
364
365 // make this thing a matrix
366 static constexpr size_type mat_rows() { return ROWS; }
367 static constexpr size_type mat_cols() { return COLS; }
368
369 row_reference mat_access ( size_type i )
370 {
371 DUNE_ASSERT_BOUNDS(i < ROWS);
372 return _data[i];
373 }
374
375 const_row_reference mat_access ( size_type i ) const
376 {
377 DUNE_ASSERT_BOUNDS(i < ROWS);
378 return _data[i];
379 }
380 };
381
382#ifndef DOXYGEN // hide specialization
385 template<class K>
386 class FieldMatrix<K,1,1> : public DenseMatrix< FieldMatrix<K,1,1> >
387 {
388 FieldVector<K,1> _data;
389 typedef DenseMatrix< FieldMatrix<K,1,1> > Base;
390 public:
391 // standard constructor and everything is sufficient ...
392
393 //===== type definitions and constants
394
396 typedef typename Base::size_type size_type;
397
400 constexpr static int blocklevel = 1;
401
402 typedef typename Base::row_type row_type;
403
404 typedef typename Base::row_reference row_reference;
405 typedef typename Base::const_row_reference const_row_reference;
406
409 constexpr static int rows = 1;
412 constexpr static int cols = 1;
413
414 //===== constructors
417 constexpr FieldMatrix() = default;
418
421 FieldMatrix(std::initializer_list<Dune::FieldVector<K, 1>> const &l)
422 {
423 std::copy_n(l.begin(), std::min(static_cast< std::size_t >( 1 ), l.size()), &_data);
424 }
425
426 template <class T,
427 typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>>
428 FieldMatrix(T const& rhs)
429 {
430 *this = rhs;
431 }
432
433 using Base::operator=;
434
436 FieldMatrix<K, 1, 1> transposed() const
437 {
438 return *this;
439 }
440
442 template <class OtherScalar>
443 friend auto operator+ ( const FieldMatrix& matrixA,
444 const FieldMatrix<OtherScalar,1,1>& matrixB)
445 {
446 return FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,1,1>{matrixA[0][0] + matrixB[0][0]};
447 }
448
450 template <class Scalar,
451 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
452 friend auto operator+ ( const FieldMatrix& matrix,
453 const Scalar& scalar)
454 {
455 return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1>{matrix[0][0] + scalar};
456 }
457
459 template <class Scalar,
460 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
461 friend auto operator+ ( const Scalar& scalar,
462 const FieldMatrix& matrix)
463 {
464 return FieldMatrix<typename PromotionTraits<Scalar,K>::PromotedType,1,1>{scalar + matrix[0][0]};
465 }
466
468 template <class OtherScalar>
469 friend auto operator- ( const FieldMatrix& matrixA,
470 const FieldMatrix<OtherScalar,1,1>& matrixB)
471 {
472 return FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,1,1>{matrixA[0][0] - matrixB[0][0]};
473 }
474
476 template <class Scalar,
477 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
478 friend auto operator- ( const FieldMatrix& matrix,
479 const Scalar& scalar)
480 {
481 return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1>{matrix[0][0] - scalar};
482 }
483
485 template <class Scalar,
486 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
487 friend auto operator- ( const Scalar& scalar,
488 const FieldMatrix& matrix)
489 {
490 return FieldMatrix<typename PromotionTraits<Scalar,K>::PromotedType,1,1>{scalar - matrix[0][0]};
491 }
492
494 template <class Scalar,
495 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
496 friend auto operator* ( const FieldMatrix& matrix, Scalar scalar)
497 {
498 return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1> {matrix[0][0] * scalar};
499 }
500
502 template <class Scalar,
503 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
504 friend auto operator* ( Scalar scalar, const FieldMatrix& matrix)
505 {
506 return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1> {scalar * matrix[0][0]};
507 }
508
510 template <class Scalar,
511 std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
512 friend auto operator/ ( const FieldMatrix& matrix, Scalar scalar)
513 {
514 return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1> {matrix[0][0] / scalar};
515 }
516
517 //===== solve
518
521 template <class OtherScalar, int otherCols>
522 friend auto operator* ( const FieldMatrix& matrixA,
523 const FieldMatrix<OtherScalar, 1, otherCols>& matrixB)
524 {
525 FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,1,otherCols> result;
526
527 for (size_type j = 0; j < matrixB.mat_cols(); ++j)
528 result[0][j] = matrixA[0][0] * matrixB[0][j];
529
530 return result;
531 }
532
539 template <class OtherMatrix, std::enable_if_t<
540 Impl::IsStaticSizeMatrix_v<OtherMatrix>
541 and not Impl::IsFieldMatrix_v<OtherMatrix>
542 and (OtherMatrix::rows==1)
543 , int> = 0>
544 friend auto operator* ( const FieldMatrix& matrixA,
545 const OtherMatrix& matrixB)
546 {
547 using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
549 for (std::size_t j=0; j<rows; ++j)
550 matrixB.mtv(matrixA[j], result[j]);
551 return result;
552 }
553
560 template <class OtherMatrix, std::enable_if_t<
561 Impl::IsStaticSizeMatrix_v<OtherMatrix>
562 and not Impl::IsFieldMatrix_v<OtherMatrix>
563 and (OtherMatrix::cols==1)
564 , int> = 0>
565 friend auto operator* ( const OtherMatrix& matrixA,
566 const FieldMatrix& matrixB)
567 {
568 using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
570 for (std::size_t j=0; j<cols; ++j)
571 {
572 auto B_j = Impl::ColumnVectorView(matrixB, j);
573 auto result_j = Impl::ColumnVectorView(result, j);
574 matrixA.mv(B_j, result_j);
575 }
576 return result;
577 }
578
580 template<int l>
581 FieldMatrix<K,l,1> leftmultiplyany (const FieldMatrix<K,l,1>& M) const
582 {
583 FieldMatrix<K,l,1> C;
584 for (size_type j=0; j<l; j++)
585 C[j][0] = M[j][0]*(*this)[0][0];
586 return C;
587 }
588
591 {
592 _data[0] *= M[0][0];
593 return *this;
594 }
595
597 template<int l>
598 FieldMatrix<K,1,l> rightmultiplyany (const FieldMatrix<K,1,l>& M) const
599 {
600 FieldMatrix<K,1,l> C;
601
602 for (size_type j=0; j<l; j++)
603 C[0][j] = M[0][j]*_data[0];
604 return C;
605 }
606
607 // make this thing a matrix
608 static constexpr size_type mat_rows() { return 1; }
609 static constexpr size_type mat_cols() { return 1; }
610
611 row_reference mat_access ([[maybe_unused]] size_type i)
612 {
613 DUNE_ASSERT_BOUNDS(i == 0);
614 return _data;
615 }
616
617 const_row_reference mat_access ([[maybe_unused]] size_type i) const
618 {
619 DUNE_ASSERT_BOUNDS(i == 0);
620 return _data;
621 }
622
624 FieldMatrix& operator+= (const K& k)
625 {
626 _data[0] += k;
627 return (*this);
628 }
629
631 FieldMatrix& operator-= (const K& k)
632 {
633 _data[0] -= k;
634 return (*this);
635 }
636
638 FieldMatrix& operator*= (const K& k)
639 {
640 _data[0] *= k;
641 return (*this);
642 }
643
645 FieldMatrix& operator/= (const K& k)
646 {
647 _data[0] /= k;
648 return (*this);
649 }
650
651 //===== conversion operator
652
653 operator const K& () const { return _data[0]; }
654
655 };
656
658 template<typename K>
659 std::ostream& operator<< (std::ostream& s, const FieldMatrix<K,1,1>& a)
660 {
661 s << a[0][0];
662 return s;
663 }
664
665#endif // DOXYGEN
666
667 namespace FMatrixHelp {
668
670 template <typename K>
671 static inline K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
672 {
673 using real_type = typename FieldTraits<K>::real_type;
674 inverse[0][0] = real_type(1.0)/matrix[0][0];
675 return matrix[0][0];
676 }
677
679 template <typename K>
680 static inline K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
681 {
682 return invertMatrix(matrix,inverse);
683 }
684
685
687 template <typename K>
688 static inline K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
689 {
690 using real_type = typename FieldTraits<K>::real_type;
691 // code generated by maple
692 K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
693 K det_1 = real_type(1.0)/det;
694 inverse[0][0] = matrix[1][1] * det_1;
695 inverse[0][1] = - matrix[0][1] * det_1;
696 inverse[1][0] = - matrix[1][0] * det_1;
697 inverse[1][1] = matrix[0][0] * det_1;
698 return det;
699 }
700
703 template <typename K>
704 static inline K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
705 {
706 using real_type = typename FieldTraits<K>::real_type;
707 // code generated by maple
708 K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
709 K det_1 = real_type(1.0)/det;
710 inverse[0][0] = matrix[1][1] * det_1;
711 inverse[1][0] = - matrix[0][1] * det_1;
712 inverse[0][1] = - matrix[1][0] * det_1;
713 inverse[1][1] = matrix[0][0] * det_1;
714 return det;
715 }
716
718 template <typename K>
719 static inline K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
720 {
721 using real_type = typename FieldTraits<K>::real_type;
722 // code generated by maple
723 K t4 = matrix[0][0] * matrix[1][1];
724 K t6 = matrix[0][0] * matrix[1][2];
725 K t8 = matrix[0][1] * matrix[1][0];
726 K t10 = matrix[0][2] * matrix[1][0];
727 K t12 = matrix[0][1] * matrix[2][0];
728 K t14 = matrix[0][2] * matrix[2][0];
729
730 K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
731 t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
732 K t17 = real_type(1.0)/det;
733
734 inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
735 inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
736 inverse[0][2] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
737 inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
738 inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
739 inverse[1][2] = -(t6-t10) * t17;
740 inverse[2][0] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
741 inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
742 inverse[2][2] = (t4-t8) * t17;
743
744 return det;
745 }
746
748 template <typename K>
749 static inline K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
750 {
751 using real_type = typename FieldTraits<K>::real_type;
752 // code generated by maple
753 K t4 = matrix[0][0] * matrix[1][1];
754 K t6 = matrix[0][0] * matrix[1][2];
755 K t8 = matrix[0][1] * matrix[1][0];
756 K t10 = matrix[0][2] * matrix[1][0];
757 K t12 = matrix[0][1] * matrix[2][0];
758 K t14 = matrix[0][2] * matrix[2][0];
759
760 K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
761 t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
762 K t17 = real_type(1.0)/det;
763
764 inverse[0][0] = (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
765 inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
766 inverse[2][0] = (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
767 inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
768 inverse[1][1] = (matrix[0][0] * matrix[2][2] - t14) * t17;
769 inverse[2][1] = -(t6-t10) * t17;
770 inverse[0][2] = (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
771 inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
772 inverse[2][2] = (t4-t8) * t17;
773
774 return det;
775 }
776
778 template< class K, int m, int n, int p >
779 static inline void multMatrix ( const FieldMatrix< K, m, n > &A,
780 const FieldMatrix< K, n, p > &B,
782 {
783 typedef typename FieldMatrix< K, m, p > :: size_type size_type;
784
785 for( size_type i = 0; i < m; ++i )
786 {
787 for( size_type j = 0; j < p; ++j )
788 {
789 ret[ i ][ j ] = K( 0 );
790 for( size_type k = 0; k < n; ++k )
791 ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
792 }
793 }
794 }
795
797 template <typename K, int rows, int cols>
799 {
800 typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
801
802 for(size_type i=0; i<cols; i++)
803 for(size_type j=0; j<cols; j++)
804 {
805 ret[i][j]=0.0;
806 for(size_type k=0; k<rows; k++)
807 ret[i][j]+=matrix[k][i]*matrix[k][j];
808 }
809 }
810
811 using Dune::DenseMatrixHelp::multAssign;
812
814 template <typename K, int rows, int cols>
816 {
817 typedef typename FieldMatrix<K,rows,cols>::size_type size_type;
818
819 for(size_type i=0; i<cols; ++i)
820 {
821 ret[i] = 0.0;
822 for(size_type j=0; j<rows; ++j)
823 ret[i] += matrix[j][i]*x[j];
824 }
825 }
826
828 template <typename K, int rows, int cols>
830 {
832 multAssign(matrix,x,ret);
833 return ret;
834 }
835
837 template <typename K, int rows, int cols>
839 {
841 multAssignTransposed( matrix, x, ret );
842 return ret;
843 }
844
845 } // end namespace FMatrixHelp
846
849} // end namespace
850
851#include "fmatrixev.hh"
852#endif
Macro for wrapping boundary checks.
A dense n x m matrix.
Definition: densematrix.hh:140
derived_type operator-() const
Matrix negation.
Definition: densematrix.hh:298
void mtv(const X &x, Y &y) const
y = A^T x
Definition: densematrix.hh:387
constexpr size_type M() const
number of columns
Definition: densematrix.hh:703
FieldMatrix< K, ROWS, COLS > & rightmultiply(const DenseMatrix< M2 > &M)
Multiplies M from the right to this matrix.
Definition: densematrix.hh:645
derived_type & operator/=(const field_type &k)
vector space division by scalar
Definition: densematrix.hh:329
derived_type & operator*=(const field_type &k)
vector space multiplication with scalar
Definition: densematrix.hh:321
derived_type & operator-=(const DenseMatrix< Other > &x)
vector space subtraction
Definition: densematrix.hh:312
static constexpr int blocklevel
The number of block levels we contain. This is the leaf, that is, 1.
Definition: densematrix.hh:178
Traits::row_type row_type
The type used to represent a row (must fulfill the Dune::DenseVector interface)
Definition: densematrix.hh:169
Traits::size_type size_type
The type used for the index access and size operation.
Definition: densematrix.hh:166
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:175
Traits::row_reference row_reference
The type used to represent a reference to a row (usually row_type &)
Definition: densematrix.hh:172
derived_type & operator+=(const DenseMatrix< Other > &x)
vector space addition
Definition: densematrix.hh:289
A dense n x m matrix.
Definition: fmatrix.hh:117
constexpr FieldMatrix()=default
Default constructor.
FieldMatrix & operator=(const FieldMatrix< T, ROWS, COLS > &x)
copy assignment from FieldMatrix over a different field
Definition: fmatrix.hh:161
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:351
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:316
FieldMatrix & rightmultiply(const FieldMatrix< K, r, c > &M)
Multiplies M from the right to this matrix.
Definition: fmatrix.hh:334
friend auto operator*(const FieldMatrix &matrix, Scalar scalar)
vector space multiplication with scalar
Definition: fmatrix.hh:212
FieldMatrix & operator=(FieldMatrix< T, rows, cols > const &)=delete
no copy assignment from FieldMatrix of different size
static constexpr int rows
The number of rows.
Definition: fmatrix.hh:123
FieldMatrix< K, COLS, ROWS > transposed() const
Return transposed of the matrix as FieldMatrix.
Definition: fmatrix.hh:172
FieldMatrix(std::initializer_list< Dune::FieldVector< K, cols > > const &l)
Constructor initializing the matrix from a list of vector.
Definition: fmatrix.hh:140
static constexpr int cols
The number of columns.
Definition: fmatrix.hh:125
friend auto operator/(const FieldMatrix &matrix, Scalar scalar)
vector space division by scalar
Definition: fmatrix.hh:240
friend auto operator+(const FieldMatrix &matrixA, const FieldMatrix< OtherScalar, ROWS, COLS > &matrixB)
vector space addition – two-argument version
Definition: fmatrix.hh:183
FieldMatrix & operator=(const FieldMatrix &)=default
copy assignment operator
Implements a matrix constructed from a given type representing a field and a compile-time given numbe...
A few common exception classes.
static FieldVector< K, cols > multTransposed(const FieldMatrix< K, rows, cols > &matrix, const FieldVector< K, rows > &x)
calculates ret = matrix^T * x
Definition: fmatrix.hh:838
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:680
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:779
static K invertMatrix(const FieldMatrix< K, 1, 1 > &matrix, FieldMatrix< K, 1, 1 > &inverse)
invert scalar without changing the original matrix
Definition: fmatrix.hh:671
static FieldVector< K, rows > mult(const FieldMatrix< K, rows, cols > &matrix, const FieldVector< K, cols > &x)
calculates ret = matrix * x
Definition: fmatrix.hh:829
static void multTransposedMatrix(const FieldMatrix< K, rows, cols > &matrix, FieldMatrix< K, cols, cols > &ret)
calculates ret= A_t*A
Definition: fmatrix.hh:798
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:815
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:30
auto min(ADLTag< 0 >, const V &v1, const V &v2)
implements binary Simd::min()
Definition: defaults.hh:89
typename Overloads::ScalarType< std::decay_t< V > >::type Scalar
Element type of some SIMD type.
Definition: interface.hh:235
Dune namespace.
Definition: alignedallocator.hh:13
Various precision settings for calculations with FieldMatrix and FieldVector.
Compute type of the result of an arithmetic operation involving two different number types.
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)