ABLATE Source Documentation  0.12.34
probes.hpp
1 #ifndef ABLATELIBRARY_PROBES_HPP
2 #define ABLATELIBRARY_PROBES_HPP
3 
4 #include <utility>
5 #include "io/interval/interval.hpp"
6 #include "monitor.hpp"
7 #include "probes/probe.hpp"
8 #include "probes/probeInitializer.hpp"
9 
10 namespace ablate::monitors {
11 
15 class Probes : public Monitor {
16  public:
20  class ProbeRecorder {
21  private:
23  const int bufferSize;
24 
26  std::filesystem::path outputPath;
27 
29  PetscReal lastOutputTime = PETSC_MIN_REAL;
30 
32  int activeIndex = -1;
33 
35  std::vector<std::vector<double>> buffer;
36 
38  std::vector<double> timeHistory;
39 
40  public:
46  ProbeRecorder(int bufferSize, const std::vector<std::string>& variables, const std::filesystem::path& outputPath);
47 
52 
57  void AdvanceTime(double time);
58 
63  void SetValue(std::size_t index, double value);
64 
68  void WriteBuffer();
69  };
70 
72  const std::shared_ptr<ablate::monitors::probes::ProbeInitializer> initializer;
73 
75  const std::vector<std::string> variableNames;
76 
78  const std::shared_ptr<io::interval::Interval> interval;
79 
81  const int bufferSize;
82 
84  std::vector<probes::Probe> localProbes;
85 
87  std::vector<domain::Field> fields;
88 
90  std::vector<int> fieldOffset;
91 
93  std::vector<DMInterpolationInfo> interpolants;
94 
96  std::vector<ProbeRecorder> recorders;
97 
98  static PetscErrorCode UpdateProbes(TS ts, PetscInt step, PetscReal crtime, Vec u, void* ctx);
99 
100  public:
108  Probes(const std::shared_ptr<ablate::monitors::probes::ProbeInitializer>&, std::vector<std::string> variableNames, const std::shared_ptr<io::interval::Interval>& interval = {},
109  const int bufferSize = 0);
110 
111  ~Probes() override;
112 
117  void Register(std::shared_ptr<solver::Solver> solver) override;
118 
123  PetscMonitorFunction GetPetscFunction() override { return UpdateProbes; }
124 };
125 
126 } // namespace ablate::monitors
127 
128 #endif // ABLATELIBRARY_PROBES_HPP
Definition: monitor.hpp:12
void SetValue(std::size_t index, double value)
Definition: probes.cpp:251
ProbeRecorder(int bufferSize, const std::vector< std::string > &variables, const std::filesystem::path &outputPath)
Definition: probes.cpp:200
void WriteBuffer()
Definition: probes.cpp:257
void AdvanceTime(double time)
Definition: probes.cpp:239
~ProbeRecorder()
Definition: probes.cpp:237
Definition: probes.hpp:15
PetscMonitorFunction GetPetscFunction() override
Definition: probes.hpp:123
void Register(std::shared_ptr< solver::Solver > solver) override
Definition: probes.cpp:15
Probes(const std::shared_ptr< ablate::monitors::probes::ProbeInitializer > &, std::vector< std::string > variableNames, const std::shared_ptr< io::interval::Interval > &interval={}, const int bufferSize=0)
Definition: probes.cpp:8
std::vector< ProbeRecorder > recorders
list of probe recorders that goe
Definition: probes.hpp:96
std::vector< int > fieldOffset
store the offset for the field in the output (needed for multiple components)
Definition: probes.hpp:90
std::vector< DMInterpolationInfo > interpolants
list of petsc intepolants
Definition: probes.hpp:93
const std::vector< std::string > variableNames
List of variables to output.
Definition: probes.hpp:75
const int bufferSize
output bufferSize
Definition: probes.hpp:81
const std::shared_ptr< ablate::monitors::probes::ProbeInitializer > initializer
Original list of all requested probe locations by name.
Definition: probes.hpp:72
std::vector< probes::Probe > localProbes
list of local probes on this rank
Definition: probes.hpp:84
const std::shared_ptr< io::interval::Interval > interval
The sampling interval.
Definition: probes.hpp:78
std::vector< domain::Field > fields
list of fields to interpolate
Definition: probes.hpp:87
Definition: boundarySolverMonitor.hpp:16