Nektar++
Functions
Python/BasicUtils/Equation.cpp File Reference
#include <LibUtilities/BasicUtils/Equation.h>
#include <LibUtilities/Interpreter/Interpreter.h>
#include <LibUtilities/Python/NekPyConfig.hpp>
#include <LibUtilities/BasicConst/NektarUnivTypeDefs.hpp>
#include <boost/python/raw_function.hpp>

Go to the source code of this file.

Functions

py::object Equation_SetConstants (py::tuple args, py::dict kwargs)
 Wrapper for Equation::SetConstants. More...
 
std::shared_ptr< EquationConstructEquation (std::shared_ptr< Interpreter > evaluator, std::string expr, std::string vlist)
 Construct an equation object from an expression string expr and a list of variables vlist. More...
 
NekDouble Equation_Evaluate1 (std::shared_ptr< Equation > equation)
 Wrapper for Equation::Evaluate (overloaded for no parameters). More...
 
NekDouble Equation_Evaluate2 (std::shared_ptr< Equation > equation, const NekDouble x, const NekDouble y=0, const NekDouble z=0, const NekDouble t=0)
 Wrapper for Equation::Evaluate (overloaded for constant parameters). More...
 
Array< OneD, NekDoubleEquation_Evaluate3 (std::shared_ptr< Equation > equation, const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z)
 Wrapper for Equation::Evaluate (overloaded for Array parameters). More...
 
Array< OneD, NekDoubleEquation_Evaluate4 (std::shared_ptr< Equation > equation, const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z, const NekDouble t)
 Wrapper for Equation::Evaluate (overloaded for Array parameters + constant time). More...
 
Array< OneD, NekDoubleEquation_Evaluate5 (std::shared_ptr< Equation > equation, const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z, const Array< OneD, const NekDouble > &t)
 Wrapper for Equation::Evaluate (overloaded for Array parameters + Array time). More...
 
void export_Equation ()
 

Function Documentation

◆ ConstructEquation()

std::shared_ptr< Equation > ConstructEquation ( std::shared_ptr< Interpreter evaluator,
std::string  expr,
std::string  vlist 
)

Construct an equation object from an expression string expr and a list of variables vlist.

Parameters
evaluatorInterpreter object
exprString contaning the expression
vlistString contining the list of variables in expr.
Returns
An #Equation object.

Definition at line 94 of file Python/BasicUtils/Equation.cpp.

96{
97 return std::make_shared<Equation>(evaluator, expr, vlist);
98}

Referenced by export_Equation().

◆ Equation_Evaluate1()

NekDouble Equation_Evaluate1 ( std::shared_ptr< Equation equation)

Wrapper for Equation::Evaluate (overloaded for no parameters).

Parameters
equationEquation object
Returns
Value of equation.

Definition at line 107 of file Python/BasicUtils/Equation.cpp.

108{
109 return equation->Evaluate();
110}

Referenced by export_Equation().

◆ Equation_Evaluate2()

NekDouble Equation_Evaluate2 ( std::shared_ptr< Equation equation,
const NekDouble  x,
const NekDouble  y = 0,
const NekDouble  z = 0,
const NekDouble  t = 0 
)

Wrapper for Equation::Evaluate (overloaded for constant parameters).

Parameters
equationEquation object from Python
xx-coordinate within the expression
yy-coordinate within the expression
zz-coordinate within the expression
tvalue of time within the expression
Returns
Value of evaluated expression.

Definition at line 123 of file Python/BasicUtils/Equation.cpp.

126{
127 return equation->Evaluate(x, y, z, t);
128}
std::vector< double > z(NPUPPER)

References Nektar::UnitTests::z().

Referenced by export_Equation().

◆ Equation_Evaluate3()

Array< OneD, NekDouble > Equation_Evaluate3 ( std::shared_ptr< Equation equation,
const Array< OneD, const NekDouble > &  x,
const Array< OneD, const NekDouble > &  y,
const Array< OneD, const NekDouble > &  z 
)

Wrapper for Equation::Evaluate (overloaded for Array parameters).

Parameters
equationEquation object from Python
xx-coordinates within the expression
yy-coordinates within the expression
zz-coordinates within the expression
Returns
Array containing values of evaluated expression.

Definition at line 140 of file Python/BasicUtils/Equation.cpp.

144{
145 Array<OneD, NekDouble> tmp(x.size());
146 equation->Evaluate(x, y, z, tmp);
147 return tmp;
148}

References Nektar::UnitTests::z().

Referenced by export_Equation().

◆ Equation_Evaluate4()

Array< OneD, NekDouble > Equation_Evaluate4 ( std::shared_ptr< Equation equation,
const Array< OneD, const NekDouble > &  x,
const Array< OneD, const NekDouble > &  y,
const Array< OneD, const NekDouble > &  z,
const NekDouble  t 
)

Wrapper for Equation::Evaluate (overloaded for Array parameters + constant time).

Parameters
equationEquation object from Python
xx-coordinates within the expression
yy-coordinates within the expression
zz-coordinates within the expression
tValue of time
Returns
Array containing values of evaluated expression.

Definition at line 162 of file Python/BasicUtils/Equation.cpp.

