Nektar++
Functions
Python/Interpreter/Interpreter.cpp File Reference
#include <LibUtilities/Interpreter/Interpreter.h>
#include <LibUtilities/BasicConst/NektarUnivTypeDefs.hpp>
#include <LibUtilities/Python/NekPyConfig.hpp>
#include <LibUtilities/BasicUtils/ErrorUtil.hpp>
#include <LibUtilities/Python/BasicUtils/SharedArray.hpp>

Go to the source code of this file.

Functions

void Interpreter_AddConstants (std::shared_ptr< Interpreter > interpreter, py::args args, const py::kwargs &kwargs)
 Wrapper for Interpreter::AddConstants. More...
 
void Interpreter_SetParameters (std::shared_ptr< Interpreter > interpreter, py::args args, const py::kwargs &kwargs)
 Wrapper for Interpreter::SetParameters. More...
 
NekDouble Interpreter_GetParameter (std::shared_ptr< Interpreter > interpreter, std::string paramName)
 Wrapper for Interpreter::GetParameter. More...
 
NekDouble Interpreter_GetConstant (std::shared_ptr< Interpreter > interpreter, std::string constantName)
 Wrapper for Interpreter::GetConstant. More...
 
NekDouble Interpreter_Evaluate (std::shared_ptr< Interpreter > interpreter, const int id)
 Wrapper for Interpreter::Evaluate (only constant or parameter). More...
 
NekDouble Interpreter_Evaluate2 (std::shared_ptr< Interpreter > interpreter, const int id, const NekDouble x, const NekDouble y, const NekDouble z, const NekDouble t)
 Wrapper for Interpreter::Evaluate (only constant parameters). More...
 
