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

Specialisation of the NodalUtil class to support nodal tetrahedral elements. More...

#include <NodalUtil.h>

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

Public Member Functions

 NodalUtilTetrahedron (size_t degree, Array< OneD, NekDouble > r, Array< OneD, NekDouble > s, Array< OneD, NekDouble > t)
 Construct the nodal utility class for a tetrahedron. More...
 
virtual ~NodalUtilTetrahedron ()
 
- Public Member Functions inherited from Nektar::LibUtilities::NodalUtil
virtual ~NodalUtil ()=default
 
NekVector< NekDoubleGetWeights ()
 Obtain the integration weights for the given nodal distribution. More...
 
SharedMatrix GetVandermonde ()
 Return the Vandermonde matrix for the nodal distribution. More...
 
SharedMatrix GetVandermondeForDeriv (size_t dir)
 Return the Vandermonde matrix of the derivative of the basis functions for the nodal distribution. More...
 
SharedMatrix GetDerivMatrix (size_t dir)
 Return the derivative matrix for the nodal distribution. More...
 
SharedMatrix GetInterpolationMatrix (Array< OneD, Array< OneD, NekDouble > > &xi)
 Construct the interpolation matrix used to evaluate the basis at the points xi inside the element. More...
 

Protected Member Functions

virtual NekVector< NekDoublev_OrthoBasis (const size_t mode) override
 Return the value of the modal functions for the tetrahedral element at the nodal points m_xi for a given mode. More...
 
virtual NekVector< NekDoublev_OrthoBasisDeriv (const size_t dir, const size_t mode) override
 Return the value of the derivative of the modal functions for the tetrahedral element at the nodal points m_xi for a given mode. More...
 
virtual std::shared_ptr< NodalUtilv_CreateUtil (Array< OneD, Array< OneD, NekDouble > > &xi) override
 Construct a NodalUtil object of the appropriate element type for a given set of points. More...
 
virtual NekDouble v_ModeZeroIntegral () override
 Return the value of the integral of the zero-th mode for this element. More...
 
virtual size_t v_NumModes () override
 Calculate the number of degrees of freedom for this element. More...
 
- Protected Member Functions inherited from Nektar::LibUtilities::NodalUtil
 NodalUtil (size_t degree, size_t dim)
 Set up the NodalUtil object. More...
 
virtual NekVector< NekDoublev_OrthoBasis (const size_t mode)=0
 Return the values of the orthogonal basis at the nodal points for a given mode. More...
 
virtual NekVector< NekDoublev_OrthoBasisDeriv (const size_t dir, const size_t mode)=0
 Return the values of the derivative of the orthogonal basis at the nodal points for a given mode. More...
 
virtual std::shared_ptr< NodalUtilv_CreateUtil (Array< OneD, Array< OneD, NekDouble > > &xi)=0
 Construct a NodalUtil object of the appropriate element type for a given set of points. More...
 
virtual NekDouble v_ModeZeroIntegral ()=0
 Return the value of the integral of the zero-th mode for this element. More...
 
virtual size_t v_NumModes ()=0
 Calculate the number of degrees of freedom for this element. More...
 

Protected Attributes

std::vector< Modem_ordering
 Mapping from the \( (i,j,k) \) indexing of the basis to a continuous ordering. More...
 
Array< OneD, Array< OneD, NekDouble > > m_eta
 Collapsed coordinates \( (\eta_1, \eta_2, \eta_3) \) of the nodal points. More...
 
- Protected Attributes inherited from Nektar::LibUtilities::NodalUtil
size_t m_dim
 Dimension of the nodal element. More...
 
size_t m_degree
 Degree of the nodal element. More...
 
size_t m_numPoints
 Total number of nodal points. More...
 
Array< OneD, Array< OneD, NekDouble > > m_xi
 Coordinates of the nodal points defining the basis. More...
 

Private Types

typedef std::tuple< int, int, int > Mode
 

Detailed Description

Specialisation of the NodalUtil class to support nodal tetrahedral elements.

Definition at line 212 of file NodalUtil.h.

Member Typedef Documentation

◆ Mode

typedef std::tuple<int, int, int> Nektar::LibUtilities::NodalUtilTetrahedron::Mode
private

Definition at line 214 of file NodalUtil.h.

Constructor & Destructor Documentation

◆ NodalUtilTetrahedron()

Nektar::LibUtilities::NodalUtilTetrahedron::NodalUtilTetrahedron ( size_t  degree,
Array< OneD, NekDouble r,
Array< OneD, NekDouble s,
Array< OneD, NekDouble t 
)

Construct the nodal utility class for a tetrahedron.

