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 
40 
41 namespace Nektar
42 {
43 namespace LibUtilities
44 {
46  PointsKey(0, ePolyEvenlySpaced), PolyEPoints::Create)};
47 
49 {
50  // Allocate the storage for points
52 
53  unsigned int npts = m_pointsKey.GetNumPoints();
54  if (npts == 1)
55  {
56  m_points[0][0] = 0.0;
57  }
58  else
59  {
60  NekDouble dx = 2.0 / (NekDouble)(npts - 1);
61  for (unsigned int i = 0; i < npts; ++i)
62  {
63  m_points[0][i] = -1.0 + i * dx;
64  }
65  }
66 }
67 
69 {
70  // Allocate the storage for points
72 
73  unsigned int npts = m_pointsKey.GetNumPoints();
74  if (npts == 1)
75  {
76  m_weights[0] = 2.0; // midpoint rule
77  }
78  else
79  {
80  PointsKey gaussKey(npts, eGaussLobattoLegendre);
81  std::shared_ptr<PointsBaseType> ptr = PointsManager()[gaussKey];
84 
85  ptr->GetZW(z, w);
86  for (unsigned int i = 0; i < npts; ++i)
87  {
88  m_weights[i] = 0.0;
89  for (unsigned j = 0; j < npts; ++j)
90  {
91  m_weights[i] += w[j] * LagrangePoly(z[j], i, npts, m_points[0]);
92  }
93  }
94  }
95 }
96 
98 {
99  // Allocate the derivative matrix.
101 
102  for (unsigned int i = 0; i < m_pointsKey.GetNumPoints(); ++i)
103  {
104  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
105  {
106  (*m_derivmatrix[0])(i, j) = LagrangePolyDeriv(
107  m_points[0][i], j, m_pointsKey.GetNumPoints(), m_points[0]);
108  }
109  }
110 }
111 
113  unsigned int npts, const Array<OneD, const NekDouble> &xpoints,
114  Array<OneD, NekDouble> &interp)
115 {
116  for (unsigned int i = 0; i < npts; ++i)
117  {
118  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
119  {
120  interp[i + j * npts] = LagrangePoly(
121  xpoints[i], j, m_pointsKey.GetNumPoints(), m_points[0]);
122  }
123  }
124 }
125 
126 std::shared_ptr<PolyEPoints::PointsBaseType> PolyEPoints::Create(
127  const PointsKey &key)
128 {
129  std::shared_ptr<PointsBaseType> returnval(
131 
132  returnval->Initialize();
133 
134  return returnval;
135 }
136 
137 const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
138  const PointsKey &pkey)
139 {
140  ASSERTL0(pkey.GetPointsDim() == 1,
141  "Gauss Points can only interp to other 1d point distributions");
142 
143  int numpoints = pkey.GetNumPoints();
145 
146  PointsManager()[pkey]->GetPoints(xpoints);
147 
148  return GetI(numpoints, xpoints);
149 }
150 
151 const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
153 {
154  int numpoints = 1;
155 
156  /// Delegate to function below.
157  return GetI(numpoints, x);
158 }
159 
160 const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
161  unsigned int numpoints, const Array<OneD, const NekDouble> &x)
162 {
163  Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
164 
165  CalculateInterpMatrix(numpoints, x, interp);
166 
167  unsigned int np = GetTotNumPoints();
168  NekDouble *d = interp.data();
169  std::shared_ptr<NekMatrix<NekDouble>> returnval(
170  MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
171  d));
172 
173  return returnval;
174 }
175 
177  NekDouble x, int npts, const Array<OneD, const NekDouble> &xpts,
178  const Array<OneD, const NekDouble> &funcvals)
179 {
180  NekDouble sum = 0.0;
181 
182  for (int i = 0; i < npts; ++i)
183  {
184  sum += funcvals[i] * LagrangePoly(x, i, npts, xpts);
185  }
186  return sum;
187 }
188 
190  const Array<OneD, const NekDouble> &xpts)
191 {
192  NekDouble h = 1.0;
193 
194  for (int i = 0; i < pt; ++i)
195  {
196  h = h * (x - xpts[i]) / (xpts[pt] - xpts[i]);
197  }
198 
199  for (int i = pt + 1; i < npts; ++i)
200  {
201  h = h * (x - xpts[i]) / (xpts[pt] - xpts[i]);
202  }
203 
204  return h;
205 }
206 
208  NekDouble x, int pt, int npts, const Array<OneD, const NekDouble> &xpts)
209 {
210  NekDouble h;
211  NekDouble y = 0.0;
212 
213  for (int j = 0; j < npts; ++j)
214  {
215  if (j != pt)
216  {
217  h = 1.0;
218  for (int i = 0; i < npts; ++i)
219  {
220  if (i != pt)
221  {
222  if (i != j)
223  {
224  h = h * (x - xpts[i]);
225  }
226  h = h / (xpts[pt] - xpts[i]);
227  }
228  }
229  y = y + h;
230  }
231  }
232  return y;
233 }
234 } // end of namespace LibUtilities
235 } // end of namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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:170
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:375
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:381
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:372
unsigned int GetNumPoints() const
Definition: Points.h:273
unsigned int GetTotNumPoints() const
Definition: Points.h:278
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:377
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:336
Defines a specification for a set of points.
Definition: Points.h:59
unsigned int GetPointsDim() const
Definition: Points.h:147
unsigned int GetNumPoints() const
Definition: Points.h:104
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
virtual void v_CalculatePoints() override
Definition: PolyEPoints.cpp:48
virtual void v_CalculateDerivMatrix() override
Definition: PolyEPoints.cpp:97
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetI(const PointsKey &pkey) override
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
virtual void v_CalculateWeights() override
Definition: PolyEPoints.cpp:68
NekDouble LagrangeInterpolant(NekDouble x, int npts, const Array< OneD, const NekDouble > &xpts, const Array< OneD, const NekDouble > &funcvals)
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:53
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble