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

Represents a basis of a given type. More...

#include <Basis.h>

Public Member Functions

virtual ~Basis ()
 Destructor. More...
 
int GetNumModes () const
 Return order of basis from the basis specification. More...
 
int GetTotNumModes () const
 Return total number of modes from the basis specification. More...
 
int GetNumPoints () const
 Return the number of points from the basis specification. More...
 
int GetTotNumPoints () const
 Return total number of points from the basis specification. More...
 
BasisType GetBasisType () const
 Return the type of expansion basis. More...
 
PointsKey GetPointsKey () const
 Return the points specification for the basis. More...
 
PointsType GetPointsType () const
 Return the type of quadrature. More...
 
const Array< OneD, const NekDouble > & GetZ () const
 
const Array< OneD, const NekDouble > & GetW () const
 
void GetZW (Array< OneD, const NekDouble > &z, Array< OneD, const NekDouble > &w) const
 
const Array< OneD, const NekDouble > & GetBaryWeights () const
 
const std::shared_ptr< NekMatrix< NekDouble > > & GetD (Direction dir=xDir) const
 
const std::shared_ptr< NekMatrix< NekDouble > > GetI (const Array< OneD, const NekDouble > &x)
 
const std::shared_ptr< NekMatrix< NekDouble > > GetI (const BasisKey &bkey)
 
bool ExactIprodInt () const
 Determine if basis has exact integration for inner product. More...
 
bool Collocation () const
 Determine if basis has collocation properties. More...
 
const Array< OneD, const NekDouble > & GetBdata () const
 Return basis definition array m_bdata. More...
 
const Array< OneD, const NekDouble > & GetDbdata () const
 Return basis definition array m_dbdata. More...
 
const BasisKey GetBasisKey () const
 Returns the specification for the Basis. More...
 
void Initialize ()
 

Static Public Member Functions

static std::shared_ptr< BasisCreate (const BasisKey &bkey)
 Returns a new instance of a Basis with given BasisKey. More...
 

Protected Attributes

BasisKey m_basisKey
 Basis specification. More...
 
PointsSharedPtr m_points
 Set of points. More...
 
Array< OneD, NekDoublem_bdata
 Basis definition. More...
 
Array< OneD, NekDoublem_dbdata
 Derivative Basis definition. More...
 
NekManager< BasisKey, NekMatrix< NekDouble >, BasisKey::opLessm_InterpManager
 

Private Member Functions

 Basis (const BasisKey &bkey)
 Private constructor with BasisKey. More...
 
 Basis ()=delete
 Private default constructor. More...
 
std::shared_ptr< NekMatrix< NekDouble > > CalculateInterpMatrix (const BasisKey &tbasis0)
 Calculate the interpolation Matrix for coefficient from one base (m_basisKey) to another (tbasis0) More...
 
void GenBasis ()
 Generate appropriate basis and their derivatives. More...
 

Static Private Attributes

static bool initBasisManager
 

Detailed Description

Represents a basis of a given type.

Definition at line 197 of file Basis.h.

Constructor & Destructor Documentation

◆ ~Basis()

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

Destructor.

Definition at line 204 of file Basis.h.

204{};

◆ Basis() [1/2]

Nektar::LibUtilities::Basis::Basis ( const BasisKey bkey)
private

Private constructor with BasisKey.

Definition at line 90 of file Foundations/Basis.cpp.

91 : m_basisKey(bkey), m_points(PointsManager()[bkey.GetPointsKey()]),
92 m_bdata(bkey.GetTotNumModes() * bkey.GetTotNumPoints()),
93 m_dbdata(bkey.GetTotNumModes() * bkey.GetTotNumPoints())
94{
95 m_InterpManager.RegisterGlobalCreator(
96 std::bind(&Basis::CalculateInterpMatrix, this, std::placeholders::_1));
97}
Array< OneD, NekDouble > m_dbdata
Derivative Basis definition.
Definition: Basis.h:324
PointsSharedPtr m_points
Set of points.
Definition: Basis.h:322
BasisKey m_basisKey
Basis specification.
Definition: Basis.h:321
std::shared_ptr< NekMatrix< NekDouble > > CalculateInterpMatrix(const BasisKey &tbasis0)
Calculate the interpolation Matrix for coefficient from one base (m_basisKey) to another (tbasis0)
NekManager< BasisKey, NekMatrix< NekDouble >, BasisKey::opLess > m_InterpManager
Definition: Basis.h:326
Array< OneD, NekDouble > m_bdata
Basis definition.
Definition: Basis.h:323
PointsManagerT & PointsManager(void)

References CalculateInterpMatrix(), and m_InterpManager.

◆ Basis() [2/2]

Nektar::LibUtilities::Basis::Basis ( )
privatedelete

Private default constructor.

Referenced by Create().

Member Function Documentation

◆ CalculateInterpMatrix()

std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::Basis::CalculateInterpMatrix ( const BasisKey tbasis0)
private

Calculate the interpolation Matrix for coefficient from one base (m_basisKey) to another (tbasis0)

Definition at line 120 of file Foundations/Basis.cpp.

