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

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

#include <NodalUtil.h>

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

Public Member Functions

 NodalUtilTriangle (int degree, Array< OneD, NekDouble > r, Array< OneD, NekDouble > s)
 Construct the nodal utility class for a triangle. More...
 
virtual ~NodalUtilTriangle ()
 
- 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 triangular 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 triangular 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...
 
Array< OneD, Array< OneD, NekDouble > > m_eta
 Collapsed coordinates \( (\eta_1, \eta_2) \) of the nodal points. 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 triangular elements.

Definition at line 169 of file NodalUtil.h.

Constructor & Destructor Documentation

◆ NodalUtilTriangle()

Nektar::LibUtilities::NodalUtilTriangle::NodalUtilTriangle ( int  degree,
Array< OneD, NekDouble r,
Array< OneD, NekDouble s 
)

Construct the nodal utility class for a triangle.

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

  • NodalUtilTriangle::m_eta is used to construct the collapsed coordinate locations of the nodal points \( (\eta_1, \eta_2) \) inside the square \([-1,1]^2\) on which the orthogonal basis functions are defined.
  • NodalUtilTriangle::m_ordering constructs a mapping from the index set \( I = \{ (i,j)\ |\ 0\leq i,j \leq P, i+j \leq P \}\) to an ordering \( 0 \leq m(ij) \leq (P+1)(P+2)/2 \) that defines the monomials \( \xi_1^i \xi_2^j \) that span the triangular space. This is then used to calculate which \( (i,j) \) pair corresponding to a column of the Vandermonde matrix when calculating the orthogonal polynomials.
Parameters
degreePolynomial order of this nodal triangle.
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 238 of file NodalUtil.cpp.

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

241  : NodalUtil(degree, 2), m_eta(2)
242 {
243  // Set up parent variables.
244  m_numPoints = r.num_elements();
245  m_xi[0] = r;
246  m_xi[1] = s;
247 
248  // Construct a mapping (i,j) -> m from the triangular tensor product space
249  // (i,j) to a single ordering m.
250  for (int i = 0; i <= m_degree; ++i)
251  {
252  for (int j = 0; j <= m_degree - i; ++j)
253  {
254  m_ordering.push_back(std::make_pair(i,j));
255  }
256  }
257 
258  // Calculate collapsed coordinates from r/s values
259  m_eta[0] = Array<OneD, NekDouble>(m_numPoints);
260  m_eta[1] = Array<OneD, NekDouble>(m_numPoints);
261 
262  for (int i = 0; i < m_numPoints; ++i)
263  {
264  if (fabs(m_xi[1][i]-1.0) < NekConstants::kNekZeroTol)
265  {
266  m_eta[0][i] = -1.0;
267  m_eta[1][i] = 1.0;
268  }
269  else
270  {
271  m_eta[0][i] = 2*(1+m_xi[0][i])/(1-m_xi[1][i])-1.0;
272  m_eta[1][i] = m_xi[1][i];
273  }
274  }
275 }
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:183
int m_degree
Degree of the nodal element.
Definition: NodalUtil.h:107
static const NekDouble kNekZeroTol
Array< OneD, Array< OneD, NekDouble > > m_eta
Collapsed coordinates of the nodal points.
Definition: NodalUtil.h:186
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

◆ ~NodalUtilTriangle()

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

Definition at line 176 of file NodalUtil.h.

177  {
178  }

Member Function Documentation

◆ v_CreateUtil()

virtual std::shared_ptr<NodalUtil> Nektar::LibUtilities::NodalUtilTriangle::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.

Reimplemented in Nektar::Utilities::NodalUtilTriMonomial.

Definition at line 192 of file NodalUtil.h.

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

194  {
196  m_degree, xi[0], xi[1]);
197  }
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::NodalUtilTriangle::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 199 of file NodalUtil.h.

200  {
201  return 2.0 * sqrt(2.0);
202  }

◆ v_NumModes()

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

Calculate the number of degrees of freedom for this element.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 204 of file NodalUtil.h.

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

205  {
206  return (m_degree + 1) * (m_degree + 2) / 2;
207  }
int m_degree
Degree of the nodal element.
Definition: NodalUtil.h:107

◆ v_OrthoBasis()

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

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

In a triangle, we use the orthogonal basis

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

where \( m(ij) \) is the mapping defined in NodalUtilTriangle::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.

Reimplemented in Nektar::Utilities::NodalUtilTriMonomial.

