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
39namespace Nektar
40{
41namespace LibUtilities
42{
62 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha0Beta2),
64 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha1Beta0),
66 PointsManager().RegisterCreator(PointsKey(0, eGaussRadauMAlpha2Beta0),
68 PointsManager().RegisterCreator(PointsKey(0, eGaussKronrodLegendre),
76
78{
79 // Allocate the storage for points and weights
82
83 size_t numpoints = m_pointsKey.GetNumPoints();
84
85 switch (m_pointsKey.GetPointsType())
86 {
88 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
89 0.0);
90 break;
91
93 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
94 0.0, 0.0);
95 break;
96
98 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
99 0.0, 0.0);
100 break;
101
103 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints, 0.0,
104 0.0);
105 break;
106
108 Polylib::zwgj(m_points[0].data(), m_weights.data(), numpoints, -0.5,
109 -0.5);
110 break;
111
113 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
114 -0.5, -0.5);
115 break;
116
118 Polylib::zwgrjp(m_points[0].data(), m_weights.data(), numpoints,
119 -0.5, -0.5);
120 break;
121
123 Polylib::zwglj(m_points[0].data(), m_weights.data(), numpoints,
124 -0.5, -0.5);
125 break;
126
128 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
129 0.0, 1.0);
130 break;
131
132 case eGaussRadauMAlpha0Beta2:
133 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
134 0.0, 2.0);
135 break;
136
137 case eGaussRadauMAlpha1Beta0:
138 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
139 1.0, 0.0);
140 break;
141
142 case eGaussRadauMAlpha2Beta0:
143 Polylib::zwgrjm(m_points[0].data(), m_weights.data(), numpoints,
144 2.0, 0.0);
145 break;
146
147 case eGaussKronrodLegendre:
148 Polylib::zwgk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
149 0.0);
150 break;
151
153 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
154 0.0);
155 break;
156
158 Polylib::zwrk(m_points[0].data(), m_weights.data(), numpoints, 1.0,
159 0.0);
160 break;
161
163 Polylib::zwlk(m_points[0].data(), m_weights.data(), numpoints, 0.0,
164 0.0);
165 break;
166
167 default:
169 "Unknown Gauss quadrature point distribution requested");
170 }
171}
172
174{
175 // For Gauss Quadrature, this is done as part of the points computation
176}
177
179{
180 // Allocate the derivative matrix
181 PointsBaseType::v_CalculateDerivMatrix();
182
183 size_t numpoints = m_pointsKey.GetNumPoints();
184 size_t totpoints = m_pointsKey.GetTotNumPoints();
185 Array<OneD, NekDouble> dmtempSharedArray =
186 Array<OneD, NekDouble>(totpoints * totpoints);
187 NekDouble *dmtemp = dmtempSharedArray.data();
188
189 switch (m_pointsKey.GetPointsType())
190 {
192 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
193 break;
194
196 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
197 break;
198
200 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
201 break;
202
204 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, 0.0, 0.0);
205 break;
206
208 Polylib::Dgj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
209 break;
210
212 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
213 break;
214
216 Polylib::Dgrjp(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
217 break;
218
220 Polylib::Dglj(dmtemp, m_points[0].data(), numpoints, -0.5, -0.5);
221 break;
222
224 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 1.0);
225 break;
226
227 case eGaussRadauMAlpha0Beta2:
228 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 0.0, 2.0);
229 break;
230
231 case eGaussRadauMAlpha1Beta0:
232 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 1.0, 0.0);
233 break;
234
235 case eGaussRadauMAlpha2Beta0:
236 Polylib::Dgrjm(dmtemp, m_points[0].data(), numpoints, 2.0, 0.0);
237 break;
238
239 case eGaussKronrodLegendre:
243 {
244 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
245 {
246 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
247 {
249 m_points[0][i], j, &m_points[0][0],
251 }
252 }
253 return;
254 }
255 break;
256
257 default:
259 "Unknown Gauss quadrature point distribution requested");
260 }
261
262 std::copy(dmtemp, dmtemp + totpoints * totpoints,
263 m_derivmatrix[0]->begin());
264}
265
267 size_t npts, const Array<OneD, const NekDouble> &xpoints,
269{
270 switch (m_pointsKey.GetPointsType())
271 {
273 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
274 GetNumPoints(), npts, 0.0, 0.0);
275 break;
276
278 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
279 GetNumPoints(), npts, 0.0, 0.0);
280 break;
281
283 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
284 GetNumPoints(), npts, 0.0, 0.0);
285 break;
286
288 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
289 GetNumPoints(), npts, 0.0, 0.0);
290 break;
291
293 Polylib::Imgj(interp.data(), m_points[0].data(), xpoints.data(),
294 GetNumPoints(), npts, -0.5, -0.5);
295 break;
296
298 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
299 GetNumPoints(), npts, -0.5, -0.5);
300 break;
301
303 Polylib::Imgrjp(interp.data(), m_points[0].data(), xpoints.data(),
304 GetNumPoints(), npts, -0.5, -0.5);
305 break;
306
308 Polylib::Imglj(interp.data(), m_points[0].data(), xpoints.data(),
309 GetNumPoints(), npts, -0.5, -0.5);
310 break;
311
313 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
314 GetNumPoints(), npts, 0.0, 1.0);
315 break;
316
317 case eGaussRadauMAlpha0Beta2:
318 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
319 GetNumPoints(), npts, 0.0, 2.0);
320 break;
321
322 case eGaussRadauMAlpha1Beta0:
323 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
324 GetNumPoints(), npts, 1.0, 0.0);
325 break;
326
327 case eGaussRadauMAlpha2Beta0:
328 Polylib::Imgrjm(interp.data(), m_points[0].data(), xpoints.data(),
329 GetNumPoints(), npts, 2.0, 0.0);
330 break;
331
332 case eGaussKronrodLegendre:
336 {
337 for (size_t i = 0; i < npts; ++i)
338 {
339 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
340 {
341 interp[i + j * npts] =
342 Polylib::laginterp(xpoints[i], j, &m_points[0][0],
344 }
345 }
346 }
347 break;
348
349 default:
351 "Unknown Gauss quadrature point distribution requested");
352 }
353}
354
355std::shared_ptr<Points<NekDouble>> GaussPoints::Create(const PointsKey &pkey)
356{
357 std::shared_ptr<Points<NekDouble>> returnval(
359
360 returnval->Initialize();
361
362 return returnval;
363}
364
365std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateMatrix(
366 const PointsKey &pkey)
367{
368 size_t numpoints = pkey.GetNumPoints();
370
371 PointsManager()[pkey]->GetPoints(xpoints);
372
373 // Delegate to function below
374 return GetI(numpoints, xpoints);
375}
376
377const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
378 const PointsKey &pkey)
379{
380 ASSERTL0(pkey.GetPointsDim() == 1,
381 "Gauss Points can only interp to other 1d point distributions");
382
383 return m_InterpManager[pkey];
384}
385
386const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
388{
389 size_t numpoints = 1;
390
391 // Delegate to function below
392 return GetI(numpoints, x);
393}
394
395const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::v_GetI(
396 size_t numpoints, const Array<OneD, const NekDouble> &x)
397{
398 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
399
400 CalculateInterpMatrix(numpoints, x, interp);
401
402 NekDouble *t = interp.data();
403 size_t np = GetNumPoints();
404 std::shared_ptr<NekMatrix<NekDouble>> returnval(
405 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
406 t));
407
408 return returnval;
409}
410
411const std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
413{
414 return m_GalerkinProjectionManager[pkey];
415}
416
417std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::CreateGPMatrix(
418 const PointsKey &pkey)
419{
420 std::shared_ptr<NekMatrix<NekDouble>> returnval =
422
423 // Delegate to function below
424 return returnval;
425}
426
427std::shared_ptr<NekMatrix<NekDouble>> GaussPoints::
429{
430 size_t numpointsfrom = pkey.GetNumPoints();
431 size_t numpointsto = GetNumPoints();
432
434
435 weightsfrom = PointsManager()[pkey]->GetW();
436
437 std::shared_ptr<NekMatrix<NekDouble>> Interp = GetI(pkey);
438
439 Array<OneD, NekDouble> GalProj(numpointsfrom * numpointsto);
440
441 // set up inner product matrix and multiply by inverse of
442 // diagaonal mass matrix
443 for (size_t i = 0; i < numpointsto; ++i)
444 {
445 Vmath::Vmul(numpointsfrom, Interp->GetPtr().get() + i * numpointsfrom,
446 1, &weightsfrom[0], 1, &GalProj[0] + i, numpointsto);
447 Vmath::Smul(numpointsfrom, 1.0 / m_weights[i], &GalProj[0] + i,
448 numpointsto, &GalProj[0] + i, numpointsto);
449 }
450
451 NekDouble *t = GalProj.data();
452 std::shared_ptr<NekMatrix<NekDouble>> returnval(
453 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(
454 numpointsto, numpointsfrom, t));
455
456 return returnval;
457}
458
459} // end of namespace LibUtilities
460} // 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)
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetGalerkinProjection(const PointsKey &pkey) override
virtual void v_CalculateWeights() override final
virtual void v_CalculateDerivMatrix() override final
virtual void v_CalculatePoints() override final
Definition: GaussPoints.cpp:77
std::shared_ptr< NekMatrix< NekDouble > > CreateMatrix(const PointsKey &pkey)
virtual 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:169
Array< OneD, DataType > m_points[3]
Storage for the point locations, allowing for up to a 3D points storage.
Definition: Points.h:361
MatrixSharedPtrType m_derivmatrix[3]
Derivative matrices.
Definition: Points.h:367
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_InterpManager
Definition: Points.h:369
NekManager< PointsKey, NekMatrix< DataType >, PointsKey::opLess > m_GalerkinProjectionManager
Definition: Points.h:371
PointsKey m_pointsKey
Points type for this points distributions.
Definition: Points.h:358
Array< OneD, DataType > m_weights
Quadrature weights for the weights.
Definition: Points.h:363
const MatrixSharedPtrType GetI(const PointsKey &key)
Definition: Points.h:322
Defines a specification for a set of points.
Definition: Points.h:55
size_t GetPointsDim() const
Definition: Points.h:133
PointsType GetPointsType() const
Definition: Points.h:95
size_t GetTotNumPoints() const
Definition: Points.h:163
size_t GetNumPoints() const
Definition: Points.h:90
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
double laginterpderiv(double z, int k, const double *zj, int np)
Definition: Polylib.cpp:100
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:647
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:223
double laginterp(double z, int j, const double *zj, int np)
Definition: Polylib.cpp:87
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:807
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:1074
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:262
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:154
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:1012
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:1043
void zwgk(double *z, double *w, const int npt, const double alpha, const double beta)
Gauss-Kronrod-Jacobi zeros and weights.
Definition: Polylib.cpp:311
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:694
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:750
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:383
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:1105
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:489
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:207
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:245