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

#include <NodalPrismElec.h>

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

Public Member Functions

virtual ~NodalPrismElec ()
 
 NodalPrismElec (const PointsKey &key)
 
- Public Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual ~Points ()
 
void Initialize (void)
 
unsigned int GetPointsDim () const
 
unsigned int GetNumPoints () const
 
unsigned int GetTotNumPoints () const
 
PointsType GetPointsType () const
 
const Array< OneD, const DataType > & GetZ () const
 
const Array< OneD, const DataType > & GetW () const
 
void GetZW (Array< OneD, const DataType > &z, Array< OneD, const DataType > &w) const
 
const Array< OneD, const NekDouble > & GetBaryWeights () const
 
void GetPoints (Array< OneD, const DataType > &x) const
 
void GetPoints (Array< OneD, const DataType > &x, Array< OneD, const DataType > &y) const
 
void GetPoints (Array< OneD, const DataType > &x, Array< OneD, const DataType > &y, Array< OneD, const DataType > &z) const
 
const MatrixSharedPtrTypeGetD (Direction dir=xDir) const
 
const MatrixSharedPtrType GetI (const PointsKey &key)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x)
 
const MatrixSharedPtrType GetI (unsigned int uint, const Array< OneD, const DataType > &x)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y)
 
const MatrixSharedPtrType GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y, const Array< OneD, const DataType > &z)
 
const MatrixSharedPtrType GetGalerkinProjection (const PointsKey &pkey)
 

Static Public Member Functions

static std::shared_ptr< PointsBaseTypeCreate (const PointsKey &key)
 

Protected Member Functions

virtual const MatrixSharedPtrType v_GetI (const PointsKey &pkey) override
 
virtual const MatrixSharedPtrType v_GetI (const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z) override
 
- Protected Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual void v_Initialize (void)
 
virtual void v_CalculateBaryWeights ()
 This function calculates the barycentric weights used for enhanced interpolation speed. More...
 
 Points (const PointsKey &key)
 
virtual const MatrixSharedPtrType v_GetI (const Array< OneD, const DataType > &x)
 
virtual const MatrixSharedPtrType v_GetI (unsigned int, const Array< OneD, const DataType > &x)
 
virtual const MatrixSharedPtrType v_GetI (const Array< OneD, const DataType > &x, const Array< OneD, const DataType > &y)
 
virtual const MatrixSharedPtrType v_GetGalerkinProjection (const PointsKey &pkey)
 

Private Member Functions

 NodalPrismElec ()
 Default constructor should not be called except by Create matrix. More...
 
void NodalPointReorder3d ()
 
virtual void v_CalculatePoints () override
 
virtual void v_CalculateWeights () override
 
virtual void v_CalculateDerivMatrix () override
 
void CalculateInterpMatrix (const Array< OneD, const NekDouble > &xi, const Array< OneD, const NekDouble > &yi, const Array< OneD, const NekDouble > &zi, Array< OneD, NekDouble > &interp)
 

Private Attributes

std::shared_ptr< NodalUtilPrismm_util
 

Static Private Attributes

static bool initPointsManager []
 

Additional Inherited Members

- Public Types inherited from Nektar::LibUtilities::Points< NekDouble >
typedef NekDouble DataType
 
typedef std::shared_ptr< NekMatrix< DataType > > MatrixSharedPtrType
 
- Protected Attributes inherited from Nektar::LibUtilities::Points< NekDouble >
PointsKey m_pointsKey
 Points type for this points distributions. More...
 
Array< OneD, DataTypem_points [3]
 Storage for the point locations, allowing for up to a 3D points storage. More...
 
Array< OneD, DataTypem_weights
 Quadrature weights for the weights. More...
 
Array< OneD, DataTypem_bcweights
 Barycentric weights. More...
 
MatrixSharedPtrType m_derivmatrix [3]
 Derivative matrices. More...
 
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLessm_InterpManager
 
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLessm_GalerkinProjectionManager
 

Detailed Description

Definition at line 51 of file NodalPrismElec.h.

