Discontinuous Galerkin Library
#include "dg/algorithm.h"
Loading...
Searching...
No Matches
operator.h
Go to the documentation of this file.
1#ifndef _DG_OPERATORS_DYN_
2#define _DG_OPERATORS_DYN_
3
4#include <vector>
5#include <algorithm>
6#include <cmath>
7#include <iterator>
8#include <stdexcept>
9#include <cassert>
10#include "dlt.h"
11#include "../blas1.h" // reproducible DOT products
12#include "../blas2.h" // Matrix tensor traits
13
14namespace dg{
15
16
29template< class T>
31{
32 public:
33 typedef T value_type;
37 SquareMatrix() = default;
43 explicit SquareMatrix( const unsigned n): n_(n), data_(n_*n_){}
50 SquareMatrix( const unsigned n, const T& value): n_(n), data_(n_*n_, value) {}
58 template< class InputIterator>
59 SquareMatrix( InputIterator first, InputIterator last, std::enable_if_t<!std::is_integral<InputIterator>::value>* = 0): data_(first, last)
60 {
61 unsigned n = std::distance( first, last);
62 n_ = (unsigned)sqrt( (double)n);
63 if( n_*n_!=n) throw Error( Message(_ping_)<<"Too few elements "<<n<<" need "<<n_*n_<<"\n");
64 }
70 SquareMatrix( const std::vector<T>& src): data_(src)
71 {
72 unsigned n = src.size();
73 n_ = (unsigned)sqrt( (double)n);
74 if( n_*n_!=n) throw Error( Message(_ping_)<<"Wrong number of elements "<<n<<" need "<<n_*n_<<"\n");
75 }
76
80 void zero() {
81 for( unsigned i=0; i<n_*n_; i++)
82 data_[i] = 0;
83 }
84
91 T& operator()(const size_t i, const size_t j)
92 {
93 return data_[ i*n_+j];
94 }
101 const T& operator()(const size_t i, const size_t j) const {
102 return data_[ i*n_+j];
103 }
104
110 unsigned size() const { return n_;}
117 void resize( unsigned m, T val = T()) {
118 n_ = m;
119 data_.resize( m*m, val);
120 }
121
127 const std::vector<value_type>& data() const {return data_;}
134 std::vector<value_type>& data() {return data_;}
135
142 void swap_lines( const size_t i, const size_t k)
143 {
144 if(!( i< n_ && k<n_)) throw Error( Message(_ping_) << "Out of range "<<i<<" "<<k<<" range is "<<n_<<"\n");
145 for( size_t j = 0; j<n_; j++)
146 {
147 std::swap( data_[i*n_+j], data_[k*n_+j]);
148 }
149 }
150
157 {
158 SquareMatrix o(*this);
159 for( unsigned i=0; i<n_; i++)
160 for( unsigned j=0; j<i; j++)
161 {
162 std::swap( o.data_[i*n_+j], o.data_[j*n_+i]);
163 }
164 return o;
165 }
166
174 template<class ContainerType1, class ContainerType2>
175 void symv( const ContainerType1& x, ContainerType2& y) const
176 {
177 for( unsigned j=0; j<n_; j++)
178 {
179 y[j] = 0;
180 for( unsigned k=0; k<n_; k++)
181 y[j] += data_[j*n_+k]*x[k];
182 }
183 }
193 template<class value_type1, class ContainerType1, class value_type2, class ContainerType2>
194 void symv( value_type1 alpha, const ContainerType1& x, value_type2 beta, ContainerType2& y) const
195 {
196 for( unsigned j=0; j<n_; j++)
197 {
198 y[j] *= beta;
199 for( unsigned k=0; k<n_; k++)
200 y[j] += alpha*data_[j*n_+k]*x[k];
201 }
202 }
203
209 bool operator!=( const SquareMatrix& rhs) const{
210 for( size_t i = 0; i < n_*n_; i++)
211 if( data_[i] != rhs.data_[i])
212 return true;
213 return false;
214 }
215
221 bool operator==( const SquareMatrix& rhs) const {return !((*this != rhs));}
222
229 {
230 SquareMatrix temp(n_, 0.);
231 for( unsigned i=0; i<n_*n_; i++)
232 temp.data_[i] = -data_[i];
233 return temp;
234 }
243 {
244 for( unsigned i=0; i<n_*n_; i++)
245 data_[i] += op.data_[i];
246 return *this;
247 }
256 {
257 for( unsigned i=0; i<n_*n_; i++)
258 data_[i] -= op.data_[i];
259 return *this;
260 }
268 SquareMatrix& operator*=( const T& value )
269 {
270 for( unsigned i=0; i<n_*n_; i++)
271 data_[i] *= value;
272 return *this;
273 }
282 friend SquareMatrix operator+( const SquareMatrix& lhs, const SquareMatrix& rhs)
283 {
284 SquareMatrix temp(lhs);
285 temp+=rhs;
286 return temp;
287 }
296 friend SquareMatrix operator-( const SquareMatrix& lhs, const SquareMatrix& rhs)
297 {
298 SquareMatrix temp(lhs);
299 temp-=rhs;
300 return temp;
301 }
310 friend SquareMatrix operator*( const T& value, const SquareMatrix& rhs )
311 {
312 SquareMatrix temp(rhs);
313 temp*=value;
314 return temp;
315 }
316
325 friend SquareMatrix operator*( const SquareMatrix& lhs, const T& value)
326 {
327 return value*lhs;
328 }
329
338 friend SquareMatrix operator*( const SquareMatrix& lhs, const SquareMatrix& rhs)
339 {
340 unsigned n_ = lhs.n_;
341 SquareMatrix temp(n_, 0.);
342 for( unsigned i=0; i< n_; i++)
343 for( unsigned j=0; j<n_; j++)
344 {
345 temp.data_[i*n_+j] = (T)0;
346 for( unsigned k=0; k<n_; k++)
347 temp.data_[i*n_+j] += lhs.data_[i*n_+k]*rhs.data_[k*n_+j];
348 }
349 return temp;
350 }
351
361 template<class ContainerType>
362 friend ContainerType operator*( const SquareMatrix& S, const ContainerType& x)
363 {
364 ContainerType out(x);
365 S.symv( x, out);
366 return out;
367 }
368
376 template< class Ostream>
377 friend Ostream& operator<<(Ostream& os, const SquareMatrix& mat)
378 {
379 unsigned n_ = mat.n_;
380 for( size_t i=0; i < n_ ; i++)
381 {
382 for( size_t j = 0;j < n_; j++)
383 os << mat(i,j) << " ";
384 os << "\n";
385 }
386 return os;
387 }
388
398 template< class Istream>
399 friend Istream& operator>> ( Istream& is, SquareMatrix<T>& mat){
400 unsigned n_ = mat.n_;
401 for( size_t i=0; i<n_; i++)
402 for( size_t j=0; j<n_; j++)
403 is >> mat(i, j);
404 return is;
405 }
406
407 private:
408 unsigned n_;
409 std::vector<T> data_;
410};
411
413template<class T>
420
421
422namespace create
423{
426
440template< class T>
441T lu_pivot( dg::SquareMatrix<T>& m, std::vector<unsigned>& p)
442{
443 //from numerical recipes
444 T pivot, determinant=(T)1;
445 unsigned pivotzeile, numberOfSwaps=0;
446 const size_t n = m.size();
447 p.resize( n);
448 for( size_t j = 0; j < n; j++) //gehe Spalten /Diagonale durch
449 {
450 //compute upper matrix except for the diagonal element (the pivot)
451 for( size_t i = 0; i< j; i++)
452 {
453 thrust::host_vector<T> mik(i), mkj(i);
454 for( size_t k=0; k<i; k++)
455 mik[k] = m(i,k), mkj[k] = m(k,j);
456 m(i,j) -= dg::blas1::dot( mik, mkj);
457 }
458 //compute canditates for pivot elements
459 for( size_t i = j; i< n; i++)
460 {
461 thrust::host_vector<T> mik(j), mkj(j);
462 for( size_t k=0; k<j; k++)
463 mik[k] = m(i,k), mkj[k] = m(k,j);
464 m(i,j) -= dg::blas1::dot( mik, mkj);
465 }
466 //search for absolute maximum of pivot candidates
467 pivot = m(j,j);
468 pivotzeile = j;
469 for( size_t i = j+1; i < n; i++)
470 if( fabs( m(i,j)) > fabs(pivot))
471 {
472 pivot = m(i,j), pivotzeile = i;
473 }
474
475 if( fabs(pivot) > 1e-15 )
476 {
477 if( pivotzeile != j)
478 {
479 m.swap_lines( pivotzeile, j);
480 numberOfSwaps++;
481 }
482 p[j] = pivotzeile;
483 //divide all elements below the diagonal by the pivot to get the lower matrix
484 for( size_t i=j+1; i<n; i++)
485 m(i,j) /= pivot;
486 determinant*=m(j,j);
487
488 }
489 else
490 throw std::runtime_error( "Matrix is singular!!");
491 }
492 if( numberOfSwaps % 2 != 0)
493 determinant*=-1.;
494 return determinant;
495
496}
497
513template<class T>
515{
516 dg::SquareMatrix<T> out(in);
517 const unsigned n = in.size();
518 std::vector<unsigned> pivot( n);
519 dg::SquareMatrix<T> lu(in);
520 T determinant = lu_pivot( lu, pivot);
521 if( fabs(determinant ) == 0)
522 throw std::runtime_error( "Determinant zero!");
523 for( unsigned i=0; i<n; i++)
524 {
525 std::vector<T> unit(n, 0);
526 unit[i] = 1;
527 lu_solve( lu, pivot, unit);
528 for( unsigned j=0; j<n; j++)
529 out(j,i) = unit[j];
530 }
531 return out;
532}
533
535
543template<class real_type>
544SquareMatrix<real_type> delta( unsigned n)
545{
547 for( unsigned i=0; i<n; i++)
548 op( i,i) = 1.;
549 return op;
550}
551
552
553// Not sure why we called this dg::create::lu_solve instead of dg::lu_solve
554template<class T>
555void lu_solve( const dg::SquareMatrix<T>& lu, const std::vector<unsigned>& p, std::vector<T>& b)
556{
557 assert(p.size() == lu.size() && p.size() == b.size());
558 const size_t n = p.size();
559 // Vorwärtseinsetzen
560 for( size_t i = 0; i<n; i++)
561 {
562 //mache Zeilentausch
563 std::swap( b[ p[i] ], b[i]);
564 thrust::host_vector<T> lui(i), bi(i);
565 for( size_t j = 0; j < i; j++)
566 lui[j] = lu(i,j), bi[j] = b[j];
567 b[i] -= dg::blas1::dot( lui, bi);
568 }
569 // Rückwärtseinsetzen
570 for( int i = n-1; i>=0; i--)
571 {
572 thrust::host_vector<T> lui(n-(i+1)), bi(n-(i+1));
573 for( size_t j = i+1; j < n; j++)
574 lui[j-(i+1)] = lu(i,j), bi[j-(i+1)] = b[j];
575 b[i] -= dg::blas1::dot( lui, bi);
576 b[i] /= lu(i,i);
577 }
578}
579
580
581// S-matrix
582template<class real_type>
583SquareMatrix<real_type> pipj( unsigned n)
584{
585 SquareMatrix<real_type> op(n, 0);
586 for( unsigned i=0; i<n; i++)
587 op( i,i) = 2./(real_type)(2*i+1);
588 return op;
589}
590// T-matrix
591template<class real_type>
592SquareMatrix<real_type> pipj_inv( unsigned n)
593{
594 SquareMatrix<real_type> op(n, 0);
595 for( unsigned i=0; i<n; i++)
596 op( i,i) = (real_type)(2*i+1)/2.;
597 return op;
598}
599// D-matrix
600template<class real_type>
601SquareMatrix<real_type> pidxpj( unsigned n)
602{
603 SquareMatrix<real_type> op(n, 0);
604 for( unsigned i=0; i<n; i++)
605 for( unsigned j=0; j<n; j++)
606 {
607 if( i < j)
608 {
609 if( (i+j)%2 != 0)
610 op( i, j) = 2;
611 }
612 }
613 return op;
614}
615// R-matrix
616template<class real_type>
617SquareMatrix<real_type> rirj( unsigned n)
618{
619 return SquareMatrix<real_type>( n, 1.);
620}
621// RL-matrix
622template<class real_type>
623SquareMatrix<real_type> rilj( unsigned n)
624{
625 SquareMatrix<real_type> op( n, -1.);
626 for( unsigned i=0; i<n; i++)
627 for( unsigned j=0; j<n; j++)
628 if( j%2 == 0)
629 op( i,j) = 1.;
630 return op;
631}
632// LR-matrix
633template<class real_type>
634SquareMatrix<real_type> lirj( unsigned n) {
635 return rilj<real_type>( n).transpose();
636}
637// L-matrix
638template<class real_type>
639SquareMatrix<real_type> lilj( unsigned n)
640{
641 SquareMatrix<real_type> op( n, -1.);
642 for( unsigned i=0; i<n; i++)
643 for( unsigned j=0; j<n; j++)
644 if( ((i+j)%2) == 0)
645 op( i,j) = 1.;
646 return op;
647}
648
649// N-matrix
650template<class real_type>
651SquareMatrix<real_type> ninj( unsigned n)
652{
653 SquareMatrix<real_type> op( n, 0.);
654 for( int i=0; i<(int)n; i++)
655 for( int j=0; j<(int)n; j++)
656 {
657 if( i == j+1)
658 op( i,j) = 2./(2*i+1)/(2*j+1);
659 if( i == j-1)
660 op( i,j) = -2./(2*i+1)/(2*j+1);
661 }
662 op(0,0) = 2;
663 return op;
664}
666
667
669}//namespace create
670
681template<class T>
682void lu_solve( const dg::SquareMatrix<T>& lu, const std::vector<unsigned>& p, std::vector<T>& b)
683{
684 dg::create::lu_solve( lu, p, b);
685}
686
690template<class T>
698template<class T>
700
701} //namespace dg
702
703#endif //_DG_OPERATORS_DYN_
class intended for the use in throw statements
Definition exceptions.h:83
small class holding a stringstream
Definition exceptions.h:29
A square nxn matrix.
Definition operator.h:31
SquareMatrix & operator-=(const SquareMatrix &op)
subtract
Definition operator.h:255
unsigned size() const
Size n of the SquareMatrix.
Definition operator.h:110
SquareMatrix & operator*=(const T &value)
scalar multiply
Definition operator.h:268
const std::vector< value_type > & data() const
access underlying data
Definition operator.h:127
SquareMatrix(const unsigned n)
allocate storage for nxn matrix
Definition operator.h:43
SquareMatrix(const unsigned n, const T &value)
Initialize elements.
Definition operator.h:50
void symv(value_type1 alpha, const ContainerType1 &x, value_type2 beta, ContainerType2 &y) const
Matrix vector multiplication .
Definition operator.h:194
const T & operator()(const size_t i, const size_t j) const
const access operator
Definition operator.h:101
void resize(unsigned m, T val=T())
Resize.
Definition operator.h:117
friend SquareMatrix operator*(const SquareMatrix &lhs, const T &value)
scalar multiplication
Definition operator.h:325
friend SquareMatrix operator*(const SquareMatrix &lhs, const SquareMatrix &rhs)
matrix multiplication
Definition operator.h:338
bool operator==(const SquareMatrix &rhs) const
two Matrices are considered equal if elements are equal
Definition operator.h:221
SquareMatrix & operator+=(const SquareMatrix &op)
add
Definition operator.h:242
void swap_lines(const size_t i, const size_t k)
Swap two lines in the square matrix.
Definition operator.h:142
void zero()
Assign zero to all elements.
Definition operator.h:80
T value_type
typically double or float
Definition operator.h:33
void symv(const ContainerType1 &x, ContainerType2 &y) const
Matrix vector multiplication .
Definition operator.h:175
friend ContainerType operator*(const SquareMatrix &S, const ContainerType &x)
matrix-vector multiplication
Definition operator.h:362
friend SquareMatrix operator*(const T &value, const SquareMatrix &rhs)
scalar multiplication
Definition operator.h:310
SquareMatrix(const std::vector< T > &src)
Copy from existing data.
Definition operator.h:70
SquareMatrix operator-() const
subtract
Definition operator.h:228
friend SquareMatrix operator+(const SquareMatrix &lhs, const SquareMatrix &rhs)
add
Definition operator.h:282
friend SquareMatrix operator-(const SquareMatrix &lhs, const SquareMatrix &rhs)
subtract
Definition operator.h:296
bool operator!=(const SquareMatrix &rhs) const
two Matrices are considered equal if elements are equal
Definition operator.h:209
friend Istream & operator>>(Istream &is, SquareMatrix< T > &mat)
Read values into a Matrix from given istream.
Definition operator.h:399
T & operator()(const size_t i, const size_t j)
access operator
Definition operator.h:91
friend Ostream & operator<<(Ostream &os, const SquareMatrix &mat)
puts a matrix linewise in output stream
Definition operator.h:377
SquareMatrix()=default
Construct empty SquareMatrix.
std::vector< value_type > & data()
Write access to underlying data.
Definition operator.h:134
SquareMatrix(InputIterator first, InputIterator last, std::enable_if_t<!std::is_integral< InputIterator >::value > *=0)
Construct from iterators.
Definition operator.h:59
SquareMatrix transpose() const
Transposition.
Definition operator.h:156
The discrete legendre trafo class.
#define _ping_
Definition exceptions.h:12
auto dot(const ContainerType1 &x, const ContainerType2 &y)
Binary reproducible Euclidean dot product between two vectors
Definition blas1.h:153
@ y
y direction
@ x
x direction
dg::SquareMatrix< T > inverse(const dg::SquareMatrix< T > &in)
Invert a square matrix.
Definition operator.h:514
void lu_solve(const dg::SquareMatrix< T > &lu, const std::vector< unsigned > &p, std::vector< T > &b)
Solve the linear system with the LU decomposition.
Definition operator.h:682
T lu_pivot(dg::SquareMatrix< T > &m, std::vector< unsigned > &p)
LU Decomposition with partial pivoting.
Definition operator.h:441
dg::SquareMatrix< T > invert(const dg::SquareMatrix< T > &in)
Compute inverse of square matrix (alias for dg::create::inverse)
Definition operator.h:691
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...
Indicates that the type has a member function with the same name and interface (up to the matrix itse...
Definition matrix_categories.h:26
Indicate sequential execution.
Definition execution_policy.h:26
T value_type
Definition operator.h:416
The vector traits.
Definition tensor_traits.h:38