File: fmatrix.hh

package info (click to toggle)
dune-common 2.11.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 6,048 kB
  • sloc: cpp: 54,403; python: 4,136; sh: 1,657; makefile: 17
file content (860 lines) | stat: -rw-r--r-- 28,765 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
#ifndef DUNE_FMATRIX_HH
#define DUNE_FMATRIX_HH

#include <cmath>
#include <cstddef>
#include <iostream>
#include <algorithm>
#include <initializer_list>

#include <dune/common/boundschecking.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/fvector.hh>
#include <dune/common/densematrix.hh>
#include <dune/common/precision.hh>
#include <dune/common/promotiontraits.hh>
#include <dune/common/typetraits.hh>
#include <dune/common/matrixconcepts.hh>

namespace Dune
{

  namespace Impl
  {

    template<class M>
    class ColumnVectorView
    {
    public:

      using value_type = typename M::value_type;
      using size_type = typename M::size_type;

      constexpr ColumnVectorView(M& matrix, size_type col) :
        matrix_(matrix),
        col_(col)
      {}

      constexpr size_type N () const {
        return matrix_.N();
      }

      template<class M_ = M,
        std::enable_if_t<std::is_same_v<M_,M> and not std::is_const_v<M_>, int> = 0>
      constexpr value_type& operator[] (size_type row) {
        return matrix_[row][col_];
      }

      constexpr const value_type& operator[] (size_type row) const {
        return matrix_[row][col_];
      }

    protected:
      M& matrix_;
      const size_type col_;
    };

  }

  template<typename M>
  struct FieldTraits< Impl::ColumnVectorView<M> >
  {
    using field_type = typename FieldTraits<M>::field_type;
    using real_type = typename FieldTraits<M>::real_type;
  };

  /**
      @addtogroup DenseMatVec
      @{
   */

  /*! \file

     \brief  Implements a matrix constructed from a given type
     representing a field and compile-time given number of rows and columns.
   */

  template< class K, int ROWS, int COLS = ROWS > class FieldMatrix;


  template< class K, int ROWS, int COLS >
  struct DenseMatVecTraits< FieldMatrix<K,ROWS,COLS> >
  {
    typedef FieldMatrix<K,ROWS,COLS> derived_type;

    // each row is implemented by a field vector
    typedef FieldVector<K,COLS> row_type;

    typedef row_type &row_reference;
    typedef const row_type &const_row_reference;

    typedef std::array<row_type,ROWS> container_type;
    typedef K value_type;
    typedef typename container_type::size_type size_type;
  };

  template< class K, int ROWS, int COLS >
  struct FieldTraits< FieldMatrix<K,ROWS,COLS> >
  {
    typedef typename FieldTraits<K>::field_type field_type;
    typedef typename FieldTraits<K>::real_type real_type;
  };

  /**
      @brief A dense n x m matrix.

     Matrices represent linear maps from a vector space V to a vector space W.
       This class represents such a linear map by storing a two-dimensional
       %array of numbers of a given field type K. The number of rows and
       columns is given at compile time.
   */
  template<class K, int ROWS, int COLS>
  class FieldMatrix : public DenseMatrix< FieldMatrix<K,ROWS,COLS> >
  {
    template<class,int,int> friend class FieldMatrix;
    std::array< FieldVector<K,COLS>, ROWS > _data;
    typedef DenseMatrix< FieldMatrix<K,ROWS,COLS> > Base;
  public:

    //! The number of rows.
    constexpr static int rows = ROWS;
    //! The number of columns.
    constexpr static int cols = COLS;

    typedef typename Base::size_type size_type;
    typedef typename Base::row_type row_type;

    typedef typename Base::row_reference row_reference;
    typedef typename Base::const_row_reference const_row_reference;

    //===== constructors
    /** \brief Default constructor
     */
    constexpr FieldMatrix() = default;

    /** \brief Constructor initializing the matrix from a list of vector
     */
    constexpr FieldMatrix(std::initializer_list<Dune::FieldVector<K, cols> > const &l) {
      assert(l.size() == rows); // Actually, this is not needed any more!
      for(std::size_t i=0; i<std::min<std::size_t>(ROWS, l.size()); ++i)
        _data[i] = std::data(l)[i];
    }

