Nektar++
Public Member Functions | Private Member Functions | Private Attributes | List of all members
Nektar::LibUtilities::PtsField Class Reference

#include <PtsField.h>

Collaboration diagram for Nektar::LibUtilities::PtsField:
Collaboration graph
[legend]

Public Member Functions

 PtsField (const int dim, const Array< OneD, Array< OneD, NekDouble > > &pts)
 
 PtsField (const int dim, const vector< std::string > fieldnames, const Array< OneD, Array< OneD, NekDouble > > &pts)
 
 PtsField (const int dim, const vector< std::string > fieldnames, const Array< OneD, Array< OneD, NekDouble > > &pts, const Array< OneD, Array< OneD, float > > &weights, const Array< OneD, Array< OneD, unsigned int > > &neighInds)
 
void CalcWeights (const Array< OneD, Array< OneD, NekDouble > > &physCoords, short int coordId=-1)
 Compute the weights for an interpolation of the field values to physical points. More...
 
void Interpolate (const Array< OneD, Array< OneD, NekDouble > > &physCoords, Array< OneD, Array< OneD, NekDouble > > &intFields, short int coordId=-1)
 Compute weights and perform the interpolate of field values to physical points. More...
 
void Interpolate (Array< OneD, Array< OneD, NekDouble > > &intFields)
 Perform the interpolate of field values to physical points. More...
 
void SetWeights (const Array< OneD, Array< OneD, float > > &weights, const Array< OneD, Array< OneD, unsigned int > > &neighbourInds)
 Set the interpolation weights for an interpolation. More...
 
void GetWeights (Array< OneD, Array< OneD, float > > &weights, Array< OneD, Array< OneD, unsigned int > > &neighbourInds) const
 Get the interpolation weights and corresponding neighbour indices. More...
 
void GetConnectivity (vector< Array< OneD, int > > &conn) const
 Set the connectivity data for ePtsTetBlock and ePtsTriBlock. More...
 
void SetConnectivity (const vector< Array< OneD, int > > &conn)
 Get the connectivity data for ePtsTetBlock and ePtsTriBlock. More...
 
void SetDim (const int ptsDim)
 
int GetDim () const
 
int GetNFields () const
 
vector< std::string > GetFieldNames () const
 
std::string GetFieldName (const int i) const
 
void SetFieldNames (const vector< std::string > fieldNames)
 
void AddField (const Array< OneD, NekDouble > &pts, const std::string fieldName)
 
int GetNpoints () const
 
NekDouble GetPointVal (const int fieldInd, const int ptInd) const
 
void GetPts (Array< OneD, Array< OneD, NekDouble > > &pts) const
 
Array< OneD, NekDoubleGetPts (const int fieldInd) const
 
void SetPts (Array< OneD, Array< OneD, NekDouble > > &pts)
 
vector< int > GetPointsPerEdge () const
 
int GetPointsPerEdge (const int i) const
 
void SetPointsPerEdge (const vector< int > nPtsPerEdge)
 Set the number of points per edge. More...
 
PtsType GetPtsType () const
 
void SetPtsType (const PtsType type)
 
template<typename FuncPointerT , typename ObjectPointerT >
void setProgressCallback (FuncPointerT func, ObjectPointerT obj)
 

Private Member Functions

void CalcW_Linear (const int physPtIdx, const NekDouble coord)
 Compute interpolation weights for a 1D physical point using linear interpolation. More...
 
void CalcW_Shepard (const int physPtIdx, const Array< OneD, NekDouble > &physPtCoords)
 Compute interpolation weights for a physical point using a modified Shepard algorithm. More...
 
void CalcW_Quadratic (const int physPtIdx, const NekDouble coord)
 Compute interpolation weights for a 1D physical point using quadratic interpolation. More...
 
NekDouble DistSq (const Array< OneD, NekDouble > &point1, const Array< OneD, NekDouble > &point2) const
 Compute the square of the euclidean distance between point1 and point2. More...
 
void FindNeighbours (const Array< OneD, NekDouble > &physPtCoords, vector< PtsPoint > &neighbourPts, const unsigned int numPts=1)
 Find nearest neighbours using a brute-force "algorithm". More...
 

Private Attributes

int m_dim
 Dimension of the pts field. More...
 
vector< std::string > m_fieldNames
 Names of the field variables. More...
 
