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 (size_t degree, Array< OneD, NekDouble > r, Array< OneD, NekDouble > s)
 Construct the nodal utility class for a triangle. More...
 
 ~NodalUtilTriangle () override
 
- 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

NekVector< NekDoublev_OrthoBasis (const size_t mode) override
 Return the value of the modal functions for the triangular element at the nodal points m_xi for a given mode. More...
 
NekVector< NekDoublev_OrthoBasisDeriv (const size_t dir, const size_t mode) override
 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...
 
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...
 
NekDouble v_ModeZeroIntegral () override
 Return the value of the integral of the zero-th mode for this element. More...
 
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...
 
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
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 triangular elements.

Definition at line 165 of file NodalUtil.h.

Constructor & Destructor Documentation

◆ NodalUtilTriangle()

Nektar::LibUtilities::NodalUtilTriangle::NodalUtilTriangle ( size_t  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 234 of file NodalUtil.cpp.

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

◆ ~NodalUtilTriangle()

Nektar::LibUtilities::NodalUtilTriangle::~NodalUtilTriangle ( )
inlineoverride

Definition at line 172 of file NodalUtil.h.

173 {
174 }

Member Function Documentation

◆ v_CreateUtil()

std::shared_ptr< NodalUtil > Nektar::LibUtilities::NodalUtilTriangle::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 188 of file NodalUtil.h.

190 {
192 m_degree, xi[0], xi[1]);
193 }
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()

NekDouble Nektar::LibUtilities::NodalUtilTriangle::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 195 of file NodalUtil.h.

196 {
197 return 2.0 * sqrt(2.0);
198 }
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References tinysimd::sqrt().

◆ v_NumModes()

size_t Nektar::LibUtilities::NodalUtilTriangle::v_NumModes ( )
inlineoverrideprotectedvirtual

Calculate the number of degrees of freedom for this element.

Implements Nektar::LibUtilities::NodalUtil.

Definition at line 200 of file NodalUtil.h.

201 {
202 return (m_degree + 1) * (m_degree + 2) / 2;
203 }

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

◆ v_OrthoBasis()

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

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.

Definition at line 289 of file NodalUtil.cpp.

290{
291 std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
292 std::pair<int, int> modes = m_ordering[mode];
293
294 // Calculate Jacobi polynomials
295 Polylib::jacobfd(m_numPoints, &m_eta[0][0], &jacobi_i[0], nullptr,
296 modes.first, 0.0, 0.0);
297 Polylib::jacobfd(m_numPoints, &m_eta[1][0], &jacobi_j[0], nullptr,
298 modes.second, 2.0 * modes.first + 1.0, 0.0);
299
300 NekVector<NekDouble> ret(m_numPoints);
301 NekDouble sqrt2 = sqrt(2.0);
302
303 for (size_t i = 0; i < m_numPoints; ++i)
304 {
305 ret[i] = sqrt2 * jacobi_i[i] * jacobi_j[i] *
306 pow(1.0 - m_eta[1][i], modes.first);
307 }
308
309 return ret;
310}
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:1248

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

◆ v_OrthoBasisDeriv()

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

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.

Definition at line 326 of file NodalUtil.cpp.

328{
329 std::vector<NekDouble> jacobi_i(m_numPoints), jacobi_j(m_numPoints);
330 std::vector<NekDouble> jacobi_di(m_numPoints), jacobi_dj(m_numPoints);
331 std::pair<int, int> modes = m_ordering[mode];
332
333 // Calculate Jacobi polynomials and their derivatives. Note that we use both
334 // jacobfd and jacobd since jacobfd is only valid for derivatives in the
335 // open interval (-1,1).
336 Polylib::jacobfd(m_numPoints, &m_eta[0][0], &jacobi_i[0], nullptr,
337 modes.first, 0.0, 0.0);
338 Polylib::jacobfd(m_numPoints, &m_eta[1][0], &jacobi_j[0], nullptr,
339 modes.second, 2.0 * modes.first + 1.0, 0.0);
340 Polylib::jacobd(m_numPoints, &m_eta[0][0], &jacobi_di[0], modes.first, 0.0,
341 0.0);
342 Polylib::jacobd(m_numPoints, &m_eta[1][0], &jacobi_dj[0], modes.second,
343 2.0 * modes.first + 1.0, 0.0);
344
345 NekVector<NekDouble> ret(m_numPoints);
346 NekDouble sqrt2 = sqrt(2.0);
347
348 if (dir == 0)
349 {
350 // d/d(\xi_1) = 2/(1-\eta_2) d/d(\eta_1)
351 for (size_t i = 0; i < m_numPoints; ++i)
352 {
353 ret[i] = 2.0 * sqrt2 * jacobi_di[i] * jacobi_j[i];
354 if (modes.first > 0)
355 {
356 ret[i] *= pow(1.0 - m_eta[1][i], modes.first - 1.0);
357 }
358 }
359 }
360 else
361 {
362 // d/d(\xi_2) = 2(1+\eta_1)/(1-\eta_2) d/d(\eta_1) + d/d(eta_2)
363 for (size_t i = 0; i < m_numPoints; ++i)
364 {
365 ret[i] = (1 + m_eta[0][i]) * sqrt2 * jacobi_di[i] * jacobi_j[i];
366 if (modes.first > 0)
367 {
368 ret[i] *= pow(1.0 - m_eta[1][i], modes.first - 1.0);
369 }
370
371 NekDouble tmp = jacobi_dj[i] * pow(1.0 - m_eta[1][i], modes.first);
372 if (modes.first > 0)
373 {
374 tmp -= modes.first * jacobi_j[i] *
375 pow(1.0 - m_eta[1][i], modes.first - 1);
376 }
377
378 ret[i] += sqrt2 * jacobi_i[i] * tmp;
379 }
380 }
381
382 return ret;
383}
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:1378

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

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

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