1#ifndef _DG_OPERATORS_DYN_
2#define _DG_OPERATORS_DYN_
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)
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");
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");
81 for(
unsigned i=0; i<n_*n_; i++)
93 return data_[ i*n_+j];
102 return data_[ i*n_+j];
110 unsigned size()
const {
return n_;}
119 data_.resize( m*m, val);
127 const std::vector<value_type>&
data()
const {
return data_;}
134 std::vector<value_type>&
data() {
return data_;}
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++)
147 std::swap( data_[i*n_+j], data_[k*n_+j]);
159 for(
unsigned i=0; i<n_; i++)
160 for(
unsigned j=0; j<i; j++)
162 std::swap( o.data_[i*n_+j], o.data_[j*n_+i]);
174 template<
class ContainerType1,
class ContainerType2>
175 void symv(
const ContainerType1& x, ContainerType2&
y)
const
177 for(
unsigned j=0; j<n_; j++)
180 for(
unsigned k=0; k<n_; k++)
181 y[j] += data_[j*n_+k]*
x[k];
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
196 for(
unsigned j=0; j<n_; j++)
199 for(
unsigned k=0; k<n_; k++)
200 y[j] += alpha*data_[j*n_+k]*
x[k];
210 for(
size_t i = 0; i < n_*n_; i++)
211 if( data_[i] != rhs.data_[i])
231 for(
unsigned i=0; i<n_*n_; i++)
232 temp.data_[i] = -data_[i];
244 for(
unsigned i=0; i<n_*n_; i++)
245 data_[i] += op.data_[i];
257 for(
unsigned i=0; i<n_*n_; i++)
258 data_[i] -= op.data_[i];
270 for(
unsigned i=0; i<n_*n_; i++)
340 unsigned n_ = lhs.n_;
342 for(
unsigned i=0; i< n_; i++)
343 for(
unsigned j=0; j<n_; j++)
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];
361 template<
class ContainerType>
364 ContainerType out(
x);
376 template<
class Ostream>
379 unsigned n_ = mat.n_;
380 for(
size_t i=0; i < n_ ; i++)
382 for(
size_t j = 0;j < n_; j++)
383 os << mat(i,j) <<
" ";
398 template<
class Istream>
400 unsigned n_ = mat.n_;
401 for(
size_t i=0; i<n_; i++)
402 for(
size_t j=0; j<n_; j++)
409 std::vector<T> data_;
444 T pivot, determinant=(T)1;
445 unsigned pivotzeile, numberOfSwaps=0;
446 const size_t n = m.
size();
448 for(
size_t j = 0; j < n; j++)
451 for(
size_t i = 0; i< j; i++)
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);
459 for(
size_t i = j; i< n; i++)
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);
469 for(
size_t i = j+1; i < n; i++)
470 if( fabs( m(i,j)) > fabs(pivot))
472 pivot = m(i,j), pivotzeile = i;
475 if( fabs(pivot) > 1e-15 )
484 for(
size_t i=j+1; i<n; i++)
490 throw std::runtime_error(
"Matrix is singular!!");
492 if( numberOfSwaps % 2 != 0)
517 const unsigned n = in.
size();
518 std::vector<unsigned> pivot( n);
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++)
525 std::vector<T> unit(n, 0);
528 for(
unsigned j=0; j<n; j++)
543template<
class real_type>
547 for(
unsigned i=0; i<n; i++)
557 assert(p.size() == lu.
size() && p.size() == b.size());
558 const size_t n = p.size();
560 for(
size_t i = 0; i<n; i++)
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];
570 for(
int i = n-1; i>=0; i--)
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];
582template<
class real_type>
583SquareMatrix<real_type> pipj(
unsigned n)
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);
591template<
class real_type>
592SquareMatrix<real_type> pipj_inv(
unsigned n)
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.;
600template<
class real_type>
601SquareMatrix<real_type> pidxpj(
unsigned n)
603 SquareMatrix<real_type> op(n, 0);
604 for(
unsigned i=0; i<n; i++)
605 for(
unsigned j=0; j<n; j++)
616template<
class real_type>
617SquareMatrix<real_type> rirj(
unsigned n)
619 return SquareMatrix<real_type>( n, 1.);
622template<
class real_type>
623SquareMatrix<real_type> rilj(
unsigned n)
625 SquareMatrix<real_type> op( n, -1.);
626 for(
unsigned i=0; i<n; i++)
627 for(
unsigned j=0; j<n; j++)
633template<
class real_type>
634SquareMatrix<real_type> lirj(
unsigned n) {
635 return rilj<real_type>( n).transpose();
638template<
class real_type>
639SquareMatrix<real_type> lilj(
unsigned n)
641 SquareMatrix<real_type> op( n, -1.);
642 for(
unsigned i=0; i<n; i++)
643 for(
unsigned j=0; j<n; j++)
650template<
class real_type>
651SquareMatrix<real_type> ninj(
unsigned n)
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++)
658 op( i,j) = 2./(2*i+1)/(2*j+1);
660 op( i,j) = -2./(2*i+1)/(2*j+1);
684 dg::create::lu_solve( lu, p, b);
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
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