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 
36 #include <LocalRegions/MatrixKey.h>
37 #include <MultiRegions/ExpList.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 
68 {
69  Array<OneD, NekDouble> out(exp->GetNcoeffs());
70  exp->FwdTrans(in, out);
71  return out;
72 }
73 
76 {
77  Array<OneD, NekDouble> out(exp->GetNpoints());
78  exp->BwdTrans(in, out);
79  return out;
80 }
81 
84 {
85  Array<OneD, NekDouble> out(exp->GetNcoeffs());
86  exp->IProductWRTBase(in, out);
87  return out;
88 }
89 
92 {
93  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
94  exp->MultiplyByInvMassMatrix(in, out);
95  return out;
96 }
97 
100  const py::object constFactorMap,
101  const py::object varCoeffMap)
102 {
103  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
104 
107 
108  if (!constFactorMap.is_none())
109  {
110  facMap = py::extract<StdRegions::ConstFactorMap>(constFactorMap);
111  }
112  if (!varCoeffMap.is_none())
113  {
114  coeffMap = py::extract<StdRegions::VarCoeffMap>(varCoeffMap);
115  }
116 
117  exp->HelmSolve(in, out, facMap, coeffMap);
118  return out;
119 }
120 
123 {
124  return exp->L2(in);
125 }
126 
129  const Array<OneD, const NekDouble> &err)
130 {
131  return exp->L2(in, err);
132 }
133 
136 {
137  return exp->Linf(in);
138 }
139 
142  const Array<OneD, const NekDouble> &err)
143 {
144  return exp->Linf(in, err);
145 }
146 
148 {
149  int nPhys = exp->GetNpoints();
150  int coordim = exp->GetCoordim(0);
151 
152  std::vector<Array<OneD, NekDouble>> coords(coordim);
153  for (int i = 0; i < coordim; ++i)
154  {
155  coords[i] = Array<OneD, NekDouble>(nPhys);
156  }
157 
158  switch (coordim)
159  {
160  case 1:
161  exp->GetCoords(coords[0]);
162  return py::make_tuple(coords[0]);
163  break;
164  case 2:
165  exp->GetCoords(coords[0], coords[1]);
166  return py::make_tuple(coords[0], coords[1]);
167  break;
168  case 3:
169  exp->GetCoords(coords[0], coords[1], coords[2]);
170  return py::make_tuple(coords[0], coords[1], coords[2]);
171  break;
172  }
173 
174  return py::tuple();
175 }
176 
178 {
179  exp->SetPhysArray(inarray);
180 }
181 
183  const Array<OneD, const NekDouble> &inarray)
184 {
185  exp->SetPhys(inarray);
186 }
187 
189 {
190  return exp->GetPhys();
191 }
192 
194 {
195  return exp->Integral();
196 }
197 
199 {
200  std::stringstream ss;
201  ss << static_cast<const void *>(&(exp->GetPhys()[0]));
202  return ss.str();
203 }
204 
206 {
207  exp->ClearGlobalLinSysManager();
209  LocalRegions::MatrixKey::opLess>::ClearManager();
211  LocalRegions::MatrixKey::opLess>::ClearManager();
212 }
213 
214 void ExpList_LoadField(ExpListSharedPtr exp, std::string filename,
215  std::string varName)
216 {
217  int nExp = exp->GetExpSize();
218  Array<OneD, int> elementGIDs(nExp);
219 
220  // Construct a map from element locations to global IDs
221  for (int i = 0; i < nExp; ++i)
222  {
223  elementGIDs[i] = exp->GetExp(i)->GetGeom()->GetGlobalID();
224  }
225 
226  std::vector<LibUtilities::FieldDefinitionsSharedPtr> def;
227  std::vector<std::vector<NekDouble>> data;
229  LibUtilities::FieldIO::CreateForFile(exp->GetSession(), filename);
230  fldIO->Import(filename, def, data, LibUtilities::NullFieldMetaDataMap,
231  elementGIDs);
232 
233  int idx = -1;
234 
235  Vmath::Zero(exp->GetNcoeffs(), exp->UpdateCoeffs(), 1);
236 
237  // Loop over all the expansions
238  for (int i = 0; i < def.size(); ++i)
239  {
240  // Find the index of the required field in the expansion segment
241  for (int j = 0; j < def[i]->m_fields.size(); ++j)
242  {
243  if (def[i]->m_fields[j] == varName)
244  {
245  idx = j;
246  }
247  }
248 
249  if (idx >= 0)
250  {
251  exp->ExtractDataToCoeffs(def[i], data[i], def[i]->m_fields[idx],
252  exp->UpdateCoeffs());
253  }
254  else
255  {
256  std::cout << "Field " + varName + " not found." << std::endl;
257  }
258  }
259 
260  exp->BwdTrans(exp->GetCoeffs(), exp->UpdatePhys());
261 }
262 
264 {
265  int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
266 
267  py::class_<ExpList, std::shared_ptr<ExpList>, boost::noncopyable>(
268  "ExpList", py::init<const LibUtilities::SessionReaderSharedPtr &,
270 
271  // Query points and offset information
272  .def("GetExp", &ExpList_GetExp)
273  .def("GetExpSize", &ExpList::GetExpSize)
274  .def("GetNpoints", &ExpList::GetNpoints)
275  .def("GetNcoeffs", GetNcoeffs)
276  .def("GetCoords", &ExpList_GetCoords)
277  .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
278  .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
279 
280  // Evaluations
281  .def("PhysEvaluate", &ExpList::PhysEvaluate)
282  .def("LoadField", &ExpList_LoadField)
283 
284  // Operators
285  .def("FwdTrans", &ExpList_FwdTrans)
286  .def("BwdTrans", &ExpList_BwdTrans)
287  .def("IProductWRTBase", &ExpList_IProductWRTBase)
288  .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
289  .def("HelmSolve", &ExpList_HelmSolve,
290  (py::arg("in"), py::arg("constFactorMap") = py::object(),
291  py::arg("varCoeffMap") = py::object()))
292 
293  // Error norms
294  .def("L2", &ExpList_L2)
295  .def("L2", &ExpList_L2_Error)
296  .def("Linf", &ExpList_Linf)
297  .def("Linf", &ExpList_Linf_Error)
298 
299  // Storage setups
300  .def("SetPhysArray", &ExpList_SetPhysArray)
301  .def("SetPhys", &ExpList_SetPhys)
302  .def("GetPhys", &ExpList_GetPhys)
303  .def("SetPhysState", &ExpList::SetPhysState)
304  .def("GetPhysState", &ExpList::GetPhysState)
305  .def("Integral", &ExpList_Integral)
306  .def("GetPhysAddress", &ExpList_GetPhysAddress)
307 
308  // Misc functions
309  .def("WriteVTK", &ExpList_WriteVTK)
310  .def("ResetManagers", &ExpList_ResetManagers);
311 }
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:227
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1641
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:2232
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2204
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1719
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:2635
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1779
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:1771
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:2240
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:301
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:172
The namespace associated with the the StdRegions library (StdRegions introduction)
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:240
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:282
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:283
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:241
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:68
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
double NekDouble
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:71