167{
168 Array<OneD, NekDouble> tmp(x.size());
169 equation->Evaluate(x, y, z, t, tmp);
170 return tmp;
171}

References Nektar::UnitTests::z().

Referenced by export_Equation().

◆ Equation_Evaluate5()

Array< OneD, NekDouble > Equation_Evaluate5 ( std::shared_ptr< Equation equation,
const Array< OneD, const NekDouble > &  x,
const Array< OneD, const NekDouble > &  y,
const Array< OneD, const NekDouble > &  z,
const Array< OneD, const NekDouble > &  t 
)

Wrapper for Equation::Evaluate (overloaded for Array parameters + Array time).

Parameters
equationEquation object from Python
xx-coordinates within the expression
yy-coordinates within the expression
zz-coordinates within the expression
tTime values within the expression.
Returns
Array containing values of evaluated expression.

Definition at line 185 of file Python/BasicUtils/Equation.cpp.

190{
191 Array<OneD, NekDouble> tmp(x.size());
192 equation->Evaluate(x, y, z, t, tmp);
193 return tmp;
194}

References Nektar::UnitTests::z().

Referenced by export_Equation().

◆ Equation_SetConstants()

py::object Equation_SetConstants ( py::tuple  args,
py::dict  kwargs 
)

Wrapper for Equation::SetConstants.

boost.python does not know how (by default) to convert from any Python datatype to a C++ map, so we add a pythonic way to set parameters through keyword arguments.

Parameters
argsPython function positional arguments.
kwargsPython function keyword arguments.
Returns
None (null py::object).

Definition at line 58 of file Python/BasicUtils/Equation.cpp.

59{
60 // To extract the 'self' object
61 Equation &equation = py::extract<Equation &>(args[0]);
62
63 // Map that will be passed to C++.
64 std::map<std::string, NekDouble> constants;
65
66 // Loop over the keys inside the kwargs dictionary
67 py::list keys = kwargs.keys();
68 for (int i = 0; i < py::len(keys); ++i)
69 {
70 py::object arg = kwargs[keys[i]];
71 if (arg)
72 {
73 constants[py::extract<std::string>(keys[i])] =
74 py::extract<NekDouble>(arg);
75 }
76 }
77
78 equation.SetConstants(constants);
79
80 // Return nothing.
81 return py::object();
82}
void SetConstants(const std::map< std::string, NekDouble > &constants)

References Nektar::LibUtilities::Equation::SetConstants().

Referenced by export_Equation().

◆ export_Equation()

void export_Equation ( )

Definition at line 196 of file Python/BasicUtils/Equation.cpp.

197{
198 py::class_<Equation, std::shared_ptr<Equation>, boost::noncopyable>(
199 "Equation", py::no_init)
200
201 .def("__init__", py::make_constructor(ConstructEquation))
202
203 .def("Evaluate", Equation_Evaluate1)
204 .def("Evaluate", Equation_Evaluate2)
205 .def("Evaluate", Equation_Evaluate3)
206 .def("Evaluate", Equation_Evaluate4)
207 .def("Evaluate", Equation_Evaluate5)
208
209 .def("SetParameter", &Equation::SetParameter)
210 .def("SetConstants", py::raw_function(Equation_SetConstants))
211 .def("GetExpression", &Equation::GetExpression)
212 .def("GetVlist", &Equation::GetVlist)
213
214 .def("GetTime", &Equation::GetTime);
215}
py::object Equation_SetConstants(py::tuple args, py::dict kwargs)
Wrapper for Equation::SetConstants.
std::shared_ptr< Equation > ConstructEquation(std::shared_ptr< Interpreter > evaluator, std::string expr, std::string vlist)
Construct an equation object from an expression string expr and a list of variables vlist.
Array< OneD, NekDouble > Equation_Evaluate5(std::shared_ptr< Equation > equation, const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z, const Array< OneD, const NekDouble > &t)
Wrapper for Equation::Evaluate (overloaded for Array parameters + Array time).
NekDouble Equation_Evaluate2(std::shared_ptr< Equation > equation, const NekDouble x, const NekDouble y=0, const NekDouble z=0, const NekDouble t=0)
Wrapper for Equation::Evaluate (overloaded for constant parameters).
NekDouble Equation_Evaluate1(std::shared_ptr< Equation > equation)
Wrapper for Equation::Evaluate (overloaded for no parameters).
Array< OneD, NekDouble > Equation_Evaluate4(std::shared_ptr< Equation > equation, const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z, const NekDouble t)
Wrapper for Equation::Evaluate (overloaded for Array parameters + constant time).
Array< OneD, NekDouble > Equation_Evaluate3(std::shared_ptr< Equation > equation, const Array< OneD, const NekDouble > &x, const Array< OneD, const NekDouble > &y, const Array< OneD, const NekDouble > &z)
Wrapper for Equation::Evaluate (overloaded for Array parameters).

References ConstructEquation(), Equation_Evaluate1(), Equation_Evaluate2(), Equation_Evaluate3(), Equation_Evaluate4(), Equation_Evaluate5(), and Equation_SetConstants().

Referenced by BOOST_PYTHON_MODULE().