The constructor of this class sets up two member variables used in the evaluation of the orthogonal basis:

  • NodalUtilTetrahedron::m_eta is used to construct the collapsed coordinate locations of the nodal points \( (\eta_1, \eta_2, \eta_3) \) inside the cube \([-1,1]^3\) on which the orthogonal basis functions are defined.
  • NodalUtilTetrahedron::m_ordering constructs a mapping from the index set \( I = \{ (i,j,k)\ |\ 0\leq i,j,k \leq P, i+j \leq P, i+j+k \leq P \}\) to an ordering \( 0 \leq m(ijk) \leq (P+1)(P+2)(P+3)/6 \) that defines the monomials \( \xi_1^i \xi_2^j \xi_3^k \) that span the tetrahedral space. This is then used to calculate which \( (i,j,k) \) triple (represented as a tuple) corresponding to a column of the Vandermonde matrix when calculating the orthogonal polynomials.
Parameters
degreePolynomial order of this nodal tetrahedron
r\( \xi_1 \)-coordinates of nodal points in the standard element.
s\( \xi_2 \)-coordinates of nodal points in the standard element.
t\( \xi_3 \)-coordinates of nodal points in the standard element.

Definition at line 412 of file NodalUtil.cpp.

416 : NodalUtil(degree, 3), m_eta(3)
417{
418 m_numPoints = r.size();
419 m_xi[0] = r;
420 m_xi[1] = s;
421 m_xi[2] = t;
422
423 for (size_t i = 0; i <= m_degree; ++i)
424 {
425 for (size_t j = 0; j <= m_degree - i; ++j)
426 {
427 for (size_t k = 0; k <= m_degree - i - j; ++k)
428 {
429 m_ordering.push_back(Mode(i, j, k));
430 }
431 }
432 }
433
434 // Calculate collapsed coordinates from r/s values
435 m_eta[0] = Array<OneD, NekDouble>(m_numPoints);
436 m_eta[1] = Array<OneD, NekDouble>(m_numPoints);
437 m_eta[2] = Array<OneD, NekDouble>(m_numPoints);
438
439 for (size_t i = 0; i < m_numPoints; ++i)
440 {
441 if (fabs(m_xi[2][i] - 1.0) < NekConstants::kNekZeroTol)
442 {
443 // Very top point of the tetrahedron
444 m_eta[0][i] = -1.0;
445 m_eta[1][i] = -1.0;
446 m_eta[2][i] = m_xi[2][i];
447 }
448 else
449 {
450 if (fabs(m_xi[1][i] - 1.0) < NekConstants::kNekZeroTol)
451 {
452 // Distant diagonal edge shared by all eta_x coordinate planes:
453 // the xi_y == -xi_z line
454 m_eta[0][i] = -1.0;
455 }
456 else if (fabs(m_xi[1][i] + m_xi[2][i]) < NekConstants::kNekZeroTol)
457 {
458 m_eta[0][i] = -1.0;
459 }
460 else
461 {
462 m_eta[0][i] =
463 2.0 * (1.0 + m_xi[0][i]) / (-m_xi[1][i] - m_xi[2][i]) - 1.0;
464 }
465 m_eta[1][i] = 2.0 * (1.0 + m_xi[1][i]) / (1.0 - m_xi[2][i]) - 1.0;
466 m_eta[2][i] = m_xi[2][i];
467 }
468 }
469}
Array< OneD, Array< OneD, NekDouble > > m_xi
Coordinates of the nodal points defining the basis.
Definition: NodalUtil.h:109
size_t m_degree
Degree of the nodal element.
Definition: NodalUtil.h:105
size_t m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:107
NodalUtil(size_t degree, size_t dim)
Set up the NodalUtil object.
Definition: NodalUtil.h:97
std::tuple< int, int, int > Mode
Definition: NodalUtil.h:214
Array< OneD, Array< OneD, NekDouble > > m_eta
Collapsed coordinates of the nodal points.
Definition: NodalUtil.h:233
std::vector< Mode > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:229
static const NekDouble kNekZeroTol

References Nektar::NekConstants::kNekZeroTol, Nektar::LibUtilities::NodalUtil::m_degree, m_eta, Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, and Nektar::LibUtilities::NodalUtil::m_xi.

◆ ~NodalUtilTetrahedron()

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

Definition at line 222 of file NodalUtil.h.

223 {
224 }

Member Function Documentation

◆ v_CreateUtil()

virtual std::shared_ptr< NodalUtil > Nektar::LibUtilities::NodalUtilTetrahedron::v_CreateUtil ( Array< OneD, Array< OneD, NekDouble > > &  xi)
inlineoverrideprotectedvirtual

Construct a NodalUtil object of the appropriate element type for a given set of points.

