ABLATE Source Documentation  0.12.33
cellSolver.hpp
1 #ifndef ABLATELIBRARY_CELLSOLVER_HPP
2 #define ABLATELIBRARY_CELLSOLVER_HPP
3 
4 #include <petsc.h>
5 #include <functional>
6 #include <vector>
7 #include "solver.hpp"
8 
9 namespace ablate::solver {
10 
11 class CellSolver : public solver::Solver {
12  public:
14  using AuxFieldUpdateFunction = PetscErrorCode (*)(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* u, const PetscInt aOff[],
15  PetscScalar* auxField, void* ctx);
16 
18  using SolutionFieldUpdateFunction = PetscErrorCode (*)(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], PetscScalar* u, void* ctx);
19 
20  private:
24  struct AuxFieldUpdateFunctionDescription {
25  AuxFieldUpdateFunction function;
26  void* context;
27  std::vector<PetscInt> inputFields;
28  std::vector<PetscInt> auxFields;
29  };
30 
32  std::vector<AuxFieldUpdateFunctionDescription> auxFieldUpdateFunctionDescriptions;
33 
37  struct SolutionFieldUpdateFunctionDescription {
39  void* context;
40  std::vector<PetscInt> inputFieldsOffsets;
41  };
42 
44  std::vector<SolutionFieldUpdateFunctionDescription> solutionFieldUpdateFunctionDescriptions;
45 
46  protected:
48  Vec cellGeomVec = nullptr;
49 
51  Vec faceGeomVec = nullptr;
52 
53  public:
59  explicit CellSolver(std::string solverId, std::shared_ptr<domain::Region> = {}, std::shared_ptr<parameters::Parameters> options = nullptr);
60  ~CellSolver() override;
61 
70  void RegisterAuxFieldUpdate(AuxFieldUpdateFunction function, void* context, const std::vector<std::string>& auxField, const std::vector<std::string>& inputFields);
71 
80  void RegisterSolutionFieldUpdate(SolutionFieldUpdateFunction function, void* context, const std::vector<std::string>& inputFields);
81 
88  void UpdateAuxFields(PetscReal time, Vec locXVec, Vec locAuxField);
89 
95  void UpdateSolutionFields(PetscReal time, Vec globXVec);
96 
100  void Setup() override;
101 
105  void Initialize() override;
106 };
107 } // namespace ablate::solver
108 
109 #endif // ABLATELIBRARY_CELLSOLVER_HPP
Definition: cellSolver.hpp:11
void RegisterAuxFieldUpdate(AuxFieldUpdateFunction function, void *context, const std::vector< std::string > &auxField, const std::vector< std::string > &inputFields)
Definition: cellSolver.cpp:16
Vec faceGeomVec
Vector used to describe the entire face geom of the dm. This is constant and does not depend upon reg...
Definition: cellSolver.hpp:51
Vec cellGeomVec
Vector used to describe the entire cell geom of the dm. This is constant and does not depend upon reg...
Definition: cellSolver.hpp:48
void UpdateAuxFields(PetscReal time, Vec locXVec, Vec locAuxField)
Definition: cellSolver.cpp:62
PetscErrorCode(*)(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *u, const PetscInt aOff[], PetscScalar *auxField, void *ctx) AuxFieldUpdateFunction
function template for updating the aux field
Definition: cellSolver.hpp:15
CellSolver(std::string solverId, std::shared_ptr< domain::Region >={}, std::shared_ptr< parameters::Parameters > options=nullptr)
Definition: cellSolver.cpp:4
PetscErrorCode(*)(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], PetscScalar *u, void *ctx) SolutionFieldUpdateFunction
function template for updating the solution field
Definition: cellSolver.hpp:18
void Initialize() override
Definition: cellSolver.cpp:198
void RegisterSolutionFieldUpdate(SolutionFieldUpdateFunction function, void *context, const std::vector< std::string > &inputFields)
Definition: cellSolver.cpp:42
void Setup() override
Definition: cellSolver.cpp:193
void UpdateSolutionFields(PetscReal time, Vec globXVec)
Definition: cellSolver.cpp:141
Definition: solver.hpp:17