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 
43 {
44  Array<OneD, NekDouble> out(exp->GetNcoeffs());
45  exp->FwdTrans(in, out);
46  return out;
47 }
48 
51 {
52  Array<OneD, NekDouble> out(exp->GetTotPoints());
53  exp->BwdTrans(in, out);
54  return out;
55 }
56 
59 {
60  Array<OneD, NekDouble> out(exp->GetNcoeffs());
61  exp->IProductWRTBase(in, out);
62  return out;
63 }
64 
67  const Array<OneD, const NekDouble> &physvals)
68 {
69  return exp->PhysEvaluate(coords, physvals);
70 }
71 
74 {
75  return exp->L2(in);
76 }
77 
81 {
82  return exp->L2(in, err);
83 }
84 
86 {
87  int nPhys = exp->GetTotPoints();
88  int coordim = exp->GetCoordim();
89 
90  std::vector<Array<OneD, NekDouble>> coords(coordim);
91  for (int i = 0; i < coordim; ++i)
92  {
93  coords[i] = Array<OneD, NekDouble>(nPhys);
94  }
95 
96  switch (coordim)
97  {
98  case 1:
99  exp->GetCoords(coords[0]);
100  return py::make_tuple(coords[0]);
101  break;
102  case 2:
103  exp->GetCoords(coords[0], coords[1]);
104  return py::make_tuple(coords[0], coords[1]);
105  break;
106  case 3:
107  exp->GetCoords(coords[0], coords[1], coords[2]);
108  return py::make_tuple(coords[0], coords[1], coords[2]);
109  break;
110  }
111 
112  return py::tuple();
113 }
114 
116  const Array<OneD, const NekDouble> &inarray)
117 {
118  int nPhys = exp->GetTotPoints();
119  int coordim = exp->GetCoordim();
120 
121  std::vector<Array<OneD, NekDouble>> derivs(coordim);
122  for (int i = 0; i < coordim; ++i)
123  {
124  derivs[i] = Array<OneD, NekDouble>(nPhys);
125  }
126 
127  switch (coordim)
128  {
129  case 1:
130  exp->PhysDeriv(inarray, derivs[0]);
131  return py::make_tuple(derivs[0]);
132  break;
133  case 2:
134  exp->PhysDeriv(inarray, derivs[0], derivs[1]);
135  return py::make_tuple(derivs[0], derivs[1]);
136  break;
137  case 3:
138  exp->PhysDeriv(inarray, derivs[0], derivs[1], derivs[2]);
139  return py::make_tuple(derivs[0], derivs[1], derivs[2]);
140  break;
141  }
142 
143  return py::tuple();
144 }
145 
147 {
148  py::class_<StdExpansion, std::shared_ptr<StdExpansion>, boost::noncopyable>(
149  "StdExpansion", py::no_init)
150 
151  .def("GetNcoeffs", &StdExpansion::GetNcoeffs)
152  .def("GetTotPoints", &StdExpansion::GetTotPoints)
153  .def("GetBasisType", &StdExpansion::GetBasisType)
154  .def("GetPointsType", &StdExpansion::GetPointsType)
155  .def("GetNverts", &StdExpansion::GetNverts)
156  .def("GetNtraces", &StdExpansion::GetNtraces)
157  .def("DetShapeType", &StdExpansion::DetShapeType)
158  .def("GetShapeDimension", &StdExpansion::GetShapeDimension)
159  .def("Integral", &StdExpansion::Integral)
160 
161  .def("GetBasis", &StdExpansion::GetBasis,
162  py::return_value_policy<py::copy_const_reference>())
163 
164  .def("GenMatrix", &StdExpansion::GenMatrix)
165  .def("GetStdMatrix", &StdExpansion::GetStdMatrix)
166 
167  .def("FwdTrans", &StdExpansion_FwdTrans)
168  .def("BwdTrans", &StdExpansion_BwdTrans)
169  .def("IProductWRTBase", &StdExpansion_IProductWRTBase)
170 
171  .def("PhysEvaluate", &StdExpansion_PhysEvaluate)
172  .def("L2", &StdExpansion_L2)
173  .def("L2", &StdExpansion_L2_Error)
174 
175  .def("GetCoords", &StdExpansion_GetCoords)
176 
177  .def("PhysDeriv", &StdExpansion_PhysDeriv);
178 }
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)
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:130
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
Definition: StdExpansion.h:140
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:118
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:162
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
int GetNtraces() const
Returns the number of trace elements connected to this element.
Definition: StdExpansion.h:357
int GetNverts() const
This function returns the number of vertices of the expansion domain.
Definition: StdExpansion.h:252
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
Definition: StdExpansion.h:373
DNekMatSharedPtr GenMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:843
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
Definition: StdExpansion.h:211
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
This function integrates the specified function over the domain.
Definition: StdExpansion.h:479
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:2
double NekDouble