1 #ifndef ABLATELIBRARY_LINEARTABLE_HPP
2 #define ABLATELIBRARY_LINEARTABLE_HPP
7 #include "mathFunction.hpp"
8 namespace ablate::mathFunctions {
19 std::vector<double> independentValues;
20 std::vector<std::vector<double>> dependentValues;
21 const std::string independentColumnName;
22 const std::vector<std::string> dependentColumnsNames;
23 const std::shared_ptr<MathFunction> independentValueFunction;
26 void ParseInputData(std::istream& inputFile);
28 void Interpolate(
double x,
size_t numInterpolations,
double* result)
const;
30 static PetscErrorCode LinearInterpolatorPetscFunction(PetscInt dim, PetscReal time,
const PetscReal x[], PetscInt Nf, PetscScalar* u,
void* ctx);
33 LinearTable(std::filesystem::path inputFile, std::string independentColumnName, std::vector<std::string> dependentColumnsNames, std::shared_ptr<MathFunction> independentValueFunction);
34 LinearTable(std::istream& inputStream, std::string independentColumnName, std::vector<std::string> dependentColumnsNames, std::shared_ptr<MathFunction> independentValueFunction);
36 double Eval(
const double& x,
const double& y,
const double& z,
const double& t)
const override;
38 double Eval(
const double* xyz,
const int& ndims,
const double& t)
const override;
40 void Eval(
const double& x,
const double& y,
const double& z,
const double& t, std::vector<double>& result)
const override;
42 void Eval(
const double* xyz,
const int& ndims,
const double& t, std::vector<double>& result)
const override;
48 const std::vector<double>& GetIndependentValues()
const {
return independentValues; }
50 const std::vector<std::vector<double>>& GetDependentValues()
const {
return dependentValues; }
Definition: linearTable.hpp:17
double Eval(const double &x, const double &y, const double &z, const double &t) const override
Definition: linearTable.cpp:106
PetscFunction GetPetscFunction() override
Definition: linearTable.hpp:46
void * GetContext() override
Definition: linearTable.hpp:44
Definition: mathFunction.hpp:13