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

 ~NodalPrismElec () override
 
 NodalPrismElec (const PointsKey &key)
 
- Public Member Functions inherited from Nektar::LibUtilities::Points< NekDouble >
virtual ~Points ()
 
void Initialize (void)
 
size_t GetPointsDim () const
 
size_t GetNumPoints () const
 
size_t 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 (size_t 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

const MatrixSharedPtrType v_GetI (const PointsKey &pkey) override
 
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_CalculatePoints ()
 
virtual void v_CalculateWeights ()
 

Private Member Functions

 NodalPrismElec ()=delete
 
 NodalPrismElec (const NodalPrismElec &points)=delete
 
void NodalPointReorder3d ()
 
void v_CalculatePoints () final
 
void v_CalculateWeights () final
 
void v_CalculateDerivMatrix () final
 
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 44 of file NodalPrismElec.h.

Constructor & Destructor Documentation

◆ ~NodalPrismElec()

Nektar::LibUtilities::NodalPrismElec::~NodalPrismElec ( )
inlineoverride

Definition at line 47 of file NodalPrismElec.h.

48 {
49 }

◆ NodalPrismElec() [1/3]

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

Definition at line 51 of file NodalPrismElec.h.

51 : PointsBaseType(key)
52 {
53 }
Points< NekDouble > PointsBaseType

◆ NodalPrismElec() [2/3]

Nektar::LibUtilities::NodalPrismElec::NodalPrismElec ( )
privatedelete

◆ NodalPrismElec() [3/3]

Nektar::LibUtilities::NodalPrismElec::NodalPrismElec ( const NodalPrismElec points)
privatedelete

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 421 of file NodalPrismElec.cpp.

425{
426 Array<OneD, Array<OneD, NekDouble>> xi(3);
427 xi[0] = xia;
428 xi[1] = yia;
429 xi[1] = zia;
430
431 std::shared_ptr<NekMatrix<NekDouble>> mat =
432 m_util->GetInterpolationMatrix(xi);
433 Vmath::Vcopy(mat->GetRows() * mat->GetColumns(), mat->GetRawPtr(), 1,
434 &interp[0], 1);
435}
std::shared_ptr< NodalUtilPrism > m_util
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825

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 449 of file NodalPrismElec.cpp.

450{
451 std::shared_ptr<PointsBaseType> returnval(
453
454 returnval->Initialize();
455
456 return returnval;
457}
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 175 of file NodalPrismElec.cpp.

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

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 ( )
finalprivate

Definition at line 439 of file NodalPrismElec.cpp.

440{
441 // Allocate the derivative matrix.
442 PointsBaseType::v_CalculateDerivMatrix();
443
444 m_derivmatrix[0] = m_util->GetDerivMatrix(0);
445 m_derivmatrix[1] = m_util->GetDerivMatrix(1);
446 m_derivmatrix[2] = m_util->GetDerivMatrix(2);
447}
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:362

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

◆ v_CalculatePoints()

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

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

Definition at line 145 of file NodalPrismElec.cpp.

146{
147 // Allocate the storage for points
149
150 // Populate m_points
151 size_t npts = GetNumPoints();
152
153 LibUtilities::PointsKey pkey1(npts, LibUtilities::eNodalTriElec);
154 Array<OneD, NekDouble> u1, v1;
155 LibUtilities::PointsManager()[pkey1]->GetPoints(u1, v1);
156 LibUtilities::PointsKey pkey2(npts, LibUtilities::eGaussLobattoLegendre);
157 Array<OneD, NekDouble> u;
158 LibUtilities::PointsManager()[pkey2]->GetPoints(u);
159
160 for (size_t y = 0, index = 0; y < npts; y++)
161 {
162 for (size_t t = 0; t < u1.size(); t++, index++)
163 {
164 m_points[0][index] = u1[t];
165 m_points[1][index] = u[y];
166 m_points[2][index] = v1[t];
167 }
168 }
169
172 npts - 1, m_points[0], m_points[1], m_points[2]);
173}
PointsManagerT & PointsManager(void)
@ eNodalTriElec
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:81
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

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 ( )
finalprivatevirtual

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

Definition at line 406 of file NodalPrismElec.cpp.

407{
408 // Allocate the storage for points
410
411 typedef DataType T;
412
413 // Solve the Vandermonde system of integrals for the weight vector
414 NekVector<T> w = m_util->GetWeights();
415
416 m_weights = Array<OneD, T>(w.GetRows(), w.GetPtr());
417}
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:358
std::vector< double > w(NPUPPER)

References m_util, Nektar::LibUtilities::Points< NekDouble >::m_weights, Nektar::LibUtilities::Points< NekDouble >::v_CalculateWeights(), and Nektar::UnitTests::w().

◆ v_GetI() [1/2]

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 
)
inlineoverrideprotected

Definition at line 69 of file NodalPrismElec.h.

73 {
74 size_t numpoints = x.size();
75 size_t np = GetTotNumPoints();
76
77 Array<OneD, NekDouble> interp(GetTotNumPoints() * numpoints);
78 CalculateInterpMatrix(x, y, z, interp);
79
80 NekDouble *d = interp.data();
81 return MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints,
82 np, d);
83 }
void CalculateInterpMatrix(const Array< OneD, const NekDouble > &xi, const Array< OneD, const NekDouble > &yi, const Array< OneD, const NekDouble > &zi, Array< OneD, NekDouble > &interp)
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble

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

◆ v_GetI() [2/2]

const MatrixSharedPtrType Nektar::LibUtilities::NodalPrismElec::v_GetI ( const PointsKey pkey)
inlineoverrideprotected

Definition at line 59 of file NodalPrismElec.h.

60 {
61 ASSERTL0(pkey.GetPointsDim() == 3,
62 "NodalPrismElec Points can only interp to "
63 "other 3d point distributions");
64 Array<OneD, const NekDouble> x, y, z;
65 PointsManager()[pkey]->GetPoints(x, y, z);
66 return GetI(x, y, z);
67 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:317

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

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:168
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
@ eNodalPrismElec
3D electrostatically spaced points on a Prism
Definition: PointsType.h:87

Definition at line 86 of file NodalPrismElec.h.

◆ m_util

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