Array< OneD, Array< OneD, NekDouble > > m_pts
 Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx]. More...
 
vector< int > m_nPtsPerEdge
 Number of points per edge. Empty if the point data has no specific shape (ePtsLine) or is a block (ePtsTetBlock, ePtsTriBlock), size=1 for ePtsLine and 2 for a ePtsPlane. More...
 
vector< Array< OneD, int > > m_ptsConn
 Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each, m_ptsConn is a vector of n arrays with 3*m (ePtsTriBlock) or 4*m (ePtsTetBlock) entries. More...
 
PtsType m_ptsType
 Type of the PtsField. More...
 
Array< OneD, Array< OneD, float > > m_weights
 Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx]. More...
 
Array< OneD, Array< OneD, unsigned int > > m_neighInds
 Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourIdx]. More...
 
boost::function< void(const int position, const int goal)> m_progressCallback
 

Detailed Description

Definition at line 90 of file PtsField.h.

Constructor & Destructor Documentation

Nektar::LibUtilities::PtsField::PtsField ( const int  dim,
const Array< OneD, Array< OneD, NekDouble > > &  pts 
)
inline

Definition at line 94 of file PtsField.h.

96  :
97  m_dim(dim),
98  m_pts(pts),
100  {
101  };
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
Nektar::LibUtilities::PtsField::PtsField ( const int  dim,
const vector< std::string >  fieldnames,
const Array< OneD, Array< OneD, NekDouble > > &  pts 
)
inline

Definition at line 103 of file PtsField.h.

106  :
107  m_dim(dim),
108  m_fieldNames(fieldnames),
109  m_pts(pts),
111  {
112  };
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
Nektar::LibUtilities::PtsField::PtsField ( const int  dim,
const vector< std::string >  fieldnames,
const Array< OneD, Array< OneD, NekDouble > > &  pts,
const Array< OneD, Array< OneD, float > > &  weights,
const Array< OneD, Array< OneD, unsigned int > > &  neighInds 
)
inline

Definition at line 114 of file PtsField.h.

119  :
120  m_dim(dim),
121  m_fieldNames(fieldnames),
122  m_pts(pts),
124  m_weights(weights),
125  m_neighInds(neighInds)
126  {
127  };
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204

Member Function Documentation

void Nektar::LibUtilities::PtsField::AddField ( const Array< OneD, NekDouble > &  pts,
const std::string  fieldName 
)

Definition at line 267 of file PtsField.cpp.

References ASSERTL1, m_fieldNames, and m_pts.

