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 (int 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
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 (int dir)
 Return the Vandermonde matrix of the derivative of the basis functions for the nodal distribution. More...
 
SharedMatrix GetDerivMatrix (int 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 int mode)
 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 int dir, const int mode)
 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)
 Construct a NodalUtil object of the appropriate element type for a given set of points. More...
 
virtual NekDouble v_ModeZeroIntegral ()
 Return the value of the integral of the zero-th mode for this element. More...
 
virtual int v_NumModes ()
 Calculate the number of degrees of freedom for this element. More...
 
- Protected Member Functions inherited from Nektar::LibUtilities::NodalUtil
 NodalUtil (int degree, int dim)
 Set up the NodalUtil object. 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
int m_dim
 Dimension of the nodal element. More...
 
int m_degree
 Degree of the nodal element. More...
 
int 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 311 of file NodalUtil.h.

Constructor & Destructor Documentation

◆ NodalUtilQuad()

Nektar::LibUtilities::NodalUtilQuad::NodalUtilQuad ( int  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 842 of file NodalUtil.cpp.

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

845  : NodalUtil(degree, 2)
846 {
847  // Set up parent variables.
848  m_numPoints = r.num_elements();
849  m_xi[0] = r;
850  m_xi[1] = s;
851 
852  // Construct a mapping (i,j) -> m from the tensor product space (i,j) to a
853  // single ordering m.
854  for (int j = 0; j <= m_degree; ++j)
855  {
856  for (int i = 0; i <= m_degree; ++i)
857  {
858  m_ordering.push_back(std::make_pair(i,j));
859  }
860  }
861 }
int m_degree
Degree of the nodal element.
Definition: NodalUtil.h:107
NodalUtil(int degree, int dim)
Set up the NodalUtil object.
Definition: NodalUtil.h:100
Array< OneD, Array< OneD, NekDouble > > m_xi
Coordinates of the nodal points defining the basis.
Definition: NodalUtil.h:111
int m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:109
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:325

◆ ~NodalUtilQuad()

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

Definition at line 318 of file NodalUtil.h.

319  {
320  }

Member Function Documentation

◆ v_CreateUtil()

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

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 331 of file NodalUtil.h.

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

333  {
335  m_degree, xi[0], xi[1]);
336  }
int m_degree
Degree of the nodal element.
Definition: NodalUtil.h:107
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ v_ModeZeroIntegral()

virtual NekDouble Nektar::LibUtilities::NodalUtilQuad::v_ModeZeroIntegral ( )
inlineprotectedvirtual

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 338 of file NodalUtil.h.

339  {
340  return 4.0;
341  }

◆ v_NumModes()

virtual int Nektar::LibUtilities::NodalUtilQuad::v_NumModes ( )
inlineprotectedvirtual

Calculate the number of degrees of freedom for this element.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 343 of file NodalUtil.h.

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

344  {
345  return (m_degree + 1) * (m_degree + 1);
346  }
int m_degree
Degree of the nodal element.
Definition: NodalUtil.h:107

◆ v_OrthoBasis()

NekVector< NekDouble > Nektar::LibUtilities::NodalUtilQuad::v_OrthoBasis ( const int  mode)
protectedvirtual

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 876 of file NodalUtil.cpp.

References Polylib::jacobfd(), Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, Nektar::LibUtilities::NodalUtil::m_xi, and CG_Iterations::modes.

877 {
878  std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
879  std::pair<int, int> modes = m_ordering[mode];
880 
881  // Calculate Jacobi polynomials
883  m_numPoints, &m_xi[0][0], &jacobi_i[0], NULL, modes.first, 0.0, 0.0);
885  m_numPoints, &m_xi[1][0], &jacobi_j[0], NULL, modes.second, 0.0, 0.0);
886 
887  NekVector<NekDouble> ret(m_numPoints);
888 
889  for (int i = 0; i < m_numPoints; ++i)
890  {
891  ret[i] = jacobi_i[i] * jacobi_j[i];
892  }
893 
894  return ret;
895 }
Array< OneD, Array< OneD, NekDouble > > m_xi
Coordinates of the nodal points defining the basis.
Definition: NodalUtil.h:111
int m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:109
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:325
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:1031

◆ v_OrthoBasisDeriv()

NekVector< NekDouble > Nektar::LibUtilities::NodalUtilQuad::v_OrthoBasisDeriv ( const int  dir,
const int  mode 
)
protectedvirtual

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 903 of file NodalUtil.cpp.

References Polylib::jacobd(), Polylib::jacobfd(), Nektar::LibUtilities::NodalUtil::m_numPoints, m_ordering, Nektar::LibUtilities::NodalUtil::m_xi, and CG_Iterations::modes.

905 {
906  std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
907  std::vector<NekDouble> jacobi_di(m_numPoints), jacobi_dj(m_numPoints);
908  std::pair<int, int> modes = m_ordering[mode];
909 
910  // Calculate Jacobi polynomials and their derivatives. Note that we use both
911  // jacobfd and jacobd since jacobfd is only valid for derivatives in the
912  // open interval (-1,1).
914  m_numPoints, &m_xi[0][0], &jacobi_i[0], NULL, modes.first, 0.0, 0.0);
916  m_numPoints, &m_xi[1][0], &jacobi_j[0], NULL, modes.second, 0.0, 0.0);
918  m_numPoints, &m_xi[0][0], &jacobi_di[0], modes.first, 0.0, 0.0);
920  m_numPoints, &m_xi[1][0], &jacobi_dj[0], modes.second, 0.0, 0.0);
921 
922  NekVector<NekDouble> ret(m_numPoints);
923 
924  if (dir == 0)
925  {
926  for (int i = 0; i < m_numPoints; ++i)
927  {
928  ret[i] = jacobi_di[i] * jacobi_j[i];
929  }
930  }
931  else
932  {
933  for (int i = 0; i < m_numPoints; ++i)
934  {
935  ret[i] = jacobi_i[i] * jacobi_dj[i];
936  }
937  }
938 
939  return ret;
940 }
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:1131
Array< OneD, Array< OneD, NekDouble > > m_xi
Coordinates of the nodal points defining the basis.
Definition: NodalUtil.h:111
int m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:109
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:325
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:1031

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 325 of file NodalUtil.h.

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