    //! copy constructor
    FieldMatrix(const FieldMatrix&) = default;

    //! copy constructor from assignable type T
    template <class T,
              typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>>
    constexpr FieldMatrix(T const& rhs)
      : _data{}
    {
      *this = rhs;
    }

    using Base::operator=;

    //! copy assignment operator
    constexpr FieldMatrix& operator=(const FieldMatrix&) = default;

    //! copy assignment from FieldMatrix over a different field
    template<typename T>
    constexpr FieldMatrix& operator=(const FieldMatrix<T, ROWS, COLS>& x)
    {
      // The copy must be done element-by-element since a std::array does not have
      // a converting assignment operator.
      for (std::size_t i = 0; i < _data.size(); ++i)
        _data[i] = x._data[i];
      return *this;
    }

    //! no copy assignment from FieldMatrix of different size
    template <typename T, int rows, int cols>
    FieldMatrix& operator=(FieldMatrix<T,rows,cols> const&) = delete;

    //! Return transposed of the matrix as FieldMatrix
    constexpr FieldMatrix<K, COLS, ROWS> transposed() const
    {
      Dune::FieldMatrix<K, COLS, ROWS> AT;
      for( int i = 0; i < ROWS; ++i )
        for( int j = 0; j < COLS; ++j )
          AT[j][i] = (*this)[i][j];
      return AT;
    }

    //! vector space addition -- two-argument version
    template <class OtherScalar>
    friend constexpr auto operator+ ( const FieldMatrix& matrixA,
                            const FieldMatrix<OtherScalar,ROWS,COLS>& matrixB)
    {
      FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,ROWS,COLS> result;

      for (size_type i = 0; i < ROWS; ++i)
        for (size_type j = 0; j < COLS; ++j)
          result[i][j] = matrixA[i][j] + matrixB[i][j];

      return result;
    }

    //! vector space subtraction -- two-argument version
    template <class OtherScalar>
    friend constexpr auto operator- ( const FieldMatrix& matrixA,
                            const FieldMatrix<OtherScalar,ROWS,COLS>& matrixB)
    {
      FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,ROWS,COLS> result;

      for (size_type i = 0; i < ROWS; ++i)
        for (size_type j = 0; j < COLS; ++j)
          result[i][j] = matrixA[i][j] - matrixB[i][j];

      return result;
    }

    //! vector space multiplication with scalar
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator* ( const FieldMatrix& matrix, Scalar scalar)
    {
      FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,ROWS,COLS> result;

      for (size_type i = 0; i < ROWS; ++i)
        for (size_type j = 0; j < COLS; ++j)
          result[i][j] = matrix[i][j] * scalar;

      return result;
    }

    //! vector space multiplication with scalar
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator* ( Scalar scalar, const FieldMatrix& matrix)
    {
      FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,ROWS,COLS> result;

      for (size_type i = 0; i < ROWS; ++i)
        for (size_type j = 0; j < COLS; ++j)
          result[i][j] = scalar * matrix[i][j];

      return result;
    }

    //! vector space division by scalar
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator/ ( const FieldMatrix& matrix, Scalar scalar)
    {
      FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,ROWS,COLS> result;

      for (size_type i = 0; i < ROWS; ++i)
        for (size_type j = 0; j < COLS; ++j)
          result[i][j] = matrix[i][j] / scalar;

      return result;
    }

    /** \brief Matrix-matrix multiplication
     */
    template <class OtherScalar, int otherCols>
    friend constexpr auto operator* ( const FieldMatrix& matrixA,
                            const FieldMatrix<OtherScalar, COLS, otherCols>& matrixB)
    {
      FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,ROWS,otherCols> result;

      for (size_type i = 0; i < matrixA.mat_rows(); ++i)
        for (size_type j = 0; j < matrixB.mat_cols(); ++j)
        {
          result[i][j] = 0;
          for (size_type k = 0; k < matrixA.mat_cols(); ++k)
            result[i][j] += matrixA[i][k] * matrixB[k][j];
        }

      return result;
    }