This function is used inside NodalUtil::GetInterpolationMatrix so that the (potentially non-square) Vandermonde matrix can be constructed to create the interpolation matrix at an arbitrary set of points in the domain.

Parameters
xiDistribution of nodal points to create utility with.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 239 of file NodalUtil.h.

241 {
243 m_degree, xi[0], xi[1], xi[2]);
244 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::LibUtilities::NodalUtil::m_degree.

◆ v_ModeZeroIntegral()

virtual NekDouble Nektar::LibUtilities::NodalUtilTetrahedron::v_ModeZeroIntegral ( )
inlineoverrideprotectedvirtual

Return the value of the integral of the zero-th mode for this element.

Note that for the orthogonal basis under consideration, all modes integrate to zero asides from the zero-th mode. This function is used in NodalUtil::GetWeights to determine integration weights.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 246 of file NodalUtil.h.

247 {
248 return 8.0 * sqrt(2.0) / 3.0;
249 }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References tinysimd::sqrt().

◆ v_NumModes()

virtual size_t Nektar::LibUtilities::NodalUtilTetrahedron::v_NumModes ( )
inlineoverrideprotectedvirtual

Calculate the number of degrees of freedom for this element.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 251 of file NodalUtil.h.

252 {
253 return (m_degree + 1) * (m_degree + 2) * (m_degree + 3) / 6;
254 }

References Nektar::LibUtilities::NodalUtil::m_degree.

◆ v_OrthoBasis()

NekVector< NekDouble > Nektar::LibUtilities::NodalUtilTetrahedron::v_OrthoBasis ( const size_t  mode)
overrideprotectedvirtual

Return the value of the modal functions for the tetrahedral element at the nodal points m_xi for a given mode.

In a tetrahedron, we use the orthogonal basis

\[ \psi_{m(ijk)} = \sqrt{8} P^{(0,0)}_i(\xi_1) P_j^{(2i+1,0)}(\xi_2) P_k^{(2i+2j+2,0)}(\xi_3) (1-\xi_2)^i (1-\xi_3)^{i+j} \]

where \( m(ijk) \) is the mapping defined in m_ordering and \( J_n^{(\alpha,\beta)}(z) \) denotes the standard Jacobi polynomial.

Parameters
modeThe mode of the orthogonal basis to evaluate.
Returns
Vector containing orthogonal basis evaluated at the points m_xi.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 487 of file NodalUtil.cpp.

488{
489 std::vector<NekDouble> jacA(m_numPoints), jacB(m_numPoints);
490 std::vector<NekDouble> jacC(m_numPoints);
491
492 size_t I, J, K;
493 std::tie(I, J, K) = m_ordering[mode];
494
495 // Calculate Jacobi polynomials
496 Polylib::jacobfd(m_numPoints, &m_eta[0][0], &jacA[0], NULL, I, 0.0, 0.0);
497 Polylib::jacobfd(m_numPoints, &m_eta[1][0], &jacB[0], NULL, J,
498 2.0 * I + 1.0, 0.0);
499 Polylib::jacobfd(m_numPoints, &m_eta[2][0], &jacC[0], NULL, K,
500 2.0 * (I + J) + 2.0, 0.0);
501
502 NekVector<NekDouble> ret(m_numPoints);
503 NekDouble sqrt8 = sqrt(8.0);
504
505 for (size_t i = 0; i < m_numPoints; ++i)
506 {
507 ret[i] = sqrt8 * jacA[i] * jacB[i] * jacC[i] *
508 pow(1.0 - m_eta[1][i], I) * pow(1.0 - m_eta[2][i], I + J);
509 }
510
511 return ret;
512}
double NekDouble
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:1206

References Polylib::jacobfd(), m_eta, Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, and tinysimd::sqrt().

◆ v_OrthoBasisDeriv()

NekVector< NekDouble > Nektar::LibUtilities::NodalUtilTetrahedron::v_OrthoBasisDeriv ( const size_t  dir,
const size_t  mode 
)
overrideprotectedvirtual

Return the value of the derivative of the modal functions for the tetrahedral element at the nodal points m_xi for a given mode.

Note that this routine must use the chain rule combined with the collapsed coordinate derivatives as described in Sherwin & Karniadakis (2nd edition), pg 152.

Parameters
dirCoordinate direction in which to evaluate the derivative.
modeThe mode of the orthogonal basis to evaluate.
Returns
Vector containing the derivative of the orthogonal basis evaluated at the points m_xi.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 528 of file NodalUtil.cpp.

