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);
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  exp->HelmSolve(in, out, facMap, coeffMap);
123  return out;
124 }
125 
127  ExpListSharedPtr exp,
129 {
130  return exp->L2(in);
131 }
132 
134  ExpListSharedPtr exp,
136  const Array<OneD, const NekDouble> &err)
137 {
138  return exp->L2(in, err);
139 }
140 
142  ExpListSharedPtr exp,
144 {
145  return exp->Linf(in);
146 }
147 
149  ExpListSharedPtr exp,
151  const Array<OneD, const NekDouble> &err)
152 {
153  return exp->Linf(in, err);
154 }
155 
157 {
158  int nPhys = exp->GetNpoints();
159  int coordim = exp->GetCoordim(0);
160 
161  std::vector<Array<OneD, NekDouble> > coords(coordim);
162  for (int i = 0; i < coordim; ++i)
163  {
164  coords[i] = Array<OneD, NekDouble>(nPhys);
165  }
166 
167  switch (coordim)
168  {
169  case 1:
170  exp->GetCoords(coords[0]);
171  return py::make_tuple(coords[0]);
172  break;
173  case 2:
174  exp->GetCoords(coords[0], coords[1]);
175  return py::make_tuple(coords[0], coords[1]);
176  break;
177  case 3:
178  exp->GetCoords(coords[0], coords[1], coords[2]);
179  return py::make_tuple(coords[0], coords[1], coords[2]);
180  break;
181  }
182 
183  return py::tuple();
184 }
185 
187  ExpListSharedPtr exp,
188  Array<OneD, NekDouble> inarray)
189 {
190  exp->SetPhysArray(inarray);
191 }
192 
194  ExpListSharedPtr exp,
195  const Array<OneD, const NekDouble> &inarray)
196 {
197  exp->SetPhys(inarray);
198 }
199 
201 {
202  return exp->GetPhys();
203 }
204 
206 {
207  return exp->Integral();
208 }
209 
211 {
212  std::stringstream ss;
213  ss << static_cast<const void*>(&(exp->GetPhys()[0]));
214  return ss.str();
215 }
216 
218 {
219  exp->ClearGlobalLinSysManager();
222  LocalRegions::MatrixKey::opLess>::ClearManager();
225  LocalRegions::MatrixKey::opLess>::ClearManager();
226 }
227 
228 void ExpList_LoadField(ExpListSharedPtr exp, std::string filename, std::string varName)
229 {
230  int nExp = exp->GetExpSize();
231  Array<OneD, int> elementGIDs(nExp);
232 
233  // Construct a map from element locations to global IDs
234  for (int i = 0; i < nExp; ++i)
235  {
236  elementGIDs[i] = exp->GetExp(i)->GetGeom()->GetGlobalID();
237  }
238 
239  std::vector<LibUtilities::FieldDefinitionsSharedPtr> def;
240  std::vector<std::vector<NekDouble> > data;
242  LibUtilities::FieldIO::CreateForFile(exp->GetSession(), filename);
243  fldIO->Import(filename, def, data, LibUtilities::NullFieldMetaDataMap,
244  elementGIDs);
245 
246  int idx = -1;
247 
248  Vmath::Zero(exp->GetNcoeffs(), exp->UpdateCoeffs(), 1);
249 
250  // Loop over all the expansions
251  for (int i = 0; i < def.size(); ++i)
252  {
253  // Find the index of the required field in the expansion segment
254  for (int j = 0; j < def[i]->m_fields.size(); ++j)
255  {
256  if (def[i]->m_fields[j] == varName)
257  {
258  idx = j;
259  }
260  }
261 
262  if (idx >= 0)
263  {
264  exp->ExtractDataToCoeffs(
265  def[i], data[i], def[i]->m_fields[idx], exp->UpdateCoeffs());
266  }
267  else
268  {
269  std::cout << "Field " + varName + " not found." << std::endl;
270  }
271  }
272 
273  exp->BwdTrans_IterPerExp(exp->GetCoeffs(), exp->UpdatePhys());
274 }
275 
277 {
278  int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
279 
280  py::class_<ExpList,
281  std::shared_ptr<ExpList>,
282  boost::noncopyable>(
283  "ExpList", py::init<
286 
287  // Query points and offset information
288  .def("GetExp", &ExpList_GetExp)
289  .def("GetExpSize", &ExpList::GetExpSize)
290  .def("GetNpoints", &ExpList::GetNpoints)
291  .def("GetNcoeffs", GetNcoeffs)
292  .def("GetCoords", &ExpList_GetCoords)
293  .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
294  .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
295 
296  // Evaluations
297  .def("PhysEvaluate", &ExpList::PhysEvaluate)
298  .def("LoadField", &ExpList_LoadField)
299 
300  // Operators
301  .def("FwdTrans", &ExpList_FwdTrans)
302  .def("BwdTrans", &ExpList_BwdTrans)
303  .def("IProductWRTBase", &ExpList_IProductWRTBase)
304  .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
305  .def("HelmSolve", &ExpList_HelmSolve, (
306  py::arg("in"),
307  py::arg("constFactorMap") = py::object(),
308  py::arg("varCoeffMap") = py::object()
309  ))
310 
311  // Error norms
312  .def("L2", &ExpList_L2)
313  .def("L2", &ExpList_L2_Error)
314  .def("Linf", &ExpList_Linf)
315  .def("Linf", &ExpList_Linf_Error)
316 
317  // Storage setups
318  .def("SetPhysArray", &ExpList_SetPhysArray)
319  .def("SetPhys", &ExpList_SetPhys)
320  .def("GetPhys", &ExpList_GetPhys)
321  .def("SetPhysState", &ExpList::SetPhysState)
322  .def("GetPhysState", &ExpList::GetPhysState)
323  .def("Integral", &ExpList_Integral)
324  .def("GetPhysAddress", &ExpList_GetPhysAddress)
325 
326  // Misc functions
327  .def("WriteVTK", &ExpList_WriteVTK)
328  .def("ResetManagers", &ExpList_ResetManagers)
329  ;
330 }
NekDouble ExpList_Linf_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
void ExpList_ResetManagers(ExpListSharedPtr exp)
NekDouble ExpList_Linf(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
Array< OneD, NekDouble > ExpList_IProductWRTBase(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
NekDouble ExpList_Integral(ExpListSharedPtr exp)
void export_ExpList()
ExpansionSharedPtr ExpList_GetExp(ExpListSharedPtr exp, int i)
void ExpList_SetPhysArray(ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
const Array< OneD, const NekDouble > ExpList_GetPhys(ExpListSharedPtr exp)
void ExpList_LoadField(ExpListSharedPtr exp, std::string filename, std::string varName)
Array< OneD, NekDouble > ExpList_MultiplyByInvMassMatrix(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
void ExpList_WriteVTK(ExpListSharedPtr exp, std::string filename)
Array< OneD, NekDouble > ExpList_FwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
void ExpList_SetPhys(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
py::tuple ExpList_GetCoords(ExpListSharedPtr exp)
NekDouble ExpList_L2(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
Array< OneD, NekDouble > ExpList_HelmSolve(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const py::object constFactorMap, const py::object varCoeffMap)
Array< OneD, NekDouble > ExpList_BwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
std::string ExpList_GetPhysAddress(ExpListSharedPtr exp)
NekDouble ExpList_L2_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
static std::shared_ptr< FieldIO > CreateForFile(const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
Construct a FieldIO object for a given input filename.
Definition: FieldIO.cpp:226
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:107
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1800
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:2431
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2401
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1882
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:2542
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1948
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:1939
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:2439
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:306
std::shared_ptr< SessionReader > SessionReaderSharedPtr
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
The namespace associated with the the StdRegions library (StdRegions introduction)
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:272
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:314
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:315
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:273
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:65
double NekDouble
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:68