Nektar++
Python/StdExpansion.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: StdExpansion.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: Python wrapper for StdExpansion.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 
38 using namespace Nektar;
39 using namespace Nektar::StdRegions;
40 
44 {
45  Array<OneD, NekDouble> out(exp->GetNcoeffs());
46  exp->FwdTrans(in, out);
47  return out;
48 }
49 
53 {
54  Array<OneD, NekDouble> out(exp->GetTotPoints());
55  exp->BwdTrans(in, out);
56  return out;
57 }
58 
62 {
63  Array<OneD, NekDouble> out(exp->GetNcoeffs());
64  exp->BwdTrans(in, out);
65  return out;
66 }
67 
70  const Array<OneD, const NekDouble> &coords,
71  const Array<OneD, const NekDouble> &physvals)
72 {
73  return exp->PhysEvaluate(coords, physvals);
74 }
75 
79 {
80  return exp->L2(in);
81 }
82 
87 {
88  return exp->L2(in, err);
89 }
90 
92 {
93  int nPhys = exp->GetTotPoints();
94  int coordim = exp->GetCoordim();
95 
96  std::vector<Array<OneD, NekDouble> > coords(coordim);
97  for (int i = 0; i < coordim; ++i)
98  {
99  coords[i] = Array<OneD, NekDouble>(nPhys);
100  }
101 
102  switch (coordim)
103  {
104  case 1:
105  exp->GetCoords(coords[0]);
106  return py::make_tuple(coords[0]);
107  break;
108  case 2:
109  exp->GetCoords(coords[0], coords[1]);
110  return py::make_tuple(coords[0], coords[1]);
111  break;
112  case 3:
113  exp->GetCoords(coords[0], coords[1], coords[2]);
114  return py::make_tuple(coords[0], coords[1], coords[2]);
115  break;
116  }
117 
118  return py::tuple();
119 }
120 
122  const Array<OneD, const NekDouble> &inarray)
123 {
124  int nPhys = exp->GetTotPoints();
125  int coordim = exp->GetCoordim();
126 
127  std::vector<Array<OneD, NekDouble> > derivs(coordim);
128  for (int i = 0; i < coordim; ++i)
129  {
130  derivs[i] = Array<OneD, NekDouble>(nPhys);
131  }
132 
133  switch (coordim)
134  {
135  case 1:
136  exp->PhysDeriv(inarray, derivs[0]);
137  return py::make_tuple(derivs[0]);
138  break;
139  case 2:
140  exp->PhysDeriv(inarray, derivs[0], derivs[1]);
141  return py::make_tuple(derivs[0], derivs[1]);
142  break;
143  case 3:
144  exp->PhysDeriv(inarray, derivs[0], derivs[1], derivs[2]);
145  return py::make_tuple(derivs[0], derivs[1], derivs[2]);
146  break;
147  }
148 
149  return py::tuple();
150 }
151 
153 {
154  py::class_<StdExpansion,
155  std::shared_ptr<StdExpansion>,
156  boost::noncopyable>(
157  "StdExpansion", py::no_init)
158 
159  .def("GetNcoeffs", &StdExpansion::GetNcoeffs)
160  .def("GetTotPoints", &StdExpansion::GetTotPoints)
161  .def("GetBasisType", &StdExpansion::GetBasisType)
162  .def("GetPointsType", &StdExpansion::GetPointsType)
163  .def("GetNverts", &StdExpansion::GetNverts)
164  .def("GetNtraces", &StdExpansion::GetNtraces)
165  .def("DetShapeType", &StdExpansion::DetShapeType)
166  .def("GetShapeDimension", &StdExpansion::GetShapeDimension)
167  .def("Integral", &StdExpansion::Integral)
168 
169  .def("GetBasis", &StdExpansion::GetBasis,
170  py::return_value_policy<py::copy_const_reference>())
171 
172  .def("GenMatrix", &StdExpansion::GenMatrix)
173  .def("GetStdMatrix", &StdExpansion::GetStdMatrix)
174 
175  .def("FwdTrans", &StdExpansion_FwdTrans)
176  .def("BwdTrans", &StdExpansion_BwdTrans)
177  .def("IProductWRTBase", &StdExpansion_IProductWRTBase)
178 
179  .def("PhysEvaluate", &StdExpansion_PhysEvaluate)
180  .def("L2", &StdExpansion_L2)
181  .def("L2", &StdExpansion_L2_Error)
182 
183  .def("GetCoords", &StdExpansion_GetCoords)
184 
185  .def("PhysDeriv", &StdExpansion_PhysDeriv)
186  ;
187 }
NekDouble StdExpansion_L2_Error(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
py::tuple StdExpansion_PhysDeriv(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
Array< OneD, NekDouble > StdExpansion_IProductWRTBase(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble StdExpansion_L2(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
py::tuple StdExpansion_GetCoords(StdExpansionSharedPtr exp)
Array< OneD, NekDouble > StdExpansion_BwdTrans(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble StdExpansion_PhysEvaluate(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
void export_StdExpansion()
Array< OneD, NekDouble > StdExpansion_FwdTrans(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
The base class for all shapes.
Definition: StdExpansion.h:63
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:134
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:111
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:158
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:617
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:360
int GetNverts() const
This function returns the number of vertices of the expansion domain.
Definition: StdExpansion.h:249
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:376
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:850
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:208
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
This function integrates the specified function over the domain.
Definition: StdExpansion.h:482
The namespace associated with the the StdRegions library (StdRegions introduction)
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
double NekDouble