Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GaussPoints.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File GaussPoints.h
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Header file of GaussPoints Distributions
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef GAUSSPOINTS_H
37 #define GAUSSPOINTS_H
38 
43 #include <boost/shared_ptr.hpp>
44 #include <boost/bind.hpp>
45 
46 namespace Nektar
47 {
48  namespace LibUtilities
49  {
50  class GaussPoints: public Points<NekDouble>
51  {
52  public:
53  virtual ~GaussPoints()
54  {
55  }
56 
57  LIB_UTILITIES_EXPORT static boost::shared_ptr< Points<NekDouble> > Create(const PointsKey &pkey);
58 
59  LIB_UTILITIES_EXPORT boost::shared_ptr< NekMatrix<NekDouble> > CreateMatrix(const PointsKey &pkey);
60 
61  LIB_UTILITIES_EXPORT const boost::shared_ptr<NekMatrix<NekDouble> > GetI(const PointsKey &pkey);
62  LIB_UTILITIES_EXPORT const boost::shared_ptr<NekMatrix<NekDouble> > GetI(const Array<OneD, const NekDouble>& x);
63  LIB_UTILITIES_EXPORT const boost::shared_ptr<NekMatrix<NekDouble> > GetI(unsigned int numpoints, const Array<OneD, const NekDouble>& x);
64 
65  LIB_UTILITIES_EXPORT boost::shared_ptr< NekMatrix<NekDouble> > CreateGPMatrix(const PointsKey &pkey);
66 
67  LIB_UTILITIES_EXPORT const boost::shared_ptr<NekMatrix<NekDouble> > GetGalerkinProjection(const PointsKey &pkey);
68 
69 
70 
72  {
74  boost::bind(&GaussPoints::CreateMatrix, this, _1));
76  boost::bind(&GaussPoints::CreateMatrix, this, _1));
78  boost::bind(&GaussPoints::CreateMatrix, this, _1));
80  boost::bind(&GaussPoints::CreateMatrix, this, _1));
82  boost::bind(&GaussPoints::CreateMatrix, this, _1));
84  boost::bind(&GaussPoints::CreateMatrix, this, _1));
86  boost::bind(&GaussPoints::CreateMatrix, this, _1));
88  boost::bind(&GaussPoints::CreateMatrix, this, _1));
90  boost::bind(&GaussPoints::CreateMatrix, this, _1));
92  boost::bind(&GaussPoints::CreateMatrix, this, _1));
94  boost::bind(&GaussPoints::CreateMatrix, this, _1));
96  boost::bind(&GaussPoints::CreateMatrix, this, _1));
98  boost::bind(&GaussPoints::CreateMatrix, this, _1));
100  boost::bind(&GaussPoints::CreateMatrix, this, _1));
102  boost::bind(&GaussPoints::CreateMatrix, this, _1));
104  boost::bind(&GaussPoints::CreateMatrix, this, _1));
106  boost::bind(&GaussPoints::CreateMatrix, this, _1));
108  boost::bind(&GaussPoints::CreateMatrix, this, _1));
110  boost::bind(&GaussPoints::CreateMatrix, this, _1));
112  boost::bind(&GaussPoints::CreateMatrix, this, _1));
114  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
116  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
118  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
120  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
122  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
124  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
126  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
128  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
130  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
132  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
134  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
136  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
138  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
140  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
142  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
144  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
146  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
148  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
150  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
152  boost::bind(&GaussPoints::CreateGPMatrix, this, _1));
153  }
154 
155 
156  private:
157  /// These should not be called. All creation is done
158  /// using the constructor requiring the key, declared
159  /// above.
160  GaussPoints();
161  GaussPoints(const GaussPoints &points);
162 
163  void CalculatePoints();
164  void CalculateWeights();
165  void CalculateDerivMatrix();
166  void CalculateInterpMatrix(unsigned int npts, const Array<OneD, const NekDouble>& xpoints, Array<OneD, NekDouble>& interp);
167 
168 
169  boost::shared_ptr<NekMatrix<NekDouble> > CalculateGalerkinProjectionMatrix(const PointsKey &pkey);
170 
171 
172  /// functions used by the Kronrod points
173  NekDouble LagrangeInterpolant(NekDouble x, int npts, const Array<OneD, const NekDouble>& xpts, const Array<OneD, const NekDouble>& funcvals);
174  NekDouble LagrangePoly(NekDouble x, int pt, int npts, const Array<OneD, const NekDouble>& xpts);
175  NekDouble LagrangePolyDeriv(NekDouble x, int pt, int npts, const Array<OneD, const NekDouble>& xpts);
176 
177  }; // class GaussPoints
178  } // end of namespace
179 } // end of namespace
180 
181 #endif //GAUSSPOINTS_H