1 #ifndef ABLATELIBRARY_SOLVER_HPP
2 #define ABLATELIBRARY_SOLVER_HPP
4 #include <domain/region.hpp>
5 #include <domain/subDomain.hpp>
7 #include <parameters/parameters.hpp>
10 #include "domain/range.hpp"
11 #include "io/serializable.hpp"
13 namespace ablate::solver {
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;
26 const std::string solverId;
29 const std::shared_ptr<domain::Region> region;
33 PetscOptions petscOptions;
36 std::shared_ptr<ablate::domain::SubDomain> subDomain;
39 explicit Solver(std::string solverId, std::shared_ptr<domain::Region> = {}, std::shared_ptr<parameters::Parameters> options =
nullptr);
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);
49 virtual void Register(std::shared_ptr<ablate::domain::SubDomain> subDomain);
55 virtual void Initialize() = 0;
58 [[nodiscard]]
inline const std::string&
GetSolverId()
const {
return solverId; }
76 [[nodiscard]]
inline std::shared_ptr<domain::Region>
GetRegion() const noexcept {
return region; }
79 void PreStage(TS ts, PetscReal stagetime);
82 void PostEvaluate(TS ts);
88 inline void RegisterPreStep(
const std::function<
void(TS ts,
Solver&)>& preStep) { this->preStepFunctions.push_back(preStep); }
94 inline void RegisterPreStage(
const std::function<
void(TS ts,
Solver&, PetscReal)>& preStage) { this->preStageFunctions.push_back(preStage); }
100 inline void RegisterPostStep(
const std::function<
void(TS ts,
Solver&)>& postStep) { this->postStepFunctions.push_back(postStep); }
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
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