ABLATE Source Documentation  0.12.33
peak.hpp
1 #ifndef ABLATELIBRARY_PEAK_FUNTION_HPP
2 #define ABLATELIBRARY_PEAK_FUNTION_HPP
3 #include <muParser.h>
4 #include "formulaBase.hpp"
5 
6 namespace ablate::mathFunctions {
11 class Peak : public MathFunction {
12  private:
13  const std::vector<double> startValue;
14  const std::vector<double> peakValue;
15  const std::vector<double> endValue;
16  const double start;
17  const double peak;
18  const double end;
19  const int dir;
20 
21  static PetscErrorCode PeakPetscFunction(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar* u, void* ctx);
22 
34  inline static double Interpolate(double x, double xS, double xP, double xE, double yS, double yP, double yE) {
35  if (x < xS) {
36  return yS;
37  } else if (x > xE) {
38  return yE;
39  } else if (x < xP) {
40  return yS + (x - xS) * (yP - yS) / (xP - xS);
41  } else {
42  return yP + (x - xP) * (yE - yP) / (xE - xP);
43  }
44  }
45 
55  inline double DetermineDirectionValue(double x, double y, double z) const {
56  switch (dir) {
57  case 0:
58  return x;
59  case 1:
60  return y;
61  case 2:
62  return z;
63  default:
64  return x;
65  }
66  }
67 
68  public:
69  Peak(const Peak&) = delete;
70  void operator=(const Peak&) = delete;
71 
80  explicit Peak(std::vector<double> startValue, std::vector<double> peakValue, std::vector<double> endValue, double start, double peak, double end, int dir);
81 
82  double Eval(const double& x, const double& y, const double& z, const double& t) const override;
83 
84  double Eval(const double* xyz, const int& ndims, const double& t) const override;
85 
86  void Eval(const double& x, const double& y, const double& z, const double& t, std::vector<double>& result) const override;
87 
88  void Eval(const double* xyz, const int& ndims, const double& t, std::vector<double>& result) const override;
89 
90  void* GetContext() override { return this; }
91 
92  PetscFunction GetPetscFunction() override { return PeakPetscFunction; }
93 };
94 } // namespace ablate::mathFunctions
95 
96 #endif // ABLATELIBRARY_SIMPLEFORMULA_HPP
Definition: mathFunction.hpp:13
Definition: peak.hpp:11
PetscFunction GetPetscFunction() override
Definition: peak.hpp:92
void * GetContext() override
Definition: peak.hpp:90
double Eval(const double &x, const double &y, const double &z, const double &t) const override
Definition: peak.cpp:18