Nektar++
Classes | Public Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
Nektar::LibUtilities::Interpolator Class Reference

A class that contains algorithms for interpolation between pts fields, expansions and different meshes. More...

#include <Interpolator.h>

Inheritance diagram for Nektar::LibUtilities::Interpolator:
[legend]

Classes

class  PtsPoint
 

Public Member Functions

 Interpolator (InterpMethod method=eNoMethod, short int coordId=-1, NekDouble filtWidth=0.0, int maxPts=1000)
 Constructor of the Interpolator class. More...
 
void CalcWeights (const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
 Compute interpolation weights without doing any interpolation. More...
 
void Interpolate (const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
 Interpolate from a pts field to a pts field. More...
 
int GetDim () const
 returns the dimension of the Interpolator. Should be higher than the dimensions of the interpolated fields More...
 
NekDouble GetFiltWidth () const
 Returns the filter width. More...
 
int GetCoordId () const
 Returns the coordinate id along which the interpolation should be performed. More...
 
InterpMethod GetInterpMethod () const
 Returns the interpolation method used by this interpolator. More...
 
LibUtilities::PtsFieldSharedPtr GetInField () const
 Returns the input field. More...
 
LibUtilities::PtsFieldSharedPtr GetOutField () const
 Returns the output field. More...
 
void PrintStatistics ()
 Returns if the weights have already been computed. More...
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetProgressCallback (FuncPointerT func, ObjectPointerT obj)
 sets a callback funtion which gets called every time the interpolation progresses More...
 

Protected Attributes

LibUtilities::PtsFieldSharedPtr m_ptsInField
 input field More...
 
LibUtilities::PtsFieldSharedPtr m_ptsOutField
 output field More...
 
std::function< void(const int position, const int goal)> m_progressCallback
 

Private Types

typedef boost::geometry::model::point< NekDouble, m_dim, boost::geometry::cs::cartesian > BPoint
 
typedef std::pair< BPoint, unsigned int > PtsPointPair
 
typedef boost::geometry::index::rtree< PtsPointPair, boost::geometry::index::rstar< 16 > > PtsRtree
 

Private Member Functions

void CalcW_Gauss (const PtsPoint &searchPt, const NekDouble sigma, const int maxPts=250)
 Computes interpolation weights using gaussian interpolation. More...
 
void CalcW_Linear (const PtsPoint &searchPt, int coordId)
 Computes interpolation weights using linear interpolation. More...
 
void CalcW_NNeighbour (const PtsPoint &searchPt)
 Computes interpolation weights using nearest neighbour interpolation. More...
 
void CalcW_Shepard (const PtsPoint &searchPt, int numPts)
 Computes interpolation weights using linear interpolation. More...
 
void CalcW_Quadratic (const PtsPoint &searchPt, int coordId)
 Computes interpolation weights using quadratic interpolation. More...
 
void SetupTree ()
 
void FindNeighbours (const PtsPoint &searchPt, std::vector< PtsPoint > &neighbourPts, const NekDouble dist, const unsigned int maxPts=1)
 Finds the neares neighbours of a point. More...
 
void FindNNeighbours (const PtsPoint &searchPt, std::vector< PtsPoint > &neighbourPts, const unsigned int numPts=1)
 Finds the neares neighbours of a point. More...
 

Private Attributes

InterpMethod m_method
 Interpolation Method. More...
 
std::shared_ptr< PtsRtreem_rtree
 A tree structure to speed up the neighbour search. Note that we fill it with an iterator, so instead of rstar, the packing algorithm is used. More...
 
Array< TwoD, NekDoublem_weights
 Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx]. More...
 
Array< TwoD, unsigned int > m_neighInds
 Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourIdx]. More...
 
NekDouble m_filtWidth
 Filter width used for some interpolation algorithms. More...
 
int m_maxPts
 Max number of interpolation points. More...
 
short int m_coordId
 coordinate id along which the interpolation should be performed More...
 

Static Private Attributes

static const int m_dim = 3
 dimension of this interpolator. Hardcoded to 3 More...
 

Detailed Description

A class that contains algorithms for interpolation between pts fields, expansions and different meshes.

Definition at line 70 of file LibUtilities/BasicUtils/Interpolator.h.

Member Typedef Documentation

◆ BPoint

Definition at line 172 of file LibUtilities/BasicUtils/Interpolator.h.

◆ PtsPointPair

typedef std::pair<BPoint, unsigned int> Nektar::LibUtilities::Interpolator::PtsPointPair
private

Definition at line 173 of file LibUtilities/BasicUtils/Interpolator.h.

◆ PtsRtree

typedef boost::geometry::index::rtree<PtsPointPair, boost::geometry::index::rstar<16> > Nektar::LibUtilities::Interpolator::PtsRtree
private

Definition at line 174 of file LibUtilities/BasicUtils/Interpolator.h.

Constructor & Destructor Documentation

◆ Interpolator()

Nektar::LibUtilities::Interpolator::Interpolator ( InterpMethod  method = eNoMethod,
short int  coordId = -1,
NekDouble  filtWidth = 0.0,
int  maxPts = 1000 
)
inline

