Nektar++
Functions
Python/ExpList.cpp File Reference
#include <LibUtilities/Python/BasicUtils/SharedArray.hpp>
#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)
 
void ExpList_SetCoeffsArray (ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
 
const Array< OneD, const NekDoubleExpList_GetCoeffs (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 77 of file Python/ExpList.cpp.

79{
80 Array<OneD, NekDouble> out(exp->GetNpoints());
81 exp->BwdTrans(in, out);
82 return out;
83}

Referenced by export_ExpList().

◆ ExpList_FwdTrans()

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

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

71{
72 Array<OneD, NekDouble> out(exp->GetNcoeffs());
73 exp->FwdTrans(in, out);
74 return out;
75}

Referenced by export_ExpList().

◆ ExpList_GetCoeffs()

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

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

203{
204 return exp->GetCoeffs();
205}

Referenced by export_ExpList().

◆ ExpList_GetCoords()

py::tuple ExpList_GetCoords ( ExpListSharedPtr  exp)

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

151{
152 size_t nPhys = exp->GetNpoints();
153 size_t coordim = exp->GetCoordim(0);
154
155 std::vector<Array<OneD, NekDouble>> coords(coordim);
156 for (size_t i = 0; i < coordim; ++i)
157 {
158 coords[i] = Array<OneD, NekDouble>(nPhys);
159 }
160
161 switch (coordim)
162 {
163 case 1:
164 exp->GetCoords(coords[0]);
165 return py::make_tuple(coords[0]);
166 break;
167 case 2:
168 exp->GetCoords(coords[0], coords[1]);
169 return py::make_tuple(coords[0], coords[1]);
170 break;
171 case 3:
172 exp->GetCoords(coords[0], coords[1], coords[2]);
173 return py::make_tuple(coords[0], coords[1], coords[2]);
174 break;
175 }
176
177 return py::tuple();
178}

Referenced by export_ExpList().

◆ ExpList_GetExp()

ExpansionSharedPtr ExpList_GetExp ( ExpListSharedPtr  exp,
int  i 
)

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

52{
53 return exp->GetExp(i);
54}

Referenced by export_ExpList().

◆ ExpList_GetPhys()

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

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

192{
193 return exp->GetPhys();
194}

Referenced by export_ExpList().

◆ ExpList_GetPhysAddress()

std::string ExpList_GetPhysAddress ( ExpListSharedPtr  exp)

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

213{
214 std::stringstream ss;
215 ss << static_cast<const void *>(&(exp->GetPhys()[0]));
216 return ss.str();
217}

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 101 of file Python/ExpList.cpp.

105{
106 Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
107
110
111 if (!constFactorMap.is_none())
112 {
113 facMap = py::extract<StdRegions::ConstFactorMap>(constFactorMap);
114 }
115 if (!varCoeffMap.is_none())
116 {
117 coeffMap = py::extract<StdRegions::VarCoeffMap>(varCoeffMap);
118 }
119
120 exp->HelmSolve(in, out, facMap, coeffMap);
121 return out;
122}
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:431
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375

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

Referenced by export_ExpList().

◆ ExpList_Integral()

NekDouble ExpList_Integral ( ExpListSharedPtr  exp)

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

208{
209 return exp->Integral();
210}

Referenced by export_ExpList().

◆ ExpList_IProductWRTBase()

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

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

87{
88 Array<OneD, NekDouble> out(exp->GetNcoeffs());
89 exp->IProductWRTBase(in, out);
90 return out;
91}

Referenced by export_ExpList().

◆ ExpList_L2()

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

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

126{
127 return exp->L2(in);
128}

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 130 of file Python/ExpList.cpp.

133{
134 return exp->L2(in, err);
135}

Referenced by export_ExpList().

◆ ExpList_Linf()

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

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

139{
140 return exp->Linf(in);
141}

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 143 of file Python/ExpList.cpp.

146{
147 return exp->Linf(in, err);
148}

Referenced by export_ExpList().

◆ ExpList_LoadField()

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

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

230{
231 size_t nExp = exp->GetExpSize();
232 Array<OneD, int> elementGIDs(nExp);
233
234 // Construct a map from element locations to global IDs
235 for (size_t i = 0; i < nExp; ++i)
236 {
237 elementGIDs[i] = exp->GetExp(i)->GetGeom()->GetGlobalID();
238 }
239
240 std::vector<LibUtilities::FieldDefinitionsSharedPtr> def;
241 std::vector<std::vector<NekDouble>> data;
243 LibUtilities::FieldIO::CreateForFile(exp->GetSession(), filename);
244 fldIO->Import(filename, def, data, LibUtilities::NullFieldMetaDataMap,
245 elementGIDs);
246
247 int idx = -1;
248
249 Vmath::Zero(exp->GetNcoeffs(), exp->UpdateCoeffs(), 1);
250
251 // Loop over all the expansions
252 for (size_t i = 0; i < def.size(); ++i)
253 {
254 // Find the index of the required field in the expansion segment
255 for (size_t j = 0; j < def[i]->m_fields.size(); ++j)
256 {
257 if (def[i]->m_fields[j] == varName)
258 {
259 idx = j;
260 }
261 }
262
263 if (idx >= 0)
264 {
265 exp->ExtractDataToCoeffs(def[i], data[i], def[i]->m_fields[idx],
266 exp->UpdateCoeffs());
267 }
268 else
269 {
270 std::cout << "Field " + varName + " not found." << std::endl;
271 }
272 }
273
274 exp->BwdTrans(exp->GetCoeffs(), exp->UpdatePhys());
275}
std::shared_ptr< FieldIO > FieldIOSharedPtr
Definition: FieldIO.h:322
static FieldMetaDataMap NullFieldMetaDataMap
Definition: FieldIO.h:51
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273

References 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 93 of file Python/ExpList.cpp.

95{
96 Array<OneD, NekDouble> out(exp->GetNcoeffs(), 0.0);
97 exp->MultiplyByInvMassMatrix(in, out);
98 return out;
99}

Referenced by export_ExpList().

◆ ExpList_ResetManagers()

void ExpList_ResetManagers ( ExpListSharedPtr  exp)

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

220{
221 exp->ClearGlobalLinSysManager();
223 LocalRegions::MatrixKey::opLess>::ClearManager();
225 LocalRegions::MatrixKey::opLess>::ClearManager();
226}
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:69

Referenced by export_ExpList().

◆ ExpList_SetCoeffsArray()

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

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

198{
199 exp->SetCoeffsArray(inarray);
200}

Referenced by export_ExpList().

◆ ExpList_SetPhys()

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

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

187{
188 exp->SetPhys(inarray);
189}

Referenced by export_ExpList().

◆ ExpList_SetPhysArray()

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

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

181{
182 exp->SetPhysArray(inarray);
183}

Referenced by export_ExpList().

◆ ExpList_WriteVTK()

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

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

57{
58 std::ofstream out(filename.c_str());
59 exp->WriteVtkHeader(out);
60 size_t nExp = exp->GetExpSize();
61 for (size_t i = 0; i < nExp; ++i)
62 {
63 exp->WriteVtkPieceHeader(out, i);
64 exp->WriteVtkPieceFooter(out, i);
65 }
66 exp->WriteVtkFooter(out);
67}

Referenced by export_ExpList().

◆ export_ExpList()

void export_ExpList ( )

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

278{
279 int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
280
281 py::class_<ExpList, std::shared_ptr<ExpList>, boost::noncopyable>(
282 "ExpList", py::init<const LibUtilities::SessionReaderSharedPtr &,
284
285 // Query points and offset information
286 .def("GetExp", &ExpList_GetExp)
287 .def("GetExpSize", &ExpList::GetExpSize)
288 .def("GetNpoints", &ExpList::GetNpoints)
289 .def("GetNcoeffs", GetNcoeffs)
290 .def("GetCoords", &ExpList_GetCoords)
291 .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
292 .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
293
294 // Evaluations
295 .def("PhysEvaluate", &ExpList::PhysEvaluate)
296 .def("LoadField", &ExpList_LoadField)
297
298 // Operators
299 .def("FwdTrans", &ExpList_FwdTrans)
300 .def("BwdTrans", &ExpList_BwdTrans)
301 .def("IProductWRTBase", &ExpList_IProductWRTBase)
302 .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
303 .def("HelmSolve", &ExpList_HelmSolve,
304 (py::arg("in"), py::arg("constFactorMap") = py::object(),
305 py::arg("varCoeffMap") = py::object()))
306
307 // Error norms
308 .def("L2", &ExpList_L2)
309 .def("L2", &ExpList_L2_Error)
310 .def("Linf", &ExpList_Linf)
311 .def("Linf", &ExpList_Linf_Error)
312
313 // Storage setups
314 .def("SetPhysArray", &ExpList_SetPhysArray)
315 .def("SetPhys", &ExpList_SetPhys)
316 .def("GetPhys", &ExpList_GetPhys)
317 .def("SetCoeffsArray", &ExpList_SetCoeffsArray)
318 .def("GetCoeffs", &ExpList_GetCoeffs)
319 .def("SetPhysState", &ExpList::SetPhysState)
320 .def("GetPhysState", &ExpList::GetPhysState)
321 .def("Integral", &ExpList_Integral)
322 .def("GetPhysAddress", &ExpList_GetPhysAddress)
323
324 .add_property("phys", &ExpList_GetPhys, &ExpList_SetPhys)
325 .add_property("coeffs", &ExpList_GetCoeffs, &ExpList_SetCoeffsArray)
326
327 // Misc functions
328 .def("WriteVTK", &ExpList_WriteVTK)
329 .def("ResetManagers", &ExpList_ResetManagers);
330
331 // Export Array<OneD, ExpListSharedPtr>
332 export_SharedArray<std::shared_ptr<ExpList>>();
333}
const Array< OneD, const NekDouble > ExpList_GetCoeffs(ExpListSharedPtr exp)
NekDouble ExpList_Linf_Error(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const Array< OneD, const NekDouble > &err)
Array< OneD, NekDouble > ExpList_MultiplyByInvMassMatrix(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
void ExpList_ResetManagers(ExpListSharedPtr exp)
NekDouble ExpList_Linf(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
const Array< OneD, const NekDouble > ExpList_GetPhys(ExpListSharedPtr exp)
Array< OneD, NekDouble > ExpList_HelmSolve(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const py::object constFactorMap, const py::object varCoeffMap)
NekDouble ExpList_Integral(ExpListSharedPtr exp)
Array< OneD, NekDouble > ExpList_IProductWRTBase(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
ExpansionSharedPtr ExpList_GetExp(ExpListSharedPtr exp, int i)
void ExpList_SetPhysArray(ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
Array< OneD, NekDouble > ExpList_BwdTrans(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
void ExpList_LoadField(ExpListSharedPtr exp, std::string filename, std::string varName)
void ExpList_WriteVTK(ExpListSharedPtr exp, std::string filename)
void ExpList_SetPhys(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &inarray)
py::tuple ExpList_GetCoords(ExpListSharedPtr exp)
void ExpList_SetCoeffsArray(ExpListSharedPtr exp, Array< OneD, NekDouble > inarray)
NekDouble ExpList_L2(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in)
std::string ExpList_GetPhysAddress(ExpListSharedPtr exp)
Array< OneD, NekDouble > ExpList_FwdTrans(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)
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:99
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174

References ExpList_BwdTrans(), ExpList_FwdTrans(), ExpList_GetCoeffs(), 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_SetCoeffsArray(), ExpList_SetPhys(), ExpList_SetPhysArray(), and ExpList_WriteVTK().

Referenced by BOOST_PYTHON_MODULE().