    /** \brief Matrix-matrix multiplication
     *
     * This implements multiplication of a FieldMatrix with another matrix
     * of type OtherMatrix. The latter has to provide
     * OtherMatrix::field_type, OtherMatrix::cols, and OtherMatrix::mtv(x,y).
     */
    template <class OtherMatrix, std::enable_if_t<
      Impl::IsStaticSizeMatrix_v<OtherMatrix>
      and not Impl::IsFieldMatrix_v<OtherMatrix>
      , int> = 0>
    friend constexpr auto operator* ( const FieldMatrix& matrixA,
                            const OtherMatrix& matrixB)
    {
      using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
      Dune::FieldMatrix<Field, rows ,OtherMatrix::cols> result;
      for (std::size_t j=0; j<rows; ++j)
        matrixB.mtv(matrixA[j], result[j]);
      return result;
    }

    /** \brief Matrix-matrix multiplication
     *
     * This implements multiplication of another matrix
     * of type OtherMatrix with a FieldMatrix. The former has to provide
     * OtherMatrix::field_type, OtherMatrix::rows, and OtherMatrix::mv(x,y).
     */
    template <class OtherMatrix, std::enable_if_t<
      Impl::IsStaticSizeMatrix_v<OtherMatrix>
      and not Impl::IsFieldMatrix_v<OtherMatrix>
      , int> = 0>
    friend constexpr auto operator* ( const OtherMatrix& matrixA,
                            const FieldMatrix& matrixB)
    {
      using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
      Dune::FieldMatrix<Field, OtherMatrix::rows, cols> result;
      for (std::size_t j=0; j<cols; ++j)
      {
        auto B_j = Impl::ColumnVectorView(matrixB, j);
        auto result_j = Impl::ColumnVectorView(result, j);
        matrixA.mv(B_j, result_j);
      }
      return result;
    }

    //! Multiplies M from the left to this matrix, this matrix is not modified
    template<int l>
    constexpr FieldMatrix<K,l,cols> leftmultiplyany (const FieldMatrix<K,l,rows>& M) const
    {
      FieldMatrix<K,l,cols> C;

      for (size_type i=0; i<l; i++) {
        for (size_type j=0; j<cols; j++) {
          C[i][j] = 0;
          for (size_type k=0; k<rows; k++)
            C[i][j] += M[i][k]*(*this)[k][j];
        }
      }
      return C;
    }

    using Base::rightmultiply;

    //! Multiplies M from the right to this matrix
    template <int r, int c>
    constexpr FieldMatrix& rightmultiply (const FieldMatrix<K,r,c>& M)
    {
      static_assert(r == c, "Cannot rightmultiply with non-square matrix");
      static_assert(r == cols, "Size mismatch");
      FieldMatrix<K,rows,cols> C(*this);

      for (size_type i=0; i<rows; i++)
        for (size_type j=0; j<cols; j++) {
          (*this)[i][j] = 0;
          for (size_type k=0; k<cols; k++)
            (*this)[i][j] += C[i][k]*M[k][j];
        }
      return *this;
    }

    //! Multiplies M from the right to this matrix, this matrix is not modified
    template<int l>
    constexpr FieldMatrix<K,rows,l> rightmultiplyany (const FieldMatrix<K,cols,l>& M) const
    {
      FieldMatrix<K,rows,l> C;

      for (size_type i=0; i<rows; i++) {
        for (size_type j=0; j<l; j++) {
          C[i][j] = 0;
          for (size_type k=0; k<cols; k++)
            C[i][j] += (*this)[i][k]*M[k][j];
        }
      }
      return C;
    }

    // make this thing a matrix
    static constexpr size_type mat_rows() { return ROWS; }
    static constexpr size_type mat_cols() { return COLS; }

    constexpr row_reference mat_access ( size_type i )
    {
      DUNE_ASSERT_BOUNDS(i < ROWS);
      return _data[i];
    }