Constructor of the Interpolator class.

Parameters
methodinterpolation method, defaults to a sensible value if not set
coordIdcoordinate id along which the interpolation should be performed
filtWidthfilter width, required by some algorithms such as eGauss
maxPtslimit number of considered points

if method is not specified, the best algorithm is chosen autpomatically.

If coordId is not specified, a full 1D/2D/3D interpolation is performed without collapsing any coordinate.

filtWidth must be specified for the eGauss algorithm only.

Definition at line 91 of file LibUtilities/BasicUtils/Interpolator.h.

References CalcWeights(), GetCoordId(), GetDim(), GetFiltWidth(), GetInField(), GetInterpMethod(), GetOutField(), Interpolate(), LIB_UTILITIES_EXPORT, and PrintStatistics().

95  : m_method(method), m_filtWidth(filtWidth), m_maxPts(maxPts),
96  m_coordId(coordId){};
short int m_coordId
coordinate id along which the interpolation should be performed
int m_maxPts
Max number of interpolation points.
NekDouble m_filtWidth
Filter width used for some interpolation algorithms.

Member Function Documentation

◆ CalcW_Gauss()

void Nektar::LibUtilities::Interpolator::CalcW_Gauss ( const PtsPoint searchPt,
const NekDouble  sigma,
const int  maxPts = 250 
)
private

Computes interpolation weights using gaussian interpolation.

Parameters
searchPtpoint for which the weights are computed
sigmastandard deviation of the gauss function

Performs an interpolation using gauss weighting. Ideal for filtering fields. The filter width should be half the FWHM (= 1.1774 sigma) and must be set in the constructor of the Interpolator class.

Definition at line 345 of file LibUtilities/BasicUtils/Interpolator.cpp.

References Nektar::LibUtilities::Interpolator::PtsPoint::idx.

348 {
349  // find nearest neighbours
350  vector<PtsPoint> neighbourPts;
351  FindNeighbours(searchPt, neighbourPts, 4 * sigma, maxPts);
352  size_t numPts = neighbourPts.size();
353 
354  // handle the cases that there was no or just one point within 4 * sigma
355  if (numPts == 0)
356  {
357  return;
358  }
359  if (numPts == 1)
360  {
361  m_neighInds[searchPt.idx][0] = neighbourPts.front().idx;
362  m_weights[searchPt.idx][0] = 1.0;
363 
364  return;
365  }
366 
367  NekDouble sigmaNew = 0.25 * neighbourPts.back().dist;
368 
369  for (size_t i = 0; i < numPts; i++)
370  {
371  m_neighInds[searchPt.idx][i] = neighbourPts.at(i).idx;
372  }
373 
374  NekDouble wSum = 0.0;
375  NekDouble ts2 = 2.0 * sigmaNew * sigmaNew;
376  for (size_t i = 0; i < numPts; ++i)
377  {
378  m_weights[searchPt.idx][i] =
379  exp(-1.0 * neighbourPts[i].dist * neighbourPts[i].dist / ts2);
380  wSum += m_weights[searchPt.idx][i];
381  }
382 
383  for (size_t i = 0; i < numPts; ++i)
384  {
385  m_weights[searchPt.idx][i] = m_weights[searchPt.idx][i] / wSum;
386  }
387 }
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
void FindNeighbours(const PtsPoint &searchPt, std::vector< PtsPoint > &neighbourPts, const NekDouble dist, const unsigned int maxPts=1)
Finds the neares neighbours of a point.
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
double NekDouble

◆ CalcW_Linear()

void Nektar::LibUtilities::Interpolator::CalcW_Linear ( const PtsPoint searchPt,
int  m_coordId 
)
private

Computes interpolation weights using linear interpolation.

Parameters
searchPtpoint for which the weights are computed
m_coordIdcoordinate id along which the interpolation should be performed

Currently, only implemented for 1D

Definition at line 398 of file LibUtilities/BasicUtils/Interpolator.cpp.

References ASSERTL0, Nektar::LibUtilities::Interpolator::PtsPoint::coords, Nektar::LibUtilities::Interpolator::PtsPoint::idx, and Nektar::NekConstants::kNekZeroTol.

