Newer
Older
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
// $Id$
#ifndef DUNE_FMATRIX_HH
#define DUNE_FMATRIX_HH
#include <math.h>
#include <complex>
#include <iostream>
#include "exceptions.hh"
//#include <dune/istl/allocator.hh>
#include "fvector.hh"
#include "precision.hh"
/*! \file
\brief This file implements a matrix constructed from a given type
representing a field and compile-time given number of rows and columns.
*/
namespace Dune {
/**
@addtogroup DenseMatVec
@ingroup Common
@{
*/
template<class K, int n, int m> class FieldMatrix;
/** @brief Error thrown if operations of a FieldMatrix fail. */
class FMatrixError : public Exception {};
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
// template meta program for assignment from scalar
template<int I>
struct fmmeta_assignscalar {
template<class T, class K>
static void assignscalar (T* x, const K& k)
{
fmmeta_assignscalar<I-1>::assignscalar(x,k);
x[I] = k;
}
};
template<>
struct fmmeta_assignscalar<0> {
template<class T, class K>
static void assignscalar (T* x, const K& k)
{
x[0] = k;
}
};
// template meta program for operator+=
template<int I>
struct fmmeta_plusequal {
template<class T>
static void plusequal (T& x, const T& y)
{
x[I] += y[I];
fmmeta_plusequal<I-1>::plusequal(x,y);
}
};
template<>
struct fmmeta_plusequal<0> {
template<class T>
static void plusequal (T& x, const T& y)
{
x[0] += y[0];
}
};
// template meta program for operator-=
template<int I>
struct fmmeta_minusequal {
template<class T>
static void minusequal (T& x, const T& y)
{
x[I] -= y[I];
fmmeta_minusequal<I-1>::minusequal(x,y);
}
};
template<>
struct fmmeta_minusequal<0> {
template<class T>
static void minusequal (T& x, const T& y)
{
x[0] -= y[0];
}
};
// template meta program for operator*=
template<int I>
struct fmmeta_multequal {
template<class T, class K>
static void multequal (T& x, const K& k)
{
x[I] *= k;
fmmeta_multequal<I-1>::multequal(x,k);
}
};
template<>
struct fmmeta_multequal<0> {
template<class T, class K>
static void multequal (T& x, const K& k)
{
x[0] *= k;
}
};
// template meta program for operator/=
template<int I>
struct fmmeta_divequal {
template<class T, class K>
static void divequal (T& x, const K& k)
{
x[I] /= k;
fmmeta_divequal<I-1>::divequal(x,k);
}
};
template<>
struct fmmeta_divequal<0> {
template<class T, class K>
static void divequal (T& x, const K& k)
{
x[0] /= k;
}
};
// template meta program for dot
template<int I>
struct fmmeta_dot {
template<class X, class Y, class K>
static K dot (const X& x, const Y& y)
{
return x[I]*y[I] + fmmeta_dot<I-1>::template dot<X,Y,K>(x,y);
}
};
template<>
struct fmmeta_dot<0> {
template<class X, class Y, class K>
static K dot (const X& x, const Y& y)
{
return x[0]*y[0];
}
};
// template meta program for umv(x,y)
template<int I>
struct fmmeta_umv {
template<class Mat, class X, class Y, int c>
static void umv (const Mat& A, const X& x, Y& y)
{
typedef typename Mat::row_type R;
typedef typename Mat::field_type K;
y[I] += fmmeta_dot<c>::template dot<R,X,K>(A[I],x);
fmmeta_umv<I-1>::template umv<Mat,X,Y,c>(A,x,y);
}
};
template<>
struct fmmeta_umv<0> {
template<class Mat, class X, class Y, int c>
static void umv (const Mat& A, const X& x, Y& y)
{
typedef typename Mat::row_type R;
typedef typename Mat::field_type K;
y[0] += fmmeta_dot<c>::template dot<R,X,K>(A[0],x);
}
};
// template meta program for mmv(x,y)
template<int I>
struct fmmeta_mmv {
template<class Mat, class X, class Y, int c>
static void mmv (const Mat& A, const X& x, Y& y)
{
typedef typename Mat::row_type R;
typedef typename Mat::field_type K;
y[I] -= fmmeta_dot<c>::template dot<R,X,K>(A[I],x);
fmmeta_mmv<I-1>::template mmv<Mat,X,Y,c>(A,x,y);
}
};
template<>
struct fmmeta_mmv<0> {
template<class Mat, class X, class Y, int c>
static void mmv (const Mat& A, const X& x, Y& y)
{
typedef typename Mat::row_type R;
typedef typename Mat::field_type K;
y[0] -= fmmeta_dot<c>::template dot<R,X,K>(A[0],x);
}
};
template<class K, int n, int m, class X, class Y>
inline void fm_mmv (const FieldMatrix<K,n,m>& A, const X& x, Y& y)
{
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[i] -= A[i][j]*x[j];
}
template<class K>
inline void fm_mmv (const FieldMatrix<K,1,1>& A, const FieldVector<K,1>& x, FieldVector<K,1>& y)
{
y[0] -= A[0][0]*x[0];
}
// template meta program for usmv(x,y)
template<int I>
struct fmmeta_usmv {
template<class Mat, class K, class X, class Y, int c>
static void usmv (const Mat& A, const K& alpha, const X& x, Y& y)
{
typedef typename Mat::row_type R;
y[I] += alpha*fmmeta_dot<c>::template dot<R,X,K>(A[I],x);
fmmeta_usmv<I-1>::template usmv<Mat,K,X,Y,c>(A,alpha,x,y);
}
};
template<>
struct fmmeta_usmv<0> {
template<class Mat, class K, class X, class Y, int c>
static void usmv (const Mat& A, const K& alpha, const X& x, Y& y)
{
typedef typename Mat::row_type R;
y[0] += alpha*fmmeta_dot<c>::template dot<R,X,K>(A[0],x);
}
};
// conjugate komplex does nothing for non-complex types
template<class K>
inline K fm_ck (const K& k)
{
return k;
}
// conjugate komplex
template<class K>
inline std::complex<K> fm_ck (const std::complex<K>& c)
{
return std::complex<K>(c.real(),-c.imag());
}
//! solve small system
template<class K, int n, class V>
void fm_solve (const FieldMatrix<K,n,n>& Ain, V& x, const V& b)
{
// make a copy of a to store factorization
FieldMatrix<K,n,n> A(Ain);
// Gaussian elimination with maximum column pivot
double norm=A.infinity_norm_real(); // for relative thresholds
double pivthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::pivoting_limit());
double singthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::singular_limit());
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
V& rhs = x; // use x to store rhs
rhs = b; // copy data
// elimination phase
for (int i=0; i<n; i++) // loop over all rows
{
double pivmax=fvmeta_absreal(A[i][i]);
// pivoting ?
if (pivmax<pivthres)
{
// compute maximum of row
int imax=i; double abs;
for (int k=i+1; k<n; k++)
if ((abs=fvmeta_absreal(A[k][i]))>pivmax)
{
pivmax = abs; imax = k;
}
// swap rows
if (imax!=i)
for (int j=i; j<n; j++)
std::swap(A[i][j],A[imax][j]);
}
// singular ?
if (pivmax<singthres)
DUNE_THROW(FMatrixError,"matrix is singular");
// eliminate
for (int k=i+1; k<n; k++)
{
K factor = -A[k][i]/A[i][i];
for (int j=i+1; j<n; j++)
A[k][j] += factor*A[i][j];
rhs[k] += factor*rhs[i];
}
}
// backsolve
for (int i=n-1; i>=0; i--)
{
for (int j=i+1; j<n; j++)
rhs[i] -= A[i][j]*x[j];
x[i] = rhs[i]/A[i][i];
}
}
//! special case for 1x1 matrix, x and b may be identical
template<class K, class V>
inline void fm_solve (const FieldMatrix<K,1,1>& A, V& x, const V& b)
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (fvmeta_absreal(A[0][0])<FMatrixPrecision<>::absolute_limit())
DUNE_THROW(FMatrixError,"matrix is singular");
#endif
x[0] = b[0]/A[0][0];
}
//! special case for 2x2 matrix, x and b may be identical
template<class K, class V>
inline void fm_solve (const FieldMatrix<K,2,2>& A, V& x, const V& b)
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (fvmeta_absreal(detinv)<FMatrixPrecision<>::absolute_limit())
DUNE_THROW(FMatrixError,"matrix is singular");
detinv = 1/detinv;
#else
K detinv = 1.0/(A[0][0]*A[1][1]-A[0][1]*A[1][0]);
#endif
K temp = b[0];
x[0] = detinv*(A[1][1]*b[0]-A[0][1]*b[1]);
x[1] = detinv*(A[0][0]*b[1]-A[1][0]*temp);
}
//! compute inverse
template<class K, int n>
void fm_invert (FieldMatrix<K,n,n>& B)
{
FieldMatrix<K,n,n> A(B);
FieldMatrix<K,n,n>& L=A;
FieldMatrix<K,n,n>& U=A;
double norm=A.infinity_norm_real(); // for relative thresholds
double pivthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::pivoting_limit());
double singthres = std::max(FMatrixPrecision<>::absolute_limit(),norm*FMatrixPrecision<>::singular_limit());
// LU decomposition of A in A
for (int i=0; i<n; i++) // loop over all rows
{
double pivmax=fvmeta_absreal(A[i][i]);
// pivoting ?
if (pivmax<pivthres)
{
// compute maximum of column
int imax=i; double abs;
for (int k=i+1; k<n; k++)
if ((abs=fvmeta_absreal(A[k][i]))>pivmax)
{
pivmax = abs; imax = k;
}
// swap rows
if (imax!=i)
for (int j=i; j<n; j++)
std::swap(A[i][j],A[imax][j]);
}
// singular ?
if (pivmax<singthres)
DUNE_THROW(FMatrixError,"matrix is singular");
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
// eliminate
for (int k=i+1; k<n; k++)
{
K factor = A[k][i]/A[i][i];
L[k][i] = factor;
for (int j=i+1; j<n; j++)
A[k][j] -= factor*A[i][j];
}
}
// initialize inverse
B = 0;
for (int i=0; i<n; i++) B[i][i] = 1;
// L Y = I; multiple right hand sides
for (int i=0; i<n; i++)
for (int j=0; j<i; j++)
for (int k=0; k<n; k++)
B[i][k] -= L[i][j]*B[j][k];
// U A^{-1} = Y
for (int i=n-1; i>=0; i--)
for (int k=0; k<n; k++)
{
for (int j=i+1; j<n; j++)
B[i][k] -= U[i][j]*B[j][k];
B[i][k] /= U[i][i];
}
}
//! compute inverse n=1
template<class K>
void fm_invert (FieldMatrix<K,1,1>& A)
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (fvmeta_absreal(A[0][0])<FMatrixPrecision<>::absolute_limit())
DUNE_THROW(FMatrixError,"matrix is singular");
#endif
A[0][0] = 1/A[0][0];
}
//! compute inverse n=2
template<class K>
void fm_invert (FieldMatrix<K,2,2>& A)
{
K detinv = A[0][0]*A[1][1]-A[0][1]*A[1][0];
#ifdef DUNE_FMatrix_WITH_CHECKING
if (fvmeta_absreal(detinv)<FMatrixPrecision<>::absolute_limit())
DUNE_THROW(FMatrixError,"matrix is singular");
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
#endif
detinv = 1/detinv;
K temp=A[0][0];
A[0][0] = A[1][1]*detinv;
A[0][1] = -A[0][1]*detinv;
A[1][0] = -A[1][0]*detinv;
A[1][1] = temp*detinv;
}
//! left multiplication with matrix
template<class K, int n, int m>
void fm_leftmultiply (const FieldMatrix<K,n,n>& M, FieldMatrix<K,n,m>& A)
{
FieldMatrix<K,n,m> C(A);
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
{
A[i][j] = 0;
for (int k=0; k<n; k++)
A[i][j] += M[i][k]*C[k][j];
}
}
//! left multiplication with matrix, n=1
template<class K>
void fm_leftmultiply (const FieldMatrix<K,1,1>& M, FieldMatrix<K,1,1>& A)
{
A[0][0] *= M[0][0];
}
//! left multiplication with matrix, n=2
template<class K>
void fm_leftmultiply (const FieldMatrix<K,2,2>& M, FieldMatrix<K,2,2>& A)
{
FieldMatrix<K,2,2> C(A);
A[0][0] = M[0][0]*C[0][0] + M[0][1]*C[1][0];
A[0][1] = M[0][0]*C[0][1] + M[0][1]*C[1][1];
A[1][0] = M[1][0]*C[0][0] + M[1][1]*C[1][0];
A[1][1] = M[1][0]*C[0][1] + M[1][1]*C[1][1];
}
//! right multiplication with matrix
template<class K, int n, int m>
void fm_rightmultiply (const FieldMatrix<K,m,m>& M, FieldMatrix<K,n,m>& A)
{
FieldMatrix<K,n,m> C(A);
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
{
A[i][j] = 0;
for (int k=0; k<m; k++)
A[i][j] += C[i][k]*M[k][j];
}
}
//! right multiplication with matrix, n=1
template<class K>
void fm_rightmultiply (const FieldMatrix<K,1,1>& M, FieldMatrix<K,1,1>& A)
{
A[0][0] *= M[0][0];
}
//! right multiplication with matrix, n=2
template<class K>
void fm_rightmultiply (const FieldMatrix<K,2,2>& M, FieldMatrix<K,2,2>& A)
{
FieldMatrix<K,2,2> C(A);
A[0][0] = C[0][0]*M[0][0] + C[0][1]*M[1][0];
A[0][1] = C[0][0]*M[0][1] + C[0][1]*M[1][1];
A[1][0] = C[1][0]*M[0][0] + C[1][1]*M[1][0];
A[1][1] = C[1][0]*M[0][1] + C[1][1]*M[1][1];
}
/** 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.
Implementation of all members uses template meta programs where appropriate
*/
template<class K, int n, int m>
class FieldMatrix
{
public:
// standard constructor and everything is sufficient ...
//===== type definitions and constants
//! export the type representing the field
typedef K field_type;
//! export the type representing the components
typedef K block_type;
//! We are at the leaf of the block recursion
enum {blocklevel = 1};
//! Each row is implemented by a field vector
typedef FieldVector<K,m> row_type;
//! export size
enum {rows = n, cols = m};
//===== constructors
/** \brief Default constructor
*/
FieldMatrix () {}
/** \brief Constructor initializing the whole matrix with a scalar
*/
FieldMatrix (const K& k)
{
for (int i=0; i<n; i++) p[i] = k;
}
//===== random access interface to rows of the matrix
//! random access to the rows
row_type& operator[] (int i)
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
#endif
return p[i];
}
//! same for read only access
const row_type& operator[] (int i) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (i<0 || i>=n) DUNE_THROW(FMatrixError,"index out of range");
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
#endif
return p[i];
}
//===== iterator interface to rows of the matrix
// forward declaration
class ConstIterator;
//! Iterator access to rows
class Iterator
{
public:
//! constructor
Iterator (row_type* _p, int _i)
{
p = _p;
i = _i;
}
//! empty constructor, use with care!
Iterator ()
{ }
//! prefix increment
Iterator& operator++()
{
++i;
return *this;
}
//! prefix decrement
Iterator& operator--()
{
--i;
return *this;
}
//! equality
bool operator== (const Iterator& it) const
{
return (p+i)==(it.p+it.i);
}
//! inequality
bool operator!= (const Iterator& it) const
{
return (p+i)!=(it.p+it.i);
}
//! dereferencing
row_type& operator* ()
{
return p[i];
}
//! arrow
row_type* operator-> ()
{
return p+i;
}
//! return index
int index ()
{
return i;
}
friend class ConstIterator;
private:
row_type* p;
int i;
};
//! begin iterator
Iterator begin ()
{
return Iterator(p,0);
}
//! end iterator
Iterator end ()
{
return Iterator(p,n);
}
//! begin iterator
Iterator rbegin ()
{
return Iterator(p,n-1);
}
//! end iterator
Iterator rend ()
{
return Iterator(p,-1);
}
//! rename the iterators for easier access
typedef Iterator RowIterator;
typedef typename row_type::Iterator ColIterator;
//! Iterator access to rows
class ConstIterator
{
public:
//! constructor
ConstIterator (const row_type* _p, int _i) : p(_p), i(_i)
{ }
//! empty constructor, use with care!
ConstIterator ()
{
p = 0;
i = 0;
}
//! prefix increment
ConstIterator& operator++()
{
++i;
return *this;
}
//! prefix decrement
ConstIterator& operator--()
{
--i;
return *this;
}
//! equality
bool operator== (const ConstIterator& it) const
{
return (p+i)==(it.p+it.i);
}
//! inequality
bool operator!= (const ConstIterator& it) const
{
return (p+i)!=(it.p+it.i);
}
//! equality
bool operator== (const Iterator& it) const
{
return (p+i)==(it.p+it.i);
}
//! inequality
bool operator!= (const Iterator& it) const
{
return (p+i)!=(it.p+it.i);
}
//! dereferencing
const row_type& operator* () const
{
return p[i];
}
//! arrow
const row_type* operator-> () const
{
return p+i;
}
//! return index
int index () const
{
return i;
}
friend class Iterator;
private:
const row_type* p;
int i;
};
//! begin iterator
ConstIterator begin () const
{
return ConstIterator(p,0);
}
//! end iterator
ConstIterator end () const
{
return ConstIterator(p,n);
}
//! begin iterator
ConstIterator rbegin () const
{
return ConstIterator(p,n-1);
}
//! end iterator
ConstIterator rend () const
{
return ConstIterator(p,-1);
}
//! rename the iterators for easier access
typedef ConstIterator ConstRowIterator;
typedef typename row_type::ConstIterator ConstColIterator;
//===== assignment from scalar
FieldMatrix& operator= (const K& k)
{
fmmeta_assignscalar<n-1>::assignscalar(p,k);
return *this;
}
//===== vector space arithmetic
//! vector space addition
FieldMatrix& operator+= (const FieldMatrix& y)
{
fmmeta_plusequal<n-1>::plusequal(*this,y);
return *this;
}
//! vector space subtraction
FieldMatrix& operator-= (const FieldMatrix& y)
{
fmmeta_minusequal<n-1>::minusequal(*this,y);
return *this;
}
//! vector space multiplication with scalar
FieldMatrix& operator*= (const K& k)
{
fmmeta_multequal<n-1>::multequal(*this,k);
return *this;
}
//! vector space division by scalar
FieldMatrix& operator/= (const K& k)
{
fmmeta_divequal<n-1>::divequal(*this,k);
return *this;
}
//===== linear maps
//! y += A x
template<class X, class Y>
void umv (const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
#endif
fmmeta_umv<n-1>::template umv<FieldMatrix,X,Y,m-1>(*this,x,y);
}
//! y += A^T x
template<class X, class Y>
void umtv (const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
#endif
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[j] += p[i][j]*x[i];
}
//! y += A^H x
template<class X, class Y>
void umhv (const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
#endif
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[j] += fm_ck(p[i][j])*x[i];
}
//! y -= A x
template<class X, class Y>
void mmv (const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
#endif
fmmeta_mmv<n-1>::template mmv<FieldMatrix,X,Y,m-1>(*this,x,y);
//fm_mmv(*this,x,y);
}
//! y -= A^T x
template<class X, class Y>
void mmtv (const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
#endif
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[j] -= p[i][j]*x[i];
}
//! y -= A^H x
template<class X, class Y>
void mmhv (const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
#endif
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[j] -= fm_ck(p[i][j])*x[i];
}
//! y += alpha A x
template<class X, class Y>
void usmv (const K& alpha, const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
#endif
fmmeta_usmv<n-1>::template usmv<FieldMatrix,K,X,Y,m-1>(*this,alpha,x,y);
}
//! y += alpha A^T x
template<class X, class Y>
void usmtv (const K& alpha, const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
#endif
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[j] += alpha*p[i][j]*x[i];
}
//! y += alpha A^H x
template<class X, class Y>
void usmhv (const K& alpha, const X& x, Y& y) const
{
#ifdef DUNE_FMatrix_WITH_CHECKING
if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
#endif
for (int i=0; i<n; i++)
for (int j=0; j<m; j++)
y[j] += alpha*fm_ck(p[i][j])*x[i];
}
//===== norms
//! frobenius norm: sqrt(sum over squared values of entries)
double frobenius_norm () const
{
double sum=0;
for (int i=0; i<n; ++i) sum += p[i].two_norm2();
return sqrt(sum);
}
//! square of frobenius norm, need for block recursion
double frobenius_norm2 () const
{
double sum=0;
for (int i=0; i<n; ++i) sum += p[i].two_norm2();
return sum;
}
//! infinity norm (row sum norm, how to generalize for blocks?)
double infinity_norm () const
{
double max=0;
for (int i=0; i<n; ++i) max = std::max(max,p[i].one_norm());
return max;
}
//! simplified infinity norm (uses Manhattan norm for complex values)
double infinity_norm_real () const
{
double max=0;
for (int i=0; i<n; ++i) max = std::max(max,p[i].one_norm_real());
return max;
}
//===== solve
/** \brief Solve system A x = b
*
* \exception FMatrixError if the matrix is singular
*/
template<class V>
void solve (V& x, const V& b) const
{
fm_solve(*this,x,b);
}
/** \brief Compute inverse
*
* \exception FMatrixError if the matrix is singular
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
*/
void invert ()
{
fm_invert(*this);
}
//! calculates the determinant of this matrix
K determinant () const;
//! Multiplies M from the left to this matrix
FieldMatrix& leftmultiply (const FieldMatrix<K,n,n>& M)
{
fm_leftmultiply(M,*this);
return *this;
}
//! Multiplies M from the right to this matrix
FieldMatrix& rightmultiply (const FieldMatrix<K,n,n>& M)
{
fm_rightmultiply(M,*this);
return *this;
}
//===== sizes
//! number of blocks in row direction
int N () const
{