    constexpr const_row_reference mat_access ( size_type i ) const
    {
      DUNE_ASSERT_BOUNDS(i < ROWS);
      return _data[i];
    }
  };

#ifndef DOXYGEN // hide specialization
  /** \brief Special type for 1x1 matrices
   */
  template<class K>
  class FieldMatrix<K,1,1> : public DenseMatrix< FieldMatrix<K,1,1> >
  {
    FieldVector<K,1> _data;
    typedef DenseMatrix< FieldMatrix<K,1,1> > Base;
  public:
    // standard constructor and everything is sufficient ...

    //===== type definitions and constants

    //! The type used for index access and size operations
    typedef typename Base::size_type size_type;

    //! The number of block levels we contain.
    //! This is always one for this type.
    constexpr static int blocklevel = 1;

    typedef typename Base::row_type row_type;

    typedef typename Base::row_reference row_reference;
    typedef typename Base::const_row_reference const_row_reference;

    //! \brief The number of rows.
    //! This is always one for this type.
    constexpr static int rows = 1;
    //! \brief The number of columns.
    //! This is always one for this type.
    constexpr static int cols = 1;

    //===== constructors
    /** \brief Default constructor
     */
    constexpr FieldMatrix() = default;

    /** \brief Constructor initializing the matrix from a list of vector
     */
    FieldMatrix(std::initializer_list<Dune::FieldVector<K, 1>> const &l)
    {
      std::copy_n(l.begin(), std::min<std::size_t>(1, l.size()), &_data);
    }

    template <class T,
              typename = std::enable_if_t<HasDenseMatrixAssigner<FieldMatrix, T>::value>>
    constexpr FieldMatrix(T const& rhs)
    {
      *this = rhs;
    }

    using Base::operator=;

    //! Return transposed of the matrix as FieldMatrix
    constexpr FieldMatrix<K, 1, 1> transposed() const
    {
      return *this;
    }

    //! vector space addition -- two-argument version
    template <class OtherScalar>
    friend constexpr auto operator+ ( const FieldMatrix& matrixA,
                            const FieldMatrix<OtherScalar,1,1>& matrixB)
    {
      return FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,1,1>{matrixA[0][0] + matrixB[0][0]};
    }

    //! Binary addition when treating FieldMatrix<K,1,1> like K
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator+ ( const FieldMatrix& matrix,
                            const Scalar& scalar)
    {
      return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1>{matrix[0][0] + scalar};
    }

    //! Binary addition when treating FieldMatrix<K,1,1> like K
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator+ ( const Scalar& scalar,
                            const FieldMatrix& matrix)
    {
      return FieldMatrix<typename PromotionTraits<Scalar,K>::PromotedType,1,1>{scalar + matrix[0][0]};
    }

    //! vector space subtraction -- two-argument version
    template <class OtherScalar>
    friend constexpr auto operator- ( const FieldMatrix& matrixA,
                            const FieldMatrix<OtherScalar,1,1>& matrixB)
    {
      return FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,1,1>{matrixA[0][0] - matrixB[0][0]};
    }

    //! Binary subtraction when treating FieldMatrix<K,1,1> like K
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator- ( const FieldMatrix& matrix,
                            const Scalar& scalar)
    {
      return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1>{matrix[0][0] - scalar};
    }

    //! Binary subtraction when treating FieldMatrix<K,1,1> like K
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator- ( const Scalar& scalar,
                            const FieldMatrix& matrix)
    {
      return FieldMatrix<typename PromotionTraits<Scalar,K>::PromotedType,1,1>{scalar - matrix[0][0]};
    }

    //! vector space multiplication with scalar
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator* ( const FieldMatrix& matrix, Scalar scalar)
    {
      return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1> {matrix[0][0] * scalar};
    }

    //! vector space multiplication with scalar
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator* ( Scalar scalar, const FieldMatrix& matrix)
    {
      return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1> {scalar * matrix[0][0]};
    }

    //! vector space division by scalar
    template <class Scalar,
              std::enable_if_t<IsNumber<Scalar>::value, int> = 0>
    friend constexpr auto operator/ ( const FieldMatrix& matrix, Scalar scalar)
    {
      return FieldMatrix<typename PromotionTraits<K,Scalar>::PromotedType,1,1> {matrix[0][0] / scalar};
    }

