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 51 of file Python/ContField.cpp.

52 {
53  py::class_<ContField, py::bases<ExpList>, std::shared_ptr<ContField>>(
54  "ContField", py::no_init)
55  .def("__init__", py::make_constructor(
57  py::default_call_policies(),
58  (py::arg("session"), py::arg("graph"), py::arg("var"),
59  py::arg("checkSingular") = true)));
60 }
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 62 of file Python/DisContField.cpp.

63 {
64  export_DisContField_Helper<DisContField, ExpList>("DisContField");
65 }

Referenced by BOOST_PYTHON_MODULE().

◆ export_ExpList()

void export_ExpList ( )

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

277 {
278  int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
279 
280  py::class_<ExpList,
281  std::shared_ptr<ExpList>,
282  boost::noncopyable>(
283  "ExpList", py::init<
286 
287  // Query points and offset information
288  .def("GetExp", &ExpList_GetExp)
289  .def("GetExpSize", &ExpList::GetExpSize)
290  .def("GetNpoints", &ExpList::GetNpoints)
291  .def("GetNcoeffs", GetNcoeffs)
292  .def("GetCoords", &ExpList_GetCoords)
293  .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
294  .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
295 
296  // Evaluations
297  .def("PhysEvaluate", &ExpList::PhysEvaluate)
298  .def("LoadField", &ExpList_LoadField)
299 
300  // Operators
301  .def("FwdTrans", &ExpList_FwdTrans)
302  .def("BwdTrans", &ExpList_BwdTrans)
303  .def("IProductWRTBase", &ExpList_IProductWRTBase)
304  .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
305  .def("HelmSolve", &ExpList_HelmSolve, (
306  py::arg("in"),
307  py::arg("constFactorMap") = py::object(),
308  py::arg("varCoeffMap") = py::object()
309  ))
310 
311  // Error norms
312  .def("L2", &ExpList_L2)
313  .def("L2", &ExpList_L2_Error)
314  .def("Linf", &ExpList_Linf)
315  .def("Linf", &ExpList_Linf_Error)
316 
317  // Storage setups
318  .def("SetPhysArray", &ExpList_SetPhysArray)
319  .def("SetPhys", &ExpList_SetPhys)
320  .def("GetPhys", &ExpList_GetPhys)
321  .def("SetPhysState", &ExpList::SetPhysState)
322  .def("GetPhysState", &ExpList::GetPhysState)
323  .def("Integral", &ExpList_Integral)
324  .def("GetPhysAddress", &ExpList_GetPhysAddress)
325 
326  // Misc functions
327  .def("WriteVTK", &ExpList_WriteVTK)
328  .def("ResetManagers", &ExpList_ResetManagers)
329  ;
330 }
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:107
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(), 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().