399 {
400  int npts = m_ptsInField->GetNpoints();
401  int i;
402 
403  NekDouble coord = searchPt.coords[m_coordId];
404 
405  for (i = 0; i < npts - 1; ++i)
406  {
407  if ((m_ptsInField->GetPointVal(0, i) <=
408  (coord + NekConstants::kNekZeroTol)) &&
409  (coord <=
410  (m_ptsInField->GetPointVal(0, i + 1) + NekConstants::kNekZeroTol)))
411  {
412  NekDouble pdiff = m_ptsInField->GetPointVal(0, i + 1) -
413  m_ptsInField->GetPointVal(0, i);
414 
415  m_neighInds[searchPt.idx][0] = i;
416  m_neighInds[searchPt.idx][1] = i + 1;
417 
418  m_weights[searchPt.idx][0] =
419  (m_ptsInField->GetPointVal(0, i + 1) - coord) / pdiff;
420  m_weights[searchPt.idx][1] =
421  (coord - m_ptsInField->GetPointVal(0, i)) / pdiff;
422 
423  break;
424  }
425  }
426  ASSERTL0(i != npts - 1, "Failed to find coordinate " +
427  boost::lexical_cast<string>(coord) +
428  " within provided input points");
429 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
short int m_coordId
coordinate id along which the interpolation should be performed
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
static const NekDouble kNekZeroTol
double NekDouble

◆ CalcW_NNeighbour()

void Nektar::LibUtilities::Interpolator::CalcW_NNeighbour ( const PtsPoint searchPt)
private

Computes interpolation weights using nearest neighbour interpolation.

Parameters
searchPtpoint for which the weights are computed
m_coordIdcoordinate id along which the interpolation should be performed

Definition at line 439 of file LibUtilities/BasicUtils/Interpolator.cpp.

References Nektar::LibUtilities::Interpolator::PtsPoint::idx.

440 {
441  // find nearest neighbours
442  vector<PtsPoint> neighbourPts;
443  // TODO: we currently dont handle the case when there are more than one
444  // most distant points (of same distance)
445  FindNNeighbours(searchPt, neighbourPts, 1);
446 
447  m_neighInds[searchPt.idx][0] = neighbourPts[0].idx;
448  m_weights[searchPt.idx][0] = 1.0;
449 }
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
void FindNNeighbours(const PtsPoint &searchPt, std::vector< PtsPoint > &neighbourPts, const unsigned int numPts=1)
Finds the neares neighbours of a point.
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].

◆ CalcW_Quadratic()

void Nektar::LibUtilities::Interpolator::CalcW_Quadratic ( const PtsPoint searchPt,
int  m_coordId 
)
private

Computes interpolation weights using quadratic interpolation.

Parameters
searchPtpoint for which the weights are computed
m_coordIdcoordinate id along which the interpolation should be performed

Currently, only implemented for 1D. Falls back to linear interpolation if only 2 values are available.

Definition at line 513 of file LibUtilities/BasicUtils/Interpolator.cpp.

References ASSERTL0, Nektar::LibUtilities::Interpolator::PtsPoint::coords, Nektar::LibUtilities::Interpolator::PtsPoint::idx, and Nektar::NekConstants::kNekZeroTol.