269 {
270  int nTotvars = m_pts.num_elements();
271  int nPts = m_pts[0].num_elements();
272 
273  ASSERTL1(pts.num_elements() == nPts, "Field size mismatch");
274 
275  // redirect existing pts
276  Array<OneD, Array<OneD, NekDouble> > newpts(nTotvars + 1);
277  for (int i = 0; i < nTotvars; ++i)
278  {
279  newpts[i] = m_pts[i];
280  }
281  newpts[nTotvars] = pts;
282 
283  m_pts = newpts;
284 
285  m_fieldNames.push_back(fieldName);
286 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LibUtilities::PtsField::CalcW_Linear ( const int  physPtIdx,
const NekDouble  coord 
)
private

Compute interpolation weights for a 1D physical point using linear interpolation.

Parameters
physPtIdxThe index of the physical point in its storage array
coordThe coordinate of the physical point

Definition at line 377 of file PtsField.cpp.

References ASSERTL0, m_neighInds, m_pts, m_weights, and npts.

Referenced by CalcWeights().

378 {
379  int npts = m_pts[0].num_elements();
380  int i;
381 
382  int numPts = 2;
383  m_neighInds[physPtIdx] = Array<OneD, unsigned int> (numPts);
384  m_weights[physPtIdx] = Array<OneD, float> (numPts, 0.0);
385 
386  for (i = 0; i < npts - 1; ++i)
387  {
388  if ((m_pts[0][i] <= coord) && (coord <= m_pts[0][i + 1]))
389  {
390  NekDouble pdiff = m_pts[0][i + 1] - m_pts[0][i];
391 
392  m_neighInds[physPtIdx][0] = i;
393  m_neighInds[physPtIdx][1] = i + 1;
394 
395  m_weights[physPtIdx][0] = (m_pts[0][i + 1] - coord) / pdiff;
396  m_weights[physPtIdx][1] = (coord - m_pts[0][i]) / pdiff;
397 
398  break;
399  }
400  }
401  ASSERTL0(i != npts - 1, "Failed to find coordinate " +
402  boost::lexical_cast<string>(coord) +
403  " within provided input points");
404 };
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
void Nektar::LibUtilities::PtsField::CalcW_Quadratic ( const int  physPtIdx,
const NekDouble  coord 
)
private

Compute interpolation weights for a 1D physical point using quadratic interpolation.

Parameters
physPtIdxThe index of the physical point in its storage array
coordThe coordinate of the physical point

Definition at line 474 of file PtsField.cpp.

References ASSERTL0, m_neighInds, m_pts, m_weights, and npts.

Referenced by CalcWeights().

475 {
476  int npts = m_pts[0].num_elements();
477  int i;
478 
479  int numPts = 3;
480  m_neighInds[physPtIdx] = Array<OneD, unsigned int> (numPts);
481  m_weights[physPtIdx] = Array<OneD, float> (numPts, 0.0);
482 
483  for (i = 0; i < npts - 1; ++i)
484  {
485  if ((m_pts[0][i] <= coord) && (coord <= m_pts[0][i + 1]))
486  {
487  NekDouble pdiff = m_pts[0][i + 1] - m_pts[0][i];
488  NekDouble h1, h2, h3;
489 
490  if (i < npts - 2)
491  {
492  // forwards stencil
493  NekDouble pdiff2 = m_pts[0][i + 2] - m_pts[0][i + 1];
494 
495  h1 = (m_pts[0][i + 1] - coord)
496  * (m_pts[0][i + 2] - coord)
497  / (pdiff * (pdiff + pdiff2));
498  h2 = (coord - m_pts[0][i])
499  * (m_pts[0][i + 2] - coord)
500  / (pdiff * pdiff2);
501  h3 = (coord - m_pts[0][i])
502  * (coord - m_pts[0][i + 1])
503  / ((pdiff + pdiff2) * pdiff2);
504 
505  m_neighInds[physPtIdx][0] = i;
506  m_neighInds[physPtIdx][1] = i + 1;
507  m_neighInds[physPtIdx][2] = i + 2;
508  }
509  else
510  {
511  // backwards stencil
512  NekDouble pdiff2 = m_pts[0][i] - m_pts[0][i - 1];
513 
514  h1 = (m_pts[0][i + 1] - coord)
515  * (coord - m_pts[0][i - 1])
516  / (pdiff * pdiff2);
517  h2 = (coord - m_pts[0][i])
518  * (coord - m_pts[0][i - 1])
519  / (pdiff * (pdiff + pdiff2));
520  h3 = (m_pts[0][i] - coord)
521  * (m_pts[0][i + 1] - coord)
522  / ((pdiff + pdiff2) * pdiff);
523 
524  m_neighInds[physPtIdx][0] = i;
525  m_neighInds[physPtIdx][1] = i + 1;
526  m_neighInds[physPtIdx][2] = i - 1;
527  }
528 
529 
530  m_weights[physPtIdx][0] = h1;
531  m_weights[physPtIdx][1] = h2;
532  m_weights[physPtIdx][2] = h3;
533 
534  break;
535  }
536  }
537  ASSERTL0(i != npts - 1, "Failed to find coordinate " +
538  boost::lexical_cast<string>(coord) +
539  " within provided input points");
540 };
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
void Nektar::LibUtilities::PtsField::CalcW_Shepard ( const int  physPtIdx,
const Array< OneD, NekDouble > &  physPt 
)
private

Compute interpolation weights for a physical point using a modified Shepard algorithm.

Parameters
physPtIdxThe index of the physical point in its storage array
physPtThe coordinates of the physical point

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 422 of file PtsField.cpp.

References ASSERTL0, FindNeighbours(), Nektar::NekConstants::kNekSqrtTol, m_dim, m_neighInds, m_pts, m_weights, and Vmath::Nnan().

Referenced by CalcWeights().

424 {
425  // find nearest neighbours
426  vector<PtsPoint> neighbourPts;
427  int numPts = pow(float(2), m_dim);
428  numPts = min(numPts, int(m_pts[0].num_elements() / 2));
429  FindNeighbours(physPt, neighbourPts, numPts);
430 
431  m_neighInds[physPtIdx] = Array<OneD, unsigned int> (numPts);
432  for (int i = 0; i < numPts; i++)
433  {
434  m_neighInds[physPtIdx][i] = neighbourPts.at(i).m_idx;
435  }
436 
437  m_weights[physPtIdx] = Array<OneD, float> (numPts, 0.0);
438 
439  // In case d < kVertexTheSameDouble ( d^2 < kNekSqrtTol), use the exact
440  // point and return
441  for (int i = 0; i < numPts; ++i)
442  {
443  if (neighbourPts[i].m_distSq <= NekConstants::kNekSqrtTol)
444  {
445  m_weights[physPtIdx][i] = 1.0;
446  return;
447  }
448  }
449 
450  NekDouble wSum = 0.0;
451  for (int i = 0; i < numPts; ++i)
452  {
453  m_weights[physPtIdx][i] = 1 / pow(double(neighbourPts[i].m_distSq),
454  double(m_dim / float(2)));
455  wSum += m_weights[physPtIdx][i];
456  }
457 
458  for (int i = 0; i < numPts; ++i)
459  {
460  m_weights[physPtIdx][i] = m_weights[physPtIdx][i] / wSum;
461  }
462 
463  ASSERTL0(Vmath::Nnan(numPts, m_weights[physPtIdx], 1) == 0, "NaN found in weights");
464 
465 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
static const NekDouble kNekSqrtTol
void FindNeighbours(const Array< OneD, NekDouble > &physPtCoords, vector< PtsPoint > &neighbourPts, const unsigned int numPts=1)
Find nearest neighbours using a brute-force "algorithm".
Definition: PtsField.cpp:576
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.cpp:869
double NekDouble
void Nektar::LibUtilities::PtsField::CalcWeights ( const Array< OneD, Array< OneD, NekDouble > > &  physCoords,
short int  coordId = -1 
)

Compute the weights for an interpolation of the field values to physical points.

Parameters
physCoordscoordinates of the physical points
coord_idid of the coordinate to use for interpolation.

Set coord_id to -1 to use n-D interpolation for an n-dimensional field. The most suitable algorithm is chosen automatically.

Definition at line 53 of file PtsField.cpp.

References ASSERTL1, CalcW_Linear(), CalcW_Quadratic(), CalcW_Shepard(), m_dim, m_neighInds, m_progressCallback, m_pts, and m_weights.

Referenced by Interpolate().

56 {
57  ASSERTL1(physCoords.num_elements() >= m_dim,
58  "physCoords is smaller than number of dimesnions");
59 
60  int nPhysPts = physCoords[0].num_elements();
61  int lastProg = 0;
62 
63  m_weights = Array<OneD, Array<OneD, float> >(nPhysPts);
64  m_neighInds = Array<OneD, Array<OneD, unsigned int> >(nPhysPts);
65 
66  // interpolate points and transform
67  for (int i = 0; i < nPhysPts; ++i)
68  {
69  Array<OneD, NekDouble> physPt(m_dim);
70  for (int j = 0; j < m_dim; ++j)
71  {
72  physPt[j] = physCoords[j][i];
73  }
74 
75  if (m_dim == 1 || coordId >= 0)
76  {
77  if (m_dim == 1)
78  {
79  coordId = 0;
80  }
81 
82  if (m_pts[0].num_elements() <= 2)
83  {
84  CalcW_Linear(i, physPt[coordId]);
85  }
86  else
87  {
88  CalcW_Quadratic(i, physPt[coordId]);
89  }
90  }
91  else
92  {
93  CalcW_Shepard(i, physPt);
94  }
95 
96  int progress = int(100 * i / nPhysPts);
97  if (m_progressCallback && progress > lastProg)
98  {
99  m_progressCallback(i, nPhysPts);
100  lastProg = progress;
101  }
102  }
103 }
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
void CalcW_Quadratic(const int physPtIdx, const NekDouble coord)
Compute interpolation weights for a 1D physical point using quadratic interpolation.
Definition: PtsField.cpp:474
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
void CalcW_Shepard(const int physPtIdx, const Array< OneD, NekDouble > &physPtCoords)
Compute interpolation weights for a physical point using a modified Shepard algorithm.
Definition: PtsField.cpp:422
void CalcW_Linear(const int physPtIdx, const NekDouble coord)
Compute interpolation weights for a 1D physical point using linear interpolation. ...
Definition: PtsField.cpp:377
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
boost::function< void(const int position, const int goal)> m_progressCallback
Definition: PtsField.h:225
NekDouble Nektar::LibUtilities::PtsField::DistSq ( const Array< OneD, NekDouble > &  point1,
const Array< OneD, NekDouble > &  point2 
) const
private

Compute the square of the euclidean distance between point1 and point2.

Parameters
point1The first point
point2The second point

Definition at line 549 of file PtsField.cpp.

Referenced by FindNeighbours().

551 {
552  NekDouble d = 0.0;
553  NekDouble tmp;
554  for (int i = 0; i < point1.num_elements(); i++)
555  {
556  tmp = point1[i] - point2[i];
557  d += tmp * tmp;
558  }
559 
560  return d;
561 }
double NekDouble
void Nektar::LibUtilities::PtsField::FindNeighbours ( const Array< OneD, NekDouble > &  physPt,
vector< PtsPoint > &  neighbourPts,
const unsigned int  numPts = 1 
)
private

Find nearest neighbours using a brute-force "algorithm".

Parameters
physPtCoordinates of the physical point its neighbours we are looking for
neighbourPtsThe points we found
numPtsThe number of points to find

This iterates over all points, computes the (squared) euclidean distance and chooses the numPts closest points. Thus, its very expensive and inefficient.

Definition at line 576 of file PtsField.cpp.

References DistSq(), m_dim, m_pts, and npts.

Referenced by CalcW_Shepard().

579 {
580  int npts = m_pts[0].num_elements();
581 
582  // generate an initial set of intPts
583  for (int i = 0; i < numPts; ++i)
584  {
585  PtsPoint intPt = PtsPoint(-1, Array<OneD, NekDouble>(m_dim), 1E30);
586  neighbourPts.push_back(intPt);
587  }
588 
589  // generate and iterate over all intPts
590  for (int i = 0; i < npts; ++i)
591  {
592  Array<OneD, NekDouble> coords(m_dim);
593  for (int j = 0; j < m_dim; ++j)
594  {
595  coords[j] = m_pts[j][i];
596  }
597  NekDouble d = DistSq(physPt, coords);
598 
599  if (d < neighbourPts.back().m_distSq)
600  {
601  // create new point struct
602  PtsPoint intPt = PtsPoint(i, coords, d);
603 
604  // add it to list, sort the list and remove last point from the sorted
605  // list
606  neighbourPts.push_back(intPt);
607  sort(neighbourPts.begin(), neighbourPts.end());
608  neighbourPts.pop_back();
609  }
610  }
611 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
NekDouble DistSq(const Array< OneD, NekDouble > &point1, const Array< OneD, NekDouble > &point2) const
Compute the square of the euclidean distance between point1 and point2.
Definition: PtsField.cpp:549
void Nektar::LibUtilities::PtsField::GetConnectivity ( vector< Array< OneD, int > > &  conn) const

Set the connectivity data for ePtsTetBlock and ePtsTriBlock.

Parameters
connConnectivity data Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each, m_ptsConn is a vector of n arrays with 3*m (ePtsTriBlock) or 4*m (ePtsTetBlock) entries.

Definition at line 206 of file PtsField.cpp.

References m_ptsConn.

207 {
208  conn = m_ptsConn;
209 }
vector< Array< OneD, int > > m_ptsConn
Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each...
Definition: PtsField.h:215
int Nektar::LibUtilities::PtsField::GetDim ( ) const

Definition at line 234 of file PtsField.cpp.

References m_dim.

235 {
236  return m_dim;
237 }
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
std::string Nektar::LibUtilities::PtsField::GetFieldName ( const int  i) const

Definition at line 252 of file PtsField.cpp.

References m_fieldNames.

253 {
254  return m_fieldNames[i];
255 }
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
vector< std::string > Nektar::LibUtilities::PtsField::GetFieldNames ( ) const

Definition at line 246 of file PtsField.cpp.

References m_fieldNames.

247 {
248  return m_fieldNames;
249 }
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
int Nektar::LibUtilities::PtsField::GetNFields ( ) const

Definition at line 240 of file PtsField.cpp.

References m_fieldNames.

241 {
242  return m_fieldNames.size();
243 }
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
int Nektar::LibUtilities::PtsField::GetNpoints ( ) const

Definition at line 289 of file PtsField.cpp.

References m_pts.

290 {
291  return m_pts[0].num_elements();
292 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
vector< int > Nektar::LibUtilities::PtsField::GetPointsPerEdge ( ) const

Definition at line 322 of file PtsField.cpp.

References m_nPtsPerEdge.

323 {
324  return m_nPtsPerEdge;
325 }
vector< int > m_nPtsPerEdge
Number of points per edge. Empty if the point data has no specific shape (ePtsLine) or is a block (eP...
Definition: PtsField.h:211
int Nektar::LibUtilities::PtsField::GetPointsPerEdge ( const int  i) const

Definition at line 328 of file PtsField.cpp.

References m_nPtsPerEdge.

329 {
330  return m_nPtsPerEdge[i];
331 }
vector< int > m_nPtsPerEdge
Number of points per edge. Empty if the point data has no specific shape (ePtsLine) or is a block (eP...
Definition: PtsField.h:211
NekDouble Nektar::LibUtilities::PtsField::GetPointVal ( const int  fieldInd,
const int  ptInd 
) const

Definition at line 295 of file PtsField.cpp.

References m_pts.

296 {
297  return m_pts[fieldInd][ptInd];
298 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
void Nektar::LibUtilities::PtsField::GetPts ( Array< OneD, Array< OneD, NekDouble > > &  pts) const

Definition at line 301 of file PtsField.cpp.

References m_pts.

302 {
303  pts = m_pts;
304 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
Array< OneD, NekDouble > Nektar::LibUtilities::PtsField::GetPts ( const int  fieldInd) const

Definition at line 307 of file PtsField.cpp.

References m_pts.

308 {
309  return m_pts[fieldInd];
310 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
PtsType Nektar::LibUtilities::PtsField::GetPtsType ( ) const

Definition at line 358 of file PtsField.cpp.

References m_ptsType.

359 {
360  return m_ptsType;
361 }
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
void Nektar::LibUtilities::PtsField::GetWeights ( Array< OneD, Array< OneD, float > > &  weights,
Array< OneD, Array< OneD, unsigned int > > &  neighbourInds 
) const

Get the interpolation weights and corresponding neighbour indices.

Parameters
weightsInterpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx]
neighbourIndsIndices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourIdx]

Definition at line 190 of file PtsField.cpp.

References m_neighInds, and m_weights.

193 {
194  weights = m_weights;
195  neighbourInds = m_neighInds;
196 }
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
void Nektar::LibUtilities::PtsField::Interpolate ( const Array< OneD, Array< OneD, NekDouble > > &  physCoords,
Array< OneD, Array< OneD, NekDouble > > &  intFields,
short int  coordId = -1 
)

Compute weights and perform the interpolate of field values to physical points.

Parameters
physCoordscoordinates of the physical points
intFieldsinterpolated field at the physical points
coord_idid of the coordinate to use for interpolation.

Set coord_id to -1 to use n-D interpolation for an n-dimensional field. The most suitable algorithm is chosen automatically.

Definition at line 116 of file PtsField.cpp.

References CalcWeights().

120 {
121  CalcWeights(physCoords, coordId);
122  Interpolate(intFields);
123 }
void Interpolate(const Array< OneD, Array< OneD, NekDouble > > &physCoords, Array< OneD, Array< OneD, NekDouble > > &intFields, short int coordId=-1)
Compute weights and perform the interpolate of field values to physical points.
Definition: PtsField.cpp:116
void CalcWeights(const Array< OneD, Array< OneD, NekDouble > > &physCoords, short int coordId=-1)
Compute the weights for an interpolation of the field values to physical points.
Definition: PtsField.cpp:53
void Nektar::LibUtilities::PtsField::Interpolate ( Array< OneD, Array< OneD, NekDouble > > &  intFields)

Perform the interpolate of field values to physical points.

Parameters
intFieldsinterpolated field at the physical points

The weights must have already been computed by or set by .

Definition at line 134 of file PtsField.cpp.

References ASSERTL1, m_dim, m_fieldNames, m_neighInds, m_pts, and m_weights.

135 {
136  ASSERTL1(m_weights[0].num_elements() == m_neighInds[0].num_elements(),
137  "weights / neighInds mismatch")
138  int nFields = m_fieldNames.size();
139  int nPhysPts = m_weights.num_elements();
140 
141  // interpolate points and transform
142  intFields = Array<OneD, Array<OneD, NekDouble> >(nFields);
143  for (int i = 0; i < nFields; ++i)
144  {
145  intFields[i] = Array<OneD, NekDouble>(nPhysPts);
146 
147  for (int j = 0; j < nPhysPts; ++j)
148  {
149  intFields[i][j] = 0.0;
150 
151  int nPts = m_weights[j].num_elements();
152  for (int k = 0; k < nPts; ++k)
153  {
154  unsigned int nIdx = m_neighInds[j][k];
155  intFields[i][j] += m_weights[j][k] * m_pts[m_dim + i][nIdx];
156  }
157  }
158  }
159 }
for(iterator iter=lhs.begin();iter!=lhs.end();++iter)
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
double NekDouble
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LibUtilities::PtsField::SetConnectivity ( const vector< Array< OneD, int > > &  conn)

Get the connectivity data for ePtsTetBlock and ePtsTriBlock.

Parameters
connConnectivity data Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each, m_ptsConn is a vector of n arrays with 3*m (ePtsTriBlock) or 4*m (ePtsTetBlock) entries.

Definition at line 219 of file PtsField.cpp.

References ASSERTL1, Nektar::LibUtilities::ePtsTetBlock, Nektar::LibUtilities::ePtsTriBlock, m_ptsConn, and m_ptsType.

220 {
222  "ptsType must be set before connectivity");
223 
224  m_ptsConn = conn;
225 }
vector< Array< OneD, int > > m_ptsConn
Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each...
Definition: PtsField.h:215
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LibUtilities::PtsField::SetDim ( const int  ptsDim)

Definition at line 228 of file PtsField.cpp.

References m_dim.

229 {
230  m_dim = ptsDim;
231 }
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
void Nektar::LibUtilities::PtsField::SetFieldNames ( const vector< std::string >  fieldNames)

Definition at line 258 of file PtsField.cpp.

References ASSERTL0, m_dim, m_fieldNames, and m_pts.

259 {
260  ASSERTL0(fieldNames.size() == m_pts.num_elements() - m_dim,
261  "Number of given fieldNames does not match the number of stored fields");
262 
263  m_fieldNames = fieldNames;
264 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
int m_dim
Dimension of the pts field.
Definition: PtsField.h:202
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:204
void Nektar::LibUtilities::PtsField::SetPointsPerEdge ( const vector< int >  nPtsPerEdge)

Set the number of points per edge.

Parameters
nPtsPerEdgeNumber of points per edge. Empty if the point data has no specific shape (ePtsLine) or is a block (ePtsTetBlock, ePtsTriBlock), size=1 for ePtsLine and 2 for a ePtsPlane

Definition at line 340 of file PtsField.cpp.

References ASSERTL0, Nektar::LibUtilities::ePtsLine, Nektar::LibUtilities::ePtsPlane, m_nPtsPerEdge, m_pts, and m_ptsType.

341 {
343  "SetPointsPerEdge only supported for ePtsLine and ePtsPlane .");
344 
345  int totPts(1);
346  for (int i = 0; i < nPtsPerEdge.size(); ++i)
347  {
348  totPts = totPts * nPtsPerEdge.at(i);
349  }
350 
351  ASSERTL0(totPts == m_pts.num_elements(),
352  "nPtsPerEdge does not match total number of points");
353 
354  m_nPtsPerEdge = nPtsPerEdge;
355 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
vector< int > m_nPtsPerEdge
Number of points per edge. Empty if the point data has no specific shape (ePtsLine) or is a block (eP...
Definition: PtsField.h:211
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::LibUtilities::PtsField::setProgressCallback ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 193 of file PtsField.h.

195  {
196  m_progressCallback = boost::bind(func, obj, _1, _2);
197  }
boost::function< void(const int position, const int goal)> m_progressCallback
Definition: PtsField.h:225
void Nektar::LibUtilities::PtsField::SetPts ( Array< OneD, Array< OneD, NekDouble > > &  pts)

Definition at line 313 of file PtsField.cpp.

References ASSERTL1, and m_pts.

314 {
315  ASSERTL1(pts.num_elements() == m_pts.num_elements(),
316  "Pts field count mismatch");
317 
318  m_pts = pts;
319 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].
Definition: PtsField.h:207
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
void Nektar::LibUtilities::PtsField::SetPtsType ( const PtsType  type)

Definition at line 364 of file PtsField.cpp.

References m_ptsType.

365 {
366  m_ptsType = type;
367 }
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:217
void Nektar::LibUtilities::PtsField::SetWeights ( const Array< OneD, Array< OneD, float > > &  weights,
const Array< OneD, Array< OneD, unsigned int > > &  neighbourInds 
)

Set the interpolation weights for an interpolation.

Parameters
weightsInterpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx]
neighbourIndsIndices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourIdx]

Definition at line 170 of file PtsField.cpp.

References ASSERTL0, m_neighInds, and m_weights.

173 {
174  ASSERTL0(weights.num_elements() == neighbourInds.num_elements(),
175  "weights and neighbourInds not of same number of physical points")
176 
177  m_weights = weights;
178  m_neighInds = neighbourInds;
179 
180 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, Array< OneD, unsigned int > > m_neighInds
Indices of the relevant neighbours for each physical point. Structure: m_neighInds[ptIdx][neighbourId...
Definition: PtsField.h:223
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:220

Member Data Documentation

int Nektar::LibUtilities::PtsField::m_dim
private

Dimension of the pts field.

Definition at line 202 of file PtsField.h.

Referenced by CalcW_Shepard(), CalcWeights(), FindNeighbours(), GetDim(), Interpolate(), SetDim(), and SetFieldNames().

vector<std::string> Nektar::LibUtilities::PtsField::m_fieldNames
private

Names of the field variables.

Definition at line 204 of file PtsField.h.

Referenced by AddField(), GetFieldName(), GetFieldNames(), GetNFields(), Interpolate(), and SetFieldNames().

Array<OneD, Array<OneD, unsigned int> > Nektar::LibUtilities::PtsField::m_neighInds
private

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

Definition at line 223 of file PtsField.h.

Referenced by CalcW_Linear(), CalcW_Quadratic(), CalcW_Shepard(), CalcWeights(), GetWeights(), Interpolate(), and SetWeights().

vector<int> Nektar::LibUtilities::PtsField::m_nPtsPerEdge
private

Number of points per edge. Empty if the point data has no specific shape (ePtsLine) or is a block (ePtsTetBlock, ePtsTriBlock), size=1 for ePtsLine and 2 for a ePtsPlane.

Definition at line 211 of file PtsField.h.

Referenced by GetPointsPerEdge(), and SetPointsPerEdge().

boost::function<void (const int position, const int goal)> Nektar::LibUtilities::PtsField::m_progressCallback
private

Definition at line 225 of file PtsField.h.

Referenced by CalcWeights().

Array<OneD, Array<OneD, NekDouble> > Nektar::LibUtilities::PtsField::m_pts
private

Point data. For a n-dimensional field, the first n fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].

Definition at line 207 of file PtsField.h.

Referenced by AddField(), CalcW_Linear(), CalcW_Quadratic(), CalcW_Shepard(), CalcWeights(), FindNeighbours(), GetNpoints(), GetPointVal(), GetPts(), Interpolate(), SetFieldNames(), SetPointsPerEdge(), and SetPts().

vector<Array<OneD, int> > Nektar::LibUtilities::PtsField::m_ptsConn
private

Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each, m_ptsConn is a vector of n arrays with 3*m (ePtsTriBlock) or 4*m (ePtsTetBlock) entries.

Definition at line 215 of file PtsField.h.

Referenced by GetConnectivity(), and SetConnectivity().

PtsType Nektar::LibUtilities::PtsField::m_ptsType
private

Type of the PtsField.

Definition at line 217 of file PtsField.h.

Referenced by GetPtsType(), SetConnectivity(), SetPointsPerEdge(), and SetPtsType().

Array<OneD, Array<OneD, float> > Nektar::LibUtilities::PtsField::m_weights
private

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

Definition at line 220 of file PtsField.h.

Referenced by CalcW_Linear(), CalcW_Quadratic(), CalcW_Shepard(), CalcWeights(), GetWeights(), Interpolate(), and SetWeights().