59 py::args args,
const py::kwargs &kwargs)
62 std::map<std::string, NekDouble> constants;
65 for (
auto &it : kwargs)
67 const std::string &key = py::cast<const std::string>(it.first);
68 const NekDouble &val = py::cast<const NekDouble>(it.second);
72 interpreter->AddConstants(constants);
88 py::args args,
const py::kwargs &kwargs)
91 std::map<std::string, NekDouble> parameters;
94 for (
auto &it : kwargs)
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;
101 interpreter->SetParameters(parameters);
113 std::string paramName)
116 return interpreter->GetParameter(paramName);
128 std::string constantName)
131 return interpreter->GetConstant(constantName);
145 return interpreter->Evaluate(
id);
165 return interpreter->Evaluate(
id, x, y,
z, t);
183 std::shared_ptr<Interpreter> interpreter,
const int id,
190 interpreter->Evaluate(
id, x, y,
z, t, tmp);
204 std::shared_ptr<Interpreter> interpreter,
const int expression_id,
208 interpreter->Evaluate(expression_id, points, tmp);
214 py::class_<Interpreter, std::shared_ptr<Interpreter>>(m,
"Interpreter")
218 .def(
"SetRandomSeed", &Interpreter::SetRandomSeed)
221 .def(
"AddConstant", &Interpreter::AddConstant)
225 .def(
"SetParameter", &Interpreter::SetParameter)
228 .def(
"GetTime", &Interpreter::GetTime)
229 .def(
"DefineFunction", &Interpreter::DefineFunction)
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.
void export_Interpreter(py::module &m)
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.
std::vector< double > z(NPUPPER)