514 {
515  int npts = m_ptsInField->GetNpoints();
516  int i;
517 
518  NekDouble coord = searchPt.coords[m_coordId];
519 
520  for (i = 0; i < npts - 1; ++i)
521  {
522  if ((m_ptsInField->GetPointVal(0, i) <=
523  (coord + NekConstants::kNekZeroTol)) &&
524  (coord <=
525  (m_ptsInField->GetPointVal(0, i + 1) + NekConstants::kNekZeroTol)))
526  {
527  NekDouble pdiff = m_ptsInField->GetPointVal(0, i + 1) -
528  m_ptsInField->GetPointVal(0, i);
529  NekDouble h1, h2, h3;
530 
531  if (i < npts - 2)
532  {
533  // forwards stencil
534  NekDouble pdiff2 = m_ptsInField->GetPointVal(0, i + 2) -
535  m_ptsInField->GetPointVal(0, i + 1);
536 
537  h1 = (m_ptsInField->GetPointVal(0, i + 1) - coord) *
538  (m_ptsInField->GetPointVal(0, i + 2) - coord) /
539  (pdiff * (pdiff + pdiff2));
540  h2 = (coord - m_ptsInField->GetPointVal(0, i)) *
541  (m_ptsInField->GetPointVal(0, i + 2) - coord) /
542  (pdiff * pdiff2);
543  h3 = (coord - m_ptsInField->GetPointVal(0, i)) *
544  (coord - m_ptsInField->GetPointVal(0, i + 1)) /
545  ((pdiff + pdiff2) * pdiff2);
546 
547  m_neighInds[searchPt.idx][0] = i;
548  m_neighInds[searchPt.idx][1] = i + 1;
549  m_neighInds[searchPt.idx][2] = i + 2;
550  }
551  else
552  {
553  // backwards stencil
554  NekDouble pdiff2 = m_ptsInField->GetPointVal(0, i) -
555  m_ptsInField->GetPointVal(0, i - 1);
556 
557  h1 = (m_ptsInField->GetPointVal(0, i + 1) - coord) *
558  (coord - m_ptsInField->GetPointVal(0, i - 1)) /
559  (pdiff * pdiff2);
560  h2 = (coord - m_ptsInField->GetPointVal(0, i)) *
561  (coord - m_ptsInField->GetPointVal(0, i - 1)) /
562  (pdiff * (pdiff + pdiff2));
563  h3 = (m_ptsInField->GetPointVal(0, i) - coord) *
564  (m_ptsInField->GetPointVal(0, i + 1) - coord) /
565  ((pdiff + pdiff2) * pdiff);
566 
567  m_neighInds[searchPt.idx][0] = i;
568  m_neighInds[searchPt.idx][1] = i + 1;
569  m_neighInds[searchPt.idx][2] = i - 1;
570  }
571 
572  m_weights[searchPt.idx][0] = h1;
573  m_weights[searchPt.idx][1] = h2;
574  m_weights[searchPt.idx][2] = h3;
575 
576  break;
577  }
578  }
579  ASSERTL0(i != npts - 1, "Failed to find coordinate " +
580  boost::lexical_cast<string>(coord) +
581  " within provided input points");
582 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
short int m_coordId
coordinate id along which the interpolation should be performed
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
static const NekDouble kNekZeroTol
double NekDouble

◆ CalcW_Shepard()

void Nektar::LibUtilities::Interpolator::CalcW_Shepard ( const PtsPoint searchPt,
int  numPts 
)
private

Computes interpolation weights using linear interpolation.

Parameters
searchPtpoint for which the weights are computed
m_coordIdcoordinate id along which the interpolation should be performed

The algorithm is based on Shepard, D. (1968). A two-dimensional interpolation function for irregularly-spaced data. Proceedings of the 1968 23rd ACM National Conference. pp. 517–524.

In order to save memory, for n dimesnions, only 2^n points are considered. Contrary to Shepard, we use a fixed number of points with fixed weighting factors 1/d^n.

Definition at line 466 of file LibUtilities/BasicUtils/Interpolator.cpp.

References Nektar::LibUtilities::Interpolator::PtsPoint::idx, and Nektar::NekConstants::kNekZeroTol.

467 {
468  // find nearest neighbours
469  vector<PtsPoint> neighbourPts;
470  FindNNeighbours(searchPt, neighbourPts, numPts);
471 
472  for (int i = 0; i < neighbourPts.size(); i++)
473  {
474  m_neighInds[searchPt.idx][i] = neighbourPts[i].idx;
475  }
476 
477  // In case d < kVertexTheSameDouble ( d^2 < kNekSqrtTol), use the exact
478  // point and return
479  for (int i = 0; i < neighbourPts.size(); ++i)
480  {
481  if (neighbourPts[i].dist <= NekConstants::kNekZeroTol)
482  {
483  m_weights[searchPt.idx][i] = 1.0;
484  return;
485  }
486  }
487 
488  NekDouble wSum = 0.0;
489  for (int i = 0; i < neighbourPts.size(); ++i)
490  {
491  m_weights[searchPt.idx][i] = 1 / pow(double(neighbourPts[i].dist),
492  double(m_ptsInField->GetDim()));
493  wSum += m_weights[searchPt.idx][i];
494  }
495 
496  for (int i = 0; i < neighbourPts.size(); ++i)
497  {
498  m_weights[searchPt.idx][i] = m_weights[searchPt.idx][i] / wSum;
499  }
500 }
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
void FindNNeighbours(const PtsPoint &searchPt, std::vector< PtsPoint > &neighbourPts, const unsigned int numPts=1)
Finds the neares neighbours of a point.
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
static const NekDouble kNekZeroTol
double NekDouble

◆ CalcWeights()

void Nektar::LibUtilities::Interpolator::CalcWeights ( const LibUtilities::PtsFieldSharedPtr  ptsInField,
LibUtilities::PtsFieldSharedPtr ptsOutField 
)

Compute interpolation weights without doing any interpolation.

Parameters
ptsInFieldinput field
ptsOutFieldoutput field

In and output fields must have the same dimension. The most suitable algorithm is chosen automatically if it wasnt set explicitly.

Definition at line 58 of file LibUtilities/BasicUtils/Interpolator.cpp.

References ASSERTL0, Nektar::LibUtilities::eGauss, Nektar::LibUtilities::eNearestNeighbour, Nektar::LibUtilities::eNoMethod, Nektar::LibUtilities::eQuadratic, Nektar::LibUtilities::eShepard, Nektar::NekConstants::kNekZeroTol, and NEKERROR.

Referenced by Interpolator().

60 {
61  ASSERTL0(ptsInField->GetDim() <= m_dim, "too many dimesions in inField");
62  ASSERTL0(ptsOutField->GetDim() <= m_dim, "too many dimesions in outField");
63 
64  m_ptsInField = ptsInField;
65  m_ptsOutField = ptsOutField;
66 
67  size_t nOutPts = m_ptsOutField->GetNpoints();
68  int lastProg = 0;
69 
70  // set a default method
71  if (m_method == eNoMethod)
72  {
73  if (m_ptsInField->GetDim() == 1 || m_coordId >= 0)
74  {
76  }
77  else
78  {
80  }
81  }
82 
83  if (m_method != eQuadratic)
84  {
85  SetupTree();
86  }
87 
88  switch (m_method)
89  {
90  case eNearestNeighbour:
91  {
92  m_weights = Array<TwoD, NekDouble>(nOutPts, 1, 0.0);
93  m_neighInds = Array<TwoD, unsigned int>(nOutPts, 1, (unsigned int) 0);
94 
95  for (size_t i = 0; i < nOutPts; ++i)
96  {
97  Array<OneD, NekDouble> tmp(m_dim, 0.0);
98  for (size_t j = 0; j < m_ptsOutField->GetDim(); ++j)
99  {
100  tmp[j] = m_ptsOutField->GetPointVal(j, i);
101  }
102  PtsPoint searchPt(i, tmp, 1E30);
103 
104  CalcW_NNeighbour(searchPt);
105 
106  int progress = int(100 * i / nOutPts);
107  if (m_progressCallback && progress > lastProg)
108  {
109  m_progressCallback(i, nOutPts);
110  lastProg = progress;
111  }
112  }
113 
114  break;
115  }
116 
117  case eQuadratic:
118  {
119  ASSERTL0(m_ptsInField->GetDim() == 1 || m_coordId >= 0,
120  "not implemented");
121 
122  m_weights = Array<TwoD, NekDouble>(nOutPts, 3, 0.0);
123  m_neighInds = Array<TwoD, unsigned int>(nOutPts, 3, (unsigned int) 0);
124 
125  for (size_t i = 0; i < nOutPts; ++i)
126  {
127  Array<OneD, NekDouble> tmp(m_dim, 0.0);
128  for (size_t j = 0; j < m_ptsOutField->GetDim(); ++j)
129  {
130  tmp[j] = m_ptsOutField->GetPointVal(j, i);
131  }
132  PtsPoint searchPt(i, tmp, 1E30);
133 
134  if (m_ptsInField->GetNpoints() <= 2)
135  {
136  CalcW_Linear(searchPt, m_coordId);
137  }
138  else
139  {
140  CalcW_Quadratic(searchPt, m_coordId);
141  }
142 
143  int progress = int(100 * i / nOutPts);
144  if (m_progressCallback && progress > lastProg)
145  {
146  m_progressCallback(i, nOutPts);
147  lastProg = progress;
148  }
149  }
150 
151  break;
152  }
153 
154  case eShepard:
155  {
156  int numPts = m_ptsInField->GetDim();
157  numPts = 2 << numPts; // 2 ^ numPts
158  numPts = min(numPts, int(m_ptsInField->GetNpoints() / 2));
159 
160  m_weights = Array<TwoD, NekDouble>(nOutPts, numPts, 0.0);
161  m_neighInds = Array<TwoD, unsigned int>(nOutPts, numPts, (unsigned int) 0);
162 
163  for (size_t i = 0; i < nOutPts; ++i)
164  {
165  Array<OneD, NekDouble> tmp(m_dim, 0.0);
166  for (size_t j = 0; j < m_ptsOutField->GetDim(); ++j)
167  {
168  tmp[j] = m_ptsOutField->GetPointVal(j, i);
169  }
170  PtsPoint searchPt(i, tmp, 1E30);
171 
172  CalcW_Shepard(searchPt, numPts);
173 
174  int progress = int(100 * i / nOutPts);
175  if (m_progressCallback && progress > lastProg)
176  {
177  m_progressCallback(i, nOutPts);
178  lastProg = progress;
179  }
180  }
181 
182  break;
183  }
184 
185  case eGauss:
186  {
188  "No filter width set");
189  // use m_filtWidth as FWHM
190  NekDouble sigma = m_filtWidth * 0.4246609001;
191 
192  m_maxPts = min(m_maxPts, int(m_ptsInField->GetNpoints() / 2));
193 
194  m_weights = Array<TwoD, NekDouble>(nOutPts, m_maxPts, 0.0);
195  m_neighInds = Array<TwoD, unsigned int>(nOutPts, m_maxPts, (unsigned int) 0);
196 
197  for (size_t i = 0; i < nOutPts; ++i)
198  {
199  Array<OneD, NekDouble> tmp(m_dim, 0.0);
200  for (size_t j = 0; j < m_ptsOutField->GetDim(); ++j)
201  {
202  tmp[j] = m_ptsOutField->GetPointVal(j, i);
203  }
204  PtsPoint searchPt(i, tmp, 1E30);
205 
206  CalcW_Gauss(searchPt, sigma, m_maxPts);
207 
208  int progress = int(100 * i / nOutPts);
209  if (m_progressCallback && progress > lastProg)
210  {
211  m_progressCallback(i, nOutPts);
212  lastProg = progress;
213  }
214  }
215 
216  break;
217  }
218 
219  default:
220  NEKERROR(ErrorUtil::efatal, "Invalid interpolation m_method");
221  break;
222  }
223 }
LibUtilities::PtsFieldSharedPtr m_ptsOutField
output field
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209
short int m_coordId
coordinate id along which the interpolation should be performed
void CalcW_Gauss(const PtsPoint &searchPt, const NekDouble sigma, const int maxPts=250)
Computes interpolation weights using gaussian interpolation.
void CalcW_Linear(const PtsPoint &searchPt, int coordId)
Computes interpolation weights using linear interpolation.
int m_maxPts
Max number of interpolation points.
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
void CalcW_Shepard(const PtsPoint &searchPt, int numPts)
Computes interpolation weights using linear interpolation.
void CalcW_NNeighbour(const PtsPoint &searchPt)
Computes interpolation weights using nearest neighbour interpolation.
static const int m_dim
dimension of this interpolator. Hardcoded to 3
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
static const NekDouble kNekZeroTol
void CalcW_Quadratic(const PtsPoint &searchPt, int coordId)
Computes interpolation weights using quadratic interpolation.
double NekDouble
std::function< void(const int position, const int goal)> m_progressCallback
NekDouble m_filtWidth
Filter width used for some interpolation algorithms.

