1 #ifndef ABLATELIBRARY_FINITEVOLUMESOLVER_HPP
2 #define ABLATELIBRARY_FINITEVOLUMESOLVER_HPP
6 #include "boundaryConditions/boundaryCondition.hpp"
7 #include "cellInterpolant.hpp"
9 #include "faceInterpolant.hpp"
10 #include "mathFunctions/fieldFunction.hpp"
11 #include "solver/cellSolver.hpp"
12 #include "solver/solver.hpp"
13 #include "solver/timeStepper.hpp"
14 #include "utilities/vectorUtilities.hpp"
16 namespace ablate::finiteVolume {
30 using PreRHSFunctionDefinition = PetscErrorCode (*)(
FiniteVolumeSolver&, TS ts, PetscReal time,
bool initialStage, Vec locX,
void* ctx);
31 using RHSArbitraryFunction = PetscErrorCode (*)(
const FiniteVolumeSolver&, DM dm, PetscReal time, Vec locXVec, Vec locFVec,
void* ctx);
41 struct ComputeTimeStepDescription {
42 ComputeTimeStepFunction
function;
48 std::vector<CellInterpolant::DiscontinuousFluxFunctionDescription> discontinuousFluxFunctionDescriptions;
49 std::vector<FaceInterpolant::ContinuousFluxFunctionDescription> continuousFluxFunctionDescriptions;
50 std::vector<CellInterpolant::PointFunctionDescription> pointFunctionDescriptions;
53 std::vector<std::pair<RHSArbitraryFunction, void*>> rhsArbitraryFunctions;
56 std::vector<std::pair<PreRHSFunctionDefinition, void*>> preRhsFunctions;
59 std::vector<ComputeTimeStepDescription> timeStepFunctions;
62 std::vector<std::shared_ptr<processes::Process>> processes;
65 const std::vector<std::shared_ptr<boundaryConditions::BoundaryCondition>> boundaryConditions;
68 std::unique_ptr<FaceInterpolant> faceInterpolant =
nullptr;
71 std::unique_ptr<CellInterpolant> cellInterpolant =
nullptr;
74 std::shared_ptr<domain::Region> solverRegionMinusGhost;
77 DM meshCharacteristicsDm =
nullptr;
80 Vec meshCharacteristicsLocalVec =
nullptr;
83 FiniteVolumeSolver(std::string solverId, std::shared_ptr<domain::Region>, std::shared_ptr<parameters::Parameters> options, std::vector<std::shared_ptr<processes::Process>> flowProcesses,
84 std::vector<std::shared_ptr<boundaryConditions::BoundaryCondition>> boundaryConditions);
90 void Setup()
override;
111 PetscErrorCode
ComputeBoundary(PetscReal time, Vec locX, Vec locX_t)
override;
122 const std::vector<std::string>& auxFields);
133 const std::vector<std::string>& auxFields);
144 const std::vector<std::string>& auxFields);
198 PetscErrorCode
Save(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time)
override;
206 PetscErrorCode
Restore(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time)
override;
224 return utilities::VectorUtilities::Find<T>(processes);
233 PetscErrorCode
PreRHSFunction(TS ts, PetscReal time,
bool initialStage, Vec locX)
override;
239 dm = meshCharacteristicsDm;
240 vec = meshCharacteristicsLocalVec;
PetscErrorCode(*)(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[], const PetscInt aOff[], const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar flux[], void *ctx) DiscontinuousFluxFunction
Definition: cellInterpolant.hpp:17
PetscErrorCode(*)(PetscInt dim, PetscReal time, const PetscFVCellGeom *cg, const PetscInt uOff[], const PetscScalar u[], const PetscInt aOff[], const PetscScalar a[], PetscScalar f[], void *ctx) PointFunction
Definition: cellInterpolant.hpp:23
PetscErrorCode(*)(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void *ctx) ContinuousFluxFunction
Definition: faceInterpolant.hpp:79
Definition: finiteVolumeSolver.hpp:28
void Setup() override
Definition: finiteVolumeSolver.cpp:28
void GetCellRangeWithoutGhost(ablate::domain::Range &faceRange) const
Definition: finiteVolumeSolver.cpp:439
std::shared_ptr< T > FindProcess()
Definition: finiteVolumeSolver.hpp:223
void RegisterRHSFunction(CellInterpolant::DiscontinuousFluxFunction function, void *context, const std::string &field, const std::vector< std::string > &inputFields, const std::vector< std::string > &auxFields)
Definition: finiteVolumeSolver.cpp:293
PetscErrorCode ComputeBoundary(PetscReal time, Vec locX, Vec locX_t) override
Definition: finiteVolumeSolver.cpp:457
const std::string & GetId() const override
Definition: finiteVolumeSolver.hpp:190
~FiniteVolumeSolver() override
cleanup
Definition: finiteVolumeSolver.cpp:19
void RegisterPreRHSFunction(PreRHSFunctionDefinition function, void *context)
Definition: finiteVolumeSolver.cpp:360
PetscErrorCode Save(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override
Definition: finiteVolumeSolver.cpp:393
PetscErrorCode ComputeRHSFunction(PetscReal time, Vec locXVec, Vec locFVec) override
Definition: finiteVolumeSolver.cpp:233
double ComputePhysicsTimeStep(TS) override
Definition: finiteVolumeSolver.cpp:366
PetscErrorCode PreRHSFunction(TS ts, PetscReal time, bool initialStage, Vec locX) override
Definition: finiteVolumeSolver.cpp:466
void GetMeshCharacteristics(DM &dm, Vec &vec)
Definition: finiteVolumeSolver.hpp:238
std::map< std::string, double > ComputePhysicsTimeSteps(TS) override
Definition: finiteVolumeSolver.cpp:376
MeshCharacteristics
store an enum for the fields in the meshCharacteristicsDm
Definition: finiteVolumeSolver.hpp:35
SerializerType Serialize() const override
Definition: finiteVolumeSolver.cpp:391
PetscErrorCode Restore(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override
Definition: finiteVolumeSolver.cpp:427
void RegisterComputeTimeStepFunction(ComputeTimeStepFunction function, void *ctx, std::string name)
Definition: finiteVolumeSolver.cpp:362
void Initialize() override
Definition: finiteVolumeSolver.cpp:47
Definition: serializable.hpp:13
SerializerType
Definition: serializable.hpp:18
Definition: boundaryFunction.hpp:8
Definition: cellSolver.hpp:11
Definition: physicsTimeStepFunction.hpp:7
Definition: rhsFunction.hpp:7
const std::string & GetSolverId() const
Definition: solver.hpp:58
Definition: loggable.hpp:9