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

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

Referenced by export_ExpList().

◆ ExpList_FwdTrans()

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

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

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

Referenced by export_ExpList().

◆ ExpList_GetCoords()

py::tuple ExpList_GetCoords ( ExpListSharedPtr  exp)

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

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}

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

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

Referenced by export_ExpList().

◆ ExpList_GetPhysAddress()

std::string ExpList_GetPhysAddress ( ExpListSharedPtr  exp)

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

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

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

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}
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:403
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:347
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:346

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

Referenced by export_ExpList().

◆ ExpList_Integral()

NekDouble ExpList_Integral ( ExpListSharedPtr  exp)

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

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

Referenced by export_ExpList().

◆ ExpList_IProductWRTBase()

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

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

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

Referenced by export_ExpList().

◆ ExpList_L2()

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

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

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

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

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

Referenced by export_ExpList().

◆ ExpList_Linf()

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

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

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

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

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

Referenced by export_ExpList().

◆ ExpList_LoadField()

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

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

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

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

Referenced by export_ExpList().

◆ ExpList_ResetManagers()

void ExpList_ResetManagers ( ExpListSharedPtr  exp)

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

207{
208 exp->ClearGlobalLinSysManager();
210 LocalRegions::MatrixKey::opLess>::ClearManager();
212 LocalRegions::MatrixKey::opLess>::ClearManager();
213}
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_SetPhys()

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

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

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

Referenced by export_ExpList().

◆ ExpList_SetPhysArray()

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

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

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

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

Referenced by export_ExpList().

◆ export_ExpList()

void export_ExpList ( )

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

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)
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)
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:98
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174

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(), and ExpList_WriteVTK().

Referenced by BOOST_PYTHON_MODULE().