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

Go to the source code of this file.

Functions

void export_ExpList (py::module &m)
 
void export_DisContField (py::module &m)
 
void export_ContField (py::module &m)
 
 PYBIND11_MODULE (_MultiRegions, m)
 

Function Documentation

◆ export_ContField()

void export_ContField ( py::module &  m)

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

51{
52 py::class_<ContField, DisContField, std::shared_ptr<ContField>>(m,
53 "ContField")
54 .def(py::init<>(&CreateContField), py::arg("session"), py::arg("graph"),
55 py::arg("var"), py::arg("checkSingular") = true);
56}
std::shared_ptr< ContField > CreateContField(const LibUtilities::SessionReaderSharedPtr &session, const SpatialDomains::MeshGraphSharedPtr &graph, const std::string &var, const bool checkSingular)

References CreateContField().

Referenced by PYBIND11_MODULE().

◆ export_DisContField()

void export_DisContField ( py::module &  m)

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

59{
60 export_DisContField_Helper<DisContField, ExpList>(m, "DisContField");
61}

Referenced by PYBIND11_MODULE().

◆ export_ExpList()

void export_ExpList ( py::module &  m)

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

268{
269 int (ExpList::*GetNcoeffs)() const = &ExpList::GetNcoeffs;
270
271 py::class_<ExpList, std::shared_ptr<ExpList>>(m, "ExpList")
272 .def(py::init<const LibUtilities::SessionReaderSharedPtr &,
274
275 // Query points and offset information
276 .def("GetExp", &ExpList_GetExp)
277 .def("GetExpSize", &ExpList::GetExpSize)
278 .def("GetNpoints", &ExpList::GetNpoints)
279 .def("GetNcoeffs", GetNcoeffs)
280 .def("GetCoords", &ExpList_GetCoords)
281 .def("GetPhys_Offset", &ExpList::GetPhys_Offset)
282 .def("GetCoeff_Offset", &ExpList::GetCoeff_Offset)
283
284 // Evaluations
285 .def("PhysEvaluate", &ExpList::PhysEvaluate)
286 .def("LoadField", &ExpList_LoadField)
287
288 // Operators
289 .def("FwdTrans", &ExpList_FwdTrans)
290 .def("BwdTrans", &ExpList_BwdTrans)
291 .def("IProductWRTBase", &ExpList_IProductWRTBase)
292 .def("MultiplyByInvMassMatrix", &ExpList_MultiplyByInvMassMatrix)
293 .def("HelmSolve", &ExpList_HelmSolve, py::arg("in"),
294 py::arg("constFactorMap") = StdRegions::NullConstFactorMap,
295 py::arg("varCoeffMap") = StdRegions::NullVarCoeffMap)
296
297 // Error norms
298 .def("L2", &ExpList_L2)
299 .def("L2", &ExpList_L2_Error)
300 .def("Linf", &ExpList_Linf)
301 .def("Linf", &ExpList_Linf_Error)
302
303 // Storage setups
304 .def("SetPhysArray", &ExpList_SetPhysArray)
305 .def("SetPhys", &ExpList_SetPhys)
306 .def("GetPhys", &ExpList_GetPhys)
307 .def("SetCoeffsArray", &ExpList_SetCoeffsArray)
308 .def("GetCoeffs", &ExpList_GetCoeffs)
309 .def("SetPhysState", &ExpList::SetPhysState)
310 .def("GetPhysState", &ExpList::GetPhysState)
311 .def("Integral", &ExpList_Integral)
312 .def("GetPhysAddress", &ExpList_GetPhysAddress)
313
314 .def_property("phys", &ExpList_GetPhys, &ExpList_SetPhys)
315 .def_property("coeffs", &ExpList_GetCoeffs, &ExpList_SetCoeffsArray)
316
317 // Misc functions
318 .def("WriteVTK", &ExpList_WriteVTK)
319 .def("ResetManagers", &ExpList_ResetManagers);
320}
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)
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_HelmSolve(ExpListSharedPtr exp, const Array< OneD, const NekDouble > &in, const StdRegions::ConstFactorMap &constFactorMap, const StdRegions::VarCoeffMap &varCoeffMap)
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
static ConstFactorMap NullConstFactorMap
Definition: StdRegions.hpp:431
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376

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(), ExpList_WriteVTK(), Nektar::StdRegions::NullConstFactorMap, and Nektar::StdRegions::NullVarCoeffMap.

Referenced by PYBIND11_MODULE().

◆ PYBIND11_MODULE()

PYBIND11_MODULE ( _MultiRegions  ,
 
)

Definition at line 41 of file MultiRegions.cpp.

42{
46}
void export_ContField(py::module &m)
void export_DisContField(py::module &m)
void export_ExpList(py::module &m)

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