    //===== solve

        /** \brief Matrix-matrix multiplication
     */
    template <class OtherScalar, int otherCols>
    friend constexpr auto operator* ( const FieldMatrix& matrixA,
                            const FieldMatrix<OtherScalar, 1, otherCols>& matrixB)
    {
      FieldMatrix<typename PromotionTraits<K,OtherScalar>::PromotedType,1,otherCols> result;

      for (size_type j = 0; j < matrixB.mat_cols(); ++j)
        result[0][j] = matrixA[0][0] * matrixB[0][j];

      return result;
    }

    /** \brief Matrix-matrix multiplication
     *
     * This implements multiplication of a FieldMatrix with another matrix
     * of type OtherMatrix. The latter has to provide
     * OtherMatrix::field_type, OtherMatrix::cols, and OtherMatrix::mtv(x,y).
     */
    template <class OtherMatrix, std::enable_if_t<
      Impl::IsStaticSizeMatrix_v<OtherMatrix>
      and not Impl::IsFieldMatrix_v<OtherMatrix>
      and (OtherMatrix::rows==1)
      , int> = 0>
    friend constexpr auto operator* ( const FieldMatrix& matrixA,
                            const OtherMatrix& matrixB)
    {
      using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
      Dune::FieldMatrix<Field, rows ,OtherMatrix::cols> result;
      for (std::size_t j=0; j<rows; ++j)
        matrixB.mtv(matrixA[j], result[j]);
      return result;
    }

    /** \brief Matrix-matrix multiplication
     *
     * This implements multiplication of another matrix
     * of type OtherMatrix with a FieldMatrix. The former has to provide
     * OtherMatrix::field_type, OtherMatrix::rows, and OtherMatrix::mv(x,y).
     */
    template <class OtherMatrix, std::enable_if_t<
      Impl::IsStaticSizeMatrix_v<OtherMatrix>
      and not Impl::IsFieldMatrix_v<OtherMatrix>
      and (OtherMatrix::cols==1)
      , int> = 0>
    friend constexpr auto operator* ( const OtherMatrix& matrixA,
                            const FieldMatrix& matrixB)
    {
      using Field = typename PromotionTraits<K, typename OtherMatrix::field_type>::PromotedType;
      Dune::FieldMatrix<Field, OtherMatrix::rows, cols> result;
      for (std::size_t j=0; j<cols; ++j)
      {
        auto B_j = Impl::ColumnVectorView(matrixB, j);
        auto result_j = Impl::ColumnVectorView(result, j);
        matrixA.mv(B_j, result_j);
      }
      return result;
    }

    //! Multiplies M from the left to this matrix, this matrix is not modified
    template<int l>
    constexpr FieldMatrix<K,l,1> leftmultiplyany (const FieldMatrix<K,l,1>& M) const
    {
      FieldMatrix<K,l,1> C;
      for (size_type j=0; j<l; j++)
        C[j][0] = M[j][0]*(*this)[0][0];
      return C;
    }

    //! left multiplication
    constexpr FieldMatrix& rightmultiply (const FieldMatrix& M)
    {
      _data[0] *= M[0][0];
      return *this;
    }

    //! Multiplies M from the right to this matrix, this matrix is not modified
    template<int l>
    constexpr FieldMatrix<K,1,l> rightmultiplyany (const FieldMatrix<K,1,l>& M) const
    {
      FieldMatrix<K,1,l> C;

      for (size_type j=0; j<l; j++)
        C[0][j] = M[0][j]*_data[0];
      return C;
    }

    // make this thing a matrix
    static constexpr size_type mat_rows() { return 1; }
    static constexpr size_type mat_cols() { return 1; }

    constexpr row_reference mat_access ([[maybe_unused]] size_type i)
    {
      DUNE_ASSERT_BOUNDS(i == 0);
      return _data;
    }

    constexpr const_row_reference mat_access ([[maybe_unused]] size_type i) const
    {
      DUNE_ASSERT_BOUNDS(i == 0);
      return _data;
    }

    //! add scalar
    constexpr FieldMatrix& operator+= (const K& k)
    {
      _data[0] += k;
      return (*this);
    }

