ABLATE Source Documentation  0.12.33
hdf5Initializer.hpp
1 #ifndef ABLATELIBRARY_DOMAIN_HDF5INITIALIZER_HPP
2 #define ABLATELIBRARY_DOMAIN_HDF5INITIALIZER_HPP
3 
4 #include <filesystem>
5 #include <memory>
6 #include <utility>
7 #include <vector>
8 #include "field.hpp"
9 #include "initializer.hpp"
10 #include "mathFunctions/fieldFunction.hpp"
11 
12 namespace ablate::domain {
16 class Hdf5Initializer : public Initializer {
17  protected:
19  std::filesystem::path hdf5Path;
20 
22  const std::shared_ptr<ablate::domain::Region> region;
23 
24  public:
28  explicit Hdf5Initializer(std::filesystem::path hdf5Path, std::shared_ptr<ablate::domain::Region> region = {});
29 
33  [[nodiscard]] std::vector<std::shared_ptr<mathFunctions::FieldFunction>> GetFieldFunctions(const std::vector<domain::Field>& fields) const override;
34 
35  protected:
39  class Hdf5Mesh {
40  private:
41  // Pointer to the viewer that can be used to load the mesh or other vectors
42  PetscViewer petscViewer = nullptr;
43 
44  // Pointer to the base dm
45  DM dm = nullptr;
46 
47  // count the number of cells in this mesh
48  PetscInt numberCells = -1;
49 
50  public:
55  explicit Hdf5Mesh(const std::filesystem::path& hdf5Path);
56 
60  ~Hdf5Mesh();
61 
62  // Allow the Hdf5Initializer to access hdf5 mesh private variables
63  friend class Hdf5Initializer;
64  };
65 
70  protected:
71  // the field used to represent this math function (useful for debugging)
72  const std::string field;
73 
74  // The dm for this specific field/math function
75  DM fieldDm = nullptr;
76 
77  // create an empty vec
78  Vec fieldVec = nullptr;
79 
80  // Hold a reference to the base mesh,viewer
81  std::shared_ptr<Hdf5Mesh> baseMesh;
82 
83  // The number of components in the field
84  PetscInt components = -1;
85 
86  // The number of dimensions in the dm
87  PetscInt dim = -1;
88 
89  // The result size, this is normally the number of components but can be overwritten by sub classes
90  PetscInt resultSize = -1;
91 
92  public:
98  Hdf5MathFunction(std::shared_ptr<Hdf5Mesh> baseMesh, std::string field);
99 
103  ~Hdf5MathFunction() override;
104 
113  virtual PetscErrorCode Eval(PetscInt dim, const PetscReal x[], PetscScalar* u) const;
114 
115  private:
126  static PetscErrorCode Hdf5PetscFunction(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar* u, void* ctx);
127 
128  public:
137  [[nodiscard]] double Eval(const double& x, const double& y, const double& z, const double& t) const override;
138 
146  [[nodiscard]] double Eval(const double* xyz, const int& ndims, const double& t) const override;
147 
156  void Eval(const double& x, const double& y, const double& z, const double& t, std::vector<double>& result) const override;
157 
165  void Eval(const double* xyz, const int& ndims, const double& t, std::vector<double>& result) const override;
166 
171  ablate::mathFunctions::PetscFunction GetPetscFunction() override { return Hdf5PetscFunction; }
172 
177  void* GetContext() override { return this; }
178  };
179 };
180 
181 } // namespace ablate::domain
182 
183 #endif // ABLATELIBRARY_DOMAIN_HDF5INITIALIZER_HPP
Definition: hdf5Initializer.hpp:69
virtual PetscErrorCode Eval(PetscInt dim, const PetscReal x[], PetscScalar *u) const
Definition: hdf5Initializer.cpp:91
void * GetContext() override
Definition: hdf5Initializer.hpp:177
ablate::mathFunctions::PetscFunction GetPetscFunction() override
Definition: hdf5Initializer.hpp:171
~Hdf5MathFunction() override
Definition: hdf5Initializer.cpp:86
Hdf5MathFunction(std::shared_ptr< Hdf5Mesh > baseMesh, std::string field)
Definition: hdf5Initializer.cpp:52
Definition: hdf5Initializer.hpp:39
Hdf5Mesh(const std::filesystem::path &hdf5Path)
Definition: hdf5Initializer.cpp:27
~Hdf5Mesh()
Definition: hdf5Initializer.cpp:44
Definition: hdf5Initializer.hpp:16
const std::shared_ptr< ablate::domain::Region > region
the region to apply this initializer
Definition: hdf5Initializer.hpp:22
std::vector< std::shared_ptr< mathFunctions::FieldFunction > > GetFieldFunctions(const std::vector< domain::Field > &fields) const override
Definition: hdf5Initializer.cpp:8
std::filesystem::path hdf5Path
the path to the hdf5 that contains the stored data
Definition: hdf5Initializer.hpp:19
Hdf5Initializer(std::filesystem::path hdf5Path, std::shared_ptr< ablate::domain::Region > region={})
Definition: hdf5Initializer.cpp:6
Definition: initializer.hpp:14
Definition: mathFunction.hpp:13