Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Public Attributes | 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, map< PtsInfo, int > ptsInfo=NullPtsInfoMap)
 
 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)
 
vector< NekDoubleGetBoxSize () const
 
void SetBoxSize (const vector< NekDouble > boxsize)
 
template<typename FuncPointerT , typename ObjectPointerT >
void setProgressCallback (FuncPointerT func, ObjectPointerT obj)
 

Public Attributes

map< PtsInfo, int > m_ptsInfo
 map for information about points that can be added through PtsInfo enum More...
 

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 m_dim 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...
 
vector< NekDoublem_boxSize
 vector of box size xmin,xmax,ymin,ymax,zmin,zmax More...
 
boost::function< void(const
int position, const int goal)> 
m_progressCallback
 

Detailed Description

Definition at line 99 of file PtsField.h.

Constructor & Destructor Documentation

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

Definition at line 198 of file PtsField.cpp.

References GetNFields(), and m_fieldNames.

198  :
199  m_dim(dim),
200  m_pts(pts),
202 {
203  for (int i = 0; i < GetNFields(); ++i)
204  {
205  m_fieldNames.push_back("NA");
206  }
207 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
Nektar::LibUtilities::PtsField::PtsField ( const int  dim,
const vector< std::string >  fieldnames,
const Array< OneD, Array< OneD, NekDouble > > &  pts,
map< PtsInfo, int >  ptsInfo = NullPtsInfoMap 
)
inline

Definition at line 106 of file PtsField.h.

110  :
111  m_ptsInfo(ptsInfo),
112  m_dim(dim),
113  m_fieldNames(fieldnames),
114  m_pts(pts),
116  {
117  };
map< PtsInfo, int > m_ptsInfo
map for information about points that can be added through PtsInfo enum
Definition: PtsField.h:209
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
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 119 of file PtsField.h.

124  :
126  m_dim(dim),
127  m_fieldNames(fieldnames),
128  m_pts(pts),
131  m_neighInds(neighInds)
132  {
133  };
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:235
map< PtsInfo, int > m_ptsInfo
map for information about points that can be added through PtsInfo enum
Definition: PtsField.h:209
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
map< pair< int, int >, NekDouble > weights(set< pair< int, int > > springs, Array< OneD, NekDouble > u, Array< OneD, NekDouble > v)
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
static map< PtsInfo, int > NullPtsInfoMap
Definition: PtsField.h:72

Member Function Documentation

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

Definition at line 280 of file PtsField.cpp.

References ASSERTL1, m_fieldNames, and m_pts.

282 {
283  int nTotvars = m_pts.num_elements();
284 
285  ASSERTL1(pts.num_elements() == m_pts[0].num_elements(),
286  "Field size mismatch");
287 
288  // redirect existing pts
289  Array<OneD, Array<OneD, NekDouble> > newpts(nTotvars + 1);
290  for (int i = 0; i < nTotvars; ++i)
291  {
292  newpts[i] = m_pts[i];
293  }
294  newpts[nTotvars] = pts;
295 
296  m_pts = newpts;
297 
298  m_fieldNames.push_back(fieldName);
299 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
#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 394 of file PtsField.cpp.

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

Referenced by CalcWeights().

395 {
396  int npts = m_pts[0].num_elements();
397  int i;
398 
399  int numPts = 2;
400  m_neighInds[physPtIdx] = Array<OneD, unsigned int> (numPts);
401  m_weights[physPtIdx] = Array<OneD, float> (numPts, 0.0);
402 
403  for (i = 0; i < npts - 1; ++i)
404  {
405  if ((m_pts[0][i] <= coord) && (coord <= m_pts[0][i + 1]))
406  {
407  NekDouble pdiff = m_pts[0][i + 1] - m_pts[0][i];
408 
409  m_neighInds[physPtIdx][0] = i;
410  m_neighInds[physPtIdx][1] = i + 1;
411 
412  m_weights[physPtIdx][0] = (m_pts[0][i + 1] - coord) / pdiff;
413  m_weights[physPtIdx][1] = (coord - m_pts[0][i]) / pdiff;
414 
415  break;
416  }
417  }
418  ASSERTL0(i != npts - 1, "Failed to find coordinate " +
419  boost::lexical_cast<string>(coord) +
420  " within provided input points");
421 };
#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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
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 491 of file PtsField.cpp.

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

Referenced by CalcWeights().

492 {
493  int npts = m_pts[0].num_elements();
494  int i;
495 
496  int numPts = 3;
497  m_neighInds[physPtIdx] = Array<OneD, unsigned int> (numPts);
498  m_weights[physPtIdx] = Array<OneD, float> (numPts, 0.0);
499 
500  for (i = 0; i < npts - 1; ++i)
501  {
502  if ((m_pts[0][i] <= coord) && (coord <= m_pts[0][i + 1]))
503  {
504  NekDouble pdiff = m_pts[0][i + 1] - m_pts[0][i];
505  NekDouble h1, h2, h3;
506 
507  if (i < npts - 2)
508  {
509  // forwards stencil
510  NekDouble pdiff2 = m_pts[0][i + 2] - m_pts[0][i + 1];
511 
512  h1 = (m_pts[0][i + 1] - coord)
513  * (m_pts[0][i + 2] - coord)
514  / (pdiff * (pdiff + pdiff2));
515  h2 = (coord - m_pts[0][i])
516  * (m_pts[0][i + 2] - coord)
517  / (pdiff * pdiff2);
518  h3 = (coord - m_pts[0][i])
519  * (coord - m_pts[0][i + 1])
520  / ((pdiff + pdiff2) * pdiff2);
521 
522  m_neighInds[physPtIdx][0] = i;
523  m_neighInds[physPtIdx][1] = i + 1;
524  m_neighInds[physPtIdx][2] = i + 2;
525  }
526  else
527  {
528  // backwards stencil
529  NekDouble pdiff2 = m_pts[0][i] - m_pts[0][i - 1];
530 
531  h1 = (m_pts[0][i + 1] - coord)
532  * (coord - m_pts[0][i - 1])
533  / (pdiff * pdiff2);
534  h2 = (coord - m_pts[0][i])
535  * (coord - m_pts[0][i - 1])
536  / (pdiff * (pdiff + pdiff2));
537  h3 = (m_pts[0][i] - coord)
538  * (m_pts[0][i + 1] - coord)
539  / ((pdiff + pdiff2) * pdiff);
540 
541  m_neighInds[physPtIdx][0] = i;
542  m_neighInds[physPtIdx][1] = i + 1;
543  m_neighInds[physPtIdx][2] = i - 1;
544  }
545 
546 
547  m_weights[physPtIdx][0] = h1;
548  m_weights[physPtIdx][1] = h2;
549  m_weights[physPtIdx][2] = h3;
550 
551  break;
552  }
553  }
554  ASSERTL0(i != npts - 1, "Failed to find coordinate " +
555  boost::lexical_cast<string>(coord) +
556  " within provided input points");
557 };
#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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
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 439 of file PtsField.cpp.

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

Referenced by CalcWeights().

441 {
442  // find nearest neighbours
443  vector<PtsPoint> neighbourPts;
444  int numPts = pow(float(2), m_dim);
445  numPts = min(numPts, int(m_pts[0].num_elements() / 2));
446  FindNeighbours(physPt, neighbourPts, numPts);
447 
448  m_neighInds[physPtIdx] = Array<OneD, unsigned int> (numPts);
449  for (int i = 0; i < numPts; i++)
450  {
451  m_neighInds[physPtIdx][i] = neighbourPts.at(i).m_idx;
452  }
453 
454  m_weights[physPtIdx] = Array<OneD, float> (numPts, 0.0);
455 
456  // In case d < kVertexTheSameDouble ( d^2 < kNekSqrtTol), use the exact
457  // point and return
458  for (int i = 0; i < numPts; ++i)
459  {
460  if (neighbourPts[i].m_distSq <= NekConstants::kNekSqrtTol)
461  {
462  m_weights[physPtIdx][i] = 1.0;
463  return;
464  }
465  }
466 
467  NekDouble wSum = 0.0;
468  for (int i = 0; i < numPts; ++i)
469  {
470  m_weights[physPtIdx][i] = 1 / pow(double(neighbourPts[i].m_distSq),
471  double(m_dim / float(2)));
472  wSum += m_weights[physPtIdx][i];
473  }
474 
475  for (int i = 0; i < numPts; ++i)
476  {
477  m_weights[physPtIdx][i] = m_weights[physPtIdx][i] / wSum;
478  }
479 
480  ASSERTL0(Vmath::Nnan(numPts, m_weights[physPtIdx], 1) == 0, "NaN found in weights");
481 
482 }
#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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
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:593
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
int Nnan(int n, const T *x, const int incx)
Return number of NaN elements of x.
Definition: Vmath.cpp:878
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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
void CalcW_Quadratic(const int physPtIdx, const NekDouble coord)
Compute interpolation weights for a 1D physical point using quadratic interpolation.
Definition: PtsField.cpp:491
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
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:439
void CalcW_Linear(const int physPtIdx, const NekDouble coord)
Compute interpolation weights for a 1D physical point using linear interpolation. ...
Definition: PtsField.cpp:394
#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:240
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 566 of file PtsField.cpp.

Referenced by FindNeighbours().

568 {
569  NekDouble d = 0.0;
570  NekDouble tmp;
571  for (int i = 0; i < point1.num_elements(); i++)
572  {
573  tmp = point1[i] - point2[i];
574  d += tmp * tmp;
575  }
576 
577  return d;
578 }
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 593 of file PtsField.cpp.

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

Referenced by CalcW_Shepard().

596 {
597  int npts = m_pts[0].num_elements();
598 
599  // generate an initial set of intPts
600  for (int i = 0; i < numPts; ++i)
601  {
602  PtsPoint intPt = PtsPoint(-1, Array<OneD, NekDouble>(m_dim), 1E30);
603  neighbourPts.push_back(intPt);
604  }
605 
606  // generate and iterate over all intPts
607  for (int i = 0; i < npts; ++i)
608  {
609  Array<OneD, NekDouble> coords(m_dim);
610  for (int j = 0; j < m_dim; ++j)
611  {
612  coords[j] = m_pts[j][i];
613  }
614  NekDouble d = DistSq(physPt, coords);
615 
616  if (d < neighbourPts.back().m_distSq)
617  {
618  // create new point struct
619  PtsPoint intPt = PtsPoint(i, coords, d);
620 
621  // add it to list, sort the list and remove last point from the sorted
622  // list
623  neighbourPts.push_back(intPt);
624  sort(neighbourPts.begin(), neighbourPts.end());
625  neighbourPts.pop_back();
626  }
627  }
628 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
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:566
vector< NekDouble > Nektar::LibUtilities::PtsField::GetBoxSize ( ) const

Definition at line 374 of file PtsField.cpp.

References m_boxSize.

375 {
376  return m_boxSize;
377 }
vector< NekDouble > m_boxSize
vector of box size xmin,xmax,ymin,ymax,zmin,zmax
Definition: PtsField.h:238
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 219 of file PtsField.cpp.

References m_ptsConn.

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

Definition at line 247 of file PtsField.cpp.

References m_dim.

248 {
249  return m_dim;
250 }
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
std::string Nektar::LibUtilities::PtsField::GetFieldName ( const int  i) const

Definition at line 265 of file PtsField.cpp.

References m_fieldNames.

266 {
267  return m_fieldNames[i];
268 }
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
vector< std::string > Nektar::LibUtilities::PtsField::GetFieldNames ( ) const

Definition at line 259 of file PtsField.cpp.

References m_fieldNames.

260 {
261  return m_fieldNames;
262 }
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
int Nektar::LibUtilities::PtsField::GetNFields ( ) const

Definition at line 253 of file PtsField.cpp.

References m_fieldNames.

Referenced by PtsField().

254 {
255  return m_fieldNames.size();
256 }
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
int Nektar::LibUtilities::PtsField::GetNpoints ( ) const

Definition at line 302 of file PtsField.cpp.

References m_pts.

303 {
304  return m_pts[0].num_elements();
305 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
vector< int > Nektar::LibUtilities::PtsField::GetPointsPerEdge ( ) const

Definition at line 335 of file PtsField.cpp.

References m_nPtsPerEdge.

336 {
337  return m_nPtsPerEdge;
338 }
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:223
int Nektar::LibUtilities::PtsField::GetPointsPerEdge ( const int  i) const

Definition at line 341 of file PtsField.cpp.

References m_nPtsPerEdge.

342 {
343  return m_nPtsPerEdge[i];
344 }
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:223
NekDouble Nektar::LibUtilities::PtsField::GetPointVal ( const int  fieldInd,
const int  ptInd 
) const

Definition at line 308 of file PtsField.cpp.

References m_pts.

309 {
310  return m_pts[fieldInd][ptInd];
311 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
void Nektar::LibUtilities::PtsField::GetPts ( Array< OneD, Array< OneD, NekDouble > > &  pts) const

Definition at line 314 of file PtsField.cpp.

References m_pts.

315 {
316  pts = m_pts;
317 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
Array< OneD, NekDouble > Nektar::LibUtilities::PtsField::GetPts ( const int  fieldInd) const

Definition at line 320 of file PtsField.cpp.

References m_pts.

321 {
322  return m_pts[fieldInd];
323 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
PtsType Nektar::LibUtilities::PtsField::GetPtsType ( ) const

Definition at line 363 of file PtsField.cpp.

References m_ptsType.

364 {
365  return m_ptsType;
366 }
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
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, m_weights, and Nektar::NekMeshUtils::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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
map< pair< int, int >, NekDouble > weights(set< pair< int, int > > springs, Array< OneD, NekDouble > u, Array< OneD, NekDouble > v)
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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
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::SetBoxSize ( const vector< NekDouble boxsize)

Definition at line 379 of file PtsField.cpp.

References m_boxSize.

380 {
381  m_boxSize = boxSize;
382 }
vector< NekDouble > m_boxSize
vector of box size xmin,xmax,ymin,ymax,zmin,zmax
Definition: PtsField.h:238
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 232 of file PtsField.cpp.

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

233 {
235  "ptsType must be set before connectivity");
236 
237  m_ptsConn = conn;
238 }
vector< Array< OneD, int > > m_ptsConn
Connectivity data needed for ePtsTetBlock and ePtsTriBlock. For n Blocks with m elements each...
Definition: PtsField.h:227
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
#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 241 of file PtsField.cpp.

References m_dim.

242 {
243  m_dim = ptsDim;
244 }
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
void Nektar::LibUtilities::PtsField::SetFieldNames ( const vector< std::string >  fieldNames)

Definition at line 271 of file PtsField.cpp.

References ASSERTL0, m_dim, m_fieldNames, and m_pts.

272 {
273  ASSERTL0(fieldNames.size() == m_pts.num_elements() - m_dim,
274  "Number of given fieldNames does not match the number of stored fields");
275 
276  m_fieldNames = fieldNames;
277 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
int m_dim
Dimension of the pts field.
Definition: PtsField.h:214
vector< std::string > m_fieldNames
Names of the field variables.
Definition: PtsField.h:216
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, 2 for ePtsPlane and 3 for ePtsBox

Definition at line 353 of file PtsField.cpp.

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

354 {
356  m_ptsType == ePtsBox,
357  "SetPointsPerEdge only supported for ePtsLine, ePtsPlane and ePtsBox.");
358 
359  m_nPtsPerEdge = nPtsPerEdge;
360 }
#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:223
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
template<typename FuncPointerT , typename ObjectPointerT >
void Nektar::LibUtilities::PtsField::setProgressCallback ( FuncPointerT  func,
ObjectPointerT  obj 
)
inline

Definition at line 202 of file PtsField.h.

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

Definition at line 326 of file PtsField.cpp.

References ASSERTL1, and m_pts.

327 {
328  ASSERTL1(pts.num_elements() == m_pts.num_elements(),
329  "Pts field count mismatch");
330 
331  m_pts = pts;
332 }
Array< OneD, Array< OneD, NekDouble > > m_pts
Point data. For a n-dimensional field, the first m_dim fields are the points spatial coordinates...
Definition: PtsField.h:219
#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 369 of file PtsField.cpp.

References m_ptsType.

370 {
371  m_ptsType = type;
372 }
PtsType m_ptsType
Type of the PtsField.
Definition: PtsField.h:229
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, m_weights, and Nektar::NekMeshUtils::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:235
Array< OneD, Array< OneD, float > > m_weights
Interpolation weights for each neighbour. Structure: m_weights[physPtIdx][neighbourIdx].
Definition: PtsField.h:232
map< pair< int, int >, NekDouble > weights(set< pair< int, int > > springs, Array< OneD, NekDouble > u, Array< OneD, NekDouble > v)

Member Data Documentation

vector<NekDouble> Nektar::LibUtilities::PtsField::m_boxSize
private

vector of box size xmin,xmax,ymin,ymax,zmin,zmax

Definition at line 238 of file PtsField.h.

Referenced by GetBoxSize(), and SetBoxSize().

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

Dimension of the pts field.

Definition at line 214 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 216 of file PtsField.h.

Referenced by AddField(), GetFieldName(), GetFieldNames(), GetNFields(), Interpolate(), PtsField(), 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 235 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 223 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 240 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 m_dim fields are the points spatial coordinates. Structure: m_pts[fieldIdx][ptIdx].

Definition at line 219 of file PtsField.h.

Referenced by AddField(), CalcW_Linear(), CalcW_Quadratic(), CalcW_Shepard(), CalcWeights(), FindNeighbours(), GetNpoints(), GetPointVal(), GetPts(), Interpolate(), SetFieldNames(), 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 227 of file PtsField.h.

Referenced by GetConnectivity(), and SetConnectivity().

map<PtsInfo,int> Nektar::LibUtilities::PtsField::m_ptsInfo

map for information about points that can be added through PtsInfo enum

Definition at line 209 of file PtsField.h.

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

Type of the PtsField.

Definition at line 229 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 232 of file PtsField.h.

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