122{
123 size_t dim = m_basisKey.GetNumModes();
124 const PointsKey pkey(dim, LibUtilities::eGaussLobattoLegendre);
125 BasisKey fbkey(m_basisKey.GetBasisType(), dim, pkey);
126 BasisKey tbkey(tbasis0.GetBasisType(), dim, pkey);
127
128 // "Constructur" of the basis
129 BasisSharedPtr fbasis = BasisManager()[fbkey];
130 BasisSharedPtr tbasis = BasisManager()[tbkey];
131
132 // Get B Matrices
133 Array<OneD, NekDouble> fB_data = fbasis->GetBdata();
134 Array<OneD, NekDouble> tB_data = tbasis->GetBdata();
135
136 // Convert to a NekMatrix
137 NekMatrix<NekDouble> fB(dim, dim, fB_data);
138 NekMatrix<NekDouble> tB(dim, dim, tB_data);
139
140 // Invert the "to" matrix: tu = tB^(-1)*fB fu = ftB fu
141 tB.Invert();
142
143 // Compute transformation matrix
144 Array<OneD, NekDouble> zero1D(dim * dim, 0.0);
145 std::shared_ptr<NekMatrix<NekDouble>> ftB(
146 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(dim, dim,
147 zero1D));
148 (*ftB) = tB * fB;
149
150 return ftB;
151}
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:131
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:74
BasisManagerT & BasisManager(void)
std::shared_ptr< Basis > BasisSharedPtr
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