◆ FindNeighbours()

void Nektar::LibUtilities::Interpolator::FindNeighbours ( const PtsPoint searchPt,
std::vector< PtsPoint > &  neighbourPts,
const NekDouble  dist,
const unsigned int  maxPts = 1 
)
private

Finds the neares neighbours of a point.

Parameters
searchPtpoint for which the neighbours are searched
neighbourPtspossible neighbour points
distlimits the distance of the neighbours

Definition at line 669 of file LibUtilities/BasicUtils/Interpolator.cpp.

References Nektar::LibUtilities::Interpolator::PtsPoint::coords.

673 {
674  BPoint searchBPoint(searchPt.coords[0], searchPt.coords[1],
675  searchPt.coords[2]);
676  BPoint bbMin(searchPt.coords[0] - dist, searchPt.coords[1] - dist,
677  searchPt.coords[2] - dist);
678  BPoint bbMax(searchPt.coords[0] + dist, searchPt.coords[1] + dist,
679  searchPt.coords[2] + dist);
680  bg::model::box<BPoint> bbox(bbMin, bbMax);
681 
682  // find points within the distance box
683  std::vector<PtsPointPair> result;
684  if (maxPts >= 1)
685  {
686  m_rtree->query(bgi::within(bbox) && bgi::nearest(searchBPoint, maxPts),
687  std::back_inserter(result));
688  }
689  else
690  {
691  m_rtree->query(bgi::within(bbox), std::back_inserter(result));
692  }
693 
694  // massage into or own format
695  for (int i = 0; i < result.size(); ++i)
696  {
697  int idx = result[i].second;
698  Array<OneD, NekDouble> coords(m_dim, 0.0);
699  for (int j = 0; j < m_ptsInField->GetDim(); ++j)
700  {
701  coords[j] = m_ptsInField->GetPointVal(j, idx);
702  }
703  NekDouble d = bg::distance(searchBPoint, result[i].first);
704 
705  // discard points beyonf dist
706  if (d <= dist)
707  {
708  neighbourPts.push_back(PtsPoint(idx, coords, d));
709  }
710  }
711 
712  sort(neighbourPts.begin(), neighbourPts.end());
713 }
boost::geometry::model::point< NekDouble, m_dim, boost::geometry::cs::cartesian > BPoint
std::shared_ptr< PtsRtree > m_rtree
A tree structure to speed up the neighbour search. Note that we fill it with an iterator, so instead of rstar, the packing algorithm is used.
static const int m_dim
dimension of this interpolator. Hardcoded to 3
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
double NekDouble

