Discontinuous Galerkin Library
#include "dg/algorithm.h"
Loading...
Searching...
No Matches
interpolation.h
Go to the documentation of this file.
1#pragma once
2//#include <iomanip>
3
5#include "dg/backend/view.h"
6#include "grid.h"
7#include "evaluation.h"
8#include "functions.h"
9#include "operator_tensor.h"
10#include "xspacelib.h"
11
17namespace dg{
18
19namespace create{
21namespace detail{
35template<class real_type>
36void shift( bool& negative, real_type & x, dg::bc bc, real_type x0, real_type x1)
37{
38 // This is now a free function because not all Topologies may have it
39 if( bc == dg::PER)
40 {
41 real_type N0 = floor((x-x0)/(x1-x0)); // ... -2[ -1[ 0[ 1[ 2[ ...
42 x = x - N0*(x1-x0); //shift
43 }
44 //mirror along boundary as often as necessary
45 while( (x<x0) || (x>x1) )
46 {
47 if( x < x0){
48 x = 2.*x0 - x;
49 //every mirror swaps the sign if Dirichlet
50 if( bc == dg::DIR || bc == dg::DIR_NEU)
51 negative = !negative;//swap sign
52 }
53 if( x > x1){
54 x = 2.*x1 - x;
55 if( bc == dg::DIR || bc == dg::NEU_DIR) //notice the different boundary NEU_DIR to the above DIR_NEU !
56 negative = !negative; //swap sign
57 }
58 }
59}
60
69template<class real_type>
70std::vector<real_type> legendre( real_type xn, unsigned n)
71{
72 assert( xn <= 1. && xn >= -1.);
73 std::vector<real_type> px(n);
74 if( xn == -1)
75 {
76 for( unsigned u=0; u<n; u++)
77 px[u] = (u%2 == 0) ? +1. : -1.;
78 }
79 else if( xn == 1)
80 {
81 for( unsigned i=0; i<n; i++)
82 px[i] = 1.;
83 }
84 else
85 {
86 px[0] = 1.;
87 if( n > 1)
88 {
89 px[1] = xn;
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);
92 }
93 }
94 return px;
95}
96
97// evaluate n base polynomials for n given abscissas
98template<class real_type>
99std::vector<real_type> lagrange( real_type x, const std::vector<real_type>& xi)
100{
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++)
105 {
106 if ( k != i)
107 l[i] *= (x-xi[k])/(xi[i]-xi[k]);
108 }
109 return l;
110}
111
112//THERE IS A BUG FOR PERIODIC BC !!
113template<class real_type>
114std::vector<real_type> choose_1d_abscissas( real_type X,
115 unsigned points_per_line, double lx,// const RealGrid1d<real_type>& g,
116 const thrust::host_vector<real_type>& abs, dg::bc bcx,
117 thrust::host_vector<unsigned>& cols)
118{
119 // Select points for nearest, linear or cubic interpolation
120 // lx is needed for PER bondary conditions
121 assert( abs.size() >= points_per_line && "There must be more points to interpolate\n");
122 //determine which cell (X) lies in
123 // abs must be sorted for std::lower_bound to work
124 auto it = std::lower_bound( abs.begin(), abs.end(), X);
125
126 std::vector<real_type> xs( points_per_line, 0);
127 cols.resize( points_per_line, 0);
128 switch( points_per_line)
129 {
130 case 1: xs[0] = 1.;
131 if( it == abs.begin())
132 cols[0] = 0;
133 else if( it == abs.end())
134 cols[0] = it - abs.begin() - 1;
135 else
136 {
137 if ( fabs(X - *it) < fabs( X - *(it-1)))
138 cols[0] = it - abs.begin();
139 else
140 cols[0] = it - abs.begin() -1;
141 }
142 break;
143 case 2: if( it == abs.begin())
144 {
145 if( bcx == dg::PER)
146 {
147 xs[0] = *it;
148 xs[1] = *(abs.end() -1)-lx;;
149 cols[0] = 0, cols[1] = abs.end()-abs.begin()-1;
150 }
151 else
152 {
153 //xs[0] = *it;
154 //xs[1] = *(it+1);
155 //cols[0] = 0, cols[1] = 1;
156 // This makes it consistent with fem_t
157 xs.resize(1);
158 xs[0] = *it;
159 cols[0] = 0;
160 }
161 }
162 else if( it == abs.end())
163 {
164 if( bcx == dg::PER)
165 {
166 xs[0] = *(abs.begin())+lx;
167 xs[1] = *(it-1);
168 cols[0] = 0, cols[1] = it-abs.begin()-1;
169 }
170 else
171 {
172 //xs[0] = *(it-2);
173 //xs[0] = *(it-1);
174 //cols[0] = it - abs.begin() - 2;
175 //cols[1] = it - abs.begin() - 1;
176 // This makes it consistent with fem_t
177 xs.resize(1);
178 xs[0] = *(it-1);
179 cols[0] = it-abs.begin()-1;
180 }
181 }
182 else
183 {
184 xs[0] = *(it-1);
185 xs[1] = *it;
186 cols[0] = it - abs.begin() - 1;
187 cols[1] = cols[0]+1;
188 }
189 break;
190 case 4: if( it <= abs.begin() +1)
191 {
192 if( bcx == dg::PER)
193 {
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;
200 }
201 else
202 {
203 it = abs.begin();
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;
208 }
209 }
210 else if( it >= abs.end() -1)
211 {
212 if( bcx == dg::PER)
213 {
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;
219 }
220 else
221 {
222 it = abs.end();
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;
226 cols[1] = cols[0]+1;
227 cols[2] = cols[1]+1;
228 cols[3] = cols[2]+1;
229 }
230 }
231 else
232 {
233 xs[0] = *(it-2), xs[1] = *(it-1);
234 xs[2] = *(it ), xs[3] = *(it+1);
235 cols[0] = it - abs.begin() - 2;
236 cols[1] = cols[0]+1;
237 cols[2] = cols[1]+1;
238 cols[3] = cols[2]+1;
239 }
240 break;
241 }
242 return xs;
243}
244
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)
250{
251 bool negative = false;
252 detail::shift( negative, X, bcx, g.x0(), g.x1());
253
254 //determine which cell (x) lies in
255 real_type xnn = (X-g.x0())/g.h();
256 unsigned nn = (unsigned)floor(xnn);
257 //determine normalized coordinates
258 real_type xn = 2.*xnn - (real_type)(2*nn+1);
259 //intervall correction
260 if (nn==g.N()) {
261 nn-=1;
262 xn = 1.;
263 }
264 //Test if the point is a Gauss point since then no interpolation is needed
265 int idxX =-1;
266 std::vector<real_type> gauss_nodesx = dg::DLT<real_type>::abscissas(g.n());
267 for( unsigned k=0; k<g.nx(); k++)
268 {
269 // Due to rounding errors xn will not perfectly be gauss_nodex[k] ( errors ~ 1e-14 possible)
270 // even if x are the abscissas of the grid
271 if( fabs( xn - gauss_nodesx[k]) < 1e-13)
272 idxX = nn*g.nx() + k; //determine which grid column it is
273 }
274 if( idxX < 0 or sp == dg::lspace) //there is no corresponding point
275 {
276 //evaluate 1d Legendre polynomials at (xn, yn)...
277 std::vector<real_type> px = detail::legendre( xn, g.n());
278 //...these are the matrix coefficients with which to multiply
279 if( sp == dg::xspace)
280 {
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];
287 px.swap( pxF);
288 }
289 for ( unsigned l=0; l<g.n(); l++)
290 {
291 cols.push_back( nn*g.n() + l);
292 vals.push_back( negative ? -px[l] : px[l]);
293 }
294 }
295 else //the point already exists
296 {
297 cols.push_back( idxX);
298 vals.push_back( negative ? -1. : 1.);
299 }
300}
301
302
303// dG interpolation
304template<class host_vector, class real_type>
306 dg::space sp,
307 const host_vector& x, // can be a view...
308 const RealGrid1d<real_type>& g,
309 dg::bc bcx )
310{
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++)
317 {
318 interpolation_row( sp, *ptr, g, bcx, column_indices, values);
319 row_offsets.push_back( values.size());
320 ptr++;
321 }
322 return {x.size(), g.size(), row_offsets, column_indices, values};
323}
324
325// nearest, linear or cubic interpolation
326template<class host_vector1, class host_vector2 >
327dg::SparseMatrix<int, dg::get_value_type<host_vector2>, thrust::host_vector> interpolation1d(
328 const host_vector1& x,
329 const host_vector2& abs, // must be sorted
331 std::string method )
332{
333 using real_type = dg::get_value_type<host_vector2>;
334 // boundary condidions for dg::Box likely won't work | Box is now removed
335 // from library ...
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")
341 points_per_line = 1;
342 else if( method == "linear")
343 points_per_line = 2;
344 else if( method == "cubic")
345 points_per_line = 4;
346 else
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++)
351 {
352 row_offsets.push_back(row_offsets[i]);
353 real_type X = *ptr;
354 ptr++;
355 bool negative = false;
356 detail::shift( negative, X, bcx, x0, x1);
357 // Test if point already exists since then no interpolation is needed
358 int idxX = -1;
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;
364 // THIS IS A VERY BAD IDEA PERFORMANCE WISE
365 //for( unsigned u=0; u<abs.size(); u++)
366 // if( fabs( X - abs[u]) <1e-13)
367 // idxX = u;
368 if( idxX < 0) //no corresponding point
369 {
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);
373
374 std::vector<real_type> px = detail::lagrange( X, xs);
375 // px may have size != points_per_line (at boundary)
376 for ( unsigned l=0; l<px.size(); l++)
377 {
378 row_offsets[i+1]++;
379 column_indices.push_back( cols[l]);
380 values.push_back(negative ? -px[l] : px[l]);
381 }
382 }
383 else //the point already exists
384 {
385 row_offsets[i+1]++;
386 column_indices.push_back(idxX);
387 values.push_back( negative ? -1. : 1.);
388 }
389 }
390 return {x.size(), abs.size(), row_offsets, column_indices, values};
391}
392
393}//namespace detail
397
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")
438{
439
440 std::array<dg::SparseMatrix<int,real_type,thrust::host_vector>,Nd> axes;
441 for( unsigned u=0; u<Nd; u++)
442 {
443 if( x[u].size() != x[0].size())
444 throw dg::Error( dg::Message(_ping_)<<"All coordinate lists must have same size "<<x[0].size());
445 if( method == "dg")
446 {
447 axes[u] = detail::interpolation1d( dg::xspace, x[u], g.grid(u), bcx[u]);
448 }
449 else
450 {
451 axes[u] = detail::interpolation1d( x[u], g.abscissas(u), bcx[u], g.p(u), g.q(u), method);
452 }
453 }
454 for( unsigned u=1; u<Nd; u++)
455 {
456 axes[0] = dg::tensorproduct_cols( axes[u], axes[0]);
457 }
458 return axes[0];
459}
460
479template<class host_vector, class real_type, typename = std::enable_if_t<dg::is_vector_v<host_vector>>>
481 const host_vector& x,
482 const RealGrid1d<real_type>& g,
483 dg::bc bcx = dg::NEU,
484 std::string method = "dg")
485{
486 dg::View<const host_vector> vx( x.data(), x.size());
487 return interpolation(
488 std::vector{vx}, g,
489 std::array<dg::bc,1>{bcx}, method);
490}
491
513template<class host_vector, class real_type>
515 const host_vector& x,
516 const host_vector& y,
518 dg::bc bcx = dg::NEU, dg::bc bcy = dg::NEU,
519 std::string method = "dg")
520{
521 dg::View<const host_vector> vx( x.data(), x.size());
522 dg::View<const host_vector> vy( y.data(), y.size());
523 return interpolation( std::vector{vx,vy}, g, {bcx,bcy}, method);
524}
525
526
527
551template<class host_vector, class real_type>
553 const host_vector& x,
554 const host_vector& y,
555 const host_vector& z,
557 dg::bc bcx = dg::NEU, dg::bc bcy = dg::NEU, dg::bc bcz = dg::PER,
558 std::string method = "dg")
559{
560 dg::View<const host_vector> vx( x.data(), x.size());
561 dg::View<const host_vector> vy( y.data(), y.size());
562 dg::View<const host_vector> vz( z.data(), z.size());
563 return interpolation( std::vector{vx,vy,vz}, g, {bcx,bcy,bcz}, method);
564}
585template<class real_type, size_t Nd>
587 const aRealTopology<real_type,Nd>& g_new,
588 const aRealTopology<real_type,Nd>& g_old, std::string method = "dg")
589{
590 //assert both grids are on the same box
591 for( unsigned u=0; u<Nd; u++)
592 {
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));
599 }
600 std::array<dg::SparseMatrix<int,real_type,thrust::host_vector>,Nd> axes;
601 for( unsigned u=0; u<Nd; u++)
602 {
603 if( method == "dg")
604 {
605 axes[u] = detail::interpolation1d( dg::xspace, g_new.abscissas(u),
606 g_old.grid(u), g_old.bc(u));
607 }
608 else
609 {
610 axes[u] = detail::interpolation1d( g_new.abscissas(u),
611 g_old.abscissas(u), g_old.bc(u), g_old.p(u), g_old.q(u),
612 method);
613 }
614 }
615 for( unsigned u=1; u<Nd; u++)
616 axes[0] = dg::tensorproduct( axes[u], axes[0]);
617 return axes[0];
618
619}
621
622
623}//namespace create
624
625
643template<class host_vector, class real_type>
644real_type interpolate(
645 dg::space sp,
646 const host_vector& v,
647 real_type x,
648 const RealGrid1d<real_type>& g,
649 dg::bc bcx = dg::NEU)
650{
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);
655 //multiply x
656 real_type value = 0;
657 for( unsigned j=0; j<vals.size(); j++)
658 value += v[cols[j]]*vals[j];
659 return value;
660}
661
681template<class host_vector, class real_type>
682real_type interpolate(
683 dg::space sp,
684 const host_vector& v,
685 real_type x, real_type y,
687 dg::bc bcx = dg::NEU, dg::bc bcy = dg::NEU )
688{
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);
696 //multiply x
697 real_type value = 0;
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];
701 return value;
702}
703
704
705} //namespace dg
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.
base topology classes
bc
Switch between boundary conditions.
Definition enums.h:15
space
Space of DG coefficients.
Definition enums.h:164
@ z
z direction
@ 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
@ y
y direction
@ x
x direction
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.
utility functions