Nektar++
Functions
Python/ExpList.cpp File Reference
#include <LibUtilities/Python/NekPyConfig.hpp>
#include <LocalRegions/MatrixKey.h>
#include <MultiRegions/ExpList.h>
#include <fstream>
#include <sstream>
#include <string>

Go to the source code of this file.

Functions

ExpansionSharedPtr ExpList_GetExp (ExpListSharedPtr exp, int i)
 
void ExpList_WriteVTK (ExpListSharedPtr exp, std::string filename)
 
Array< OneD, NekDoubleExpList_FwdTrans (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
Array< OneD, NekDoubleExpList_BwdTrans (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
Array< OneD, NekDoubleExpList_IProductWRTBase (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
Array< OneD, NekDoubleExpList_MultiplyByInvMassMatrix (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
Array< OneD, NekDoubleExpList_HelmSolve (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const py::object constFactorMap, const py::object varCoeffMap)
 
NekDouble ExpList_L2 (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
 
NekDouble ExpList_L2_Error (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
 
NekDouble ExpList_Linf (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)
 
py::tuple ExpList_GetCoords (ExpListSharedPtr exp)
 
void ExpList_SetPhysArray (ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
 
void ExpList_SetPhys (ExpListSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
 
const Array< OneD, const NekDoubleExpList_GetPhys (ExpListSharedPtr exp)
 
NekDouble ExpList_Integral (ExpListSharedPtr exp)
 
std::string ExpList_GetPhysAddress (ExpListSharedPtr exp)
 
void ExpList_ResetManagers (ExpListSharedPtr exp)
 
void ExpList_LoadField (ExpListSharedPtr exp, std::string filename, std::string varName)
 
void export_ExpList ()
 

Function Documentation

◆ ExpList_BwdTrans()

Array<OneD, NekDouble> ExpList_BwdTrans ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 74 of file Python/ExpList.cpp.

76 {
77  Array<OneD, NekDouble> out(exp->GetNpoints());
78  exp->BwdTrans(in, out);
79  return out;
80 }

Referenced by export_ExpList().

◆ ExpList_FwdTrans()

Array<OneD, NekDouble> ExpList_FwdTrans ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 66 of file Python/ExpList.cpp.

68 {
69  Array<OneD, NekDouble> out(exp->GetNcoeffs());
70  exp->FwdTrans(in, out);
71  return out;
72 }

Referenced by export_ExpList().

◆ ExpList_GetCoords()

py::tuple ExpList_GetCoords ( ExpListSharedPtr  exp)

Definition at line 147 of file Python/ExpList.cpp.

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 }

Referenced by export_ExpList().

◆ ExpList_GetExp()

ExpansionSharedPtr ExpList_GetExp ( ExpListSharedPtr  exp,
int  i 
)

Definition at line 49 of file Python/ExpList.cpp.

50 {
51  return exp->GetExp(i);
52 }

Referenced by export_ExpList().

◆ ExpList_GetPhys()

const Array<OneD, const NekDouble> ExpList_GetPhys ( ExpListSharedPtr  exp)

Definition at line 188 of file Python/ExpList.cpp.

189 {
190  return exp->GetPhys();
191 }

Referenced by export_ExpList().

◆ ExpList_GetPhysAddress()

std::string ExpList_GetPhysAddress ( ExpListSharedPtr  exp)

Definition at line 198 of file Python/ExpList.cpp.

199 {
200  std::stringstream ss;
201  ss << static_cast<const void *>(&(exp->GetPhys()[0]));
202  return ss.str();
203 }

Referenced by export_ExpList().

◆ ExpList_HelmSolve()

Array<OneD, NekDouble> ExpList_HelmSolve ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in,
const py::object  constFactorMap,
const py::object  varCoeffMap 
)

Definition at line 98 of file Python/ExpList.cpp.

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 }
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

References Nektar::StdRegions::NullConstFactorMap, and Nektar::StdRegions::NullVarCoeffMap.

Referenced by export_ExpList().

◆ ExpList_Integral()

NekDouble ExpList_Integral ( ExpListSharedPtr  exp)

Definition at line 193 of file Python/ExpList.cpp.

194 {
195  return exp->Integral();
196 }

Referenced by export_ExpList().

◆ ExpList_IProductWRTBase()

Array<OneD, NekDouble> ExpList_IProductWRTBase ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 82 of file Python/ExpList.cpp.

84 {
85  Array<OneD, NekDouble> out(exp->GetNcoeffs());
86  exp->IProductWRTBase(in, out);
87  return out;
88 }

Referenced by export_ExpList().

◆ ExpList_L2()

NekDouble ExpList_L2 ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 121 of file Python/ExpList.cpp.

123 {
124  return exp->L2(in);
125 }

Referenced by export_ExpList().

◆ ExpList_L2_Error()

NekDouble ExpList_L2_Error ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in,
const Array< OneD, const NekDouble > &  err 
)

Definition at line 127 of file Python/ExpList.cpp.

130 {
131  return exp->L2(in, err);
132 }

Referenced by export_ExpList().

◆ ExpList_Linf()

NekDouble ExpList_Linf ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 134 of file Python/ExpList.cpp.

136 {
137  return exp->Linf(in);
138 }

Referenced by export_ExpList().

◆ ExpList_Linf_Error()

NekDouble ExpList_Linf_Error ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in,
const Array< OneD, const NekDouble > &  err 
)

Definition at line 140 of file Python/ExpList.cpp.

143 {
144  return exp->Linf(in, err);
145 }

Referenced by export_ExpList().

◆ ExpList_LoadField()

void ExpList_LoadField ( ExpListSharedPtr  exp,
std::string  filename,
std::string  varName 
)

Definition at line 214 of file Python/ExpList.cpp.

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 }
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:301
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:53
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492

References Nektar::LibUtilities::FieldIO::CreateForFile(), Nektar::LibUtilities::NullFieldMetaDataMap, and Vmath::Zero().

Referenced by export_ExpList().

◆ ExpList_MultiplyByInvMassMatrix()

Array<OneD, NekDouble> ExpList_MultiplyByInvMassMatrix ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  in 
)

Definition at line 90 of file Python/ExpList.cpp.

92 {
93  Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
94  exp->MultiplyByInvMassMatrix(in, out);
95  return out;
96 }

Referenced by export_ExpList().

◆ ExpList_ResetManagers()

void ExpList_ResetManagers ( ExpListSharedPtr  exp)

Definition at line 205 of file Python/ExpList.cpp.

206 {
207  exp->ClearGlobalLinSysManager();
209  LocalRegions::MatrixKey::opLess>::ClearManager();
211  LocalRegions::MatrixKey::opLess>::ClearManager();
212 }
NekMatrix< NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag >, BlockMatrixTag > DNekScalBlkMat
Definition: NekTypeDefs.hpp:68
NekMatrix< NekMatrix< NekDouble, StandardMatrixTag >, ScaledMatrixTag > DNekScalMat
Used to lookup the create function in NekManager.
Definition: MatrixKey.h:71

Referenced by export_ExpList().

◆ ExpList_SetPhys()

void ExpList_SetPhys ( ExpListSharedPtr  exp,
const Array< OneD, const NekDouble > &  inarray 
)

Definition at line 182 of file Python/ExpList.cpp.

184 {
185  exp->SetPhys(inarray);
186 }

Referenced by export_ExpList().

◆ ExpList_SetPhysArray()

void ExpList_SetPhysArray ( ExpListSharedPtr  exp,
Array< OneD, NekDouble inarray 
)

Definition at line 177 of file Python/ExpList.cpp.

178 {
179  exp->SetPhysArray(inarray);
180 }

Referenced by export_ExpList().

◆ ExpList_WriteVTK()

void ExpList_WriteVTK ( ExpListSharedPtr  exp,
std::string  filename 
)

Definition at line 54 of file Python/ExpList.cpp.

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 }

