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
38
40{
60 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta2),
62 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha1Beta0),
64 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha2Beta0),
66 PointsManager().RegisterCreator(PointsKey(0, eGaussKronrodLegendre),
74
76{
77 // Allocate the storage for points and weights
80
81 size_t numpoints = m_pointsKey.GetNumPoints();
82
83 switch (m_pointsKey.GetPointsType())
84 {
86 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
87 0.0);
88 break;
89
91 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
92 0.0, 0.0);
93 break;
94
96 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
97 0.0, 0.0);
98 break;
99
101 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
102 0.0);
103 break;
104
106 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
107 -0.5);
108 break;
109
111 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
112 -0.5, -0.5);
113 break;
114
116 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
117 -0.5, -0.5);
118 break;
119
121 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
122 -0.5, -0.5);
123 break;
124
126 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
127 0.0, 1.0);
128 break;
129
130 case eGaussRadauMAlpha0Beta2:
131 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
132 0.0, 2.0);
133 break;
134
135 case eGaussRadauMAlpha1Beta0:
136 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
137 1.0, 0.0);
138 break;
139
140 case eGaussRadauMAlpha2Beta0:
141 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
142 2.0, 0.0);
143 break;
144
145 case eGaussKronrodLegendre:
146 Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
147 0.0);
148 break;
149
151 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
152 0.0);
153 break;
154
156 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
157 0.0);
158 break;
159
161 Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
162 0.0);
163 break;
164
165 default:
167 "Unknown Gauss quadrature point distribution requested");
168 }
169}
170
172{
173 // For Gauss Quadrature, this is done as part of the points computation
174}
175
177{
178 // Allocate the derivative matrix
179 PointsBaseType::v_CalculateDerivMatrix();
180
181 size_t numpoints = m_pointsKey.GetNumPoints();
182 size_t totpoints = m_pointsKey.GetTotNumPoints();
183 Array<OneD, NekDouble> dmtempSharedArray =
184 Array<OneD, NekDouble>(totpoints * totpoints);
185 NekDouble *dmtemp = dmtempSharedArray.data();
186
187 switch (m_pointsKey.GetPointsType())
188 {
190 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
191 break;
192
194 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
195 break;
196
198 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
199 break;
200
202 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
203 break;
204
206 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
207 break;
208
210 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
211 break;
212
214 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
215 break;
216
218 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
219 break;
220
222 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
223 break;
224
225 case eGaussRadauMAlpha0Beta2:
226 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
227 break;
228
229 case eGaussRadauMAlpha1Beta0:
230 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
231 break;
232
233 case eGaussRadauMAlpha2Beta0:
234 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
235 break;
236
237 case eGaussKronrodLegendre:
241 {
242 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
243 {
244 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
245 {
247 m_points[0][i], j, &m_points[0][0],
249 }
250 }
251 return;
252 }
253 break;
254
255 default:
257 "Unknown Gauss quadrature point distribution requested");
258 }
259
260 std::copy(dmtemp, dmtemp + totpoints * totpoints,
261 m_derivmatrix[0]->begin());
262}
263
265 size_t npts, const Array<OneD, const NekDouble> &xpoints,
267{
268 switch (m_pointsKey.GetPointsType())
269 {
271 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
272 GetNumPoints(), npts, 0.0, 0.0);
273 break;
274
276 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
277 GetNumPoints(), npts, 0.0, 0.0);
278 break;
279
281 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
282 GetNumPoints(), npts, 0.0, 0.0);
283 break;
284
286 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
287 GetNumPoints(), npts, 0.0, 0.0);
288 break;
289
291 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
292 GetNumPoints(), npts, -0.5, -0.5);
293 break;
294
296 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
297 GetNumPoints(), npts, -0.5, -0.5);
298 break;
299
301 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
302 GetNumPoints(), npts, -0.5, -0.5);
303 break;
304
306 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
307 GetNumPoints(), npts, -0.5, -0.5);
308 break;
309
311 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
312 GetNumPoints(), npts, 0.0, 1.0);
313 break;
314
315 case eGaussRadauMAlpha0Beta2:
316 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
317 GetNumPoints(), npts, 0.0, 2.0);
318 break;
319
320 case eGaussRadauMAlpha1Beta0:
321 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
322 GetNumPoints(), npts, 1.0, 0.0);
323 break;
324
325 case eGaussRadauMAlpha2Beta0:
326 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
327 GetNumPoints(), npts, 2.0, 0.0);
328 break;
329
330 case eGaussKronrodLegendre:
334 {
335 for (size_t i = 0; i < npts; ++i)
336 {
337 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
338 {
339 interp[i + j * npts] =
340 Polylib::laginterp(xpoints[i], j, &m_points[0][0],
342 }
343 }
344 }
345 break;
346
347 default:
349 "Unknown Gauss quadrature point distribution requested");
350 }
351}
352
353std::shared_ptr<Points<NekDouble>> GaussPoints::Create(const PointsKey &pkey)
354{
355 std::shared_ptr<Points<NekDouble>> returnval(
357
358 returnval->Initialize();
359
360 return returnval;
361}
362
363std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateMatrix(
364 const PointsKey &pkey)
365{
366 size_t numpoints = pkey.GetNumPoints();
368
369 PointsManager()[pkey]->GetPoints(xpoints);
370
371 // Delegate to function below
372 return GetI(numpoints, xpoints);
373}
374
375const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
376 const PointsKey &pkey)
377{
378 ASSERTL0(pkey.GetPointsDim() == 1,
379 "Gauss Points can only interp to other 1d point distributions");
380
381 return m_InterpManager[pkey];
382}
383
384const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
386{
387 size_t numpoints = 1;
388
389 // Delegate to function below
390 return GetI(numpoints, x);
391}
392
393const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
394 size_t numpoints, const Array<OneD, const NekDouble> &x)
395{
396 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
397
398 CalculateInterpMatrix(numpoints, x, interp);
399
400 NekDouble *t = interp.data();
401 size_t np = GetNumPoints();
402 std::shared_ptr<NekMatrix<NekDouble>> returnval(
403 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
404 t));
405
406 return returnval;
407}
408
409const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
411{
412 return m_GalerkinProjectionManager[pkey];
413}
414
415std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateGPMatrix(
416 const PointsKey &pkey)
417{
418 std::shared_ptr<NekMatrix<NekDouble>> returnval =
420
421 // Delegate to function below
422 return returnval;
423}
424
425std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
427{
428 size_t numpointsfrom = pkey.GetNumPoints();
429 size_t numpointsto = GetNumPoints();
430
432
433 weightsfrom = PointsManager()[pkey]->GetW();
434
435 std::shared_ptr<NekMatrix<NekDouble>> Interp = GetI(pkey);
436
437 Array<OneD, NekDouble> GalProj(numpointsfrom * numpointsto);
438
439 // set up inner product matrix and multiply by inverse of
440 // diagaonal mass matrix
441 for (size_t i = 0; i < numpointsto; ++i)
442 {
443 Vmath::Vmul(numpointsfrom, Interp->GetPtr().get() + i * numpointsfrom,
444 1, &weightsfrom[0], 1, &GalProj[0] + i, numpointsto);
445 Vmath::Smul(numpointsfrom, 1.0 / m_weights[i], &GalProj[0] + i,
446 numpointsto, &GalProj[0] + i, numpointsto);
447 }
448
449 NekDouble *t = GalProj.data();
450 std::shared_ptr<NekMatrix<NekDouble>> returnval(
451 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(
452 numpointsto, numpointsfrom, t));
453
454 return returnval;
455}
456
457} // namespace Nektar::LibUtilities
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
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)
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:168
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...
def copy(self)
Definition: pycml.py:2663
PointsManagerT & PointsManager(void)
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:47
@ 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
@ 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 NekDouble
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