Nektar++
Python/ExpList.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: ExpList.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 ExpList.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <MultiRegions/ExpList.h>
36 #include <LocalRegions/MatrixKey.h>
38 
39 #include <fstream>
40 #include <sstream>
41 #include <string>
42 
43 using namespace Nektar;
44 using namespace Nektar::StdRegions;
45 using namespace Nektar::LocalRegions;
46 using namespace Nektar::SpatialDomains;
47 using namespace Nektar::MultiRegions;
48 
50 {
51  return exp->GetExp(i);
52 }
53 
54 void ExpList_WriteVTK(ExpListSharedPtr exp, std::string filename)
55 {
56  std::ofstream out(filename.c_str());
57  exp->WriteVtkHeader(out);
58  for (int i = 0; i < exp->GetExpSize(); ++i)
59  {
60  exp->WriteVtkPieceHeader(out, i);
61  exp->WriteVtkPieceFooter(out, i);
62  }
63  exp->WriteVtkFooter(out);
64 }
65 
67  ExpListSharedPtr exp,
69 {
70  Array<OneD, NekDouble> out(exp->GetNcoeffs());
71  exp->FwdTrans(in, out);
72  return out;
73 }
74 
76  ExpListSharedPtr exp,
78 {
79  Array<OneD, NekDouble> out(exp->GetNpoints());
80  exp->BwdTrans(in, out);
81  return out;
82 }
83 
85  ExpListSharedPtr exp,
87 {
88  Array<OneD, NekDouble> out(exp->GetNcoeffs());
89  exp->IProductWRTBase(in, out);
90  return out;
91 }
92 
94  ExpListSharedPtr exp,
96 {
97  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
98  exp->MultiplyByInvMassMatrix(in, out, eLocal);
99  return out;
100 }
101 
103  ExpListSharedPtr exp,
105  const py::object constFactorMap,
106  const py::object varCoeffMap)
107 {
108  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
109 
112 
113  if (!constFactorMap.is_none())
114  {
115  facMap = py::extract<StdRegions::ConstFactorMap>(constFactorMap);
116  }
117  if (!varCoeffMap.is_none())
118  {
119  coeffMap = py::extract<StdRegions::VarCoeffMap>(varCoeffMap);
120  }
121 
122  FlagList notUsed;
123 
124  exp->HelmSolve(in, out, notUsed, facMap, coeffMap);
125  return out;
126 }
127 
129  ExpListSharedPtr exp,
131 {
132  return exp->L2(in);
133 }
134 
136  ExpListSharedPtr exp,
138  const Array<OneD, const NekDouble> &err)
139 {
140  return exp->L2(in, err);
141 }
142 
144  ExpListSharedPtr exp,
146 {
147  return exp->Linf(in);
148 }
149 
151  ExpListSharedPtr exp,
153  const Array<OneD, const NekDouble> &err)
154 {
155  return exp->Linf(in, err);
156 }
157 
159 {
160  int nPhys = exp->GetNpoints();
161  int coordim = exp->GetCoordim(0);
162 
163  std::vector<Array<OneD, NekDouble> > coords(coordim);
164  for (int i = 0; i < coordim; ++i)
165  {
166  coords[i] = Array<OneD, NekDouble>(nPhys);
167  }
168 
169  switch (coordim)
170  {
171  case 1:
172  exp->GetCoords(coords[0]);
173  return py::make_tuple(coords[0]);
174  break;
175  case 2:
176  exp->GetCoords(coords[0], coords[1]);
177  return py::make_tuple(coords[0], coords[1]);
178  break;
179  case 3:
180  exp->GetCoords(coords[0], coords[1], coords[2]);
181  return py::make_tuple(coords[0], coords[1], coords[2]);
182  break;
183  }
184 
185  return py::tuple();
186 }
187 
189  ExpListSharedPtr exp,
190  Array<OneD, NekDouble> inarray)
191 {
192  exp->SetPhysArray(inarray);
193 }
194 
196  ExpListSharedPtr exp,
197  const Array<OneD, const NekDouble> &inarray)
198 {
199  exp->SetPhys(inarray);
200 }
201 
203 {
204  return exp->GetPhys();
205 }
206 
208 {
209  return exp->PhysIntegral();
210 }
211 
213 {
214  std::stringstream ss;
215  ss << static_cast<const void*>(&(exp->GetPhys()[0]));
216  return ss.str();
217 }
218 
220 {
221  exp->ClearGlobalLinSysManager();
224  LocalRegions::MatrixKey::opLess>::ClearManager();
226  LocalRegions::MatrixKey, DNekScalBlkMat,
227  LocalRegions::MatrixKey::opLess>::ClearManager();
228 }
229 
231 {
232  int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
233 
234  py::class_<ExpList,
235  std::shared_ptr<ExpList>,
236  boost::noncopyable>(
237  "ExpList", py::no_init)
238 
239  // Expansions
240  .def("GetExp", ExpList_GetExp)
241  .def("GetExpSize", &ExpList::GetExpSize)
242 
243  // Query points and offset information
244  .def("GetNpoints", &ExpList::GetNpoints)
245  .def("GetNcoeffs", GetNcoeffs)
246  .def("GetCoords", &ExpList_GetCoords)
247  .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
248  .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
249 
250  // Evaluations
251  .def("PhysEvaluate", &ExpList::PhysEvaluate)
252 
253  // Operators
254  .def("FwdTrans", &ExpList_FwdTrans)
255  .def("BwdTrans", &ExpList_BwdTrans)
256  .def("IProductWRTBase", &ExpList_IProductWRTBase)
257  .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
258  .def("HelmSolve", &ExpList_HelmSolve, (
259  py::arg("in"),
260  py::arg("constFactorMap") = py::object(),
261  py::arg("varCoeffMap") = py::object()
262  ))
263 
264  // Error norms
265  .def("L2", &ExpList_L2)
266  .def("L2", &ExpList_L2_Error)
267  .def("Linf", &ExpList_Linf)
268  .def("Linf", &ExpList_Linf_Error)
269 
270  // Storage setups
271  .def("SetPhysArray", &ExpList_SetPhysArray)
272  .def("SetPhys", &ExpList_SetPhys)
273  .def("GetPhys", &ExpList_GetPhys)
274  .def("SetPhysState", &ExpList::SetPhysState)
275  .def("GetPhysState", &ExpList::GetPhysState)
276  .def("PhysIntegral", &ExpList_PhysIntegral)
277  .def("GetPhysAddress", &ExpList_GetPhysAddress)
278 
279  // Misc functions
280  .def("WriteVTK", &ExpList_WriteVTK)
281  .def("ResetManagers", &ExpList_ResetManagers)
282  ;
283 }
Array< OneD, NekDouble > ExpList_MultiplyByInvMassMatrix(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
Array< OneD, NekDouble > ExpList_FwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
Local coefficients.
Array< OneD, NekDouble > ExpList_HelmSolve(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const py::object constFactorMap, const py::object varCoeffMap)
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
NekDouble ExpList_L2(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble ExpList_Linf_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1640
ExpansionSharedPtr ExpList_GetExp(ExpListSharedPtr exp, int i)
void ExpList_ResetManagers(ExpListSharedPtr exp)
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
std::string ExpList_GetPhysAddress(ExpListSharedPtr exp)
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2170
void ExpList_SetPhys(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
NekDouble ExpList_L2_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:103
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1706
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:65
NekDouble ExpList_PhysIntegral(ExpListSharedPtr exp)
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:264
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1558
const Array< OneD, const NekDouble > ExpList_GetPhys(ExpListSharedPtr exp)
py::tuple ExpList_GetCoords(ExpListSharedPtr exp)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:1545
double NekDouble
Defines a list of flags.
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
Array< OneD, NekDouble > ExpList_IProductWRTBase(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
The namespace associated with the the StdRegions library (StdRegions introduction) ...
Definition: IndexMapKey.cpp:40
void SetPhysState(const bool physState)
This function manually sets whether the array of physical values (implemented as m_phys) is filled o...
Definition: ExpList.h:1697
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:67
Array< OneD, NekDouble > ExpList_BwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
int GetCoeff_Offset(int n) const
Get the start offset position for a global list of m_coeffs correspoinding to element n...
Definition: ExpList.h:2200
void ExpList_WriteVTK(ExpListSharedPtr exp, std::string filename)
void export_ExpList()
int GetPhys_Offset(int n) const
Get the start offset position for a global list of m_phys correspoinding to element n...
Definition: ExpList.h:2208
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
void ExpList_SetPhysArray(ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:265
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:295
NekDouble ExpList_Linf(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)