◆ FindNNeighbours()

void Nektar::LibUtilities::Interpolator::FindNNeighbours ( const PtsPoint searchPt,
std::vector< PtsPoint > &  neighbourPts,
const unsigned int  numPts = 1 
)
private

Finds the neares neighbours of a point.

Parameters
searchPtpoint for which the neighbours are searched
neighbourPtspossible neighbour points
numPtslimits the number of neighbours found to the numPts nearest ones

Definition at line 635 of file LibUtilities/BasicUtils/Interpolator.cpp.

References Nektar::LibUtilities::Interpolator::PtsPoint::coords.

638 {
639  std::vector<PtsPointPair> result;
640  BPoint searchBPoint(searchPt.coords[0], searchPt.coords[1],
641  searchPt.coords[2]);
642  m_rtree->query(bgi::nearest(searchBPoint, numPts),
643  std::back_inserter(result));
644 
645  // massage into or own format
646  for (int i = 0; i < result.size(); ++i)
647  {
648  int idx = result[i].second;
649  Array<OneD, NekDouble> coords(m_dim, 0.0);
650  for (int j = 0; j < m_ptsInField->GetDim(); ++j)
651  {
652  coords[j] = m_ptsInField->GetPointVal(j, idx);
653  }
654  NekDouble d = bg::distance(searchBPoint, result[i].first);
655  neighbourPts.push_back(PtsPoint(idx, coords, d));
656  }
657 
658  sort(neighbourPts.begin(), neighbourPts.end());
659 }
boost::geometry::model::point< NekDouble, m_dim, boost::geometry::cs::cartesian > BPoint
std::shared_ptr< PtsRtree > m_rtree
A tree structure to speed up the neighbour search. Note that we fill it with an iterator, so instead of rstar, the packing algorithm is used.
static const int m_dim
dimension of this interpolator. Hardcoded to 3
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
double NekDouble

◆ GetCoordId()

