Nektar++
PolyEPoints.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PolyEPoints.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: 1D Evenly-Spaced Points
32//
33///////////////////////////////////////////////////////////////////////////////
34
38
40{
41
44
46{
47 // Allocate the storage for points
49
50 size_t npts = m_pointsKey.GetNumPoints();
51 if (npts == 1)
52 {
53 m_points[0][0] = 0.0;
54 }
55 else
56 {
57 NekDouble dx = 2.0 / (NekDouble)(npts - 1);
58 for (size_t i = 0; i < npts; ++i)
59 {
60 m_points[0][i] = -1.0 + i * dx;
61 }
62 }
63}
64
66{
67 // Allocate the storage for points
69
70 size_t npts = m_pointsKey.GetNumPoints();
71 if (npts == 1)
72 {
73 m_weights[0] = 2.0; // midpoint rule
74 }
75 else
76 {
77 PointsKey gaussKey(npts, eGaussLobattoLegendre);
78 std::shared_ptr<PointsBaseType> ptr = PointsManager()[gaussKey];
81
82 ptr->GetZW(z, w);
83 for (size_t i = 0; i < npts; ++i)
84 {
85 m_weights[i] = 0.0;
86 for (unsigned j = 0; j < npts; ++j)
87 {
88 m_weights[i] +=
89 w[j] * Polylib::laginterp(z[j], i, &m_points[0][0], npts);
90 }
91 }
92 }
93}
94
96{
97 // Allocate the derivative matrix.
98 PointsBaseType::v_CalculateDerivMatrix();
99
100 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
101 {
102 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
103 {
105 m_points[0][i], j, &m_points[0][0], m_pointsKey.GetNumPoints());
106 }
107 }
108}
109
111 size_t npts, const Array<OneD, const NekDouble> &xpoints,
113{
114 for (size_t i = 0; i < npts; ++i)
115 {
116 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
117 {
118 interp[i + j * npts] = Polylib::laginterp(
119 xpoints[i], j, &m_points[0][0], m_pointsKey.GetNumPoints());
120 }
121 }
122}
123
124std::shared_ptr<PolyEPoints::PointsBaseType> PolyEPoints::Create(
125 const PointsKey &key)
126{
127 std::shared_ptr<PointsBaseType> returnval(
129
130 returnval->Initialize();
131
132 return returnval;
133}
134
135const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
136 const PointsKey &pkey)
137{
138 ASSERTL0(pkey.GetPointsDim() == 1,
139 "Gauss Points can only interp to other 1d point distributions");
140
141 size_t numpoints = pkey.GetNumPoints();
143
144 PointsManager()[pkey]->GetPoints(xpoints);
145
146 return GetI(numpoints, xpoints);
147}
148
149const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
151{
152 size_t numpoints = 1;
153
154 /// Delegate to function below.
155 return GetI(numpoints, x);
156}
157
158const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
159 size_t numpoints, const Array<OneD, const NekDouble> &x)
160{
161 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
162
163 CalculateInterpMatrix(numpoints, x, interp);
164
165 size_t np = GetTotNumPoints();
166 NekDouble *d = interp.data();
167 std::shared_ptr<NekMatrix<NekDouble>> returnval(
168 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
169 d));
170
171 return returnval;
172}
173
174} // namespace Nektar::LibUtilities
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
bool RegisterCreator(const KeyType &key, const CreateFuncType &createFunc)
Register the given function and associate it with the key. The return value is just to facilitate cal...
Definition: NekManager.hpp:168
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:356
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:362
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:353
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:358
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:317
Defines a specification for a set of points.
Definition: Points.h:50
size_t GetPointsDim() const
Definition: Points.h:128
size_t GetNumPoints() const
Definition: Points.h:85
void CalculateInterpMatrix(size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
const std::shared_ptr< NekMatrix< NekDouble > > v_GetI(const PointsKey &pkey) override
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
PointsManagerT & PointsManager(void)
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:73
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
double NekDouble
double laginterpderiv(double z, int k, const double *zj, int np)
Definition: Polylib.cpp:98
double laginterp(double z, int j, const double *zj, int np)
Definition: Polylib.cpp:85