530{
531 std::vector<NekDouble> jacA(m_numPoints), jacB(m_numPoints);
532 std::vector<NekDouble> jacC(m_numPoints);
533 std::vector<NekDouble> jacDerivA(m_numPoints), jacDerivB(m_numPoints);
534 std::vector<NekDouble> jacDerivC(m_numPoints);
535
536 size_t I, J, K;
537 std::tie(I, J, K) = m_ordering[mode];
538
539 // Calculate Jacobi polynomials
540 Polylib::jacobfd(m_numPoints, &m_eta[0][0], &jacA[0], NULL, I, 0.0, 0.0);
541 Polylib::jacobfd(m_numPoints, &m_eta[1][0], &jacB[0], NULL, J,
542 2.0 * I + 1.0, 0.0);
543 Polylib::jacobfd(m_numPoints, &m_eta[2][0], &jacC[0], NULL, K,
544 2.0 * (I + J) + 2.0, 0.0);
545 Polylib::jacobd(m_numPoints, &m_eta[0][0], &jacDerivA[0], I, 0.0, 0.0);
546 Polylib::jacobd(m_numPoints, &m_eta[1][0], &jacDerivB[0], J, 2.0 * I + 1.0,
547 0.0);
548 Polylib::jacobd(m_numPoints, &m_eta[2][0], &jacDerivC[0], K,
549 2.0 * (I + J) + 2.0, 0.0);
550
551 NekVector<NekDouble> ret(m_numPoints);
552 NekDouble sqrt8 = sqrt(8.0);
553
554 // Always compute x-derivative since this term appears in the latter two
555 // terms.
556 for (size_t i = 0; i < m_numPoints; ++i)
557 {
558 ret[i] = 4.0 * sqrt8 * jacDerivA[i] * jacB[i] * jacC[i];
559
560 if (I > 0)
561 {
562 ret[i] *= pow(1 - m_eta[1][i], I - 1);
563 }
564
565 if (I + J > 0)
566 {
567 ret[i] *= pow(1 - m_eta[2][i], I + J - 1);
568 }
569 }
570
571 if (dir >= 1)
572 {
573 // Multiply by (1+a)/2
574 NekVector<NekDouble> tmp(m_numPoints);
575
576 for (size_t i = 0; i < m_numPoints; ++i)
577 {
578 ret[i] *= 0.5 * (m_eta[0][i] + 1.0);
579
580 tmp[i] = 2.0 * sqrt8 * jacA[i] * jacC[i];
581 if (I + J > 0)
582 {
583 tmp[i] *= pow(1.0 - m_eta[2][i], I + J - 1);
584 }
585
586 NekDouble tmp2 = jacDerivB[i] * pow(1.0 - m_eta[1][i], I);
587 if (I > 0)
588 {
589 tmp2 -= I * jacB[i] * pow(1.0 - m_eta[1][i], I - 1);
590 }
591
592 tmp[i] *= tmp2;
593 }
594
595 if (dir == 1)
596 {
597 ret += tmp;
598 return ret;
599 }
600
601 for (size_t i = 0; i < m_numPoints; ++i)
602 {
603 ret[i] += 0.5 * (1.0 + m_eta[1][i]) * tmp[i];
604
605 NekDouble tmp2 = jacDerivC[i] * pow(1.0 - m_eta[2][i], I + J);
606 if (I + J > 0)
607 {
608 tmp2 -= jacC[i] * (I + J) * pow(1.0 - m_eta[2][i], I + J - 1);
609 }
610
611 ret[i] +=
612 sqrt8 * jacA[i] * jacB[i] * pow(1.0 - m_eta[1][i], I) * tmp2;
613 }
614 }
615
616 return ret;
617}
void jacobd(const int np, const double *z, double *polyd, const int n, const double alpha, const double beta)
Calculate the derivative of Jacobi polynomials.
Definition: Polylib.cpp:1318

References Polylib::jacobd(), Polylib::jacobfd(), m_eta, Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, and tinysimd::sqrt().

Member Data Documentation

◆ m_eta

Array<OneD, Array<OneD, NekDouble> > Nektar::LibUtilities::NodalUtilTetrahedron::m_eta
protected

Collapsed coordinates \( (\eta_1, \eta_2, \eta_3) \) of the nodal points.

Definition at line 233 of file NodalUtil.h.

Referenced by NodalUtilTetrahedron(), v_OrthoBasis(), and v_OrthoBasisDeriv().

◆ m_ordering

std::vector<Mode> Nektar::LibUtilities::NodalUtilTetrahedron::m_ordering
protected

Mapping from the \( (i,j,k) \) indexing of the basis to a continuous ordering.

Definition at line 229 of file NodalUtil.h.

Referenced by NodalUtilTetrahedron(), v_OrthoBasis(), and v_OrthoBasisDeriv().