Nektar++
Loading...
Searching...
No Matches
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
38
40{
42 PointsManager().RegisterCreator(PointsKey(0, eGaussGaussLegendre),
44 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMLegendre),
46 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauPLegendre),
48 PointsManager().RegisterCreator(PointsKey(0, eGaussLobattoLegendre),
50 PointsManager().RegisterCreator(PointsKey(0, eGaussGaussChebyshev),
52 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMChebyshev),
54 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauPChebyshev),
56 PointsManager().RegisterCreator(PointsKey(0, eGaussLobattoChebyshev),
58 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta1),
60 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta2),
62 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha1Beta0),
64 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha2Beta0),
66 PointsManager().RegisterCreator(PointsKey(0, eGaussKronrodLegendre),
68 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauKronrodMLegendre),
70 PointsManager().RegisterCreator(
72 PointsManager().RegisterCreator(PointsKey(0, eGaussLobattoKronrodLegendre),
74 PointsManager().RegisterCreator(PointsKey(0, eGaussLegendreWithMP),
76 PointsManager().RegisterCreator(PointsKey(0, eGaussLegendreWithM),
78
80{
81 // Allocate the storage for points and weights
84
85 size_t numpoints = m_pointsKey.GetNumPoints();
86
87 switch (m_pointsKey.GetPointsType())
88 {
90 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
91 0.0);
92 break;
94 Polylib::zwgj(m_points[0].data() + 1, m_weights.data() + 1,
95 numpoints - 2, 0.0, 0.0);
96 m_points[0][0] = -1.0;
97 m_points[0][numpoints - 1] = 1.0;
98 m_weights[0] = 0.0;
99 m_weights[numpoints - 1] = 0.0;
100 break;
102 Polylib::zwgj(m_points[0].data() + 1, m_weights.data() + 1,
103 numpoints - 1, 0.0, 0.0);
104 m_points[0][0] = -1.0;
105 m_weights[0] = 0.0;
106 break;
108 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
109 0.0, 0.0);
110 break;
111
113 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
114 0.0, 0.0);
115 break;
116
118 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
119 0.0);
120 break;
121
123 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
124 -0.5);
125 break;
126
128 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
129 -0.5, -0.5);
130 break;
131
133 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
134 -0.5, -0.5);
135 break;
136
138 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
139 -0.5, -0.5);
140 break;
141
143 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
144 0.0, 1.0);
145 break;
146
147 case eGaussRadauMAlpha0Beta2:
148 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
149 0.0, 2.0);
150 break;
151
152 case eGaussRadauMAlpha1Beta0:
153 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
154 1.0, 0.0);
155 break;
156
157 case eGaussRadauMAlpha2Beta0:
158 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
159 2.0, 0.0);
160 break;
161
162 case eGaussKronrodLegendre:
163 Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
164 0.0);
165 break;
166
168 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
169 0.0);
170 break;
171
173 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
174 0.0);
175 break;
176
178 Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
179 0.0);
180 break;
181
182 default:
184 "Unknown Gauss quadrature point distribution requested");
185 }
186}
187
189{
190 // For Gauss Quadrature, this is done as part of the points computation
191}
192
194{
195 // Allocate the derivative matrix
196 PointsBaseType::v_CalculateDerivMatrix();
197
198 size_t numpoints = m_pointsKey.GetNumPoints();
199 size_t totpoints = m_pointsKey.GetTotNumPoints();
200 Array<OneD, NekDouble> dmtempSharedArray =
201 Array<OneD, NekDouble>(totpoints * totpoints);
202 NekDouble *dmtemp = dmtempSharedArray.data();
203
204 switch (m_pointsKey.GetPointsType())
205 {
207 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
208 break;
209
211 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
212 break;
213
215 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
216 break;
217
219 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
220 break;
221
223 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
224 break;
225
227 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
228 break;
229
231 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
232 break;
233
235 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
236 break;
237
239 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
240 break;
241
242 case eGaussRadauMAlpha0Beta2:
243 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
244 break;
245
246 case eGaussRadauMAlpha1Beta0:
247 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
248 break;
249
250 case eGaussRadauMAlpha2Beta0:
251 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
252 break;
253
256 case eGaussKronrodLegendre:
260 {
261 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
262 {
263 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
264 {
266 m_points[0][i], j, &m_points[0][0],
268 }
269 }
270 return;
271 }
272 break;
273
274 default:
276 "Unknown Gauss quadrature point distribution requested");
277 }
278
279 std::copy(dmtemp, dmtemp + totpoints * totpoints,
280 m_derivmatrix[0]->begin());
281}
282
284 size_t npts, const Array<OneD, const NekDouble> &xpoints,
286{
287 switch (m_pointsKey.GetPointsType())
288 {
290 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
291 GetNumPoints(), npts, 0.0, 0.0);
292 break;
293
295 // first / last line of interp are zeros:
296 for (int i = 0; i < interp.size(); ++i)
297 {
298 interp[i] = 0.0;
299 }
300 // only interpolate from interior points
301 Polylib::Imgj(interp.data() + npts, m_points[0].data() + 1,
302 xpoints.data(), GetNumPoints() - 2, npts, 0.0, 0.0);
303 break;
304
306 // first / last line of interp are zeros:
307 for (int i = 0; i < interp.size(); ++i)
308 {
309 interp[i] = 0.0;
310 }
311 // only interpolate from interior points
312 Polylib::Imgj(interp.data() + npts, m_points[0].data() + 1,
313 xpoints.data(), GetNumPoints() - 1, npts, 0.0, 0.0);
314 break;
315
317 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
318 GetNumPoints(), npts, 0.0, 0.0);
319 break;
320
322 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
323 GetNumPoints(), npts, 0.0, 0.0);
324 break;
325
327 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
328 GetNumPoints(), npts, 0.0, 0.0);
329 break;
330
332 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
333 GetNumPoints(), npts, -0.5, -0.5);
334 break;
335
337 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
338 GetNumPoints(), npts, -0.5, -0.5);
339 break;
340
342 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
343 GetNumPoints(), npts, -0.5, -0.5);
344 break;
345
347 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
348 GetNumPoints(), npts, -0.5, -0.5);
349 break;
350
352 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
353 GetNumPoints(), npts, 0.0, 1.0);
354 break;
355
356 case eGaussRadauMAlpha0Beta2:
357 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
358 GetNumPoints(), npts, 0.0, 2.0);
359 break;
360
361 case eGaussRadauMAlpha1Beta0:
362 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
363 GetNumPoints(), npts, 1.0, 0.0);
364 break;
365
366 case eGaussRadauMAlpha2Beta0:
367 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
368 GetNumPoints(), npts, 2.0, 0.0);
369 break;
370
371 case eGaussKronrodLegendre:
375 {
376 for (size_t i = 0; i < npts; ++i)
377 {
378 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
379 {
380 interp[i + j * npts] =
381 Polylib::laginterp(xpoints[i], j, &m_points[0][0],
383 }
384 }
385 }
386 break;
387
388 default:
390 "Unknown Gauss quadrature point distribution requested");
391 }
392}
393
394std::shared_ptr<Points<NekDouble>> GaussPoints::Create(const PointsKey &pkey)
395{
396 std::shared_ptr<Points<NekDouble>> returnval(
398
399 returnval->Initialize();
400
401 return returnval;
402}
403
404std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateMatrix(
405 const PointsKey &pkey)
406{
407 size_t numpoints = pkey.GetNumPoints();
409
410 PointsManager()[pkey]->GetPoints(xpoints);
411
412 // Delegate to function below
413 return GetI(numpoints, xpoints);
414}
415
416const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
417 const PointsKey &pkey)
418{
419 ASSERTL0(pkey.GetPointsDim() == 1,
420 "Gauss Points can only interp to other 1d point distributions");
421
422 return m_InterpManager[pkey];
423}
424
425const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
427{
428 size_t numpoints = 1;
429
430 // Delegate to function below
431 return GetI(numpoints, x);
432}
433
434const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
435 size_t numpoints, const Array<OneD, const NekDouble> &x)
436{
437 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
438
439 CalculateInterpMatrix(numpoints, x, interp);
440
441 NekDouble *t = interp.data();
442 size_t np = GetNumPoints();
443 std::shared_ptr<NekMatrix<NekDouble>> returnval(
444 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
445 t));
446
447 return returnval;
448}
449
450const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
452{
453 return m_GalerkinProjectionManager[pkey];
454}
455
456std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateGPMatrix(
457 const PointsKey &pkey)
458{
459 std::shared_ptr<NekMatrix<NekDouble>> returnval =
461
462 // Delegate to function below
463 return returnval;
464}
465
466std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
468{
469 size_t numpointsfrom = pkey.GetNumPoints();
470 size_t numpointsto = GetNumPoints();
471
473
474 weightsfrom = PointsManager()[pkey]->GetW();
475
476 std::shared_ptr<NekMatrix<NekDouble>> Interp = GetI(pkey);
477
478 Array<OneD, NekDouble> GalProj(numpointsfrom * numpointsto);
479
480 // set up inner product matrix and multiply by inverse of
481 // diagaonal mass matrix
482 for (size_t i = 0; i < numpointsto; ++i)
483 {
484 Vmath::Vmul(numpointsfrom, Interp->GetPtr().data() + i * numpointsfrom,
485 1, &weightsfrom[0], 1, &GalProj[0] + i, numpointsto);
486 if (m_weights[i] == 0.0)
487 {
488 Vmath::Zero(numpointsto, &GalProj[0] + i, numpointsto);
489 }
490 else
491 {
492 Vmath::Smul(numpointsfrom, 1.0 / m_weights[i], &GalProj[0] + i,
493 numpointsto, &GalProj[0] + i, numpointsto);
494 }
495 }
496
497 NekDouble *t = GalProj.data();
498 std::shared_ptr<NekMatrix<NekDouble>> returnval(
499 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(
500 numpointsto, numpointsfrom, t));
501
502 return returnval;
503}
504
505} // namespace Nektar::LibUtilities
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
static std::shared_ptr< Points< NekDouble > > Create(const PointsKey &pkey)
std::shared_ptr< NekMatrix< NekDouble > > CreateGPMatrix(const PointsKey &pkey)
const std::shared_ptr< NekMatrix< NekDouble > > v_GetGalerkinProjection(const PointsKey &pkey) override
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
const std::shared_ptr< NekMatrix< NekDouble > > v_GetI(const PointsKey &pkey) override
void CalculateInterpMatrix(size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
std::shared_ptr< NekMatrix< NekDouble > > CalculateGalerkinProjectionMatrix(const PointsKey &pkey)
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
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition Points.h:364
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition Points.h:366
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
PointsType GetPointsType() const
Definition Points.h:90
size_t GetTotNumPoints() const
Definition Points.h:158
size_t GetNumPoints() const
Definition Points.h:85
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
PointsManagerT & PointsManager(void)
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition PointsType.h:47
@ eGaussLegendreWithMP
1D Gauss-Legendre quadrature points with additional x=-1 and x=1 end points
Definition PointsType.h:95
@ eGaussRadauMAlpha0Beta1
Gauss Radau pinned at x=-1,.
Definition PointsType.h:58
@ eGaussRadauKronrodMLegendre
1D Gauss-Radau-Kronrod-Legendre quadrature points, pinned at x=-1
Definition PointsType.h:67
@ eGaussLobattoChebyshev
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:57
@ eGaussRadauPChebyshev
1D Gauss-Radau-Chebyshev quadrature points, pinned at x=1
Definition PointsType.h:55
@ 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:69
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition PointsType.h:51
@ eGaussLegendreWithM
1D Gauss-Legendre quadrature points with additional x=-1 point
Definition PointsType.h:97
@ eGaussGaussChebyshev
1D Gauss-Gauss-Chebyshev quadrature points
Definition PointsType.h:52
@ eGaussLobattoKronrodLegendre
1D Lobatto Kronrod quadrature points
Definition PointsType.h:72
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition PointsType.h:46
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition PointsType.h:49
double laginterpderiv(double z, int k, const double *zj, int np)
Definition Polylib.cpp:98
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:653
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:182
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:225
double laginterp(double z, int j, const double *zj, int np)
Definition Polylib.cpp:85
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:833
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:1116
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:266
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition Polylib.cpp:152
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:1054
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:1085
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition Polylib.cpp:317
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:704
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:768
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:389
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:1147
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:495
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.hpp:72
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.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition Vmath.hpp:273