Nektar++
Loading...
Searching...
No Matches
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
36
40
41using namespace Nektar;
42using namespace Nektar::StdRegions;
43
46{
47 Array<OneD, NekDouble> out(exp->GetNcoeffs());
48 exp->FwdTrans(in, out);
49 return out;
50}
51
54{
55 Array<OneD, NekDouble> out(exp->GetTotPoints());
56 exp->BwdTrans(in, out);
57 return out;
58}
59
62{
63 Array<OneD, NekDouble> out(exp->GetNcoeffs());
64 exp->IProductWRTBase(in, out);
65 return out;
66}
67
70 const Array<OneD, const NekDouble> &physvals)
71{
72 return exp->PhysEvaluate(coords, physvals);
73}
74
77{
78 return exp->L2(in);
79}
80
84{
85 return exp->L2(in, err);
86}
87
89{
90 int nPhys = exp->GetTotPoints();
91 int coordim = exp->GetCoordim();
92
93 std::vector<Array<OneD, NekDouble>> coords(coordim);
94 for (int i = 0; i < coordim; ++i)
95 {
96 coords[i] = Array<OneD, NekDouble>(nPhys);
97 }
98
99 switch (coordim)
100 {
101 case 1:
102 exp->GetCoords(coords[0]);
103 return py::make_tuple(coords[0]);
104 break;
105 case 2:
106 exp->GetCoords(coords[0], coords[1]);
107 return py::make_tuple(coords[0], coords[1]);
108 break;
109 case 3:
110 exp->GetCoords(coords[0], coords[1], coords[2]);
111 return py::make_tuple(coords[0], coords[1], coords[2]);
112 break;
113 }
114
115 return py::tuple();
116}
117
119 const Array<OneD, const NekDouble> &inarray)
120{
121 int nPhys = exp->GetTotPoints();
122 int coordim = exp->GetCoordim();
123
124 std::vector<Array<OneD, NekDouble>> derivs(coordim);
125 for (int i = 0; i < coordim; ++i)
126 {
127 derivs[i] = Array<OneD, NekDouble>(nPhys);
128 }
129
130 switch (coordim)
131 {
132 case 1:
133 exp->PhysDeriv(inarray, derivs[0]);
134 return py::make_tuple(derivs[0]);
135 break;
136 case 2:
137 exp->PhysDeriv(inarray, derivs[0], derivs[1]);
138 return py::make_tuple(derivs[0], derivs[1]);
139 break;
140 case 3:
141 exp->PhysDeriv(inarray, derivs[0], derivs[1], derivs[2]);
142 return py::make_tuple(derivs[0], derivs[1], derivs[2]);
143 break;
144 }
145
146 return py::tuple();
147}
148
149py::array StdExpansion_GenMatrix(std::shared_ptr<StdExpansion> &stdExp,
150 const StdMatrixKey &mkey)
151{
152 auto mat = stdExp->GenMatrix(mkey);
153 // A bit wasteful but at the moment we can't easily keep hold of the
154 // NekMatrix returned as a std::shared_ptr, because pybind11 doesn't like
155 // std::shared_ptr around types that have custom casters (see issue #787).
156 return py::array({mat->GetRows(), mat->GetColumns()},
157 {sizeof(NekDouble), mat->GetRows() * sizeof(NekDouble)},
158 mat->GetRawPtr());
159}
160
161py::array StdExpansion_GetStdMatrix(std::shared_ptr<StdExpansion> &stdExp,
162 const StdMatrixKey &mkey)
163{
164 auto mat = stdExp->GetStdMatrix(mkey);
165 return py::array({mat->GetRows(), mat->GetColumns()},
166 {sizeof(NekDouble), mat->GetRows() * sizeof(NekDouble)},
167 mat->GetRawPtr());
168}
169
170void export_StdExpansion(py::module &m)
171{
172 py::class_<StdExpansion, std::shared_ptr<StdExpansion>>(m, "StdExpansion")
173
174 .def("GetNcoeffs", &StdExpansion::GetNcoeffs)
175 .def("GetTotPoints", &StdExpansion::GetTotPoints)
176 .def("GetBasisType", &StdExpansion::GetBasisType)
177 .def("GetPointsType", &StdExpansion::GetPointsType)
178 .def("GetNverts", &StdExpansion::GetNverts)
179 .def("GetNtraces", &StdExpansion::GetNtraces)
180 .def("DetShapeType", &StdExpansion::DetShapeType)
181 .def("GetShapeDimension", &StdExpansion::GetShapeDimension)
182 .def("Integral", &StdExpansion::Integral)
183
184 .def("GetBasis", &StdExpansion::GetBasis)
185
186 .def("GenMatrix", &StdExpansion_GenMatrix)
187 .def("GetStdMatrix", &StdExpansion_GetStdMatrix)
188
189 .def("FwdTrans", &StdExpansion_FwdTrans)
190 .def("BwdTrans", &StdExpansion_BwdTrans)
191 .def("IProductWRTBase", &StdExpansion_IProductWRTBase)
192
193 .def("PhysEvaluate", &StdExpansion_PhysEvaluate)
194 .def("L2", &StdExpansion_L2)
195 .def("L2", &StdExpansion_L2_Error)
196
197 .def("GetCoords", &StdExpansion_GetCoords)
198
199 .def("PhysDeriv", &StdExpansion_PhysDeriv);
200}
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)
py::array StdExpansion_GetStdMatrix(std::shared_ptr< StdExpansion > &stdExp, const StdMatrixKey &mkey)
Array< OneD, NekDouble > StdExpansion_FwdTrans(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
void export_StdExpansion(py::module &m)
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)
py::array StdExpansion_GenMatrix(std::shared_ptr< StdExpansion > &stdExp, const StdMatrixKey &mkey)
Array< OneD, NekDouble > StdExpansion_IProductWRTBase(StdExpansionSharedPtr exp, const Array< OneD, const NekDouble > &in)
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
int GetTotPoints() const
This function returns the total number of quadrature points used in the element.
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
int GetNtraces() const
Returns the number of trace elements connected to this element.
int GetNverts() const
This function returns the number of vertices of the expansion domain.
LibUtilities::ShapeType DetShapeType() const
This function returns the shape of the expansion domain.
LibUtilities::PointsType GetPointsType(const int dir) const
This function returns the type of quadrature points used in the dir direction.
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
This function integrates the specified function over the domain.
std::shared_ptr< StdExpansion > StdExpansionSharedPtr