Array< OneD, NekDoubleInterpreter_Evaluate3 (std::shared_ptr< Interpreter > interpreter, const int id, 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 Interpreter::Evaluate (vectorised version of the evaluation method that will allow the same function to be evaluated many times). More...
 
Array< OneD, NekDoubleInterpreter_Evaluate4 (std::shared_ptr< Interpreter > interpreter, const int expression_id, const std::vector< Array< OneD, const NekDouble > > &points)
 Wrapper for Interpreter::Evaluate (zero or more arrays). More...
 
void export_Interpreter (py::module &m)
 

Function Documentation

◆ export_Interpreter()

void export_Interpreter ( py::module &  m)

Definition at line 212 of file Python/Interpreter/Interpreter.cpp.

213{
214 py::class_<Interpreter, std::shared_ptr<Interpreter>>(m, "Interpreter")
215
216 .def(py::init<>())
217
218 .def("SetRandomSeed", &Interpreter::SetRandomSeed)
219
220 .def("AddConstants", &Interpreter_AddConstants)
221 .def("AddConstant", &Interpreter::AddConstant)
222 .def("GetConstant", Interpreter_GetConstant)
223
224 .def("SetParameters", &Interpreter_SetParameters)
225 .def("SetParameter", &Interpreter::SetParameter)
226 .def("GetParameter", Interpreter_GetParameter)
227
228 .def("GetTime", &Interpreter::GetTime)
229 .def("DefineFunction", &Interpreter::DefineFunction)
230
231 .def("Evaluate", Interpreter_Evaluate)
232 .def("Evaluate", Interpreter_Evaluate2)
233 .def("Evaluate", Interpreter_Evaluate3)
234 .def("Evaluate", Interpreter_Evaluate4)
235
236 ;
237}
NekDouble Interpreter_Evaluate2(std::shared_ptr< Interpreter > interpreter, const int id, const NekDouble x, const NekDouble y, const NekDouble z, const NekDouble t)
Wrapper for Interpreter::Evaluate (only constant parameters).
Array< OneD, NekDouble > Interpreter_Evaluate4(std::shared_ptr< Interpreter > interpreter, const int expression_id, const std::vector< Array< OneD, const NekDouble > > &points)
Wrapper for Interpreter::Evaluate (zero or more arrays).
NekDouble Interpreter_Evaluate(std::shared_ptr< Interpreter > interpreter, const int id)
Wrapper for Interpreter::Evaluate (only constant or parameter).
NekDouble Interpreter_GetParameter(std::shared_ptr< Interpreter > interpreter, std::string paramName)
Wrapper for Interpreter::GetParameter.
Array< OneD, NekDouble > Interpreter_Evaluate3(std::shared_ptr< Interpreter > interpreter, const int id, 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 Interpreter::Evaluate (vectorised version of the evaluation method that will allow the sa...
NekDouble Interpreter_GetConstant(std::shared_ptr< Interpreter > interpreter, std::string constantName)
Wrapper for Interpreter::GetConstant.
void Interpreter_SetParameters(std::shared_ptr< Interpreter > interpreter, py::args args, const py::kwargs &kwargs)
Wrapper for Interpreter::SetParameters.
void Interpreter_AddConstants(std::shared_ptr< Interpreter > interpreter, py::args args, const py::kwargs &kwargs)
Wrapper for Interpreter::AddConstants.

References Interpreter_AddConstants(), Interpreter_Evaluate(), Interpreter_Evaluate2(), Interpreter_Evaluate3(), Interpreter_Evaluate4(), Interpreter_GetConstant(), Interpreter_GetParameter(), and Interpreter_SetParameters().

Referenced by PYBIND11_MODULE().

◆ Interpreter_AddConstants()

void Interpreter_AddConstants ( std::shared_ptr< Interpreter interpreter,
py::args  args,
const py::kwargs &  kwargs 
)

Wrapper for Interpreter::AddConstants.

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/Interpreter/Interpreter.cpp.

60{
61 // Map that will be passed to C++
62 std::map<std::string, NekDouble> constants;
63
64 // Loop over the keys inside the kwargs dictionary
65 for (auto &it : kwargs)
66 {
67 const std::string &key = py::cast<const std::string>(it.first);
68 const NekDouble &val = py::cast<const NekDouble>(it.second);
69 constants[key] = val;
70 }
71
72 interpreter->AddConstants(constants);
73}
double NekDouble

Referenced by export_Interpreter().

◆ Interpreter_Evaluate()

NekDouble Interpreter_Evaluate ( std::shared_ptr< Interpreter interpreter,
const int  id 
)

Wrapper for Interpreter::Evaluate (only constant or parameter).

Parameters
interpreterInterpreter object
idid of the expression
Returns
Value corresponding to id

Definition at line 142 of file Python/Interpreter/Interpreter.cpp.

144{
145 return interpreter->Evaluate(id);
146}

Referenced by export_Interpreter().

◆ Interpreter_Evaluate2()

NekDouble Interpreter_Evaluate2 ( std::shared_ptr< Interpreter interpreter,
const int  id,
const NekDouble  x,
const NekDouble  y,
const NekDouble  z,
const NekDouble  t 
)

Wrapper for Interpreter::Evaluate (only constant parameters).

Parameters
interpreterInterpreter object
idid of the expression
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 160 of file Python/Interpreter/Interpreter.cpp.

164{
165 return interpreter->Evaluate(id, x, y, z, t);
166}
std::vector< double > z(NPUPPER)

References Nektar::UnitTests::z().

Referenced by export_Interpreter().

◆ Interpreter_Evaluate3()

Array< OneD, NekDouble > Interpreter_Evaluate3 ( std::shared_ptr< Interpreter interpreter,
const int  id,
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 Interpreter::Evaluate (vectorised version of the evaluation method that will allow the same function to be evaluated many times).

Parameters
interpreterInterpreter object
idid of the expression
xx-coordinates within the expression
yy-coordinates within the expression
zz-coordinates within the expression
tvalues of time within the expression
Returns
Values of evaluated expression.

Definition at line 182 of file Python/Interpreter/Interpreter.cpp.

188{
189 Array<OneD, NekDouble> tmp(x.size());
190 interpreter->Evaluate(id, x, y, z, t, tmp);
191 return tmp;
192}

References Nektar::UnitTests::z().

Referenced by export_Interpreter().

◆ Interpreter_Evaluate4()

Array< OneD, NekDouble > Interpreter_Evaluate4 ( std::shared_ptr< Interpreter interpreter,
const int  expression_id,
const std::vector< Array< OneD, const NekDouble > > &  points 
)

Wrapper for Interpreter::Evaluate (zero or more arrays).

Parameters
expression_idid of the expression
pointsvector containing arrays of values required for the expression.
Returns
Values of evaluated expression.

Definition at line 203 of file Python/Interpreter/Interpreter.cpp.

206{
207 Array<OneD, NekDouble> tmp(points.size());
208 interpreter->Evaluate(expression_id, points, tmp);
209 return tmp;
210}

Referenced by export_Interpreter().

◆ Interpreter_GetConstant()

NekDouble Interpreter_GetConstant ( std::shared_ptr< Interpreter interpreter,
std::string  constantName 
)

Wrapper for Interpreter::GetConstant.

Parameters
interpreterInterpreter object
constantNameName of constant
Returns
Constant defined by constantName.

Definition at line 127 of file Python/Interpreter/Interpreter.cpp.

129{
130
131 return interpreter->GetConstant(constantName);
132}

Referenced by export_Interpreter().

◆ Interpreter_GetParameter()

NekDouble Interpreter_GetParameter ( std::shared_ptr< Interpreter interpreter,
std::string  paramName 
)

Wrapper for Interpreter::GetParameter.

Parameters
interpreterInterpreter object
paramNameName of parameter
Returns
Parameter defined by paramName.

Definition at line 112 of file Python/Interpreter/Interpreter.cpp.

114{
115
116 return interpreter->GetParameter(paramName);
117}

Referenced by export_Interpreter().

◆ Interpreter_SetParameters()

void Interpreter_SetParameters ( std::shared_ptr< Interpreter interpreter,
py::args  args,
const py::kwargs &  kwargs 
)

Wrapper for Interpreter::SetParameters.

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 87 of file Python/Interpreter/Interpreter.cpp.

89{
90 // Map that will be passed to C++
91 std::map<std::string, NekDouble> parameters;
92
93 // Loop over the keys inside the kwargs dictionary
94 for (auto &it : kwargs)
95 {
96 const std::string &key = py::cast<const std::string>(it.first);
97 const NekDouble &val = py::cast<const NekDouble>(it.second);
98 parameters[key] = val;
99 }
100
101 interpreter->SetParameters(parameters);
102}

Referenced by export_Interpreter().