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 namespace Nektar
46 {
47 namespace LibUtilities
48 {
68  PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta2),
70  PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha1Beta0),
72  PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha2Beta0),
74  PointsManager().RegisterCreator(PointsKey(0, eGaussKronrodLegendre),
82 
84 {
85  // Allocate the storage for points and weights
88 
89  int numpoints = m_pointsKey.GetNumPoints();
90 
91  switch (m_pointsKey.GetPointsType())
92  {
94  Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
95  0.0);
96  break;
97 
99  Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
100  0.0, 0.0);
101  break;
102 
104  Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
105  0.0, 0.0);
106  break;
107 
109  Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
110  0.0);
111  break;
112 
114  Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
115  -0.5);
116  break;
117 
119  Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
120  -0.5, -0.5);
121  break;
122 
124  Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
125  -0.5, -0.5);
126  break;
127 
129  Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
130  -0.5, -0.5);
131  break;
132 
134  Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
135  0.0, 1.0);
136  break;
137 
138  case eGaussRadauMAlpha0Beta2:
139  Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
140  0.0, 2.0);
141  break;
142 
143  case eGaussRadauMAlpha1Beta0:
144  Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
145  1.0, 0.0);
146  break;
147 
148  case eGaussRadauMAlpha2Beta0:
149  Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
150  2.0, 0.0);
151  break;
152 
153  case eGaussKronrodLegendre:
154  Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
155  0.0);
156  break;
157 
159  Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
160  0.0);
161  break;
162 
164  Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
165  0.0);
166  break;
167 
169  Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
170  0.0);
171  break;
172 
173  default:
175  "Unknown Gauss quadrature point distribution requested");
176  }
177 }
178 
180 {
181  // For Gauss Quadrature, this is done as part of the points computation
182 }
183 
185 {
186  // Allocate the derivative matrix
188 
189  int numpoints = m_pointsKey.GetNumPoints();
190  int totpoints = m_pointsKey.GetTotNumPoints();
191  Array<OneD, NekDouble> dmtempSharedArray =
192  Array<OneD, NekDouble>(totpoints * totpoints);
193  NekDouble *dmtemp = dmtempSharedArray.data();
194 
195  switch (m_pointsKey.GetPointsType())
196  {
197  case eGaussGaussLegendre:
198  Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
199  break;
200 
202  Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
203  break;
204 
206  Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
207  break;
208 
210  Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
211  break;
212 
214  Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
215  break;
216 
218  Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
219  break;
220 
222  Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
223  break;
224 
226  Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
227  break;
228 
230  Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
231  break;
232 
233  case eGaussRadauMAlpha0Beta2:
234  Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
235  break;
236 
237  case eGaussRadauMAlpha1Beta0:
238  Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
239  break;
240 
241  case eGaussRadauMAlpha2Beta0:
242  Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
243  break;
244 
245  case eGaussKronrodLegendre:
249  {
250  for (unsigned int i = 0; i < m_pointsKey.GetNumPoints(); ++i)
251  {
252  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
253  {
254  (*m_derivmatrix[0])(i, j) = LagrangePolyDeriv(
255  m_points[0][i], j, m_pointsKey.GetNumPoints(),
256  m_points[0]);
257  }
258  }
259  return;
260  }
261  break;
262 
263  default:
265  "Unknown Gauss quadrature point distribution requested");
266  }
267 
268  std::copy(dmtemp, dmtemp + totpoints * totpoints,
269  m_derivmatrix[0]->begin());
270 }
271 
273  unsigned int npts, const Array<OneD, const NekDouble> &xpoints,
274  Array<OneD, NekDouble> &interp)
275 {
276  switch (m_pointsKey.GetPointsType())
277  {
278  case eGaussGaussLegendre:
279  Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
280  GetNumPoints(), npts, 0.0, 0.0);
281  break;
282 
284  Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
285  GetNumPoints(), npts, 0.0, 0.0);
286  break;
287 
289  Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
290  GetNumPoints(), npts, 0.0, 0.0);
291  break;
292 
294  Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
295  GetNumPoints(), npts, 0.0, 0.0);
296  break;
297 
299  Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
300  GetNumPoints(), npts, -0.5, -0.5);
301  break;
302 
304  Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
305  GetNumPoints(), npts, -0.5, -0.5);
306  break;
307 
309  Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
310  GetNumPoints(), npts, -0.5, -0.5);
311  break;
312 
314  Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
315  GetNumPoints(), npts, -0.5, -0.5);
316  break;
317 
319  Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
320  GetNumPoints(), npts, 0.0, 1.0);
321  break;
322 
323  case eGaussRadauMAlpha0Beta2:
324  Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
325  GetNumPoints(), npts, 0.0, 2.0);
326  break;
327 
328  case eGaussRadauMAlpha1Beta0:
329  Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
330  GetNumPoints(), npts, 1.0, 0.0);
331  break;
332 
333  case eGaussRadauMAlpha2Beta0:
334  Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
335  GetNumPoints(), npts, 2.0, 0.0);
336  break;
337 
338  case eGaussKronrodLegendre:
342  {
343  for (unsigned int i = 0; i < npts; ++i)
344  {
345  for (unsigned int j = 0; j < m_pointsKey.GetNumPoints(); ++j)
346  {
347  interp[i + j * npts] = LagrangePoly(
348  xpoints[i], j, m_pointsKey.GetNumPoints(), m_points[0]);
349  }
350  }
351  }
352  break;
353 
354  default:
356  "Unknown Gauss quadrature point distribution requested");
357  }
358 }
359 
360 std::shared_ptr<Points<NekDouble>> GaussPoints::Create(const PointsKey &pkey)
361 {
362  std::shared_ptr<Points<NekDouble>> returnval(
364 
365  returnval->Initialize();
366 
367  return returnval;
368 }
369 
370 std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateMatrix(
371  const PointsKey &pkey)
372 {
373  int numpoints = pkey.GetNumPoints();
375 
376  PointsManager()[pkey]->GetPoints(xpoints);
377 
378  // Delegate to function below
379  return GetI(numpoints, xpoints);
380 }
381 
382 const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
383  const PointsKey &pkey)
384 {
385  ASSERTL0(pkey.GetPointsDim() == 1,
386  "Gauss Points can only interp to other 1d point distributions");
387 
388  return m_InterpManager[pkey];
389 }
390 
391 const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
393 {
394  int numpoints = 1;
395 
396  // Delegate to function below
397  return GetI(numpoints, x);
398 }
399 
400 const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
401  unsigned int numpoints, const Array<OneD, const NekDouble> &x)
402 {
403  Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
404 
405  CalculateInterpMatrix(numpoints, x, interp);
406 
407  NekDouble *t = interp.data();
408  unsigned int np = GetNumPoints();
409  std::shared_ptr<NekMatrix<NekDouble>> returnval(
410  MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
411  t));
412 
413  return returnval;
414 }
415 
417  NekDouble x, int npts, const Array<OneD, const NekDouble> &xpts,
418  const Array<OneD, const NekDouble> &funcvals)
419 {
420  NekDouble sum = 0.0;
421 
422  for (int i = 0; i < npts; ++i)
423  {
424  sum += funcvals[i] * LagrangePoly(x, i, npts, xpts);
425  }
426  return sum;
427 }
428 
430  const Array<OneD, const NekDouble> &xpts)
431 {
432  NekDouble h = 1.0;
433 
434  for (int i = 0; i < pt; ++i)
435  {
436  h = h * (x - xpts[i]) / (xpts[pt] - xpts[i]);
437  }
438 
439  for (int i = pt + 1; i < npts; ++i)
440  {
441  h = h * (x - xpts[i]) / (xpts[pt] - xpts[i]);
442  }
443 
444  return h;
445 }
446 
448  NekDouble x, int pt, int npts, const Array<OneD, const NekDouble> &xpts)
449 {
450  NekDouble h;
451  NekDouble y = 0.0;
452 
453  for (int j = 0; j < npts; ++j)
454  {
455  if (j != pt)
456  {
457  h = 1.0;
458  for (int i = 0; i < npts; ++i)
459  {
460  if (i != pt)
461  {
462  if (i != j)
463  {
464  h = h * (x - xpts[i]);
465  }
466  h = h / (xpts[pt] - xpts[i]);
467  }
468  }
469  y = y + h;
470  }
471  }
472  return y;
473 }
474 
475 const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
477 {
478  return m_GalerkinProjectionManager[pkey];
479 }
480 
481 std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateGPMatrix(
482  const PointsKey &pkey)
483 {
484  std::shared_ptr<NekMatrix<NekDouble>> returnval =
486 
487  // Delegate to function below
488  return returnval;
489 }
490 
491 std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
493 {
494  int numpointsfrom = pkey.GetNumPoints();
495  int numpointsto = GetNumPoints();
496 
497  Array<OneD, const NekDouble> weightsfrom;
498 
499  weightsfrom = PointsManager()[pkey]->GetW();
500 
501  std::shared_ptr<NekMatrix<NekDouble>> Interp = GetI(pkey);
502 
503  Array<OneD, NekDouble> GalProj(numpointsfrom * numpointsto);
504 
505  // set up inner product matrix and multiply by inverse of
506  // diagaonal mass matrix
507  for (int i = 0; i < numpointsto; ++i)
508  {
509  Vmath::Vmul(numpointsfrom, Interp->GetPtr().get() + i * numpointsfrom,
510  1, &weightsfrom[0], 1, &GalProj[0] + i, numpointsto);
511  Vmath::Smul(numpointsfrom, 1.0 / m_weights[i], &GalProj[0] + i,
512  numpointsto, &GalProj[0] + i, numpointsto);
513  }
514 
515  NekDouble *t = GalProj.data();
516  std::shared_ptr<NekMatrix<NekDouble>> returnval(
517  MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(
518  numpointsto, numpointsfrom, t));
519 
520  return returnval;
521 }
522 
523 } // end of namespace LibUtilities
524 } // end of namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#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)
std::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix(const PointsKey &pkey)
NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetGalerkinProjection(const PointsKey &pkey) override
NekDouble LagrangeInterpolant(NekDouble x, int npts, const Array< OneD, const NekDouble > &xpts, const Array< OneD, const NekDouble > &funcvals)
functions used by the Kronrod points
virtual void v_CalculateWeights() override
virtual void v_CalculatePoints() override
Definition: GaussPoints.cpp:83
void CalculateInterpMatrix(unsigned int npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetI(const PointsKey &pkey) override
NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array< OneD, const NekDouble > &xpts)
std::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix(const PointsKey &pkey)
virtual void v_CalculateDerivMatrix() override
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
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:383
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition: Points.h:385
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:372
unsigned int GetNumPoints() const
Definition: Points.h:273
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 GetTotNumPoints() const
Definition: Points.h:177
PointsType GetPointsType() const
Definition: Points.h:109
unsigned int GetPointsDim() const
Definition: Points.h:147
unsigned int GetNumPoints() const
Definition: Points.h:104
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
@ eGaussRadauMAlpha0Beta1
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:60
@ eGaussRadauKronrodMLegendre
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:69
@ eGaussLobattoChebyshev
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:59
@ eGaussRadauPChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition: PointsType.h:57
@ eGaussRadauMChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=-1
Definition: PointsType.h:55
@ eGaussRadauKronrodMAlpha1Beta0
1D Gauss-Radau-Kronrod-Legendre pinned at x=-1,
Definition: PointsType.h:71
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ eGaussGaussChebyshev
1D Gauss-Gauss-Chebyshev quadrature points
Definition: PointsType.h:54
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
Definition: PointsType.h:74
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:51
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:627
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:161
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:202
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:787
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:1054
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:241
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:133
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:992
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:1023
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:290
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:674
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:730
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:362
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:1085
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:468
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:209
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:248