Constructor & Destructor Documentation

◆ ~NodalPrismElec()

virtual Nektar::LibUtilities::NodalPrismElec::~NodalPrismElec ( )
inlinevirtual

Definition at line 54 of file NodalPrismElec.h.

55  {
56  }

◆ NodalPrismElec() [1/2]

Nektar::LibUtilities::NodalPrismElec::NodalPrismElec ( const PointsKey key)
inline

Definition at line 58 of file NodalPrismElec.h.

58  : PointsBaseType(key)
59  {
60  }
Points< NekDouble > PointsBaseType

◆ NodalPrismElec() [2/2]

Nektar::LibUtilities::NodalPrismElec::NodalPrismElec ( )
inlineprivate

Default constructor should not be called except by Create matrix.

Definition at line 98 of file NodalPrismElec.h.

99  {
100  }
static const PointsKey NullPointsKey(0, eNoPointsType)

Member Function Documentation

◆ CalculateInterpMatrix()

void Nektar::LibUtilities::NodalPrismElec::CalculateInterpMatrix ( const Array< OneD, const NekDouble > &  xi,
const Array< OneD, const NekDouble > &  yi,
const Array< OneD, const NekDouble > &  zi,
Array< OneD, NekDouble > &  interp 
)
private

Definition at line 432 of file NodalPrismElec.cpp.

436 {
437  Array<OneD, Array<OneD, NekDouble>> xi(3);
438  xi[0] = xia;
439  xi[1] = yia;
440  xi[1] = zia;
441 
442  std::shared_ptr<NekMatrix<NekDouble>> mat =
443  m_util->GetInterpolationMatrix(xi);
444  Vmath::Vcopy(mat->GetRows() * mat->GetColumns(), mat->GetRawPtr(), 1,
445  &interp[0], 1);
446 }
std::shared_ptr< NodalUtilPrism > m_util
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255

References m_util, and Vmath::Vcopy().

Referenced by v_GetI().

◆ Create()

std::shared_ptr< PointsBaseType > Nektar::LibUtilities::NodalPrismElec::Create ( const PointsKey key)
static

Definition at line 460 of file NodalPrismElec.cpp.