References Nektar::LibUtilities::BasisManager(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::LibUtilities::BasisKey::GetNumModes(), and m_basisKey.

Referenced by Basis().

◆ Collocation()

bool Nektar::LibUtilities::Basis::Collocation ( ) const
inline

Determine if basis has collocation properties.

Definition at line 295 of file Basis.h.

296 {
297 return m_basisKey.Collocation();
298 }
bool Collocation() const
Determine if basis has collocation properties.

References Nektar::LibUtilities::BasisKey::Collocation(), and m_basisKey.

◆ Create()

std::shared_ptr< Basis > Nektar::LibUtilities::Basis::Create ( const BasisKey bkey)
static

Returns a new instance of a Basis with given BasisKey.

Definition at line 99 of file Foundations/Basis.cpp.

100{
101 std::shared_ptr<Basis> returnval(new Basis(bkey));
102 returnval->Initialize();
103
104 return returnval;
105}
Basis()=delete
Private default constructor.

References Basis().

◆ ExactIprodInt()

bool Nektar::LibUtilities::Basis::ExactIprodInt ( void  ) const
inline

Determine if basis has exact integration for inner product.

Definition at line 289 of file Basis.h.

290 {
291 return m_basisKey.ExactIprodInt();
292 }
bool ExactIprodInt() const
Determine if basis has exact integration for inner product.

References Nektar::LibUtilities::BasisKey::ExactIprodInt(), and m_basisKey.

◆ GenBasis()

void Nektar::LibUtilities::Basis::GenBasis ( )
private

Generate appropriate basis and their derivatives.

The following expansions are generated depending on the enum type defined in m_basisKey.m_basistype:

NOTE: This definition does not follow the order in the Karniadakis & Sherwin book since this leads to a more compact hierarchical pattern for implementation purposes. The order of these modes dictates the ordering of the expansion coefficients.

In the following m_numModes = P

eModified_A:

m_bdata[i + j*m_numpoints] = \( \phi^a_i(z_j) = \left \{ \begin{array}{ll} \left ( \frac{1-z_j}{2}\right ) & i = 0 \\ \\ \left ( \frac{1+z_j}{2}\right ) & i = 1 \\ \\ \left ( \frac{1-z_j}{2}\right )\left ( \frac{1+z_j}{2}\right ) P^{1,1}_{i-2}(z_j) & 2\leq i < P\\ \end{array} \right . \)

eModified_B:

m_bdata[n(i,j) + k*m_numpoints] = \( \phi^b_{ij}(z_k) = \left \{ \begin{array}{lll} \phi^a_j(z_k) & i = 0, & 0\leq j < P \\ \\ \left ( \frac{1-z_k}{2}\right )^{i} & 1 \leq i < P,& j = 0 \\ \\ \left ( \frac{1-z_k}{2}\right )^{i} \left ( \frac{1+z_k}{2}\right ) P^{2i-1,1}_{j-1}(z_k) & 1 \leq i < P,\ & 1\leq j < P-i\ \\ \end{array} \right . , \)

where \( n(i,j) \) is a consecutive ordering of the triangular indices \( 0 \leq i, i+j < P \) where j runs fastest.

eModified_C:

m_bdata[n(i,j,k) + l*m_numpoints] = \( \phi^c_{ij,k}(z_l) = \phi^b_{i+j,k}(z_l) = \left \{ \begin{array}{llll} \phi^b_{j,k}(z_l) & i = 0, & 0\leq j < P & 0\leq k < P-j\\ \\ \left ( \frac{1-z_l}{2}\right )^{i+j} & 1\leq i < P,\ & 0\leq j < P-i,\ & k = 0 \\ \\ \left ( \frac{1-z_l}{2}\right )^{i+j} \left ( \frac{1+z_l}{2}\right ) P^{2i+2j-1,1}_{k-1}(z_k) & 1\leq i < P,& 0\leq j < P-i& 1\leq k < P-i-j \\ \\ \end{array} \right . , \)

where \( n(i,j,k) \) is a consecutive ordering of the triangular indices \( 0 \leq i, i+j, i+j+k < P \) where k runs fastest, then j and finally i.

Orthogonal basis A

\(\tilde \psi_p^a (\eta_1) = L_p(\eta_1) = P_p^{0,0}(\eta_1)\)

Orthogonal basis B

\(\tilde \psi_{pq}^b(\eta_2) = \left ( {1 - \eta_2} \over 2 \right)^p P_q^{2p+1,0}(\eta_2)\) \

Orthogonal basis C

\(\tilde \psi_{pqr}^c = \left ( {1 - \eta_3} \over 2 \right)^{p+q} P_r^{2p+2q+2, 0}(\eta_3)\) \ \

Orthogonal basis C for Pyramid expansion (which is richer than tets)

\(\tilde \psi_{pqr}^c = \left ( {1 - \eta_3} \over 2\right)^{pq} P_r^{2pq+2, 0}(\eta_3)\) \( \mbox{where }pq = max(p+q,0) \)

This orthogonal expansion has modes that are always in the Cartesian space, however the equivalent ModifiedPyr_C has vertex modes that do not lie in this space. If one chooses \(pq = max(p+q-1,0)\) then the expansion will space the same space as the vertices but the order of the expanion in 'r' is reduced by one.

1) Eta_z values are the changing the fastest, then r, q, and finally p. 2) r index increases by the stride of numPoints.

Definition at line 216 of file Foundations/Basis.cpp.

217{
218 size_t i, p, q;
219 NekDouble scal;
220 Array<OneD, NekDouble> modeSharedArray;
221 NekDouble *mode;
222 Array<OneD, const NekDouble> z;
223 Array<OneD, const NekDouble> w;
224
225 m_points->GetZW(z, w);
226
227 const NekDouble *D = &(m_points->GetD()->GetPtr())[0];
228 size_t numModes = GetNumModes();
229 size_t numPoints = GetNumPoints();
230
231 switch (GetBasisType())
232 {
233
234 /** \brief Orthogonal basis A
235
236 \f$\tilde \psi_p^a (\eta_1) = L_p(\eta_1) = P_p^{0,0}(\eta_1)\f$
237
238 */
239 case eOrtho_A:
240 case eLegendre:
241 {
242 mode = m_bdata.data();
243
244 for (p = 0; p < numModes; ++p, mode += numPoints)
245 {
246 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, p, 0.0,
247 0.0);
248 // normalise
249 scal = sqrt(0.5 * (2.0 * p + 1.0));
250 for (i = 0; i < numPoints; ++i)
251 {
252 mode[i] *= scal;
253 }
254 }
255
256 // Define derivative basis
257 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
258 numPoints, m_bdata.data(), numPoints, 0.0,
259 m_dbdata.data(), numPoints);
260 } // end scope
261 break;
262
263 /** \brief Orthogonal basis B
264
265 \f$\tilde \psi_{pq}^b(\eta_2) = \left ( {1 - \eta_2} \over 2
266 \right)^p
267 P_q^{2p+1,0}(\eta_2)\f$ \\
268
269 */
270
271 // This is tilde psi_pq in Spencer's book, page 105
272 // The 3-dimensional array is laid out in memory such that
273 // 1) Eta_y values are the changing the fastest, then q and p.
274 // 2) q index increases by the stride of numPoints.
275 case eOrtho_B:
276 {
277 mode = m_bdata.data();
278
279 for (size_t p = 0; p < numModes; ++p)
280 {
281 for (size_t q = 0; q < numModes - p; ++q, mode += numPoints)
282 {
283 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q,
284 2 * p + 1.0, 0.0);
285 for (size_t j = 0; j < numPoints; ++j)
286 {
287 mode[j] *=
288 sqrt(p + q + 1.0) * pow(0.5 * (1.0 - z[j]), p);
289 }
290 }
291 }
292
293 // Define derivative basis
294 Blas::Dgemm('n', 'n', numPoints, numModes * (numModes + 1) / 2,
295 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
296 0.0, m_dbdata.data(), numPoints);
297 }
298 break;
299
300 /** \brief Orthogonal basis C
301
302 \f$\tilde \psi_{pqr}^c = \left ( {1 - \eta_3} \over 2 \right)^{p+q}
303 P_r^{2p+2q+2, 0}(\eta_3)\f$ \ \
304
305 */
306
307 // This is tilde psi_pqr in Spencer's book, page 105
308 // The 4-dimensional array is laid out in memory such that
309 // 1) Eta_z values are the changing the fastest, then r, q, and finally
310 // p. 2) r index increases by the stride of numPoints.
311 case eOrtho_C:
312 {
313 size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
314 mode = m_bdata.data();
315
316 for (size_t p = 0; p <= P; ++p)
317 {
318 for (size_t q = 0; q <= Q - p; ++q)
319 {
320 for (size_t r = 0; r <= R - p - q; ++r, mode += numPoints)
321 {
322 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, r,
323 2 * p + 2 * q + 2.0, 0.0);
324 for (size_t k = 0; k < numPoints; ++k)
325 {
326 // Note factor of 0.5 is part of normalisation
327 mode[k] *= pow(0.5 * (1.0 - z[k]), p + q);
328
329 // finish normalisation
330 mode[k] *= sqrt(r + p + q + 1.5);
331 }
332 }
333 }
334 }
335
336 // Define derivative basis
337 Blas::Dgemm('n', 'n', numPoints,
338 numModes * (numModes + 1) * (numModes + 2) / 6,
339 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
340 0.0, m_dbdata.data(), numPoints);
341 }
342 break;
343
344 /** \brief Orthogonal basis C for Pyramid expansion
345 (which is richer than tets)
346
347 \f$\tilde \psi_{pqr}^c = \left ( {1 - \eta_3} \over
348 2\right)^{pq} P_r^{2pq+2, 0}(\eta_3)\f$ \f$ \mbox{where
349 }pq = max(p+q,0) \f$
350
351 This orthogonal expansion has modes that are
352 always in the Cartesian space, however the
353 equivalent ModifiedPyr_C has vertex modes that do
354 not lie in this space. If one chooses \f$pq =
355 max(p+q-1,0)\f$ then the expansion will space the
356 same space as the vertices but the order of the
357 expanion in 'r' is reduced by one.
358
359 1) Eta_z values are the changing the fastest, then
360 r, q, and finally p. 2) r index increases by the
361 stride of numPoints.
362 */
363 case eOrthoPyr_C:
364 {
365 size_t P = numModes - 1, Q = numModes - 1, R = numModes - 1;
366 mode = m_bdata.data();
367
368 for (size_t p = 0; p <= P; ++p)
369 {
370 for (size_t q = 0; q <= Q; ++q)
371 {
372 for (size_t r = 0; r <= R - std::max(p, q);
373 ++r, mode += numPoints)
374 {
375 // this offset allows for orthogonal
376 // expansion to span linear FE space
377 // of modified basis but means that
378 // the cartesian polynomial space
379 // spanned by the expansion is one
380 // order lower.
381 // size_t pq = max(p + q -1,0);
382 size_t pq = std::max(p + q, size_t(0));
383
384 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, r,
385 2 * pq + 2.0, 0.0);
386 for (size_t k = 0; k < numPoints; ++k)
387 {
388 // Note factor of 0.5 is part of normalisation
389 mode[k] *= pow(0.5 * (1.0 - z[k]), pq);
390
391 // finish normalisation
392 mode[k] *= sqrt(r + pq + 1.5);
393 }
394 }
395 }
396 }
397
398 // Define derivative basis
399 Blas::Dgemm('n', 'n', numPoints,
400 numModes * (numModes + 1) * (numModes + 2) / 6,
401 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
402 0.0, m_dbdata.data(), numPoints);
403 }
404 break;
405
406 case eModified_A:
407 {
408 // Note the following packing deviates from the
409 // definition in the Book by Karniadakis in that we
410 // put the vertex degrees of freedom at the lower
411 // index range to follow a more hierarchic structure.
412
413 for (i = 0; i < numPoints; ++i)
414 {
415 m_bdata[i] = 0.5 * (1 - z[i]);
416 m_bdata[numPoints + i] = 0.5 * (1 + z[i]);
417 }
418
419 mode = m_bdata.data() + 2 * numPoints;
420
421 for (p = 2; p < numModes; ++p, mode += numPoints)
422 {
423 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, p - 2, 1.0,
424 1.0);
425
426 for (i = 0; i < numPoints; ++i)
427 {
428 mode[i] *= m_bdata[i] * m_bdata[numPoints + i];
429 }
430 }
431
432 // Define derivative basis
433 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
434 numPoints, m_bdata.data(), numPoints, 0.0,
435 m_dbdata.data(), numPoints);
436 }
437 break;
438
439 case eModified_B:
440 {
441 // Note the following packing deviates from the
442 // definition in the Book by Karniadakis in two
443 // ways. 1) We put the vertex degrees of freedom
444 // at the lower index range to follow a more
445 // hierarchic structure. 2) We do not duplicate
446 // the singular vertex definition so that only a
447 // triangular number (i.e. (modes)*(modes+1)/2) of
448 // modes are required consistent with the
449 // orthogonal basis.
450
451 // In the current structure the q index runs
452 // faster than the p index so that the matrix has
453 // a more compact structure
454
455 const NekDouble *one_m_z_pow, *one_p_z;
456
457 // bdata should be of size order*(order+1)/2*zorder
458
459 // first fow
460 for (i = 0; i < numPoints; ++i)
461 {
462 m_bdata[0 * numPoints + i] = 0.5 * (1 - z[i]);
463 m_bdata[1 * numPoints + i] = 0.5 * (1 + z[i]);
464 }
465
466 mode = m_bdata.data() + 2 * numPoints;
467
468 for (q = 2; q < numModes; ++q, mode += numPoints)
469 {
470 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q - 2, 1.0,
471 1.0);
472
473 for (i = 0; i < numPoints; ++i)
474 {
475 mode[i] *= m_bdata[i] * m_bdata[numPoints + i];
476 }
477 }
478
479 // second row
480 for (i = 0; i < numPoints; ++i)
481 {
482 mode[i] = 0.5 * (1 - z[i]);
483 }
484
485 mode += numPoints;
486
487 for (q = 2; q < numModes; ++q, mode += numPoints)
488 {
489 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q - 2, 1.0,
490 1.0);
491
492 for (i = 0; i < numPoints; ++i)
493 {
494 mode[i] *= m_bdata[i] * m_bdata[numPoints + i];
495 }
496 }
497
498 // third and higher rows
499 one_m_z_pow = m_bdata.data();
500 one_p_z = m_bdata.data() + numPoints;
501
502 for (p = 2; p < numModes; ++p)
503 {
504 for (i = 0; i < numPoints; ++i)
505 {
506 mode[i] = m_bdata[i] * one_m_z_pow[i];
507 }
508
509 one_m_z_pow = mode;
510 mode += numPoints;
511
512 for (q = 1; q < numModes - p; ++q, mode += numPoints)
513 {
514 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, q - 1,
515 2 * p - 1, 1.0);
516
517 for (i = 0; i < numPoints; ++i)
518 {
519 mode[i] *= one_m_z_pow[i] * one_p_z[i];
520 }
521 }
522 }
523
524 Blas::Dgemm('n', 'n', numPoints, numModes * (numModes + 1) / 2,
525 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
526 0.0, m_dbdata.data(), numPoints);
527 }
528 break;
529
530 case eModified_C:
531 {
532 // Note the following packing deviates from the
533 // definition in the Book by Karniadakis &
534 // Sherwin in two ways. 1) We put the vertex
535 // degrees of freedom at the lower index range to
536 // follow a more hierarchic structure. 2) We do
537 // not duplicate the singular vertex definition
538 // (or the duplicated face information in the book
539 // ) so that only a tetrahedral number
540 // (i.e. (modes)*(modes+1)*(modes+2)/6) of modes
541 // are required consistent with the orthogonal
542 // basis.
543
544 // In the current structure the r index runs
545 // fastest rollowed by q and than the p index so
546 // that the matrix has a more compact structure
547
548 // Note that eModified_C is a re-organisation/
549 // duplication of eModified_B so will get a
550 // temporary Modified_B expansion and copy the
551 // correct components.
552
553 // Generate Modified_B basis;
554 BasisKey ModBKey(eModified_B, m_basisKey.GetNumModes(),
556 BasisSharedPtr ModB = BasisManager()[ModBKey];
557
558 Array<OneD, const NekDouble> ModB_data = ModB->GetBdata();
559
560 // Copy Modified_B basis into first
561 // (numModes*(numModes+1)/2)*numPoints entires of
562 // bdata. This fills in the complete (r,p) face.
563
564 // Set up \phi^c_{p,q,r} = \phi^b_{p+q,r}
565
566 size_t N;
567 size_t B_offset = 0;
568 size_t offset = 0;
569 for (p = 0; p < numModes; ++p)
570 {
571 N = numPoints * (numModes - p) * (numModes - p + 1) / 2;
572 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1,
573 &m_bdata[0] + offset, 1);
574 B_offset += numPoints * (numModes - p);
575 offset += N;
576 }
577
578 // set up derivative of basis.
579 Blas::Dgemm('n', 'n', numPoints,
580 numModes * (numModes + 1) * (numModes + 2) / 6,
581 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
582 0.0, m_dbdata.data(), numPoints);
583 }
584 break;
585
586 case eModifiedPyr_C:
587 {
588 // Note the following packing deviates from the
589 // definition in the Book by Karniadakis &
590 // Sherwin in two ways. 1) We put the vertex
591 // degrees of freedom at the lower index range to
592 // follow a more hierarchic structure. 2) We do
593 // not duplicate the singular vertex definition
594 // so that only a pyramidic number
595 // (i.e. (modes)*(modes+1)*(2*modes+1)/6) of modes
596 // are required consistent with the orthogonal
597 // basis.
598
599 // In the current structure the r index runs
600 // fastest rollowed by q and than the p index so
601 // that the matrix has a more compact structure
602
603 // Generate Modified_B basis;
604 BasisKey ModBKey(eModified_B, m_basisKey.GetNumModes(),
606 BasisSharedPtr ModB = BasisManager()[ModBKey];
607
608 Array<OneD, const NekDouble> ModB_data = ModB->GetBdata();
609
610 // Copy Modified_B basis into first
611 // (numModes*(numModes+1)/2)*numPoints entires of
612 // bdata.
613
614 size_t N;
615 size_t B_offset = 0;
616 size_t offset = 0;
617
618 // Vertex 0,3,4, edges 3,4,7, face 4
619 N = numPoints * (numModes) * (numModes + 1) / 2;
620 Vmath::Vcopy(N, &ModB_data[0], 1, &m_bdata[0], 1);
621 offset += N;
622
623 B_offset += numPoints * (numModes);
624 // Vertex 1 edges 5
625 N = numPoints * (numModes - 1);
626 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, &m_bdata[0] + offset,
627 1);
628 offset += N;
629
630 // Vertex 2 edges 1,6, face 2
631 N = numPoints * (numModes - 1) * (numModes) / 2;
632 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, &m_bdata[0] + offset,
633 1);
634 offset += N;
635
636 B_offset += numPoints * (numModes - 1);
637
638 NekDouble *one_m_z_pow, *one_p_z;
639
640 mode = m_bdata.data() + offset;
641
642 for (p = 2; p < numModes; ++p)
643 {
644 // edges 0 2, faces 1 3
645 N = numPoints * (numModes - p);
646 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, mode, 1);
647 mode += N;
648 Vmath::Vcopy(N, &ModB_data[0] + B_offset, 1, mode, 1);
649 mode += N;
650 B_offset += N;
651
652 one_p_z = m_bdata.data() + numPoints;
653
654 for (q = 2; q < numModes; ++q)
655 {
656 // face 0
657 for (i = 0; i < numPoints; ++i)
658 {
659 // [(1-z)/2]^{p+q-2} Note in book it
660 // seems to suggest p+q-1 but that
661 // does not seem to give complete
662 // polynomial space for pyramids
663 mode[i] = pow(m_bdata[i], p + q - 2);
664 }
665
666 one_m_z_pow = mode;
667 mode += numPoints;
668
669 // interior
670 for (size_t r = 1; r < numModes - std::max(p, q); ++r)
671 {
672 Polylib::jacobfd(numPoints, z.data(), mode, nullptr,
673 r - 1, 2 * p + 2 * q - 3, 1.0);
674
675 for (i = 0; i < numPoints; ++i)
676 {
677 mode[i] *= one_m_z_pow[i] * one_p_z[i];
678 }
679 mode += numPoints;
680 }
681 }
682 }
683
684 // set up derivative of basis.
685 Blas::Dgemm('n', 'n', numPoints,
686 numModes * (numModes + 1) * (2 * numModes + 1) / 6,
687 numPoints, 1.0, D, numPoints, m_bdata.data(), numPoints,
688 0.0, m_dbdata.data(), numPoints);
689 } // end scope
690 break;
691
692 case eGLL_Lagrange:
693 {
694 mode = m_bdata.data();
695 std::shared_ptr<Points<NekDouble>> m_points =
696 PointsManager()[PointsKey(numModes, eGaussLobattoLegendre)];
697 const Array<OneD, const NekDouble> &zp(m_points->GetZ());
698
699 for (p = 0; p < numModes; ++p, mode += numPoints)
700 {
701 for (q = 0; q < numPoints; ++q)
702 {
703 mode[q] =
704 Polylib::hglj(p, z[q], zp.data(), numModes, 0.0, 0.0);
705 }
706 }
707
708 // Define derivative basis
709 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
710 numPoints, m_bdata.data(), numPoints, 0.0,
711 m_dbdata.data(), numPoints);
712
713 } // end scope
714 break;
715
716 case eGauss_Lagrange:
717 {
718 mode = m_bdata.data();
719 std::shared_ptr<Points<NekDouble>> m_points =
720 PointsManager()[PointsKey(numModes, eGaussGaussLegendre)];
721 const Array<OneD, const NekDouble> &zp(m_points->GetZ());
722
723 for (p = 0; p < numModes; ++p, mode += numPoints)
724 {
725 for (q = 0; q < numPoints; ++q)
726 {
727 mode[q] =
728 Polylib::hgj(p, z[q], zp.data(), numModes, 0.0, 0.0);
729 }
730 }
731
732 // Define derivative basis
733 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
734 numPoints, m_bdata.data(), numPoints, 0.0,
735 m_dbdata.data(), numPoints);
736
737 } // end scope
738 break;
739
740 case eFourier:
741 {
742 ASSERTL0(numModes % 2 == 0,
743 "Fourier modes should be a factor of 2");
744
745 for (i = 0; i < numPoints; ++i)
746 {
747 m_bdata[i] = 1.0;
748 m_bdata[numPoints + i] = 0.0;
749
750 m_dbdata[i] = m_dbdata[numPoints + i] = 0.0;
751 }
752
753 for (p = 1; p < numModes / 2; ++p)
754 {
755 for (i = 0; i < numPoints; ++i)
756 {
757 m_bdata[2 * p * numPoints + i] = cos(p * M_PI * (z[i] + 1));
758 m_bdata[(2 * p + 1) * numPoints + i] =
759 -sin(p * M_PI * (z[i] + 1));
760
761 m_dbdata[2 * p * numPoints + i] =
762 -p * M_PI * sin(p * M_PI * (z[i] + 1));
763 m_dbdata[(2 * p + 1) * numPoints + i] =
764 -p * M_PI * cos(p * M_PI * (z[i] + 1));
765 }
766 }
767 } // end scope
768 break;
769
770 // Fourier Single Mode (1st mode)
772 {
773 for (i = 0; i < numPoints; ++i)
774 {
775 m_bdata[i] = cos(M_PI * (z[i] + 1));
776 m_bdata[numPoints + i] = -sin(M_PI * (z[i] + 1));
777
778 m_dbdata[i] = -M_PI * sin(M_PI * (z[i] + 1));
779 m_dbdata[numPoints + i] = -M_PI * cos(M_PI * (z[i] + 1));
780 }
781
782 for (p = 1; p < numModes / 2; ++p)
783 {
784 for (i = 0; i < numPoints; ++i)
785 {
786 m_bdata[2 * p * numPoints + i] = 0.;
787 m_bdata[(2 * p + 1) * numPoints + i] = 0.;
788
789 m_dbdata[2 * p * numPoints + i] = 0.;
790 m_dbdata[(2 * p + 1) * numPoints + i] = 0.;
791 }
792 }
793 } // end scope
794 break;
795
796 // Fourier Real Half Mode
798 {
799 m_bdata[0] = cos(M_PI * z[0]);
800 m_dbdata[0] = -M_PI * sin(M_PI * z[0]);
801 } // end scope
802 break;
803
804 // Fourier Imaginary Half Mode
806 {
807 m_bdata[0] = -sin(M_PI * z[0]);
808 m_dbdata[0] = -M_PI * cos(M_PI * z[0]);
809 } // end scope
810 break;
811
812 case eChebyshev:
813 {
814 mode = m_bdata.data();
815
816 for (p = 0, scal = 1; p < numModes; ++p, mode += numPoints)
817 {
818 Polylib::jacobfd(numPoints, z.data(), mode, nullptr, p, -0.5,
819 -0.5);
820
821 for (i = 0; i < numPoints; ++i)
822 {
823 mode[i] *= scal;
824 }
825
826 scal *= 4 * (p + 1) * (p + 1) / (2 * p + 2) / (2 * p + 1);
827 }
828
829 // Define derivative basis
830 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
831 numPoints, m_bdata.data(), numPoints, 0.0,
832 m_dbdata.data(), numPoints);
833 } // end scope
834 break;
835
836 case eMonomial:
837 {
838 mode = m_bdata.data();
839
840 for (size_t p = 0; p < numModes; ++p, mode += numPoints)
841 {
842 for (size_t i = 0; i < numPoints; ++i)
843 {
844 mode[i] = pow(z[i], p);
845 }
846 }
847
848 // Define derivative basis
849 Blas::Dgemm('n', 'n', numPoints, numModes, numPoints, 1.0, D,
850 numPoints, m_bdata.data(), numPoints, 0.0,
851 m_dbdata.data(), numPoints);
852 } // end scope
853 break;
854
855 default:
856 NEKERROR(ErrorUtil::efatal, "Basis Type not known or "
857 "not implemented at this time.");
858 }
859}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
int GetNumModes() const
Return order of basis from the basis specification.
Definition: Basis.h:207
BasisType GetBasisType() const
Return the type of expansion basis.
Definition: Basis.h:231
int GetNumPoints() const
Return the number of points from the basis specification.
Definition: Basis.h:219
PointsKey GetPointsKey() const
Return distribution of points.
Definition: Basis.h:137
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
Definition: Blas.hpp:383
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:46
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:49
@ P
Monomial polynomials .
Definition: BasisType.h:62
@ eGauss_Lagrange
Lagrange Polynomials using the Gauss points.
Definition: BasisType.h:57
@ eOrtho_A
Principle Orthogonal Functions .
Definition: BasisType.h:42
@ eModified_C
Principle Modified Functions .
Definition: BasisType.h:50
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:56
@ eFourierSingleMode
Fourier ModifiedExpansion with just the first mode .
Definition: BasisType.h:64
@ eChebyshev
Chebyshev Polynomials.
Definition: BasisType.h:61
@ eOrtho_C
Principle Orthogonal Functions .
Definition: BasisType.h:46
@ eModifiedPyr_C
Principle Modified Functions.
Definition: BasisType.h:53
@ eOrtho_B
Principle Orthogonal Functions .
Definition: BasisType.h:44
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:48
@ eFourierHalfModeIm
Fourier Modified expansions with just the imaginary part of the first mode .
Definition: BasisType.h:68
@ eFourierHalfModeRe
Fourier Modified expansions with just the real part of the first mode .
Definition: BasisType.h:66
@ eOrthoPyr_C
Principle Orthogonal Functions .
Definition: BasisType.h:51
@ eFourier
Fourier Expansion .
Definition: BasisType.h:55
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)
double NekDouble
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the ar...
Definition: Polylib.cpp:914
double hglj(const int i, const double z, const double *zglj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zgrj ...
Definition: Polylib.cpp:1025
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition: Polylib.cpp:1248
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References ASSERTL0, Nektar::LibUtilities::BasisManager(), Blas::Dgemm(), Nektar::LibUtilities::eChebyshev, Nektar::ErrorUtil::efatal, Nektar::LibUtilities::eFourier, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::LibUtilities::eFourierSingleMode, Nektar::LibUtilities::eGauss_Lagrange, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eModifiedPyr_C, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::LibUtilities::eOrthoPyr_C, GetBasisType(), Nektar::LibUtilities::BasisKey::GetNumModes(), GetNumModes(), GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Polylib::hgj(), Polylib::hglj(), Polylib::jacobfd(), m_basisKey, m_bdata, m_dbdata, m_points, NEKERROR, Nektar::LibUtilities::P, CellMLToNektar.cellml_metadata::p, Nektar::LibUtilities::PointsManager(), Nektar::UnitTests::q(), tinysimd::sqrt(), Vmath::Vcopy(), Nektar::UnitTests::w(), and Nektar::UnitTests::z().

Referenced by Initialize().

◆ GetBaryWeights()

const Array< OneD, const NekDouble > & Nektar::LibUtilities::Basis::GetBaryWeights ( ) const
inline

Definition at line 264 of file Basis.h.

265 {
266 return m_points->GetBaryWeights();
267 }

References m_points.

◆ GetBasisKey()

const BasisKey Nektar::LibUtilities::Basis::GetBasisKey ( ) const
inline

Returns the specification for the Basis.

Definition at line 313 of file Basis.h.

314 {
315 return m_basisKey;
316 }

References m_basisKey.

Referenced by Nektar::MultiRegions::ExpListHomogeneous2D::GenHomogeneous2DBlockMatrix().

◆ GetBasisType()

BasisType Nektar::LibUtilities::Basis::GetBasisType ( ) const
inline

Return the type of expansion basis.

Definition at line 231 of file Basis.h.

232 {
233 return m_basisKey.GetBasisType();
234 }

References Nektar::LibUtilities::BasisKey::GetBasisType(), and m_basisKey.

Referenced by GenBasis().

◆ GetBdata()

const Array< OneD, const NekDouble > & Nektar::LibUtilities::Basis::GetBdata ( ) const
inline

Return basis definition array m_bdata.

Definition at line 301 of file Basis.h.

302 {
303 return m_bdata;
304 }

References m_bdata.

◆ GetD()

const std::shared_ptr< NekMatrix< NekDouble > > & Nektar::LibUtilities::Basis::GetD ( Direction  dir = xDir) const
inline

Definition at line 269 of file Basis.h.

271 {
272 return m_points->GetD(dir);
273 }

References m_points.

◆ GetDbdata()

const Array< OneD, const NekDouble > & Nektar::LibUtilities::Basis::GetDbdata ( ) const
inline

Return basis definition array m_dbdata.

Definition at line 307 of file Basis.h.

308 {
309 return m_dbdata;
310 }

References m_dbdata.

◆ GetI() [1/2]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::Basis::GetI ( const Array< OneD, const NekDouble > &  x)
inline

Definition at line 275 of file Basis.h.

277 {
278 return m_points->GetI(x);
279 }

References m_points.

◆ GetI() [2/2]

const std::shared_ptr< NekMatrix< NekDouble > > Nektar::LibUtilities::Basis::GetI ( const BasisKey bkey)
inline

Definition at line 281 of file Basis.h.

282 {
283 ASSERTL0(bkey.GetPointsKey().GetPointsDim() == 1,
284 "Interpolation only to other 1d basis");
285 return m_InterpManager[bkey];
286 }

References ASSERTL0, Nektar::LibUtilities::PointsKey::GetPointsDim(), Nektar::LibUtilities::BasisKey::GetPointsKey(), and m_InterpManager.

◆ GetNumModes()

int Nektar::LibUtilities::Basis::GetNumModes ( ) const
inline

Return order of basis from the basis specification.

Definition at line 207 of file Basis.h.

208 {
209 return m_basisKey.GetNumModes();
210 }

References Nektar::LibUtilities::BasisKey::GetNumModes(), and m_basisKey.

Referenced by GenBasis(), Initialize(), and main().

◆ GetNumPoints()

int Nektar::LibUtilities::Basis::GetNumPoints ( ) const
inline

Return the number of points from the basis specification.

Definition at line 219 of file Basis.h.

220 {
221 return m_basisKey.GetNumPoints();
222 }
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:120

References Nektar::LibUtilities::BasisKey::GetNumPoints(), and m_basisKey.

Referenced by GenBasis().

◆ GetPointsKey()

PointsKey Nektar::LibUtilities::Basis::GetPointsKey ( ) const
inline

Return the points specification for the basis.

Definition at line 237 of file Basis.h.

238 {
239 return m_basisKey.GetPointsKey();
240 }

References Nektar::LibUtilities::BasisKey::GetPointsKey(), and m_basisKey.

◆ GetPointsType()

PointsType Nektar::LibUtilities::Basis::GetPointsType ( ) const
inline

Return the type of quadrature.

Definition at line 243 of file Basis.h.

244 {
245 return m_basisKey.GetPointsType();
246 }
PointsType GetPointsType() const
Return type of quadrature.
Definition: Basis.h:143

References Nektar::LibUtilities::BasisKey::GetPointsType(), and m_basisKey.

◆ GetTotNumModes()

int Nektar::LibUtilities::Basis::GetTotNumModes ( ) const
inline

Return total number of modes from the basis specification.

Definition at line 213 of file Basis.h.

214 {
215 return m_basisKey.GetTotNumModes();
216 }
int GetTotNumModes() const
Definition: Basis.h:80

References Nektar::LibUtilities::BasisKey::GetTotNumModes(), and m_basisKey.

◆ GetTotNumPoints()

int Nektar::LibUtilities::Basis::GetTotNumPoints ( ) const
inline

Return total number of points from the basis specification.

Definition at line 225 of file Basis.h.

226 {
228 }
int GetTotNumPoints() const
Definition: Basis.h:125

References Nektar::LibUtilities::BasisKey::GetTotNumPoints(), and m_basisKey.

Referenced by Initialize().

◆ GetW()

const Array< OneD, const NekDouble > & Nektar::LibUtilities::Basis::GetW ( ) const
inline

Definition at line 253 of file Basis.h.

254 {
255 return m_points->GetW();
256 }

References m_points.

◆ GetZ()

const Array< OneD, const NekDouble > & Nektar::LibUtilities::Basis::GetZ ( ) const
inline

Definition at line 248 of file Basis.h.

249 {
250 return m_points->GetZ();
251 }

References m_points.

◆ GetZW()

void Nektar::LibUtilities::Basis::GetZW ( Array< OneD, const NekDouble > &  z,
Array< OneD, const NekDouble > &  w 
) const
inline

Definition at line 258 of file Basis.h.

260 {
261 m_points->GetZW(z, w);
262 }

References m_points, Nektar::UnitTests::w(), and Nektar::UnitTests::z().

◆ Initialize()

void Nektar::LibUtilities::Basis::Initialize ( )

Definition at line 107 of file Foundations/Basis.cpp.

108{
109 ASSERTL0(GetNumModes() > 0,
110 "Cannot call Basis initialisation with zero or negative order");
111 ASSERTL0(GetTotNumPoints() > 0, "Cannot call Basis initialisation with "
112 "zero or negative numbers of points");
113
114 GenBasis();
115}
int GetTotNumPoints() const
Return total number of points from the basis specification.
Definition: Basis.h:225
void GenBasis()
Generate appropriate basis and their derivatives.

References ASSERTL0, GenBasis(), GetNumModes(), and GetTotNumPoints().

Member Data Documentation

◆ initBasisManager

bool Nektar::LibUtilities::Basis::initBasisManager
staticprivate
Initial value:
= {
static std::shared_ptr< Basis > Create(const BasisKey &bkey)
Returns a new instance of a Basis with given BasisKey.
bool RegisterGlobalCreator(const CreateFuncType &createFunc)
Register the Global Create Function. The return value is just to facilitate calling statically.
Definition: NekManager.hpp:177

Definition at line 329 of file Basis.h.

◆ m_basisKey

BasisKey Nektar::LibUtilities::Basis::m_basisKey
protected

◆ m_bdata

Array<OneD, NekDouble> Nektar::LibUtilities::Basis::m_bdata
protected

Basis definition.

Definition at line 323 of file Basis.h.

Referenced by GenBasis(), and GetBdata().

◆ m_dbdata

Array<OneD, NekDouble> Nektar::LibUtilities::Basis::m_dbdata
protected

Derivative Basis definition.

Definition at line 324 of file Basis.h.

Referenced by GenBasis(), and GetDbdata().

◆ m_InterpManager

NekManager<BasisKey, NekMatrix<NekDouble>, BasisKey::opLess> Nektar::LibUtilities::Basis::m_InterpManager
protected

Definition at line 326 of file Basis.h.

Referenced by Basis(), and GetI().

◆ m_points

PointsSharedPtr Nektar::LibUtilities::Basis::m_points
protected

Set of points.

Definition at line 322 of file Basis.h.

Referenced by GenBasis(), GetBaryWeights(), GetD(), GetI(), GetW(), GetZ(), and GetZW().