3 #ifndef DUNE_SCALED_IDENTITY_MATRIX_HH
4 #define DUNE_SCALED_IDENTITY_MATRIX_HH
16 #include <dune/common/exceptions.hh>
17 #include <dune/common/fmatrix.hh>
25 template<
class K,
int n>
83 return (
this==&other);
206 return p_==other.
scalar();
212 return p_!=other.
scalar();
218 template<
class X,
class Y>
219 void mv (
const X& x, Y& y)
const
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");
230 template<
class X,
class Y>
231 void mtv (
const X& x, Y& y)
const
237 template<
class X,
class Y>
238 void umv (
const X& x, Y& y)
const
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");
249 template<
class X,
class Y>
250 void umtv (
const X& x, Y& y)
const
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");
261 template<
class X,
class Y>
262 void umhv (
const X& x, Y& y)
const
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");
269 y[i] += conjugateComplex(p_)*x[i];
273 template<
class X,
class Y>
274 void mmv (
const X& x, Y& y)
const
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");
285 template<
class X,
class Y>
286 void mmtv (
const X& x, Y& y)
const
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");
297 template<
class X,
class Y>
298 void mmhv (
const X& x, Y& y)
const
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");
305 y[i] -= conjugateComplex(p_)*x[i];
309 template<
class X,
class Y>
310 void usmv (
const K& alpha,
const X& x, Y& y)
const
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");
317 y[i] += alpha * p_ * x[i];
321 template<
class X,
class Y>
322 void usmtv (
const K& alpha,
const X& x, Y& y)
const
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");
329 y[i] += alpha * p_ * x[i];
333 template<
class X,
class Y>
334 void usmhv (
const K& alpha,
const X& x, Y& y)
const
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");
341 y[i] += alpha * conjugateComplex(p_) * x[i];
349 return fvmeta::sqrt(n*p_*p_);
367 return fvmeta::absreal(p_);
377 for (
int i=0; i<n; i++)
390 return std::pow(p_,n);
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");
422 friend std::ostream& operator<< (std::ostream& s, const ScaledIdentityMatrix<K,n>& a)
426 s << ((i==j) ? a.p_ : 0) <<
" ";
435 return reference(const_cast<K*>(&p_), i);
476 template<
class M,
class K,
int n>
480 for(
int i=0; i<n; ++i)