461 {
462  std::shared_ptr<PointsBaseType> returnval(
464 
465  returnval->Initialize();
466 
467  return returnval;
468 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ NodalPointReorder3d()

void Nektar::LibUtilities::NodalPrismElec::NodalPointReorder3d ( )
private

Definition at line 186 of file NodalPrismElec.cpp.

187 {
188  unsigned int npts = GetNumPoints();
189  using std::vector;
190  vector<int> vertex;
191  vector<int> iEdge_01; // interior edge 0
192  vector<int> iEdge_12; // interior edge 1
193  vector<int> iEdge_23; // interior edge 2
194  vector<int> iEdge_30; // interior edge 3
195  vector<int> iEdge_04; // interior edge 4
196  vector<int> iEdge_14; // interior edge 5
197  vector<int> iEdge_25; // interior edge 6
198  vector<int> iEdge_35; // interior edge 7
199  vector<int> iEdge_45; // interior edge 8
200  vector<int> iFace_0123; // interior face 0
201  vector<int> iFace_014; // interior face 1
202  vector<int> iFace_1254; // interior face 2
203  vector<int> iFace_325; // interior face 3
204  vector<int> iFace_0354; // interior face 4
205  vector<int> interiorVolumePoints; // interior volume points
206  vector<int> map;
207 
208  // Build the lattice prism left to right - bottom to top
209  for (unsigned int y = 0, index = 0; y < npts; y++)
210  {
211  for (unsigned int t = 0; t < npts * (npts + 1) / 2; t++, index++)
212  {
213  if (isVertex(t, y, npts))
214  {
215  vertex.push_back(index);
216  }
217  else if (isEdge(t, y, npts))
218  {
219  if (isEdge_01(t, y, npts))
220  {
221  iEdge_01.push_back(index);
222  }
223  else if (isEdge_12(t, y, npts))
224  {
225  iEdge_12.push_back(index);
226  }
227  else if (isEdge_23(t, y, npts))
228  {
229  iEdge_23.push_back(index);
230  }
231  else if (isEdge_30(t, y, npts))
232  {
233  iEdge_30.push_back(index);
234  }
235  else if (isEdge_04(t, y, npts))
236  {
237  iEdge_04.push_back(index);
238  }
239  else if (isEdge_14(t, y, npts))
240  {
241  iEdge_14.push_back(index);
242  }
243  else if (isEdge_25(t, y, npts))
244  {
245  iEdge_25.push_back(index);
246  }
247  else if (isEdge_35(t, y, npts))
248  {
249  iEdge_35.push_back(index);
250  }
251  else if (isEdge_45(t, y, npts))
252  {
253  iEdge_45.push_back(index);
254  }
255  }
256  else if (isFace(t, y, npts))
257  {
258  if (isFace_0123(t, y, npts))
259  {
260  iFace_0123.push_back(index);
261  }
262  else if (isFace_014(t, y, npts))
263  {
264  iFace_014.push_back(index);
265  }
266  else if (isFace_1254(t, y, npts))
267  {
268  iFace_1254.push_back(index);
269  }
270  else if (isFace_325(t, y, npts))
271  {
272  iFace_325.push_back(index);
273  }
274  else if (isFace_0354(t, y, npts))
275  {
276  iFace_0354.push_back(index);
277  }
278  }
279  else
280  {
281  interiorVolumePoints.push_back(index);
282  }
283  }
284  }
285 
286  // sort vertices
287  std::swap(vertex[2], vertex[4]);
288  // sort edges
289  std::reverse(iEdge_23.begin(), iEdge_23.end());
290  std::reverse(iEdge_30.begin(), iEdge_30.end());
291  std::reverse(iEdge_04.begin(), iEdge_04.end());
292  std::reverse(iEdge_35.begin(), iEdge_35.end());
293 
294  // faces
295  for (unsigned int i = 0; i < npts - 2; i++)
296  {
297  for (unsigned int j = i + 1; j < npts - 2; j++)
298  {
299  std::swap(iFace_1254[i * (npts - 2) + j],
300  iFace_1254[j * (npts - 2) + i]);
301  }
302  }
303  for (int i = 0; i < npts - 2; i++)
304  {
305  std::reverse(iFace_0354.begin() + (i * (npts - 2)),
306  iFace_0354.begin() + (i * (npts - 2) + npts - 2));
307  }
308  for (unsigned int i = 0; i < npts - 2; i++)
309  {
310  for (unsigned int j = i + 1; j < npts - 2; j++)
311  {
312  std::swap(iFace_0354[i * (npts - 2) + j],
313  iFace_0354[j * (npts - 2) + i]);
314  }
315  }
316 
317  for (unsigned int n = 0; n < vertex.size(); ++n)
318  {
319  map.push_back(vertex[n]);
320  }
321 
322  for (unsigned int n = 0; n < iEdge_01.size(); ++n)
323  {
324  map.push_back(iEdge_01[n]);
325  }
326 
327  for (unsigned int n = 0; n < iEdge_12.size(); ++n)
328  {
329  map.push_back(iEdge_12[n]);
330  }
331 
332  for (unsigned int n = 0; n < iEdge_23.size(); ++n)
333  {
334  map.push_back(iEdge_23[n]);
335  }
336 
337  for (unsigned int n = 0; n < iEdge_30.size(); ++n)
338  {
339  map.push_back(iEdge_30[n]);
340  }
341 
342  for (unsigned int n = 0; n < iEdge_04.size(); ++n)
343  {
344  map.push_back(iEdge_04[n]);
345  }
346 
347  for (unsigned int n = 0; n < iEdge_14.size(); ++n)
348  {
349  map.push_back(iEdge_14[n]);
350  }
351 
352  for (unsigned int n = 0; n < iEdge_25.size(); ++n)
353  {
354  map.push_back(iEdge_25[n]);
355  }
356 
357  for (unsigned int n = 0; n < iEdge_35.size(); ++n)
358  {
359  map.push_back(iEdge_35[n]);
360  }
361 
362  for (unsigned int n = 0; n < iEdge_45.size(); ++n)
363  {
364  map.push_back(iEdge_45[n]);
365  }
366 
367  for (unsigned int n = 0; n < iFace_0123.size(); ++n)
368  {
369  map.push_back(iFace_0123[n]);
370  }
371 
372  for (unsigned int n = 0; n < iFace_014.size(); ++n)
373  {
374  map.push_back(iFace_014[n]);
375  }
376 
377  for (unsigned int n = 0; n < iFace_1254.size(); ++n)
378  {
379  map.push_back(iFace_1254[n]);
380  }
381 
382  for (unsigned int n = 0; n < iFace_325.size(); ++n)
383  {
384  map.push_back(iFace_325[n]);
385  }
386 
387  for (unsigned int n = 0; n < iFace_0354.size(); ++n)
388  {
389  map.push_back(iFace_0354[n]);
390  }
391 
392  for (unsigned int n = 0; n < interiorVolumePoints.size(); ++n)
393  {
394  map.push_back(interiorVolumePoints[n]);
395  }
396 
397  Array<OneD, NekDouble> points[3];
398  points[0] = Array<OneD, NekDouble>(GetTotNumPoints());
399  points[1] = Array<OneD, NekDouble>(GetTotNumPoints());
400  points[2] = Array<OneD, NekDouble>(GetTotNumPoints());
401 
402  for (unsigned int index = 0; index < map.size(); ++index)
403  {
404  points[0][index] = m_points[0][index];
405  points[1][index] = m_points[1][index];
406  points[2][index] = m_points[2][index];
407  }
408 
409  for (unsigned int index = 0; index < map.size(); ++index)
410  {
411  m_points[0][index] = points[0][map[index]];
412  m_points[1][index] = points[1][map[index]];
413  m_points[2][index] = points[2][map[index]];
414  }
415 }
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:375
unsigned int GetNumPoints() const
Definition: Points.h:273
unsigned int GetTotNumPoints() const
Definition: Points.h:278

References Nektar::LibUtilities::Points< NekDouble >::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::GetTotNumPoints(), and Nektar::LibUtilities::Points< NekDouble >::m_points.

Referenced by v_CalculatePoints().

◆ v_CalculateDerivMatrix()

void Nektar::LibUtilities::NodalPrismElec::v_CalculateDerivMatrix ( )
overrideprivatevirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 450 of file NodalPrismElec.cpp.

451 {
452  // Allocate the derivative matrix.
454 
455  m_derivmatrix[0] = m_util->GetDerivMatrix(0);
456  m_derivmatrix[1] = m_util->GetDerivMatrix(1);
457  m_derivmatrix[2] = m_util->GetDerivMatrix(2);
458 }
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:381

References Nektar::LibUtilities::Points< NekDouble >::m_derivmatrix, m_util, and Nektar::LibUtilities::Points< NekDouble >::v_CalculateDerivMatrix().

◆ v_CalculatePoints()

void Nektar::LibUtilities::NodalPrismElec::v_CalculatePoints ( )
overrideprivatevirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 156 of file NodalPrismElec.cpp.

157 {
158  // Allocate the storage for points
160 
161  // Populate m_points
162  unsigned int npts = GetNumPoints();
163 
164  LibUtilities::PointsKey pkey1(npts, LibUtilities::eNodalTriElec);
165  Array<OneD, NekDouble> u1, v1;
166  LibUtilities::PointsManager()[pkey1]->GetPoints(u1, v1);
167  LibUtilities::PointsKey pkey2(npts, LibUtilities::eGaussLobattoLegendre);
168  Array<OneD, NekDouble> u;
169  LibUtilities::PointsManager()[pkey2]->GetPoints(u);
170 
171  for (unsigned int y = 0, index = 0; y < npts; y++)
172  {
173  for (size_t t = 0; t < u1.size(); t++, index++)
174  {
175  m_points[0][index] = u1[t];
176  m_points[1][index] = u[y];
177  m_points[2][index] = v1[t];
178  }
179  }
180 
183  npts - 1, m_points[0], m_points[1], m_points[2]);
184 }
PointsManagerT & PointsManager(void)
@ eNodalTriElec
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:83
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eNodalTriElec, Nektar::LibUtilities::Points< NekDouble >::GetNumPoints(), Nektar::LibUtilities::Points< NekDouble >::m_points, m_util, NodalPointReorder3d(), Nektar::LibUtilities::PointsManager(), and Nektar::LibUtilities::Points< NekDouble >::v_CalculatePoints().

◆ v_CalculateWeights()

void Nektar::LibUtilities::NodalPrismElec::v_CalculateWeights ( )
overrideprivatevirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 417 of file NodalPrismElec.cpp.

418 {
419  // Allocate the storage for points
421 
422  typedef DataType T;
423 
424  // Solve the Vandermonde system of integrals for the weight vector
425  NekVector<T> w = m_util->GetWeights();
426 
427  m_weights = Array<OneD, T>(w.GetRows(), w.GetPtr());
428 }
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:377

References Nektar::NekVector< DataType >::GetPtr(), Nektar::NekVector< DataType >::GetRows(), m_util, Nektar::LibUtilities::Points< NekDouble >::m_weights, and Nektar::LibUtilities::Points< NekDouble >::v_CalculateWeights().

◆ v_GetI() [1/2]

virtual const MatrixSharedPtrType Nektar::LibUtilities::NodalPrismElec::v_GetI ( const Array< OneD, const NekDouble > &  x,
const Array< OneD, const NekDouble > &  y,
const Array< OneD, const NekDouble > &  z 
)
inlineoverrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 76 of file NodalPrismElec.h.

80  {
81  size_t numpoints = x.size();
82  unsigned int np = GetTotNumPoints();
83 
84  Array<OneD, NekDouble> interp(GetTotNumPoints() * numpoints);
85  CalculateInterpMatrix(x, y, z, interp);
86 
87  NekDouble *d = interp.data();
88  return MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints,
89  np, d);
90  }
void CalculateInterpMatrix(const Array< OneD, const NekDouble > &xi, const Array< OneD, const NekDouble > &yi, const Array< OneD, const NekDouble > &zi, Array< OneD, NekDouble > &interp)
double NekDouble

