1#ifndef _DG_OPERATORS_DYN_
2#define _DG_OPERATORS_DYN_
40 explicit Operator(
const unsigned n): n_(n), data_(n_*n_){}
47 Operator(
const unsigned n,
const T& value): n_(n), data_(n_*n_, value) {}
55 template<
class InputIterator>
56 Operator( InputIterator first, InputIterator last, std::enable_if_t<!std::is_integral<InputIterator>::value>* = 0): data_(first, last)
58 unsigned n = std::distance( first, last);
59 n_ = (unsigned)sqrt( (
double)n);
60 if( n_*n_!=n)
throw Error(
Message(
_ping_)<<
"Too few elements "<<n<<
" need "<<n_*n_<<
"\n");
67 Operator(
const std::vector<T>& src): data_(src)
69 unsigned n = src.size();
70 n_ = (unsigned)sqrt( (
double)n);
71 if( n_*n_!=n)
throw Error(
Message(
_ping_)<<
"Wrong number of elements "<<n<<
" need "<<n_*n_<<
"\n");
78 for(
unsigned i=0; i<n_*n_; i++)
91 if(!(i<n_&&j<n_))
throw Error(
Message(
_ping_) <<
"You tried to access out of range "<<i<<
" "<<j<<
" size is "<<n_<<
"\n");
93 return data_[ i*n_+j];
103 if(!(i<n_&&j<n_))
throw Error(
Message(
_ping_) <<
"You tried to access out of range "<<i<<
" "<<j<<
" size is "<<n_<<
"\n");
105 return data_[ i*n_+j];
113 unsigned size()
const {
return n_;}
122 data_.resize( m*m, val);
130 const std::vector<value_type>&
data()
const {
return data_;}
140 if(!( i< n_ && k<n_))
throw Error(
Message(
_ping_) <<
"Out of range "<<i<<
" "<<k<<
" range is "<<n_<<
"\n");
141 for(
size_t j = 0; j<n_; j++)
143 std::swap( data_[i*n_+j], data_[k*n_+j]);
155 for(
unsigned i=0; i<n_; i++)
156 for(
unsigned j=0; j<i; j++)
158 std::swap( o.data_[i*n_+j], o.data_[j*n_+i]);
169 for(
size_t i = 0; i < n_*n_; i++)
170 if( data_[i] != rhs.data_[i])
190 for(
unsigned i=0; i<n_*n_; i++)
191 temp.data_[i] = -data_[i];
204 assert( op.
size() == this->size());
206 for(
unsigned i=0; i<n_*n_; i++)
207 data_[i] += op.data_[i];
220 assert( op.
size() == this->size());
222 for(
unsigned i=0; i<n_*n_; i++)
223 data_[i] -= op.data_[i];
235 for(
unsigned i=0; i<n_*n_; i++)
305 unsigned n_ = lhs.n_;
310 for(
unsigned i=0; i< n_; i++)
311 for(
unsigned j=0; j<n_; j++)
313 temp.data_[i*n_+j] = (T)0;
314 for(
unsigned k=0; k<n_; k++)
315 temp.data_[i*n_+j] += lhs.data_[i*n_+k]*rhs.data_[k*n_+j];
327 template<
class Ostream>
330 unsigned n_ = mat.n_;
331 for(
size_t i=0; i < n_ ; i++)
333 for(
size_t j = 0;j < n_; j++)
334 os << mat(i,j) <<
" ";
349 template<
class Istream>
351 unsigned n_ = mat.n_;
352 for(
size_t i=0; i<n_; i++)
353 for(
size_t j=0; j<n_; j++)
360 std::vector<T> data_;
382 unsigned pivotzeile, numberOfSwaps=0;
383 const size_t n = m.
size();
385 for(
size_t j = 0; j < n; j++)
388 for(
size_t i = 0; i< j; i++)
390 thrust::host_vector<T> mik(i), mkj(i);
391 for(
size_t k=0; k<i; k++)
392 mik[k] = m(i,k), mkj[k] = m(k,j);
396 for(
size_t i = j; i< n; i++)
398 thrust::host_vector<T> mik(j), mkj(j);
399 for(
size_t k=0; k<j; k++)
400 mik[k] = m(i,k), mkj[k] = m(k,j);
406 for(
size_t i = j+1; i < n; i++)
407 if( fabs( m(i,j)) > fabs(pivot))
409 pivot = m(i,j), pivotzeile = i;
412 if( fabs(pivot) > 1e-15 )
421 for(
size_t i=j+1; i<n; i++)
427 throw std::runtime_error(
"Matrix is singular!!");
429 if( numberOfSwaps % 2 != 0)
447 assert(p.size() == lu.
size() && p.size() == b.size());
448 const size_t n = p.size();
450 for(
size_t i = 0; i<n; i++)
453 std::swap( b[ p[i] ], b[i]);
454 thrust::host_vector<T> lui(i), bi(i);
455 for(
size_t j = 0; j < i; j++)
456 lui[j] = lu(i,j), bi[j] = b[j];
460 for(
int i = n-1; i>=0; i--)
462 thrust::host_vector<T> lui(n-(i+1)), bi(n-(i+1));
463 for(
size_t j = i+1; j < n; j++)
464 lui[j-(i+1)] = lu(i,j), bi[j-(i+1)] = b[j];
486 const unsigned n = in.
size();
487 std::vector<unsigned> pivot( n);
491 throw std::runtime_error(
"Determinant zero!");
492 for(
unsigned i=0; i<n; i++)
494 std::vector<T> unit(n, 0);
497 for(
unsigned j=0; j<n; j++)
513template<
class real_type>
517 for(
unsigned i=0; i<n; i++)
528template<
class real_type>
532 for(
unsigned i=0; i<n; i++)
533 op( i,i) = 2./(real_type)(2*i+1);
543template<
class real_type>
547 for(
unsigned i=0; i<n; i++)
548 op( i,i) = (real_type)(2*i+1)/2.;
558template<
class real_type>
562 for(
unsigned i=0; i<n; i++)
563 for(
unsigned j=0; j<n; j++)
580template<
class real_type>
592template<
class real_type>
596 for(
unsigned i=0; i<n; i++)
597 for(
unsigned j=0; j<n; j++)
609template<
class real_type>
620template<
class real_type>
624 for(
unsigned i=0; i<n; i++)
625 for(
unsigned j=0; j<n; j++)
638template<
class real_type>
642 for(
int i=0; i<(int)n; i++)
643 for(
int j=0; j<(int)n; j++)
646 op( i,j) = 2./(2*i+1)/(2*j+1);
648 op( i,j) = -2./(2*i+1)/(2*j+1);
662template<
class real_type>
665 unsigned n = dlt.
weights().size();
667 for(
unsigned i=0; i<n; i++)
678template<
class real_type>
681 unsigned n = dlt.
weights().size();
683 for(
unsigned i=0; i<n; i++)
Struct holding coefficients for Discrete Legendre Transformation (DLT) related operations.
Definition: dlt.h:23
const std::vector< real_type > & weights() const
Return Gauss-Legendre weights.
Definition: dlt.h:42
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:28
void zero()
Assign zero to all elements.
Definition: operator.h:77
T & operator()(const size_t i, const size_t j)
access operator
Definition: operator.h:89
Operator & operator*=(const T &value)
scalar multiply
Definition: operator.h:233
friend Operator operator-(const Operator &lhs, const Operator &rhs)
subtract
Definition: operator.h:261
Operator(const std::vector< T > &src)
Copy from existing data.
Definition: operator.h:67
friend Operator operator*(const Operator &lhs, const Operator &rhs)
matrix multiplication
Definition: operator.h:303
unsigned size() const
Size n of the Operator.
Definition: operator.h:113
T value_type
typically double or float
Definition: operator.h:30
friend Istream & operator>>(Istream &is, Operator< T > &mat)
Read values into a Matrix from given istream.
Definition: operator.h:350
void resize(unsigned m, T val=T())
Resize.
Definition: operator.h:120
Operator & operator-=(const Operator &op)
subtract
Definition: operator.h:217
bool operator==(const Operator &rhs) const
two Matrices are considered equal if elements are equal
Definition: operator.h:180
Operator(const unsigned n)
allocate storage for nxn matrix
Definition: operator.h:40
Operator()=default
Construct empty Operator.
friend Operator operator*(const Operator &lhs, const T &value)
scalar multiplication
Definition: operator.h:290
Operator operator-() const
subtract
Definition: operator.h:187
void swap_lines(const size_t i, const size_t k)
Swap two lines in the square matrix.
Definition: operator.h:138
const T & operator()(const size_t i, const size_t j) const
const access operator
Definition: operator.h:101
friend Operator operator*(const T &value, const Operator &rhs)
scalar multiplication
Definition: operator.h:275
Operator & operator+=(const Operator &op)
add
Definition: operator.h:201
friend Operator operator+(const Operator &lhs, const Operator &rhs)
add
Definition: operator.h:247
bool operator!=(const Operator &rhs) const
two Matrices are considered equal if elements are equal
Definition: operator.h:168
friend Ostream & operator<<(Ostream &os, const Operator &mat)
puts a matrix linewise in output stream
Definition: operator.h:328
const std::vector< value_type > & data() const
access underlying data
Definition: operator.h:130
Operator(InputIterator first, InputIterator last, std::enable_if_t<!std::is_integral< InputIterator >::value > *=0)
Construct from iterators.
Definition: operator.h:56
Operator transpose() const
Transposition.
Definition: operator.h:152
Operator(const unsigned n, const T &value)
Initialize elements.
Definition: operator.h:47
The discrete legendre trafo class.
#define _ping_
Definition: exceptions.h:12
get_value_type< ContainerType1 > dot(const ContainerType1 &x, const ContainerType2 &y)
Binary reproducible Euclidean dot product between two vectors
Definition: blas1.h:87
dg::Operator< T > inverse(const dg::Operator< T > &in)
Compute the inverse of a square matrix.
Definition: operator.h:483
T lu_pivot(dg::Operator< T > &m, std::vector< unsigned > &p)
LU Decomposition with partial pivoting.
Definition: operator.h:378
void lu_solve(const dg::Operator< T > &lu, const std::vector< unsigned > &p, std::vector< T > &b)
Solve the linear system with the LU decomposition.
Definition: operator.h:445
dg::Operator< T > invert(const dg::Operator< T > &in)
Alias for dg::create::inverse. Compute inverse of square matrix.
Definition: operator.h:695
MPI_Vector< thrust::host_vector< real_type > > weights(const aRealMPITopology2d< real_type > &g)
Nodal weight coefficients.
Definition: mpi_weights.h:22
MPI_Vector< thrust::host_vector< real_type > > inv_weights(const aRealMPITopology2d< real_type > &g)
inverse nodal weight coefficients
Definition: mpi_weights.h:29
Operator< real_type > lirj(unsigned n)
Create the LR-matrix.
Definition: operator.h:610
Operator< real_type > pipj(unsigned n)
Create the S-matrix.
Definition: operator.h:529
Operator< real_type > rirj(unsigned n)
Create the R-matrix.
Definition: operator.h:581
Operator< real_type > rilj(unsigned n)
Create the RL-matrix.
Definition: operator.h:593
Operator< real_type > delta(unsigned n)
Create the unit matrix.
Definition: operator.h:514
Operator< real_type > pidxpj(unsigned n)
Create the D-matrix.
Definition: operator.h:559
Operator< real_type > pipj_inv(unsigned n)
Create the T-matrix.
Definition: operator.h:544
Operator< real_type > lilj(unsigned n)
Create the L-matrix.
Definition: operator.h:621
Operator< real_type > ninj(unsigned n)
Create the N-matrix.
Definition: operator.h:639
ContainerType determinant(const SparseTensor< ContainerType > &t)
Definition: multiply.h:349
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...