ABLATE Source Documentation  0.12.34
solver.hpp
1 #ifndef ABLATELIBRARY_SOLVER_HPP
2 #define ABLATELIBRARY_SOLVER_HPP
3 
4 #include <domain/region.hpp>
5 #include <domain/subDomain.hpp>
6 #include <functional>
7 #include <parameters/parameters.hpp>
8 #include <string>
9 #include <vector>
10 #include "domain/range.hpp"
11 #include "io/serializable.hpp"
12 
13 namespace ablate::solver {
14 
15 class TimeStepper;
16 
17 class Solver {
18  private:
19  // pre and post step functions for the flow
20  std::vector<std::function<void(TS ts, Solver&)>> preStepFunctions;
21  std::vector<std::function<void(TS ts, Solver&, PetscReal)>> preStageFunctions;
22  std::vector<std::function<void(TS ts, Solver&)>> postStepFunctions;
23  std::vector<std::function<void(TS ts, Solver&)>> postEvaluateFunctions;
24 
25  // The name of this solver
26  const std::string solverId;
27 
28  // The region of this solver.
29  const std::shared_ptr<domain::Region> region;
30 
31  protected:
32  // an optional petscOptions that is used for this solver
33  PetscOptions petscOptions;
34 
35  // use the subDomain to setup the problem
36  std::shared_ptr<ablate::domain::SubDomain> subDomain;
37 
38  // The constructor to be call by any Solve implementation
39  explicit Solver(std::string solverId, std::shared_ptr<domain::Region> = {}, std::shared_ptr<parameters::Parameters> options = nullptr);
40 
41  // Replacement calls for PETSC versions allowing multiple DS
42  static PetscErrorCode DMPlexInsertBoundaryValues_Plex(DM dm, PetscDS ds, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM);
43  static PetscErrorCode DMPlexInsertTimeDerivativeBoundaryValues_Plex(DM dm, PetscDS ds, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM);
44 
45  public:
46  virtual ~Solver();
47 
49  virtual void Register(std::shared_ptr<ablate::domain::SubDomain> subDomain);
50 
52  virtual void Setup() = 0;
53 
54  /*** Set up mesh dependent initialization, this may be called multiple times if the mesh changes **/
55  virtual void Initialize() = 0;
56 
58  [[nodiscard]] inline const std::string& GetSolverId() const { return solverId; }
59 
64  inline ablate::domain::SubDomain& GetSubDomain() noexcept { return *subDomain; }
65 
70  [[nodiscard]] inline const ablate::domain::SubDomain& GetSubDomain() const noexcept { return *subDomain; }
71 
76  [[nodiscard]] inline std::shared_ptr<domain::Region> GetRegion() const noexcept { return region; }
77 
78  // Support for timestepping calls
79  void PreStage(TS ts, PetscReal stagetime);
80  void PreStep(TS ts);
81  void PostStep(TS ts);
82  void PostEvaluate(TS ts);
83 
88  inline void RegisterPreStep(const std::function<void(TS ts, Solver&)>& preStep) { this->preStepFunctions.push_back(preStep); }
89 
94  inline void RegisterPreStage(const std::function<void(TS ts, Solver&, PetscReal)>& preStage) { this->preStageFunctions.push_back(preStage); }
95 
100  inline void RegisterPostStep(const std::function<void(TS ts, Solver&)>& postStep) { this->postStepFunctions.push_back(postStep); }
101 
106  inline void RegisterPostEvaluate(const std::function<void(TS ts, Solver&)>& postEval) { this->postEvaluateFunctions.push_back(postEval); }
107 
112  void GetCellRange(ablate::domain::Range& cellRange) const { this->subDomain->GetCellRange(this->GetRegion(), cellRange); };
113 
118  void GetFaceRange(ablate::domain::Range& faceRange) const { this->subDomain->GetFaceRange(this->GetRegion(), faceRange); }
119 
125  void GetRange(PetscInt depth, ablate::domain::Range& range) const { this->subDomain->GetRange(this->GetRegion(), depth, range); }
126 
134  void RestoreRange(ablate::domain::Range& range) const { this->subDomain->RestoreRange(range); }
135 };
136 
137 } // namespace ablate::solver
138 #endif // ABLATELIBRARY_SOLVER_HPP
Definition: subDomain.hpp:19
Definition: solver.hpp:17
virtual void Register(std::shared_ptr< ablate::domain::SubDomain > subDomain)
Definition: solver.cpp:21
const std::string & GetSolverId() const
Definition: solver.hpp:58
virtual void Setup()=0
void GetFaceRange(ablate::domain::Range &faceRange) const
Definition: solver.hpp:118
void RegisterPostStep(const std::function< void(TS ts, Solver &)> &postStep)
Definition: solver.hpp:100
void RegisterPreStep(const std::function< void(TS ts, Solver &)> &preStep)
Definition: solver.hpp:88
void GetRange(PetscInt depth, ablate::domain::Range &range) const
Definition: solver.hpp:125
void RegisterPostEvaluate(const std::function< void(TS ts, Solver &)> &postEval)
Definition: solver.hpp:106
ablate::domain::SubDomain & GetSubDomain() noexcept
Definition: solver.hpp:64
std::shared_ptr< domain::Region > GetRegion() const noexcept
Definition: solver.hpp:76
const ablate::domain::SubDomain & GetSubDomain() const noexcept
Definition: solver.hpp:70
void GetCellRange(ablate::domain::Range &cellRange) const
Definition: solver.hpp:112
void RegisterPreStage(const std::function< void(TS ts, Solver &, PetscReal)> &preStage)
Definition: solver.hpp:94
void RestoreRange(ablate::domain::Range &range) const
Definition: solver.hpp:134
Definition: range.hpp:11