Nektar++
Functions
MultiRegions.cpp File Reference
#include <LibUtilities/Python/NekPyConfig.hpp>

Go to the source code of this file.

Functions

void export_ExpList ()
 
void export_DisContField ()
 
void export_ContField ()
 
 BOOST_PYTHON_MODULE (_MultiRegions)
 

Function Documentation

◆ BOOST_PYTHON_MODULE()

BOOST_PYTHON_MODULE ( _MultiRegions  )

Definition at line 41 of file MultiRegions.cpp.

42{
43 np::initialize();
44
48}
void export_DisContField()
void export_ExpList()
void export_ContField()

References export_ContField(), export_DisContField(), and export_ExpList().

◆ export_ContField()

void export_ContField ( )

Definition at line 50 of file Python/ContField.cpp.

51{
52 py::class_<ContField, py::bases<ExpList>, std::shared_ptr<ContField>>(
53 "ContField", py::no_init)
54 .def("__init__",
55 py::make_constructor(&CreateContField, py::default_call_policies(),
56 (py::arg("session"), py::arg("graph"),
57 py::arg("var"),
58 py::arg("checkSingular") = true)));
59}
std::shared_ptr< ContField > CreateContField(const LibUtilities::SessionReaderSharedPtr &session, const SpatialDomains::MeshGraphSharedPtr &graph, const std::string &var, const bool checkSingular)

References CreateContField().

Referenced by BOOST_PYTHON_MODULE().

◆ export_DisContField()

void export_DisContField ( )

Definition at line 60 of file Python/DisContField.cpp.

61{
62 export_DisContField_Helper<DisContField, ExpList>("DisContField");
63}

Referenced by BOOST_PYTHON_MODULE().

◆ 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().