1 #ifndef ABLATELIBRARY_DOMAIN_HDF5INITIALIZER_HPP
2 #define ABLATELIBRARY_DOMAIN_HDF5INITIALIZER_HPP
9 #include "initializer.hpp"
10 #include "mathFunctions/fieldFunction.hpp"
12 namespace ablate::domain {
22 const std::shared_ptr<ablate::domain::Region>
region;
33 [[nodiscard]] std::vector<std::shared_ptr<mathFunctions::FieldFunction>>
GetFieldFunctions(
const std::vector<domain::Field>& fields)
const override;
42 PetscViewer petscViewer =
nullptr;
48 PetscInt numberCells = -1;
72 const std::string field;
78 Vec fieldVec =
nullptr;
81 std::shared_ptr<Hdf5Mesh> baseMesh;
84 PetscInt components = -1;
90 PetscInt resultSize = -1;
113 virtual PetscErrorCode
Eval(PetscInt dim,
const PetscReal x[], PetscScalar* u)
const;
126 static PetscErrorCode Hdf5PetscFunction(PetscInt dim, PetscReal time,
const PetscReal x[], PetscInt Nf, PetscScalar* u,
void* ctx);
137 [[nodiscard]]
double Eval(
const double& x,
const double& y,
const double& z,
const double& t)
const override;
146 [[nodiscard]]
double Eval(
const double* xyz,
const int& ndims,
const double& t)
const override;
156 void Eval(
const double& x,
const double& y,
const double& z,
const double& t, std::vector<double>& result)
const override;
165 void Eval(
const double* xyz,
const int& ndims,
const double& t, std::vector<double>& result)
const override;
171 ablate::mathFunctions::PetscFunction
GetPetscFunction()
override {
return Hdf5PetscFunction; }
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