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  size_t nExp = exp->GetExpSize();
59  for (size_t i = 0; i < nExp; ++i)
60  {
61  exp->WriteVtkPieceHeader(out, i);
62  exp->WriteVtkPieceFooter(out, i);
63  }
64  exp->WriteVtkFooter(out);
65 }
66 
69 {
70  Array<OneD, NekDouble> out(exp->GetNcoeffs());
71  exp->FwdTrans(in, out);
72  return out;
73 }
74 
77 {
78  Array<OneD, NekDouble> out(exp->GetNpoints());
79  exp->BwdTrans(in, out);
80  return out;
81 }
82 
85 {
86  Array<OneD, NekDouble> out(exp->GetNcoeffs());
87  exp->IProductWRTBase(in, out);
88  return out;
89 }
90 
93 {
94  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
95  exp->MultiplyByInvMassMatrix(in, out);
96  return out;
97 }
98 
101  const py::object constFactorMap,
102  const py::object varCoeffMap)
103 {
104  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
105 
108 
109  if (!constFactorMap.is_none())
110  {
111  facMap = py::extract<StdRegions::ConstFactorMap>(constFactorMap);
112  }
113  if (!varCoeffMap.is_none())
114  {
115  coeffMap = py::extract<StdRegions::VarCoeffMap>(varCoeffMap);
116  }
117 
118  exp->HelmSolve(in, out, facMap, coeffMap);
119  return out;
120 }
121 
124 {
125  return exp->L2(in);
126 }
127 
130  const Array<OneD, const NekDouble> &err)
131 {
132  return exp->L2(in, err);
133 }
134 
137 {
138  return exp->Linf(in);
139 }
140 
143  const Array<OneD, const NekDouble> &err)
144 {
145  return exp->Linf(in, err);
146 }
147 
149 {
150  size_t nPhys = exp->GetNpoints();
151  size_t coordim = exp->GetCoordim(0);
152 
153  std::vector<Array<OneD, NekDouble>> coords(coordim);
154  for (size_t i = 0; i < coordim; ++i)
155  {
156  coords[i] = Array<OneD, NekDouble>(nPhys);
157  }
158 
159  switch (coordim)
160  {
161  case 1:
162  exp->GetCoords(coords[0]);
163  return py::make_tuple(coords[0]);
164  break;
165  case 2:
166  exp->GetCoords(coords[0], coords[1]);
167  return py::make_tuple(coords[0], coords[1]);
168  break;
169  case 3:
170  exp->GetCoords(coords[0], coords[1], coords[2]);
171  return py::make_tuple(coords[0], coords[1], coords[2]);
172  break;
173  }
174 
175  return py::tuple();
176 }
177 
179 {
180  exp->SetPhysArray(inarray);
181 }
182 
184  const Array<OneD, const NekDouble> &inarray)
185 {
186  exp->SetPhys(inarray);
187 }
188 
190 {
191  return exp->GetPhys();
192 }
193 
195 {
196  return exp->Integral();
197 }
198 
200 {
201  std::stringstream ss;
202  ss << static_cast<const void *>(&(exp->GetPhys()[0]));
203  return ss.str();
204 }
205 
207 {
208  exp->ClearGlobalLinSysManager();
210  LocalRegions::MatrixKey::opLess>::ClearManager();
212  LocalRegions::MatrixKey::opLess>::ClearManager();
213 }
214 
215 void ExpList_LoadField(ExpListSharedPtr exp, std::string filename,
216  std::string varName)
217 {
218  size_t nExp = exp->GetExpSize();
219  Array<OneD, int> elementGIDs(nExp);
220 
221  // Construct a map from element locations to global IDs
222  for (size_t i = 0; i < nExp; ++i)
223  {
224  elementGIDs[i] = exp->GetExp(i)->GetGeom()->GetGlobalID();
225  }
226 
227  std::vector<LibUtilities::FieldDefinitionsSharedPtr> def;
228  std::vector<std::vector<NekDouble>> data;
230  LibUtilities::FieldIO::CreateForFile(exp->GetSession(), filename);
231  fldIO->Import(filename, def, data, LibUtilities::NullFieldMetaDataMap,
232  elementGIDs);
233 
234  int idx = -1;
235 
236  Vmath::Zero(exp->GetNcoeffs(), exp->UpdateCoeffs(), 1);
237 
238  // Loop over all the expansions
239  for (size_t i = 0; i < def.size(); ++i)
240  {
241  // Find the index of the required field in the expansion segment
242  for (size_t j = 0; j < def[i]->m_fields.size(); ++j)
243  {
244  if (def[i]->m_fields[j] == varName)
245  {
246  idx = j;
247  }
248  }
249 
250  if (idx >= 0)
251  {
252  exp->ExtractDataToCoeffs(def[i], data[i], def[i]->m_fields[idx],
253  exp->UpdateCoeffs());
254  }
255  else
256  {
257  std::cout << "Field " + varName + " not found." << std::endl;
258  }
259  }
260 
261  exp->BwdTrans(exp->GetCoeffs(), exp->UpdatePhys());
262 }
263 
265 {
266  int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
267 
268  py::class_<ExpList, std::shared_ptr<ExpList>, boost::noncopyable>(
269  "ExpList", py::init<const LibUtilities::SessionReaderSharedPtr &,
271 
272  // Query points and offset information
273  .def("GetExp", &ExpList_GetExp)
274  .def("GetExpSize", &ExpList::GetExpSize)
275  .def("GetNpoints", &ExpList::GetNpoints)
276  .def("GetNcoeffs", GetNcoeffs)
277  .def("GetCoords", &ExpList_GetCoords)
278  .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
279  .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
280 
281  // Evaluations
282  .def("PhysEvaluate", &ExpList::PhysEvaluate)
283  .def("LoadField", &ExpList_LoadField)
284 
285  // Operators
286  .def("FwdTrans", &ExpList_FwdTrans)
287  .def("BwdTrans", &ExpList_BwdTrans)
288  .def("IProductWRTBase", &ExpList_IProductWRTBase)
289  .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
290  .def("HelmSolve", &ExpList_HelmSolve,
291  (py::arg("in"), py::arg("constFactorMap") = py::object(),
292  py::arg("varCoeffMap") = py::object()))
293 
294  // Error norms
295  .def("L2", &ExpList_L2)
296  .def("L2", &ExpList_L2_Error)
297  .def("Linf", &ExpList_Linf)
298  .def("Linf", &ExpList_Linf_Error)
299 
300  // Storage setups
301  .def("SetPhysArray", &ExpList_SetPhysArray)
302  .def("SetPhys", &ExpList_SetPhys)
303  .def("GetPhys", &ExpList_GetPhys)
304  .def("SetPhysState", &ExpList::SetPhysState)
305  .def("GetPhysState", &ExpList::GetPhysState)
306  .def("Integral", &ExpList_Integral)
307  .def("GetPhysAddress", &ExpList_GetPhysAddress)
308 
309  // Misc functions
310  .def("WriteVTK", &ExpList_WriteVTK)
311  .def("ResetManagers", &ExpList_ResetManagers);
312 }
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:103
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1492
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
Definition: ExpList.h:2037
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:1996
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1568
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:2638
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1623
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:1616
int GetPhys_Offset(int n) const
Get the start offset position for a local contiguous list of quadrature points in a full array corres...
Definition: ExpList.h:2044
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:327
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< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:400
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:344
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:343
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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