35template<
class real_type>
36void shift(
bool& negative, real_type & x,
dg::bc bc, real_type x0, real_type x1)
41 real_type N0 = floor((x-x0)/(x1-x0));
45 while( (x<x0) || (x>x1) )
69template<
class real_type>
70std::vector<real_type> legendre( real_type xn,
unsigned n)
72 assert( xn <= 1. && xn >= -1.);
73 std::vector<real_type> px(n);
76 for(
unsigned u=0; u<n; u++)
77 px[u] = (u%2 == 0) ? +1. : -1.;
81 for(
unsigned i=0; i<n; i++)
90 for(
unsigned i=1; i<n-1; i++)
91 px[i+1] = ((real_type)(2*i+1)*xn*px[i]-(real_type)i*px[i-1])/(real_type)(i+1);
98template<
class real_type>
99std::vector<real_type> lagrange( real_type x,
const std::vector<real_type>& xi)
101 unsigned n = xi.size();
102 std::vector<real_type> l( n , 1.);
103 for(
unsigned i=0; i<n; i++)
104 for(
unsigned k=0; k<n; k++)
107 l[i] *= (
x-xi[k])/(xi[i]-xi[k]);
113template<
class real_type>
114std::vector<real_type> choose_1d_abscissas( real_type X,
115 unsigned points_per_line,
double lx,
116 const thrust::host_vector<real_type>& abs,
dg::bc bcx,
117 thrust::host_vector<unsigned>& cols)
121 assert( abs.size() >= points_per_line &&
"There must be more points to interpolate\n");
124 auto it = std::lower_bound( abs.begin(), abs.end(), X);
126 std::vector<real_type> xs( points_per_line, 0);
127 cols.resize( points_per_line, 0);
128 switch( points_per_line)
131 if( it == abs.begin())
133 else if( it == abs.end())
134 cols[0] = it - abs.begin() - 1;
137 if ( fabs(X - *it) < fabs( X - *(it-1)))
138 cols[0] = it - abs.begin();
140 cols[0] = it - abs.begin() -1;
143 case 2:
if( it == abs.begin())
148 xs[1] = *(abs.end() -1)-lx;;
149 cols[0] = 0, cols[1] = abs.end()-abs.begin()-1;
162 else if( it == abs.end())
166 xs[0] = *(abs.begin())+lx;
168 cols[0] = 0, cols[1] = it-abs.begin()-1;
179 cols[0] = it-abs.begin()-1;
186 cols[0] = it - abs.begin() - 1;
190 case 4:
if( it <= abs.begin() +1)
194 xs[0] = *abs.begin(), cols[0] = 0;
195 xs[1] = *(abs.begin()+1), cols[1] = 1;
196 xs[2] = it == abs.begin() ? *(abs.end() -2) : *(abs.begin()+2);
197 cols[2] = it == abs.begin() ? abs.end()-abs.begin() -2 : 2;
198 xs[3] = *(abs.end() -1);
199 cols[3] = abs.end()-abs.begin() -1;
204 xs[0] = *it, xs[1] = *(it+1);
205 xs[2] = *(it+2), xs[3] = *(it+3);
206 cols[0] = 0, cols[1] = 1;
207 cols[2] = 2, cols[3] = 3;
210 else if( it >= abs.end() -1)
214 xs[0] = *abs.begin(), cols[0] = 0;
215 xs[1] = it == abs.end() ? *(abs.begin()+1) : *(abs.end() -3) ;
216 cols[1] = it == abs.end() ? 1 : abs.end()-abs.begin()-3 ;
217 xs[2] = *(abs.end() - 2), cols[2] = abs.end()-abs.begin()-2;
218 xs[3] = *(abs.end() - 1), cols[3] = abs.end()-abs.begin()-1;
223 xs[0] = *(it-4), xs[1] = *(it-3);
224 xs[2] = *(it-2), xs[3] = *(it-1);
225 cols[0] = it - abs.begin() - 4;
233 xs[0] = *(it-2), xs[1] = *(it-1);
234 xs[2] = *(it ), xs[3] = *(it+1);
235 cols[0] = it - abs.begin() - 2;
245template<
class real_type>
246void interpolation_row(
dg::space sp, real_type X,
247 const RealGrid1d<real_type>& g,
dg::bc bcx,
248 thrust::host_vector<int>& cols,
249 thrust::host_vector<real_type>& vals)
251 bool negative =
false;
252 detail::shift( negative, X, bcx, g.x0(), g.x1());
255 real_type xnn = (X-g.x0())/g.h();
256 unsigned nn = (unsigned)floor(xnn);
258 real_type xn = 2.*xnn - (real_type)(2*nn+1);
267 for(
unsigned k=0; k<g.nx(); k++)
271 if( fabs( xn - gauss_nodesx[k]) < 1e-13)
272 idxX = nn*g.nx() + k;
277 std::vector<real_type> px = detail::legendre( xn, g.n());
281 std::vector<real_type> pxF(px.size(), 0);
282 std::vector<real_type>
forward =
284 for(
unsigned l=0; l<g.n(); l++)
285 for(
unsigned k=0; k<g.n(); k++)
286 pxF[l]+= px[k]*forward[k*g.n() + l];
289 for (
unsigned l=0; l<g.n(); l++)
291 cols.push_back( nn*g.n() + l);
292 vals.push_back( negative ? -px[l] : px[l]);
297 cols.push_back( idxX);
298 vals.push_back( negative ? -1. : 1.);
304template<
class host_vector,
class real_type>
307 const host_vector& x,
308 const RealGrid1d<real_type>& g,
311 thrust::host_vector<real_type> values;
312 thrust::host_vector<int> row_offsets;
313 thrust::host_vector<int> column_indices;
314 auto ptr =
x.begin();
315 row_offsets.push_back(0);
316 for(
unsigned i=0; i<
x.size(); i++)
318 interpolation_row( sp, *ptr, g, bcx, column_indices, values);
319 row_offsets.push_back( values.size());
322 return {
x.size(), g.size(), row_offsets, column_indices, values};
326template<
class host_vector1,
class host_vector2 >
328 const host_vector1& x,
329 const host_vector2& abs,
336 thrust::host_vector<real_type> values;
337 thrust::host_vector<int> row_offsets;
338 thrust::host_vector<int> column_indices;
339 unsigned points_per_line = 1;
340 if( method ==
"nearest")
342 else if( method ==
"linear")
344 else if( method ==
"cubic")
347 throw std::runtime_error(
"Interpolation method "+method+
" not recognized!\n");
348 auto ptr =
x.begin();
349 row_offsets.push_back(0);
350 for(
unsigned i=0; i<
x.size(); i++)
352 row_offsets.push_back(row_offsets[i]);
355 bool negative =
false;
356 detail::shift( negative, X, bcx, x0, x1);
359 auto it = std::lower_bound( abs.begin(), abs.end(), X);
360 if( fabs( X - *it) < 1e-13)
361 idxX = it - abs.begin();
362 if( it != abs.begin() and fabs( X - *(it-1)) < 1e-13)
363 idxX = (it - abs.begin())-1;
370 thrust::host_vector<unsigned> cols;
371 std::vector<real_type> xs = detail::choose_1d_abscissas( X,
372 points_per_line, x1-x0, abs, bcx, cols);
374 std::vector<real_type> px = detail::lagrange( X, xs);
376 for (
unsigned l=0; l<px.size(); l++)
379 column_indices.push_back( cols[l]);
380 values.push_back(negative ? -px[l] : px[l]);
386 column_indices.push_back(idxX);
387 values.push_back( negative ? -1. : 1.);
390 return {
x.size(), abs.size(), row_offsets, column_indices, values};
432template<
class RecursiveHostVector,
class real_type,
size_t Nd>
434 const RecursiveHostVector& x,
436 std::array<dg::bc, Nd> bcx,
437 std::string method =
"dg")
440 std::array<dg::SparseMatrix<int,real_type,thrust::host_vector>,Nd> axes;
441 for(
unsigned u=0; u<Nd; u++)
443 if(
x[u].size() !=
x[0].size())
447 axes[u] = detail::interpolation1d(
dg::xspace,
x[u], g.
grid(u), bcx[u]);
451 axes[u] = detail::interpolation1d(
x[u], g.
abscissas(u), bcx[u], g.
p(u), g.
q(u), method);
454 for(
unsigned u=1; u<Nd; u++)
479template<
class host_vector,
class real_type,
typename = std::enable_if_t<dg::is_vector_v<host_vector>>>
481 const host_vector& x,
484 std::string method =
"dg")
489 std::array<dg::bc,1>{bcx}, method);
513template<
class host_vector,
class real_type>
515 const host_vector& x,
516 const host_vector&
y,
519 std::string method =
"dg")
523 return interpolation( std::vector{vx,vy}, g, {bcx,bcy}, method);
551template<
class host_vector,
class real_type>
553 const host_vector& x,
554 const host_vector&
y,
555 const host_vector&
z,
558 std::string method =
"dg")
563 return interpolation( std::vector{vx,vy,vz}, g, {bcx,bcy,bcz}, method);
585template<
class real_type,
size_t Nd>
591 for(
unsigned u=0; u<Nd; u++)
593 if( g_new.
p(u) < g_old.
p(u))
594 std::cerr <<
"ERROR: New grid boundary number "<<u<<
" with value "<<g_new.
p(u)<<
" lies outside old grid "<<g_old.
p(u)<<
" "<<g_old.
p(u)-g_new.
p(u)<<
"\n";
595 assert( g_new.
p(u) >= g_old.
p(u));
596 if( g_new.
q(u) > g_old.
q(u))
597 std::cerr <<
"ERROR: New grid boundary number "<<u<<
" with value "<<g_new.
q(u)<<
" lies outside old grid "<<g_old.
q(u)<<
" "<<g_old.
q(u)-g_new.
q(u)<<
"\n";
598 assert( g_new.
q(u) <= g_old.
q(u));
600 std::array<dg::SparseMatrix<int,real_type,thrust::host_vector>,Nd> axes;
601 for(
unsigned u=0; u<Nd; u++)
606 g_old.
grid(u), g_old.
bc(u));
610 axes[u] = detail::interpolation1d( g_new.
abscissas(u),
615 for(
unsigned u=1; u<Nd; u++)
643template<
class host_vector,
class real_type>
646 const host_vector& v,
651 assert( v.size() == g.
size());
652 thrust::host_vector<real_type> vals;
653 thrust::host_vector<int> cols;
654 create::detail::interpolation_row( sp,
x, g, bcx, cols, vals);
657 for(
unsigned j=0; j<vals.size(); j++)
658 value += v[cols[j]]*vals[j];
681template<
class host_vector,
class real_type>
684 const host_vector& v,
685 real_type x, real_type
y,
689 assert( v.size() == g.
size());
690 thrust::host_vector<real_type> valsx;
691 thrust::host_vector<int> colsx;
692 create::detail::interpolation_row( sp,
x, g.
gx(), bcx, colsx, valsx);
693 thrust::host_vector<real_type> valsy;
694 thrust::host_vector<int> colsy;
695 create::detail::interpolation_row( sp,
y, g.
gy(), bcy, colsy, valsy);
698 for(
unsigned i=0; i<valsy.size(); i++)
699 for(
unsigned j=0; j<valsx.size(); j++)
700 value += v[colsy[i]*g.
shape(0) + colsx[j]]*valsx[j]*valsy[i];
class intended for the use in throw statements
Definition exceptions.h:83
small class holding a stringstream
Definition exceptions.h:29
Function discretization routines.
#define _ping_
Definition exceptions.h:12
Some utility functions for the dg::evaluate routines.
bc
Switch between boundary conditions.
Definition enums.h:15
space
Space of DG coefficients.
Definition enums.h:164
@ NEU_DIR
Neumann on left, Dirichlet on right boundary.
Definition enums.h:19
@ PER
periodic boundaries
Definition enums.h:16
@ NEU
Neumann on both boundaries.
Definition enums.h:20
@ DIR
homogeneous dirichlet boundaries
Definition enums.h:17
@ DIR_NEU
Dirichlet on left, Neumann on right boundary.
Definition enums.h:18
@ xspace
Configuration space "nodal values".
Definition enums.h:166
@ lspace
DG Polynomial space "modal values".
Definition enums.h:165
@ forward
forward derivative (cell to the right and current cell)
Definition enums.h:98
typename TensorTraits< std::decay_t< Vector > >::value_type get_value_type
Definition tensor_traits.h:45
real_type interpolate(dg::space sp, const host_vector &v, real_type x, const RealGrid1d< real_type > &g, dg::bc bcx=dg::NEU)
Interpolate a vector on a single point on a 1d Grid.
Definition interpolation.h:644
dg::SparseMatrix< int, real_type, thrust::host_vector > interpolation(const RecursiveHostVector &x, const aRealTopology< real_type, Nd > &g, std::array< dg::bc, Nd > bcx, std::string method="dg")
Create interpolation matrix of a list of points in given grid.
Definition interpolation.h:433
dg::SparseMatrix< int, T, thrust::host_vector > tensorproduct_cols(const dg::SparseMatrix< int, T, thrust::host_vector > &lhs, const dg::SparseMatrix< int, T, thrust::host_vector > &rhs)
Form the tensor (Kronecker) product between two matrices in the column index
Definition xspacelib.h:89
SquareMatrix< T > tensorproduct(const SquareMatrix< T > &op1, const SquareMatrix< T > &op2)
Form the tensor product between two operators.
Definition operator_tensor.h:25
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...
static std::vector< real_type > abscissas(unsigned n)
Return Gauss-Legendre nodes on the interval [-1,1].
Definition dlt.h:27
static std::vector< real_type > forward(unsigned n)
Return forward DLT trafo matrix.
Definition dlt.h:195
The simplest implementation of aRealTopology.
Definition grid.h:710
A CSR formatted sparse matrix.
Definition sparsematrix.h:96
A vector view class, usable in dg functions.
Definition view.h:45
An abstract base class for Nd-dimensional dG grids.
Definition grid.h:95
RealGrid< real_type, 1 > gy() const
Equivalent to grid(1)
Definition grid.h:360
real_type q(unsigned u=0) const
Get right boundary point for axis u.
Definition grid.h:253
RealGrid< real_type, 1 > gx() const
Equivalent to grid(0)
Definition grid.h:354
unsigned size() const
The total number of points.
Definition grid.h:532
host_vector abscissas(unsigned u=0) const
Construct grid abscissas of the u axis.
Definition grid.h:128
RealGrid< real_type, 1 > grid(unsigned u) const
Get axis u as a 1d grid.
Definition grid.h:274
real_type p(unsigned u=0) const
Get left boundary point for axis u.
Definition grid.h:250
dg::bc bc(unsigned u=0) const
Get boundary condition for axis u.
Definition grid.h:268
unsigned shape(unsigned u=0) const
the total number of points of an axis
Definition grid.h:114
Useful typedefs of commonly used types.