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

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

#include <NodalUtil.h>

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

Public Member Functions

 NodalUtilQuad (size_t degree, Array< OneD, NekDouble > r, Array< OneD, NekDouble > s)
 Construct the nodal utility class for a quadrilateral. More...
 
virtual ~NodalUtilQuad ()
 
- 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 quad 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 quadrilateral 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< std::pair< int, int > > m_ordering
 Mapping from the \( (i,j) \) indexing of the basis to a continuous ordering. 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...
 

Detailed Description

Specialisation of the NodalUtil class to support nodal quad elements.

Definition at line 308 of file NodalUtil.h.

Constructor & Destructor Documentation

◆ NodalUtilQuad()

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

Construct the nodal utility class for a quadrilateral.

The constructor of this class sets up the m_ordering member variable used in the evaluation of the orthogonal basis.

Parameters
degreePolynomial order of this nodal quad.
r\( \xi_1 \)-coordinates of nodal points in the standard element.
s\( \xi_2 \)-coordinates of nodal points in the standard element.

Definition at line 825 of file NodalUtil.cpp.

827 : NodalUtil(degree, 2)
828{
829 // Set up parent variables.
830 m_numPoints = r.size();
831 m_xi[0] = r;
832 m_xi[1] = s;
833
834 // Construct a mapping (i,j) -> m from the tensor product space (i,j) to a
835 // single ordering m.
836 for (size_t j = 0; j <= m_degree; ++j)
837 {
838 for (size_t i = 0; i <= m_degree; ++i)
839 {
840 m_ordering.push_back(std::make_pair(i, j));
841 }
842 }
843}
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::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:321

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

◆ ~NodalUtilQuad()

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

Definition at line 314 of file NodalUtil.h.

315 {
316 }

Member Function Documentation

◆ v_CreateUtil()

virtual std::shared_ptr< NodalUtil > Nektar::LibUtilities::NodalUtilQuad::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 327 of file NodalUtil.h.

329 {
331 xi[1]);
332 }
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::NodalUtilQuad::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 334 of file NodalUtil.h.

335 {
336 return 4.0;
337 }

◆ v_NumModes()

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

Calculate the number of degrees of freedom for this element.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 339 of file NodalUtil.h.

340 {
341 return (m_degree + 1) * (m_degree + 1);
342 }

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

◆ v_OrthoBasis()

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

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

In a quad, we use the orthogonal basis

\[ \psi_{m(ij)} = P^{(0,0)}_i(\xi_1) P_j^{(0,0)}(\xi_2) \]

Parameters
modeThe mode of the orthogonal basis to evaluate.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 858 of file NodalUtil.cpp.

859{
860 std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
861 std::pair<int, int> modes = m_ordering[mode];
862
863 // Calculate Jacobi polynomials
864 Polylib::jacobfd(m_numPoints, &m_xi[0][0], &jacobi_i[0], NULL, modes.first,
865 0.0, 0.0);
866 Polylib::jacobfd(m_numPoints, &m_xi[1][0], &jacobi_j[0], NULL, modes.second,
867 0.0, 0.0);
868
869 NekVector<NekDouble> ret(m_numPoints);
870
871 for (size_t i = 0; i < m_numPoints; ++i)
872 {
873 ret[i] = jacobi_i[i] * jacobi_j[i];
874 }
875
876 return ret;
877}
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(), Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, Nektar::LibUtilities::NodalUtil::m_xi, and CG_Iterations::modes.

◆ v_OrthoBasisDeriv()

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

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

Parameters
modeThe mode of the orthogonal basis to evaluate.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 885 of file NodalUtil.cpp.

887{
888 std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
889 std::vector<NekDouble> jacobi_di(m_numPoints), jacobi_dj(m_numPoints);
890 std::pair<int, int> modes = m_ordering[mode];
891
892 // Calculate Jacobi polynomials and their derivatives. Note that we use both
893 // jacobfd and jacobd since jacobfd is only valid for derivatives in the
894 // open interval (-1,1).
895 Polylib::jacobfd(m_numPoints, &m_xi[0][0], &jacobi_i[0], NULL, modes.first,
896 0.0, 0.0);
897 Polylib::jacobfd(m_numPoints, &m_xi[1][0], &jacobi_j[0], NULL, modes.second,
898 0.0, 0.0);
899 Polylib::jacobd(m_numPoints, &m_xi[0][0], &jacobi_di[0], modes.first, 0.0,
900 0.0);
901 Polylib::jacobd(m_numPoints, &m_xi[1][0], &jacobi_dj[0], modes.second, 0.0,
902 0.0);
903
904 NekVector<NekDouble> ret(m_numPoints);
905
906 if (dir == 0)
907 {
908 for (size_t i = 0; i < m_numPoints; ++i)
909 {
910 ret[i] = jacobi_di[i] * jacobi_j[i];
911 }
912 }
913 else
914 {
915 for (size_t i = 0; i < m_numPoints; ++i)
916 {
917 ret[i] = jacobi_i[i] * jacobi_dj[i];
918 }
919 }
920
921 return ret;
922}
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(), Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, Nektar::LibUtilities::NodalUtil::m_xi, and CG_Iterations::modes.

Member Data Documentation

◆ m_ordering

std::vector<std::pair<int, int> > Nektar::LibUtilities::NodalUtilQuad::m_ordering
protected

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

Definition at line 321 of file NodalUtil.h.

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