Nektar++
PolyEPoints.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: PolyEPoints.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: 1D Evenly-Spaced Points
32//
33///////////////////////////////////////////////////////////////////////////////
34
38
39namespace Nektar
40{
41namespace LibUtilities
42{
45
47{
48 // Allocate the storage for points
50
51 size_t npts = m_pointsKey.GetNumPoints();
52 if (npts == 1)
53 {
54 m_points[0][0] = 0.0;
55 }
56 else
57 {
58 NekDouble dx = 2.0 / (NekDouble)(npts - 1);
59 for (size_t i = 0; i < npts; ++i)
60 {
61 m_points[0][i] = -1.0 + i * dx;
62 }
63 }
64}
65
67{
68 // Allocate the storage for points
70
71 size_t npts = m_pointsKey.GetNumPoints();
72 if (npts == 1)
73 {
74 m_weights[0] = 2.0; // midpoint rule
75 }
76 else
77 {
78 PointsKey gaussKey(npts, eGaussLobattoLegendre);
79 std::shared_ptr<PointsBaseType> ptr = PointsManager()[gaussKey];
82
83 ptr->GetZW(z, w);
84 for (size_t i = 0; i < npts; ++i)
85 {
86 m_weights[i] = 0.0;
87 for (unsigned j = 0; j < npts; ++j)
88 {
89 m_weights[i] +=
90 w[j] * Polylib::laginterp(z[j], i, &m_points[0][0], npts);
91 }
92 }
93 }
94}
95
97{
98 // Allocate the derivative matrix.
99 PointsBaseType::v_CalculateDerivMatrix();
100
101 for (size_t i = 0; i < m_pointsKey.GetNumPoints(); ++i)
102 {
103 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
104 {
106 m_points[0][i], j, &m_points[0][0], m_pointsKey.GetNumPoints());
107 }
108 }
109}
110
112 size_t npts, const Array<OneD, const NekDouble> &xpoints,
114{
115 for (size_t i = 0; i < npts; ++i)
116 {
117 for (size_t j = 0; j < m_pointsKey.GetNumPoints(); ++j)
118 {
119 interp[i + j * npts] = Polylib::laginterp(
120 xpoints[i], j, &m_points[0][0], m_pointsKey.GetNumPoints());
121 }
122 }
123}
124
125std::shared_ptr<PolyEPoints::PointsBaseType> PolyEPoints::Create(
126 const PointsKey &key)
127{
128 std::shared_ptr<PointsBaseType> returnval(
130
131 returnval->Initialize();
132
133 return returnval;
134}
135
136const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
137 const PointsKey &pkey)
138{
139 ASSERTL0(pkey.GetPointsDim() == 1,
140 "Gauss Points can only interp to other 1d point distributions");
141
142 size_t numpoints = pkey.GetNumPoints();
144
145 PointsManager()[pkey]->GetPoints(xpoints);
146
147 return GetI(numpoints, xpoints);
148}
149
150const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
152{
153 size_t numpoints = 1;
154
155 /// Delegate to function below.
156 return GetI(numpoints, x);
157}
158
159const std::shared_ptr<NekMatrix<NekDouble>> PolyEPoints::v_GetI(
160 size_t numpoints, const Array<OneD, const NekDouble> &x)
161{
162 Array<OneD, NekDouble> interp(GetNumPoints() * numpoints);
163
164 CalculateInterpMatrix(numpoints, x, interp);
165
166 size_t np = GetTotNumPoints();
167 NekDouble *d = interp.data();
168 std::shared_ptr<NekMatrix<NekDouble>> returnval(
169 MemoryManager<NekMatrix<NekDouble>>::AllocateSharedPtr(numpoints, np,
170 d));
171
172 return returnval;
173}
174
175} // end of namespace LibUtilities
176} // end of namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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
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
size_t GetNumPoints() const
Definition: Points.h:90
virtual void v_CalculateWeights() override final
Definition: PolyEPoints.cpp:66
void CalculateInterpMatrix(size_t npts, const Array< OneD, const NekDouble > &xpoints, Array< OneD, NekDouble > &interp)
virtual void v_CalculatePoints() override final
Definition: PolyEPoints.cpp:46
virtual const std::shared_ptr< NekMatrix< NekDouble > > v_GetI(const PointsKey &pkey) override
static std::shared_ptr< PointsBaseType > Create(const PointsKey &key)
virtual void v_CalculateDerivMatrix() override final
Definition: PolyEPoints.cpp:96
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
PointsManagerT & PointsManager(void)
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
std::vector< double > d(NPUPPER *NPUPPER)
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
double laginterp(double z, int j, const double *zj, int np)
Definition: Polylib.cpp:87