Nektar++
GaussPoints.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File GaussPoints.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: GaussPoints Definitions
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 
44 
45 
46 namespace Nektar
47 {
48  namespace LibUtilities
49  {
67  };
68 
70  {
71  // Allocate the storage for points and weights
74 
75  int numpoints = m_pointsKey.GetNumPoints();
76 
77  switch(m_pointsKey.GetPointsType())
78  {
80  Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
81  break;
82 
84  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
85  break;
86 
88  Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
89  break;
90 
92  Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
93  break;
94 
96  Polylib::zwgj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
97  break;
98 
100  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
101  break;
102 
104  Polylib::zwgrjp(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
105  break;
106 
108  Polylib::zwglj(m_points[0].data(),m_weights.data(),numpoints,-0.5,-0.5);
109  break;
110 
112  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,1.0);
113  break;
114 
116  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,0.0,2.0);
117  break;
118 
120  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
121  break;
122 
124  Polylib::zwgrjm(m_points[0].data(),m_weights.data(),numpoints,2.0,0.0);
125  break;
126 
128  Polylib::zwgk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
129  break;
130 
132  Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
133  break;
134 
136  Polylib::zwrk(m_points[0].data(),m_weights.data(),numpoints,1.0,0.0);
137  break;
138 
140  Polylib::zwlk(m_points[0].data(),m_weights.data(),numpoints,0.0,0.0);
141  break;
142 
143  default:
145  "Unknown Gauss quadrature point distribution requested");
146  }
147  }
148 
150  {
151  // For Gauss Quadrature, this is done as part of the points computation
152  }
153 
155  {
156  // Allocate the derivative matrix
158 
159  int numpoints = m_pointsKey.GetNumPoints();
160  int totpoints = m_pointsKey.GetTotNumPoints();
161  Array<OneD, NekDouble> dmtempSharedArray = Array<OneD, NekDouble>(totpoints*totpoints);
162  NekDouble *dmtemp = dmtempSharedArray.data();
163 
164  switch(m_pointsKey.GetPointsType())
165  {
166  case eGaussGaussLegendre:
167  Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
168  break;
169 
171  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
172  break;
173 
175  Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
176  break;
177 
179  Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,0.0,0.0);
180  break;
181 
183  Polylib::Dgj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
184  break;
185 
187  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
188  break;
189 
191  Polylib::Dgrjp(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
192  break;
193 
195  Polylib::Dglj(dmtemp,m_points[0].data(),numpoints,-0.5,-0.5);
196  break;
197 
199  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,1.0);
200  break;
201 
203  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,0.0,2.0);
204  break;
205 
207  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,1.0,0.0);
208  break;
209 
211  Polylib::Dgrjm(dmtemp,m_points[0].data(),numpoints,2.0,0.0);
212  break;
213 
218  {
219  for(unsigned int i=0;i<m_pointsKey.GetNumPoints();++i)
220  {
221  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
222  {
224  }
225  }
226  return;
227  }
228  break;
229 
230  default:
232  "Unknown Gauss quadrature point distribution requested");
233  }
234 
235  std::copy(dmtemp,dmtemp+totpoints*totpoints,m_derivmatrix[0]->begin());
236  }
237 
239  {
240  switch(m_pointsKey.GetPointsType())
241  {
242  case eGaussGaussLegendre:
243  Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
244  break;
245 
247  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
248  break;
249 
251  Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
252  break;
253 
255  Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,0.0);
256  break;
257 
259  Polylib::Imgj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
260  break;
261 
263  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
264  break;
265 
267  Polylib::Imgrjp(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
268  break;
269 
271  Polylib::Imglj(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,-0.5,-0.5);
272  break;
273 
275  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,1.0);
276  break;
277 
279  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,0.0,2.0);
280  break;
281 
283  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,1.0,0.0);
284  break;
285 
287  Polylib::Imgrjm(interp.data(),m_points[0].data(),xpoints.data(),GetNumPoints(),npts,2.0,0.0);
288  break;
289 
294  {
295  for(unsigned int i=0;i<npts;++i)
296  {
297  for(unsigned int j=0;j<m_pointsKey.GetNumPoints();++j)
298  {
299  interp[i + j*npts] = LagrangePoly(xpoints[i],j,m_pointsKey.GetNumPoints(),m_points[0]);
300  }
301  }
302  }
303  break;
304 
305  default:
307  "Unknown Gauss quadrature point distribution requested");
308  }
309  }
310 
311  std::shared_ptr<Points<NekDouble> > GaussPoints::Create(const PointsKey &pkey)
312  {
313  std::shared_ptr< Points<NekDouble> > returnval(MemoryManager< GaussPoints >::AllocateSharedPtr(pkey));
314 
315  returnval->Initialize();
316 
317  return returnval;
318  }
319 
320  std::shared_ptr< NekMatrix<NekDouble> > GaussPoints::CreateMatrix(const PointsKey &pkey)
321  {
322  int numpoints = pkey.GetNumPoints();
324 
325  PointsManager()[pkey]->GetPoints(xpoints);
326 
327  // Delegate to function below
328  return GetI(numpoints, xpoints);
329  }
330 
331 
332  const std::shared_ptr<NekMatrix<NekDouble> > GaussPoints::GetI(const PointsKey &pkey)
333  {
334  ASSERTL0(pkey.GetPointsDim()==1, "Gauss Points can only interp to other 1d point distributions");
335 
336  return m_InterpManager[pkey];
337  }
338 
339  const std::shared_ptr<NekMatrix<NekDouble> > GaussPoints::GetI(const Array<OneD, const NekDouble>& x)
340  {
341  int numpoints = 1;
342 
343  // Delegate to function below
344  return GetI(numpoints, x);
345  }
346 
347  const std::shared_ptr<NekMatrix<NekDouble> > GaussPoints::GetI(unsigned int numpoints, const Array<OneD, const NekDouble>& x)
348  {
349  Array<OneD, NekDouble> interp(GetNumPoints()*numpoints);
350 
351  CalculateInterpMatrix(numpoints, x, interp);
352 
353  NekDouble* t = interp.data();
354  unsigned int np = GetNumPoints();
355  std::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpoints,np,t));
356 
357  return returnval;
358  }
359 
361  const Array<OneD, const NekDouble>& funcvals)
362  {
363  NekDouble sum = 0.0;
364 
365  for(int i=0;i<npts;++i)
366  {
367  sum += funcvals[i]*LagrangePoly(x,i,npts,xpts);
368  }
369  return sum;
370  }
371 
372 
374  {
375  NekDouble h=1.0;
376 
377  for(int i=0;i<pt; ++i)
378  {
379  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
380  }
381 
382  for(int i=pt+1;i<npts;++i)
383  {
384  h = h * (x - xpts[i])/(xpts[pt]-xpts[i]);
385  }
386 
387  return h;
388  }
389 
391  {
392  NekDouble h;
393  NekDouble y=0.0;
394 
395  for(int j=0;j<npts;++j)
396  {
397  if(j!=pt)
398  {
399  h=1.0;
400  for(int i=0;i<npts;++i)
401  {
402  if(i!=pt)
403  {
404  if(i!=j)
405  {
406  h = h*(x-xpts[i]);
407  }
408  h = h/(xpts[pt]-xpts[i]);
409  }
410  }
411  y = y + h;
412  }
413  }
414  return y;
415  }
416 
417  const std::shared_ptr<NekMatrix<NekDouble> > GaussPoints::GetGalerkinProjection(const PointsKey &pkey)
418  {
419  return m_GalerkinProjectionManager[pkey];
420  }
421 
422  std::shared_ptr< NekMatrix<NekDouble> > GaussPoints::CreateGPMatrix(const PointsKey &pkey)
423  {
424  std::shared_ptr< NekMatrix<NekDouble> > returnval = CalculateGalerkinProjectionMatrix(pkey);
425 
426  // Delegate to function below
427  return returnval;
428  }
429 
430  std::shared_ptr<NekMatrix<NekDouble> > GaussPoints::CalculateGalerkinProjectionMatrix(const PointsKey &pkey)
431  {
432  int numpointsfrom = pkey.GetNumPoints();
433  int numpointsto = GetNumPoints();
434 
435  Array<OneD, const NekDouble> weightsfrom;
436 
437  weightsfrom = PointsManager()[pkey]->GetW();
438 
439  std::shared_ptr< NekMatrix<NekDouble> > Interp = GetI(pkey);
440 
441  Array<OneD, NekDouble> GalProj(numpointsfrom*numpointsto);
442 
443  // set up inner product matrix and multiply by inverse of
444  // diagaonal mass matrix
445  for(int i = 0; i < numpointsto; ++i)
446  {
447  Vmath::Vmul(numpointsfrom,Interp->GetPtr().get() +i*numpointsfrom,1,
448  &weightsfrom[0],1,&GalProj[0] +i,numpointsto);
449  Vmath::Smul(numpointsfrom,1.0/m_weights[i],&GalProj[0]+i,numpointsto,
450  &GalProj[0]+i,numpointsto);
451  }
452 
453 
454  NekDouble* t = GalProj.data();
455  std::shared_ptr< NekMatrix<NekDouble> > returnval(MemoryManager<NekMatrix<NekDouble> >::AllocateSharedPtr(numpointsto,numpointsfrom,t));
456 
457  return returnval;
458  }
459 
460  } // end of namespace LibUtilities
461 } // end of namespace Nektar
462 
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
static std::shared_ptr< Points< NekDouble > > Create(const PointsKey &pkey)
const std::shared_ptr< NekMatrix< NekDouble > > GetI(const PointsKey &pkey)
std::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix(const PointsKey &pkey)
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
NekDouble LagrangeInterpolant(NekDouble x, int npts, const Array< OneD, const NekDouble > &xpts, const Array< OneD, const NekDouble > &funcvals)
functions used by the Kronrod points
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
const std::shared_ptr< NekMatrix< NekDouble > > GetGalerkinProjection(const PointsKey &pkey)
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
std::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix(const PointsKey &pkey)
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
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:397
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition: Points.h:398
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:387
unsigned int GetNumPoints() const
Definition: Points.h:273
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 GetTotNumPoints() const
Definition: Points.h:180
PointsType GetPointsType() const
Definition: Points.h:112
unsigned int GetPointsDim() const
Definition: Points.h:150
unsigned int GetNumPoints() const
Definition: Points.h:107
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
def copy(self)
Definition: pycml.py:2663
PointsManagerT & PointsManager(void)
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:49
@ eGaussRadauMAlpha0Beta2
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
@ eGaussKronrodLegendre
1D Gauss-Kronrod-Legendre quadrature points
Definition: PointsType.h:60
@ eGaussRadauMAlpha0Beta1
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:56
@ eGaussRadauKronrodMLegendre
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:61
@ eGaussLobattoChebyshev
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:55
@ eGaussRadauPChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:54
@ eGaussRadauMAlpha2Beta0
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:59
@ eGaussRadauMChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:53
@ eGaussRadauKronrodMAlpha1Beta0
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition: PointsType.h:62
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
@ eGaussGaussChebyshev
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:52
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:63
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
@ eGaussRadauMAlpha1Beta0
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:50
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble
void Dgj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.
Definition: Polylib.cpp:546
void zwgrjm(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition: Polylib.cpp:119
void zwgrjp(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition: Polylib.cpp:158
void Dglj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros.
Definition: Polylib.cpp:690
void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at ...
Definition: Polylib.cpp:946
void zwglj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition: Polylib.cpp:195
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:91
void Imgj(double *im, const double *zgj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Jacobi points to an arbitrary distribution at points zm.
Definition: Polylib.cpp:887
void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at...
Definition: Polylib.cpp:916
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:240
void Dgrjm(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a z...
Definition: Polylib.cpp:589
void Dgrjp(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the Gauss-Radau-Jacobi zeros with a zero at z=1.
Definition: Polylib.cpp:639
void zwrk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Radau-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:313
void Imglj(double *im, const double *zglj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm.
Definition: Polylib.cpp:977
void zwlk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Lobatto-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:420
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:192
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225