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