References CalculateInterpMatrix(), and Nektar::LibUtilities::Points< NekDouble >::GetTotNumPoints().

◆ v_GetI() [2/2]

virtual const MatrixSharedPtrType Nektar::LibUtilities::NodalPrismElec::v_GetI ( const PointsKey pkey)
inlineoverrideprotectedvirtual

Reimplemented from Nektar::LibUtilities::Points< NekDouble >.

Definition at line 66 of file NodalPrismElec.h.

67  {
68  ASSERTL0(pkey.GetPointsDim() == 3,
69  "NodalPrismElec Points can only interp to "
70  "other 3d point distributions");
71  Array<OneD, const NekDouble> x, y, z;
72  PointsManager()[pkey]->GetPoints(x, y, z);
73  return GetI(x, y, z);
74  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:336

References ASSERTL0, Nektar::LibUtilities::Points< NekDouble >::GetI(), Nektar::LibUtilities::PointsKey::GetPointsDim(), and Nektar::LibUtilities::PointsManager().

Member Data Documentation

◆ initPointsManager

bool Nektar::LibUtilities::NodalPrismElec::initPointsManager
staticprivate
Initial value:
bool RegisterCreator(const KeyType &key, const CreateFuncType &createFunc)
Register the given function and associate it with the key. The return value is just to facilitate cal...
Definition: NekManager.hpp:170
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
@ eNodalPrismElec
3D electrostatically spaced points on a Prism
Definition: PointsType.h:89

Definition at line 93 of file NodalPrismElec.h.

◆ m_util

std::shared_ptr<NodalUtilPrism> Nektar::LibUtilities::NodalPrismElec::m_util
private