int Nektar::LibUtilities::Interpolator::GetCoordId ( ) const

Returns the coordinate id along which the interpolation should be performed.

Definition at line 288 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

289 {
290  return m_coordId;
291 }
short int m_coordId
coordinate id along which the interpolation should be performed

◆ GetDim()

int Nektar::LibUtilities::Interpolator::GetDim ( ) const

returns the dimension of the Interpolator. Should be higher than the dimensions of the interpolated fields

Definition at line 283 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

284 {
285  return m_dim;
286 }
static const int m_dim
dimension of this interpolator. Hardcoded to 3

◆ GetFiltWidth()

NekDouble Nektar::LibUtilities::Interpolator::GetFiltWidth ( ) const

Returns the filter width.

Definition at line 293 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

294 {
295  return m_filtWidth;
296 }
NekDouble m_filtWidth
Filter width used for some interpolation algorithms.

◆ GetInField()

LibUtilities::PtsFieldSharedPtr Nektar::LibUtilities::Interpolator::GetInField ( ) const

Returns the input field.

Definition at line 303 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

304 {
305  return m_ptsInField;
306 }
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field

◆ GetInterpMethod()

InterpMethod Nektar::LibUtilities::Interpolator::GetInterpMethod ( ) const

Returns the interpolation method used by this interpolator.

Definition at line 298 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

299 {
300  return m_method;
301 }

◆ GetOutField()

LibUtilities::PtsFieldSharedPtr Nektar::LibUtilities::Interpolator::GetOutField ( ) const

Returns the output field.

Definition at line 308 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

309 {
310  return m_ptsOutField;
311 }
LibUtilities::PtsFieldSharedPtr m_ptsOutField
output field

◆ Interpolate()

void Nektar::LibUtilities::Interpolator::Interpolate ( const LibUtilities::PtsFieldSharedPtr  ptsInField,
LibUtilities::PtsFieldSharedPtr ptsOutField 
)

Interpolate from a pts field to a pts field.

Parameters
ptsInFieldinput field
ptsOutFieldoutput field

In and output fields must have the same dimension and number of fields. The most suitable algorithm is chosen automatically if it wasnt set explicitly.

Definition at line 235 of file LibUtilities/BasicUtils/Interpolator.cpp.

References ASSERTL0.

Referenced by Nektar::Utilities::ProcessCurve::EvaluateCoordinate(), and Interpolator().