Referenced by export_ExpList().

◆ export_ExpList()

void export_ExpList ( )

Definition at line 263 of file Python/ExpList.cpp.

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)
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)
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:101
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172

References ExpList_BwdTrans(), ExpList_FwdTrans(), ExpList_GetCoords(), ExpList_GetExp(), ExpList_GetPhys(), ExpList_GetPhysAddress(), ExpList_HelmSolve(), ExpList_Integral(), ExpList_IProductWRTBase(), ExpList_L2(), ExpList_L2_Error(), ExpList_Linf(), ExpList_Linf_Error(), ExpList_LoadField(), ExpList_MultiplyByInvMassMatrix(), ExpList_ResetManagers(), ExpList_SetPhys(), ExpList_SetPhysArray(), ExpList_WriteVTK(), Nektar::MultiRegions::ExpList::GetCoeff_Offset(), Nektar::MultiRegions::ExpList::GetExpSize(), Nektar::MultiRegions::ExpList::GetNcoeffs(), Nektar::MultiRegions::ExpList::GetNpoints(), Nektar::MultiRegions::ExpList::GetPhys_Offset(), Nektar::MultiRegions::ExpList::GetPhysState(), Nektar::MultiRegions::ExpList::PhysEvaluate(), and Nektar::MultiRegions::ExpList::SetPhysState().

Referenced by BOOST_PYTHON_MODULE().