Nektar++
LibUtilities/BasicUtils/Interpolator.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File Interpolator.h
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2016 Kilian Lackhove
10 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11 // Department of Aeronautics, Imperial College London (UK), and Scientific
12 // Computing and Imaging Institute, University of Utah (USA).
13 //
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Interpolator
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef LIBUTILITIES_BASICUTILS_INTERPOLATOR_H
37 #define LIBUTILITIES_BASICUTILS_INTERPOLATOR_H
38 
39 #include <vector>
40 #include <iostream>
41 #include <functional>
42 #include <memory>
43 
44 #include <boost/geometry/geometries/box.hpp>
45 #include <boost/geometry/geometries/point.hpp>
46 #include <boost/geometry/index/rtree.hpp>
47 
53 
54 namespace Nektar
55 {
56 namespace LibUtilities
57 {
58 
60 {
66 };
67 
68 /// A class that contains algorithms for interpolation between pts fields,
69 /// expansions and different meshes
71 {
72 public:
73  /**
74  * @brief Constructor of the Interpolator class
75  *
76  * @param method interpolation method, defaults to a sensible value if not
77  * set
78  * @param coordId coordinate id along which the interpolation should be
79  * performed
80  * @param filtWidth filter width, required by some algorithms such as eGauss
81  * @param maxPts limit number of considered points
82  *
83  * if method is not specified, the best algorithm is chosen autpomatically.
84  *
85  * If coordId is not specified, a full 1D/2D/3D interpolation is performed
86  * without
87  * collapsing any coordinate.
88  *
89  * filtWidth must be specified for the eGauss algorithm only.
90  */
92  short int coordId = -1,
93  NekDouble filtWidth = 0.0,
94  int maxPts = 1000)
95  : m_method(method), m_filtWidth(filtWidth), m_maxPts(maxPts),
96  m_coordId(coordId){};
97 
98  /// Compute interpolation weights without doing any interpolation
100  const LibUtilities::PtsFieldSharedPtr ptsInField,
101  LibUtilities::PtsFieldSharedPtr &ptsOutField);
102 
103  /// Interpolate from a pts field to a pts field
105  const LibUtilities::PtsFieldSharedPtr ptsInField,
106  LibUtilities::PtsFieldSharedPtr &ptsOutField);
107 
108  /// returns the dimension of the Interpolator.
109  /// Should be higher than the dimensions of the interpolated fields
110  LIB_UTILITIES_EXPORT int GetDim() const;
111 
112  /// Returns the filter width
114 
115  /// Returns the coordinate id along which the interpolation should be
116  /// performed
117  LIB_UTILITIES_EXPORT int GetCoordId() const;
118 
119  /// Returns the interpolation method used by this interpolator
121 
122  /// Returns the input field
124 
125  /// Returns the output field
127 
128  /// Returns if the weights have already been computed
130 
131  /// sets a callback funtion which gets called every time the interpolation
132  /// progresses
133  template <typename FuncPointerT, typename ObjectPointerT>
134  void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
135  {
136  m_progressCallback = std::bind(
137  func, obj, std::placeholders::_1, std::placeholders::_2);
138  }
139 
140 protected:
141 
142  /// input field
144  /// output field
146 
147  std::function<void(const int position, const int goal)>
149 
150 private:
151 
152  class PtsPoint
153  {
154  public:
155  int idx;
158 
159  PtsPoint() : idx(-1), coords(Array<OneD, NekDouble>(3)), dist(1E30){};
160 
162  : idx(idx), coords(coords), dist(dist){};
163 
164  bool operator<(const PtsPoint &comp) const
165  {
166  return (dist < comp.dist);
167  };
168  };
169 
170  /// dimension of this interpolator. Hardcoded to 3
171  static const int m_dim = 3;
172  typedef boost::geometry::model::point<NekDouble, m_dim, boost::geometry::cs::cartesian> BPoint;
173  typedef std::pair<BPoint, unsigned int> PtsPointPair;
174  typedef boost::geometry::index::rtree<PtsPointPair, boost::geometry::index::rstar<16> > PtsRtree;
175 
176 
177 
178  /// Interpolation Method
180  /// A tree structure to speed up the neighbour search.
181  /// Note that we fill it with an iterator, so instead of rstar, the
182  /// packing algorithm is used.
183  std::shared_ptr<PtsRtree> m_rtree;
184  /// Interpolation weights for each neighbour.
185  /// Structure: m_weights[physPtIdx][neighbourIdx]
187  /// Indices of the relevant neighbours for each physical point.
188  /// Structure: m_neighInds[ptIdx][neighbourIdx]
190  /// Filter width used for some interpolation algorithms
192  /// Max number of interpolation points
193  int m_maxPts;
194  /// coordinate id along which the interpolation should be performed
195  short int m_coordId;
196 
197  LIB_UTILITIES_EXPORT void CalcW_Gauss(const PtsPoint &searchPt,
198  const NekDouble sigma,
199  const int maxPts = 250);
200 
201  LIB_UTILITIES_EXPORT void CalcW_Linear(const PtsPoint &searchPt, int coordId);
202 
203  LIB_UTILITIES_EXPORT void CalcW_NNeighbour(const PtsPoint &searchPt);
204 
205  LIB_UTILITIES_EXPORT void CalcW_Shepard(const PtsPoint &searchPt, int numPts);
206 
207  LIB_UTILITIES_EXPORT void CalcW_Quadratic(const PtsPoint &searchPt,
208  int coordId);
209 
211 
212  LIB_UTILITIES_EXPORT void FindNeighbours(const PtsPoint &searchPt,
213  std::vector<PtsPoint> &neighbourPts,
214  const NekDouble dist,
215  const unsigned int maxPts = 1);
216 
217  LIB_UTILITIES_EXPORT void FindNNeighbours(const PtsPoint &searchPt,
218  std::vector<PtsPoint> &neighbourPts,
219  const unsigned int numPts = 1);
220 };
221 
222 typedef std::shared_ptr<Interpolator> InterpolatorSharedPtr;
223 
224 }
225 }
226 
227 #endif
LibUtilities::PtsFieldSharedPtr m_ptsOutField
output field
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
Interpolator(InterpMethod method=eNoMethod, short int coordId=-1, NekDouble filtWidth=0.0, int maxPts=1000)
Constructor of the Interpolator class.
void CalcW_Gauss(const PtsPoint &searchPt, const NekDouble sigma, const int maxPts=250)
Computes interpolation weights using gaussian interpolation.
boost::geometry::index::rtree< PtsPointPair, boost::geometry::index::rstar< 16 > > PtsRtree
NekDouble GetFiltWidth() const
Returns the filter width.
A class that contains algorithms for interpolation between pts fields, expansions and different meshe...
void SetProgressCallback(FuncPointerT func, ObjectPointerT obj)
sets a callback funtion which gets called every time the interpolation progresses ...
void FindNeighbours(const PtsPoint &searchPt, std::vector< PtsPoint > &neighbourPts, const NekDouble dist, const unsigned int maxPts=1)
Finds the neares neighbours of a point.
void CalcWeights(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
Compute interpolation weights without doing any interpolation.
PtsPoint(int idx, Array< OneD, NekDouble > coords, NekDouble dist)
LibUtilities::PtsFieldSharedPtr GetInField() const
Returns the input field.
void CalcW_Linear(const PtsPoint &searchPt, int coordId)
Computes interpolation weights using linear interpolation.
int m_maxPts
Max number of interpolation points.
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].
boost::geometry::model::point< NekDouble, m_dim, boost::geometry::cs::cartesian > BPoint
void CalcW_Shepard(const PtsPoint &searchPt, int numPts)
Computes interpolation weights using linear interpolation.
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.
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
void Interpolate(const LibUtilities::PtsFieldSharedPtr ptsInField, LibUtilities::PtsFieldSharedPtr &ptsOutField)
Interpolate from a pts field to a pts field.
std::shared_ptr< PtsField > PtsFieldSharedPtr
Definition: PtsField.h:179
#define LIB_UTILITIES_EXPORT
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
std::shared_ptr< Interpolator > InterpolatorSharedPtr
InterpMethod GetInterpMethod() const
Returns the interpolation method used by this interpolator.
LibUtilities::PtsFieldSharedPtr GetOutField() const
Returns the output field.
NekDouble m_filtWidth
Filter width used for some interpolation algorithms.
int GetDim() const
returns the dimension of the Interpolator. Should be higher than the dimensions of the interpolated f...
int GetCoordId() const
Returns the coordinate id along which the interpolation should be performed.
void PrintStatistics()
Returns if the weights have already been computed.