    //! subtract scalar
    constexpr FieldMatrix& operator-= (const K& k)
    {
      _data[0] -= k;
      return (*this);
    }

    //! multiplication with scalar
    constexpr FieldMatrix& operator*= (const K& k)
    {
      _data[0] *= k;
      return (*this);
    }

    //! division by scalar
    constexpr FieldMatrix& operator/= (const K& k)
    {
      _data[0] /= k;
      return (*this);
    }

    //===== conversion operator

    constexpr operator const K& () const { return _data[0]; }

  };

  /** \brief Sends the matrix to an output stream */
  template<typename K>
  std::ostream& operator<< (std::ostream& s, const FieldMatrix<K,1,1>& a)
  {
    s << a[0][0];
    return s;
  }

#endif // DOXYGEN

  namespace FMatrixHelp {

    //! invert scalar without changing the original matrix
    template <typename K>
    static constexpr K invertMatrix (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
    {
      using real_type = typename FieldTraits<K>::real_type;
      inverse[0][0] = real_type(1.0)/matrix[0][0];
      return matrix[0][0];
    }

    //! invert scalar without changing the original matrix
    template <typename K>
    static constexpr K invertMatrix_retTransposed (const FieldMatrix<K,1,1> &matrix, FieldMatrix<K,1,1> &inverse)
    {
      return invertMatrix(matrix,inverse);
    }


    //! invert 2x2 Matrix without changing the original matrix
    template <typename K>
    static constexpr K invertMatrix (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
    {
      using real_type = typename FieldTraits<K>::real_type;
      // code generated by maple
      K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
      K det_1 = real_type(1.0)/det;
      inverse[0][0] =   matrix[1][1] * det_1;
      inverse[0][1] = - matrix[0][1] * det_1;
      inverse[1][0] = - matrix[1][0] * det_1;
      inverse[1][1] =   matrix[0][0] * det_1;
      return det;
    }

    //! invert 2x2 Matrix without changing the original matrix
    //! return transposed matrix
    template <typename K>
    static constexpr K invertMatrix_retTransposed (const FieldMatrix<K,2,2> &matrix, FieldMatrix<K,2,2> &inverse)
    {
      using real_type = typename FieldTraits<K>::real_type;
      // code generated by maple
      K det = (matrix[0][0]*matrix[1][1] - matrix[0][1]*matrix[1][0]);
      K det_1 = real_type(1.0)/det;
      inverse[0][0] =   matrix[1][1] * det_1;
      inverse[1][0] = - matrix[0][1] * det_1;
      inverse[0][1] = - matrix[1][0] * det_1;
      inverse[1][1] =   matrix[0][0] * det_1;
      return det;
    }

    //! invert 3x3 Matrix without changing the original matrix
    template <typename K>
    static constexpr K invertMatrix (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
    {
      using real_type = typename FieldTraits<K>::real_type;
      // code generated by maple
      K t4  = matrix[0][0] * matrix[1][1];
      K t6  = matrix[0][0] * matrix[1][2];
      K t8  = matrix[0][1] * matrix[1][0];
      K t10 = matrix[0][2] * matrix[1][0];
      K t12 = matrix[0][1] * matrix[2][0];
      K t14 = matrix[0][2] * matrix[2][0];

      K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
               t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
      K t17 = real_type(1.0)/det;

      inverse[0][0] =  (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
      inverse[0][1] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
      inverse[0][2] =  (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
      inverse[1][0] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
      inverse[1][1] =  (matrix[0][0] * matrix[2][2] - t14) * t17;
      inverse[1][2] = -(t6-t10) * t17;
      inverse[2][0] =  (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
      inverse[2][1] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
      inverse[2][2] =  (t4-t8) * t17;

      return det;
    }

    //! invert 3x3 Matrix without changing the original matrix
    template <typename K>
    static constexpr K invertMatrix_retTransposed (const FieldMatrix<K,3,3> &matrix, FieldMatrix<K,3,3> &inverse)
    {
      using real_type = typename FieldTraits<K>::real_type;
      // code generated by maple
      K t4  = matrix[0][0] * matrix[1][1];
      K t6  = matrix[0][0] * matrix[1][2];
      K t8  = matrix[0][1] * matrix[1][0];
      K t10 = matrix[0][2] * matrix[1][0];
      K t12 = matrix[0][1] * matrix[2][0];
      K t14 = matrix[0][2] * matrix[2][0];

      K det = (t4*matrix[2][2]-t6*matrix[2][1]-t8*matrix[2][2]+
               t10*matrix[2][1]+t12*matrix[1][2]-t14*matrix[1][1]);
      K t17 = real_type(1.0)/det;

      inverse[0][0] =  (matrix[1][1] * matrix[2][2] - matrix[1][2] * matrix[2][1])*t17;
      inverse[1][0] = -(matrix[0][1] * matrix[2][2] - matrix[0][2] * matrix[2][1])*t17;
      inverse[2][0] =  (matrix[0][1] * matrix[1][2] - matrix[0][2] * matrix[1][1])*t17;
      inverse[0][1] = -(matrix[1][0] * matrix[2][2] - matrix[1][2] * matrix[2][0])*t17;
      inverse[1][1] =  (matrix[0][0] * matrix[2][2] - t14) * t17;
      inverse[2][1] = -(t6-t10) * t17;
      inverse[0][2] =  (matrix[1][0] * matrix[2][1] - matrix[1][1] * matrix[2][0]) * t17;
      inverse[1][2] = -(matrix[0][0] * matrix[2][1] - t12) * t17;
      inverse[2][2] =  (t4-t8) * t17;

      return det;
    }

    //! calculates ret = A * B
    template< class K, int m, int n, int p >
    static constexpr void multMatrix ( const FieldMatrix< K, m, n > &A,
                                    const FieldMatrix< K, n, p > &B,
                                    FieldMatrix< K, m, p > &ret )
    {
      typedef typename FieldMatrix< K, m, p > :: size_type size_type;

      for( size_type i = 0; i < m; ++i )
      {
        for( size_type j = 0; j < p; ++j )
        {
          ret[ i ][ j ] = K( 0 );
          for( size_type k = 0; k < n; ++k )
            ret[ i ][ j ] += A[ i ][ k ] * B[ k ][ j ];
        }
      }
    }

    //! calculates ret= A_t*A
    template <typename K, int rows, int cols>
    static constexpr void multTransposedMatrix(const FieldMatrix<K,rows,cols> &matrix, FieldMatrix<K,cols,cols>& ret)
    {
      typedef typename FieldMatrix<K,rows,cols>::size_type size_type;

      for(size_type i=0; i<cols; i++)
        for(size_type j=0; j<cols; j++)
        {
          ret[i][j]=0.0;
          for(size_type k=0; k<rows; k++)
            ret[i][j]+=matrix[k][i]*matrix[k][j];
        }
    }

    using Dune::DenseMatrixHelp::multAssign;

    //! calculates ret = matrix^T * x
    template <typename K, int rows, int cols>
    static constexpr void multAssignTransposed( const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x, FieldVector<K,cols> & ret)
    {
      typedef typename FieldMatrix<K,rows,cols>::size_type size_type;

      for(size_type i=0; i<cols; ++i)
      {
        ret[i] = 0.0;
        for(size_type j=0; j<rows; ++j)
          ret[i] += matrix[j][i]*x[j];
      }
    }

    //! calculates ret = matrix * x
    template <typename K, int rows, int cols>
    static constexpr FieldVector<K,rows> mult(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,cols> & x)
    {
      FieldVector<K,rows> ret;
      multAssign(matrix,x,ret);
      return ret;
    }

    //! calculates ret = matrix^T * x
    template <typename K, int rows, int cols>
    static constexpr FieldVector<K,cols> multTransposed(const FieldMatrix<K,rows,cols> &matrix, const FieldVector<K,rows> & x)
    {
      FieldVector<K,cols> ret;
      multAssignTransposed( matrix, x, ret );
      return ret;
    }

  } // end namespace FMatrixHelp

  /** @} end documentation */

} // end namespace

#include "fmatrixev.hh"
#endif