Definition at line 294 of file NodalUtil.cpp.

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

295 {
296  std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
297  std::pair<int, int> modes = m_ordering[mode];
298 
299  // Calculate Jacobi polynomials
301  m_numPoints, &m_eta[0][0], &jacobi_i[0], NULL, modes.first, 0.0, 0.0);
303  m_numPoints, &m_eta[1][0], &jacobi_j[0], NULL, modes.second,
304  2.0 * modes.first + 1.0, 0.0);
305 
306  NekVector<NekDouble> ret(m_numPoints);
307  NekDouble sqrt2 = sqrt(2.0);
308 
309  for (int i = 0; i < m_numPoints; ++i)
310  {
311  ret[i] = sqrt2 * jacobi_i[i] * jacobi_j[i] *
312  pow(1.0 - m_eta[1][i], modes.first);
313  }
314 
315  return ret;
316 }
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:183
Array< OneD, Array< OneD, NekDouble > > m_eta
Collapsed coordinates of the nodal points.
Definition: NodalUtil.h:186
double NekDouble
int m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:109
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::NodalUtilTriangle::v_OrthoBasisDeriv ( const int  dir,
const int  mode 
)
protectedvirtual

Return the value of the derivative of the modal functions for the triangular 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 150.

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.

Reimplemented in Nektar::Utilities::NodalUtilTriMonomial.

Definition at line 332 of file NodalUtil.cpp.

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

334 {
335  std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
336  std::vector<NekDouble> jacobi_di(m_numPoints), jacobi_dj(m_numPoints);
337  std::pair<int, int> modes = m_ordering[mode];
338 
339  // Calculate Jacobi polynomials and their derivatives. Note that we use both
340  // jacobfd and jacobd since jacobfd is only valid for derivatives in the
341  // open interval (-1,1).
343  m_numPoints, &m_eta[0][0], &jacobi_i[0], NULL, modes.first, 0.0,
344  0.0);
346  m_numPoints, &m_eta[1][0], &jacobi_j[0], NULL, modes.second,
347  2.0*modes.first + 1.0, 0.0);
349  m_numPoints, &m_eta[0][0], &jacobi_di[0], modes.first, 0.0, 0.0);
351  m_numPoints, &m_eta[1][0], &jacobi_dj[0], modes.second,
352  2.0*modes.first + 1.0, 0.0);
353 
354  NekVector<NekDouble> ret(m_numPoints);
355  NekDouble sqrt2 = sqrt(2.0);
356 
357  if (dir == 0)
358  {
359  // d/d(\xi_1) = 2/(1-\eta_2) d/d(\eta_1)
360  for (int i = 0; i < m_numPoints; ++i)
361  {
362  ret[i] = 2.0 * sqrt2 * jacobi_di[i] * jacobi_j[i];
363  if (modes.first > 0)
364  {
365  ret[i] *= pow(1.0 - m_eta[1][i], modes.first - 1.0);
366  }
367  }
368  }
369  else
370  {
371  // d/d(\xi_2) = 2(1+\eta_1)/(1-\eta_2) d/d(\eta_1) + d/d(eta_2)
372  for (int i = 0; i < m_numPoints; ++i)
373  {
374  ret[i] = (1 + m_eta[0][i]) * sqrt2 * jacobi_di[i] * jacobi_j[i];
375  if (modes.first > 0)
376  {
377  ret[i] *= pow(1.0 - m_eta[1][i], modes.first - 1.0);
378  }
379 
380  NekDouble tmp = jacobi_dj[i] * pow(1.0 - m_eta[1][i], modes.first);
381  if (modes.first > 0)
382  {
383  tmp -= modes.first * jacobi_j[i] *
384  pow(1.0 - m_eta[1][i], modes.first-1);
385  }
386 
387  ret[i] += sqrt2 * jacobi_i[i] * tmp;
388  }
389  }
390 
391  return ret;
392 }
std::vector< std::pair< int, int > > m_ordering
Mapping from the indexing of the basis to a continuous ordering.
Definition: NodalUtil.h:183
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_eta
Collapsed coordinates of the nodal points.
Definition: NodalUtil.h:186
double NekDouble
int m_numPoints
Total number of nodal points.
Definition: NodalUtil.h:109
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_eta

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

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

Definition at line 186 of file NodalUtil.h.

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

◆ m_ordering

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

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

Definition at line 183 of file NodalUtil.h.

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