237 {
238  ASSERTL0(ptsInField->GetNFields() == ptsOutField->GetNFields(),
239  "number of fields does not match");
240  ASSERTL0(ptsInField->GetDim() <= m_dim, "too many dimesions in inField");
241  ASSERTL0(ptsOutField->GetDim() <= m_dim, "too many dimesions in outField");
242 
243  m_ptsInField = ptsInField;
244  m_ptsOutField = ptsOutField;
245 
246  if (m_weights.GetRows() == 0)
247  {
248  CalcWeights(m_ptsInField, m_ptsOutField);
249  }
250 
251  ASSERTL0(m_weights.GetRows() == m_ptsOutField->GetNpoints(),
252  "weights dimension mismatch");
253 
254  size_t nFields = m_ptsOutField->GetNFields();
255  size_t nOutPts = m_ptsOutField->GetNpoints();
256  size_t inDim = m_ptsInField->GetDim();
257 
258  // interpolate points and transform
259  for (size_t i = 0; i < nFields; ++i)
260  {
261  for (size_t j = 0; j < nOutPts; ++j)
262  {
263  size_t nPts = m_weights.GetColumns();
264 
265  // skip if there were no neighbours found for this point
266  if (nPts == 0)
267  {
268  continue;
269  }
270 
271  NekDouble val = 0.0;
272  for (size_t k = 0; k < nPts; ++k)
273  {
274  size_t nIdx = m_neighInds[j][k];
275  val += m_weights[j][k] *
276  m_ptsInField->GetPointVal(inDim + i, nIdx);
277  }
278  m_ptsOutField->SetPointVal(m_ptsOutField->GetDim() + i, j, val);
279  }
280  }
281 }
LibUtilities::PtsFieldSharedPtr m_ptsOutField
output field
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
void CalcWeights(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
Compute interpolation weights without doing any interpolation.
Array< TwoD, NekDouble > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
static const int m_dim
dimension of this interpolator. Hardcoded to 3
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
double NekDouble

◆ PrintStatistics()

void Nektar::LibUtilities::Interpolator::PrintStatistics ( )

Returns if the weights have already been computed.

Definition at line 313 of file LibUtilities/BasicUtils/Interpolator.cpp.

Referenced by Interpolator().

314 {
315  int meanN = 0;
316  for (int i = 0; i < m_neighInds.GetRows(); ++i)
317  {
318  for (int j = 0; j < m_neighInds.GetColumns(); ++j)
319  {
320  if (m_neighInds[i][j] > 0)
321  {
322  meanN +=1;
323  }
324  }
325  }
326 
327  cout << "Number of points: " << m_neighInds.GetRows() << endl;
328  if (m_neighInds.GetRows() > 0)
329  {
330  cout << "mean Number of Neighbours per point: "
331  << meanN / m_neighInds.GetRows() << endl;
332  }
333 }
Array< TwoD, unsigned int > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...

◆ SetProgressCallback()

template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::LibUtilities::Interpolator::SetProgressCallback ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

sets a callback funtion which gets called every time the interpolation progresses

Definition at line 134 of file LibUtilities/BasicUtils/Interpolator.h.

References m_progressCallback.

Referenced by Nektar::FieldUtils::ProcessInterpPoints::InterpolateFieldToPts(), Nektar::FieldUtils::ProcessInterpPtsToPts::InterpolatePtsToPts(), Nektar::FieldUtils::ProcessInterpField::Process(), and Nektar::FieldUtils::ProcessInterpPointDataToFld::Process().

135  {
136  m_progressCallback = std::bind(
137  func, obj, std::placeholders::_1, std::placeholders::_2);
138  }
std::function< void(const int position, const int goal)> m_progressCallback

◆ SetupTree()

void Nektar::LibUtilities::Interpolator::SetupTree ( )
private

Definition at line 584 of file LibUtilities/BasicUtils/Interpolator.cpp.

References Nektar::NekConstants::kNekZeroTol.

585 {
586  std::vector<PtsPointPair> inPoints;
587  for (int i = 0; i < m_ptsInField->GetNpoints(); ++i)
588  {
589  Array<OneD, NekDouble> coords(3, 0.0);
590  for (int j = 0; j < m_ptsInField->GetDim(); ++j)
591  {
592  coords[j] = m_ptsInField->GetPointVal(j, i);
593  }
594  inPoints.push_back(
595  PtsPointPair(BPoint(coords[0], coords[1], coords[2]), i));
596  }
598  m_rtree->insert(inPoints.begin(), inPoints.end());
599 
600  // remove duplicates from tree
601  for (auto &it : inPoints)
602  {
603  std::vector<PtsPointPair> result;
604 
605  // find nearest 2 points (2 because one of these might be the one we
606  // are
607  // checking)
608  m_rtree->query(bgi::nearest(it.first, 2),
609  std::back_inserter(result));
610 
611  // in case any of these 2 points is too close, remove the current
612  // point
613  // from the tree
614  for (auto &it2 : result)
615  {
616  if (it.second != it2.second &&
617  bg::distance(it.first, it2.first) <= NekConstants::kNekZeroTol)
618  {
619  m_rtree->remove(it);
620  break;
621  }
622  }
623  }
624 }
boost::geometry::model::point< NekDouble, m_dim, boost::geometry::cs::cartesian > BPoint
std::shared_ptr< PtsRtree > m_rtree
A tree structure to speed up the neighbour search. Note that we fill it with an iterator, so instead of rstar, the packing algorithm is used.
LibUtilities::PtsFieldSharedPtr m_ptsInField
input field
static const NekDouble kNekZeroTol
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

Member Data Documentation

◆ m_coordId

short int Nektar::LibUtilities::Interpolator::m_coordId
private

coordinate id along which the interpolation should be performed

Definition at line 195 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_dim

const int Nektar::LibUtilities::Interpolator::m_dim = 3
staticprivate

dimension of this interpolator. Hardcoded to 3

Definition at line 171 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_filtWidth

NekDouble Nektar::LibUtilities::Interpolator::m_filtWidth
private

Filter width used for some interpolation algorithms.

Definition at line 191 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_maxPts

int Nektar::LibUtilities::Interpolator::m_maxPts
private

Max number of interpolation points.

Definition at line 193 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_method

InterpMethod Nektar::LibUtilities::Interpolator::m_method
private

Interpolation Method.

Definition at line 179 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_neighInds

Array<TwoD, unsigned int> Nektar::LibUtilities::Interpolator::m_neighInds
private

Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourIdx].

Definition at line 189 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_progressCallback

std::function<void(const int position, const int goal)> Nektar::LibUtilities::Interpolator::m_progressCallback
protected

Definition at line 148 of file LibUtilities/BasicUtils/Interpolator.h.

Referenced by SetProgressCallback().

◆ m_ptsInField

LibUtilities::PtsFieldSharedPtr Nektar::LibUtilities::Interpolator::m_ptsInField
protected

input field

Definition at line 143 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_ptsOutField

LibUtilities::PtsFieldSharedPtr Nektar::LibUtilities::Interpolator::m_ptsOutField
protected

output field

Definition at line 145 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_rtree

std::shared_ptr<PtsRtree> Nektar::LibUtilities::Interpolator::m_rtree
private

A tree structure to speed up the neighbour search. Note that we fill it with an iterator, so instead of rstar, the packing algorithm is used.

Definition at line 183 of file LibUtilities/BasicUtils/Interpolator.h.

◆ m_weights

Array<TwoD, NekDouble> Nektar::LibUtilities::Interpolator::m_weights
private

Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].

Definition at line 186 of file LibUtilities/BasicUtils/Interpolator.h.