dune-istl  2.2.1
scaledidmatrix.hh
Go to the documentation of this file.
1 // -*- tab-width: 8; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2 // vi: set et ts=8 sw=4 sts=4:
3 #ifndef DUNE_SCALED_IDENTITY_MATRIX_HH
4 #define DUNE_SCALED_IDENTITY_MATRIX_HH
5 
12 #include<cmath>
13 #include<cstddef>
14 #include<complex>
15 #include<iostream>
16 #include <dune/common/exceptions.hh>
17 #include <dune/common/fmatrix.hh>
19 
20 namespace Dune {
21 
25 template<class K, int n>
27 {
29 
30 public:
31  //===== type definitions and constants
32 
34  typedef K field_type;
35 
37  typedef K block_type;
38 
40  typedef std::size_t size_type;
41 
43  enum {
46  };
47 
53 
55  enum {
57  rows = n,
59  cols = n
60  };
61 
62  //===== constructors
66 
69  ScaledIdentityMatrix (const K& k)
70  : p_(k)
71  {}
72 
73  //===== assignment from scalar
75  {
76  p_ = k;
77  return *this;
78  }
79 
80  // check if matrix is identical to other matrix (not only identical values)
81  bool identical(const ScaledIdentityMatrix<K,n>& other) const
82  {
83  return (this==&other);
84  }
85 
86  //===== iterator interface to rows of the matrix
90  typedef Iterator iterator;
94  typedef typename row_type::Iterator ColIterator;
95 
98  {
99  return Iterator(WrapperType(this),0);
100  }
101 
104  {
105  return Iterator(WrapperType(this),n);
106  }
107 
111  {
112  return Iterator(WrapperType(this),n-1);
113  }
114 
118  {
119  return Iterator(WrapperType(this),-1);
120  }
121 
122 
131 
134  {
135  return ConstIterator(WrapperType(this),0);
136  }
137 
140  {
141  return ConstIterator(WrapperType(this),n);
142  }
143 
147  {
148  return ConstIterator(WrapperType(this),n-1);
149  }
150 
154  {
155  return ConstIterator(WrapperType(this),-1);
156  }
157 
158  //===== vector space arithmetic
159 
162  {
163  p_ += y.p_;
164  return *this;
165  }
166 
169  {
170  p_ -= y.p_;
171  return *this;
172  }
173 
176  {
177  p_ += k;
178  return *this;
179  }
180 
183  {
184  p_ -= k;
185  return *this;
186  }
189  {
190  p_ *= k;
191  return *this;
192  }
193 
196  {
197  p_ /= k;
198  return *this;
199  }
200 
201  //===== comparison ops
202 
204  bool operator==(const ScaledIdentityMatrix& other) const
205  {
206  return p_==other.scalar();
207  }
208 
210  bool operator!=(const ScaledIdentityMatrix& other) const
211  {
212  return p_!=other.scalar();
213  }
214 
215  //===== linear maps
216 
218  template<class X, class Y>
219  void mv (const X& x, Y& y) const
220  {
221 #ifdef DUNE_FMatrix_WITH_CHECKING
222  if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
223  if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
224 #endif
225  for (size_type i=0; i<n; ++i)
226  y[i] = p_ * x[i];
227  }
228 
230  template<class X, class Y>
231  void mtv (const X& x, Y& y) const
232  {
233  mv(x, y);
234  }
235 
237  template<class X, class Y>
238  void umv (const X& x, Y& y) const
239  {
240 #ifdef DUNE_FMatrix_WITH_CHECKING
241  if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
242  if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
243 #endif
244  for (size_type i=0; i<n; ++i)
245  y[i] += p_ * x[i];
246  }
247 
249  template<class X, class Y>
250  void umtv (const X& x, Y& y) const
251  {
252 #ifdef DUNE_FMatrix_WITH_CHECKING
253  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
254  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
255 #endif
256  for (size_type i=0; i<n; ++i)
257  y[i] += p_ * x[i];
258  }
259 
261  template<class X, class Y>
262  void umhv (const X& x, Y& y) const
263  {
264 #ifdef DUNE_FMatrix_WITH_CHECKING
265  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
266  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
267 #endif
268  for (size_type i=0; i<n; i++)
269  y[i] += conjugateComplex(p_)*x[i];
270  }
271 
273  template<class X, class Y>
274  void mmv (const X& x, Y& y) const
275  {
276 #ifdef DUNE_FMatrix_WITH_CHECKING
277  if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
278  if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
279 #endif
280  for (size_type i=0; i<n; ++i)
281  y[i] -= p_ * x[i];
282  }
283 
285  template<class X, class Y>
286  void mmtv (const X& x, Y& y) const
287  {
288 #ifdef DUNE_FMatrix_WITH_CHECKING
289  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
290  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
291 #endif
292  for (size_type i=0; i<n; ++i)
293  y[i] -= p_ * x[i];
294  }
295 
297  template<class X, class Y>
298  void mmhv (const X& x, Y& y) const
299  {
300 #ifdef DUNE_FMatrix_WITH_CHECKING
301  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
302  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
303 #endif
304  for (size_type i=0; i<n; i++)
305  y[i] -= conjugateComplex(p_)*x[i];
306  }
307 
309  template<class X, class Y>
310  void usmv (const K& alpha, const X& x, Y& y) const
311  {
312 #ifdef DUNE_FMatrix_WITH_CHECKING
313  if (x.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
314  if (y.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
315 #endif
316  for (size_type i=0; i<n; i++)
317  y[i] += alpha * p_ * x[i];
318  }
319 
321  template<class X, class Y>
322  void usmtv (const K& alpha, const X& x, Y& y) const
323  {
324 #ifdef DUNE_FMatrix_WITH_CHECKING
325  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
326  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
327 #endif
328  for (size_type i=0; i<n; i++)
329  y[i] += alpha * p_ * x[i];
330  }
331 
333  template<class X, class Y>
334  void usmhv (const K& alpha, const X& x, Y& y) const
335  {
336 #ifdef DUNE_FMatrix_WITH_CHECKING
337  if (x.N()!=N()) DUNE_THROW(FMatrixError,"index out of range");
338  if (y.N()!=M()) DUNE_THROW(FMatrixError,"index out of range");
339 #endif
340  for (size_type i=0; i<n; i++)
341  y[i] += alpha * conjugateComplex(p_) * x[i];
342  }
343 
344  //===== norms
345 
347  double frobenius_norm () const
348  {
349  return fvmeta::sqrt(n*p_*p_);
350  }
351 
353  double frobenius_norm2 () const
354  {
355  return n*p_*p_;
356  }
357 
359  double infinity_norm () const
360  {
361  return std::abs(p_);
362  }
363 
365  double infinity_norm_real () const
366  {
367  return fvmeta::absreal(p_);
368  }
369 
370  //===== solve
371 
374  template<class V>
375  void solve (V& x, const V& b) const
376  {
377  for (int i=0; i<n; i++)
378  x[i] = b[i]/p_;
379  }
380 
383  void invert()
384  {
385  p_ = 1/p_;
386  }
387 
389  K determinant () const {
390  return std::pow(p_,n);
391  }
392 
393  //===== sizes
394 
396  size_type N () const
397  {
398  return n;
399  }
400 
402  size_type M () const
403  {
404  return n;
405  }
406 
407  //===== query
408 
410  bool exists (size_type i, size_type j) const
411  {
412 #ifdef DUNE_FMatrix_WITH_CHECKING
413  if (i<0 || i>=n) DUNE_THROW(FMatrixError,"row index out of range");
414  if (j<0 || j>=n) DUNE_THROW(FMatrixError,"column index out of range");
415 #endif
416  return i==j;
417  }
418 
419  //===== conversion operator
420 
422  friend std::ostream& operator<< (std::ostream& s, const ScaledIdentityMatrix<K,n>& a)
423  {
424  for (size_type i=0; i<n; i++) {
425  for (size_type j=0; j<n; j++)
426  s << ((i==j) ? a.p_ : 0) << " ";
427  s << std::endl;
428  }
429  return s;
430  }
431 
434  {
435  return reference(const_cast<K*>(&p_), i);
436  }
437 
440  {
441  return const_reference(const_cast<K*>(&p_), i);
442  }
443 
445  const K& diagonal(size_type i) const
446  {
447  return p_;
448  }
449 
452  {
453  return p_;
454  }
455 
458  const K& scalar() const
459  {
460  return p_;
461  }
462 
465  K& scalar()
466  {
467  return p_;
468  }
469 
470 private:
471  // the data, very simply a single number
472  K p_;
473 
474 };
475 
476 template<class M, class K, int n>
477 void istl_assign_to_fmatrix(DenseMatrix<M>& fm, const ScaledIdentityMatrix<K,n>& s)
478 {
479  fm = K();
480  for(int i=0; i<n; ++i)
481  fm[i][i] = s.scalar();
482 }
483 
484 } // end namespace
485 
486 #endif