ABLATE Source Documentation  0.12.33
sublimation.hpp
1 #ifndef ABLATELIBRARY_SUBLIMATION_HPP
2 #define ABLATELIBRARY_SUBLIMATION_HPP
3 
4 #include "boundarySolver/boundaryProcess.hpp"
5 #include "boundarySolver/physics/subModels/sublimationModel.hpp"
6 #include "eos/radiationProperties/zimmer.hpp"
7 #include "eos/transport/transportModel.hpp"
8 #include "finiteVolume/processes/navierStokesTransport.hpp"
9 #include "finiteVolume/processes/pressureGradientScaling.hpp"
10 #include "io/interval/interval.hpp"
11 #include "radiation/surfaceRadiation.hpp"
12 
13 namespace ablate::boundarySolver::physics {
14 
19  private:
20  // static name of this model
21  inline const static std::string sublimationId = "Sublimation";
22 
24  const std::shared_ptr<ablate::eos::transport::TransportModel> transportModel = nullptr;
25  const std::shared_ptr<ablate::eos::EOS> eos;
26  const std::shared_ptr<mathFunctions::MathFunction> additionalHeatFlux;
27  PetscReal currentTime = 0.0;
28 
30  const std::shared_ptr<ablate::mathFunctions::FieldFunction> massFractions;
31  const mathFunctions::PetscFunction massFractionsFunction;
32  void *massFractionsContext;
33  PetscInt numberSpecies = 0;
34 
38  const bool diffusionFlame;
39 
43  const std::shared_ptr<finiteVolume::processes::PressureGradientScaling> pressureGradientScaling;
44 
48  std::shared_ptr<ablate::radiation::SurfaceRadiation> radiation;
49 
53  eos::ThermodynamicTemperatureFunction effectiveConductivity;
54 
58  eos::ThermodynamicTemperatureFunction viscosityFunction;
59 
63  eos::ThermodynamicTemperatureFunction computeTemperatureFunction;
64 
68  eos::ThermodynamicTemperatureFunction computeSensibleEnthalpy;
69 
74 
78  const std::shared_ptr<io::interval::Interval> radiationInterval;
79 
83  const double emissivity;
84 
88  const double absorptivity;
89 
93  void UpdateSpecies(TS ts, ablate::solver::Solver &);
94 
98  std::shared_ptr<subModels::SublimationModel> sublimationModel = nullptr;
99 
100  // Hold onto a BoundaryPreRHSPointFunctionDefinition to precompute fe heat transfer if returned from the sublimation model
101  BoundarySolver::BoundaryPreRHSPointFunctionDefinition solidHeatTransferUpdateFunctionDefinition{};
102 
103  public:
104  explicit Sublimation(std::shared_ptr<subModels::SublimationModel> sublimationModel, std::shared_ptr<ablate::eos::transport::TransportModel> transportModel, std::shared_ptr<ablate::eos::EOS> eos,
105  const std::shared_ptr<ablate::mathFunctions::FieldFunction> & = {}, std::shared_ptr<mathFunctions::MathFunction> additionalHeatFlux = {},
106  std::shared_ptr<finiteVolume::processes::PressureGradientScaling> pressureGradientScaling = {}, bool diffusionFlame = false,
107  std::shared_ptr<ablate::radiation::SurfaceRadiation> radiationIn = {}, const std::shared_ptr<io::interval::Interval> &intervalIn = {},
108  const std::shared_ptr<ablate::parameters::Parameters> & = {});
109 
110  void Setup(ablate::boundarySolver::BoundarySolver &bSolver) override;
111  void Initialize(ablate::boundarySolver::BoundarySolver &bSolver) override;
112 
117  void Setup(PetscInt numberSpecies);
118 
123  [[nodiscard]] const std::string &GetId() const override { return sublimationId; }
124 
129  [[nodiscard]] SerializerType Serialize() const override { return sublimationModel ? sublimationModel->Serialize() : io::Serializable::SerializerType::none; }
130 
137  PetscErrorCode Save(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override {
138  PetscFunctionBegin;
139  PetscCall(sublimationModel->Save(viewer, sequenceNumber, time));
140  PetscFunctionReturn(PETSC_SUCCESS);
141  };
142 
149  PetscErrorCode Restore(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override {
150  PetscFunctionBegin;
151  PetscCall(sublimationModel->Restore(viewer, sequenceNumber, time));
152  PetscFunctionReturn(PETSC_SUCCESS);
153  }
154 
161  static PetscErrorCode SublimationPreRHS(BoundarySolver &, TS ts, PetscReal time, bool initialStage, Vec locX, void *ctx);
162 
182  static PetscErrorCode SublimationFunction(PetscInt dim, const ablate::boundarySolver::BoundarySolver::BoundaryFVFaceGeom *fg, const PetscFVCellGeom *boundaryCell, const PetscInt uOff[],
183  const PetscScalar *boundaryValues, const PetscScalar *stencilValues[], const PetscInt aOff[], const PetscScalar *auxValues,
184  const PetscScalar *stencilAuxValues[], PetscInt stencilSize, const PetscInt stencil[], const PetscScalar stencilWeights[], const PetscInt sOff[],
185  PetscScalar source[], void *ctx);
186 
206  static PetscErrorCode SublimationOutputFunction(PetscInt dim, const ablate::boundarySolver::BoundarySolver::BoundaryFVFaceGeom *fg, const PetscFVCellGeom *boundaryCell, const PetscInt uOff[],
207  const PetscScalar *boundaryValues, const PetscScalar *stencilValues[], const PetscInt aOff[], const PetscScalar *auxValues,
208  const PetscScalar *stencilAuxValues[], PetscInt stencilSize, const PetscInt stencil[], const PetscScalar stencilWeights[], const PetscInt sOff[],
209  PetscScalar source[], void *ctx);
210 
230  static PetscErrorCode UpdateBoundaryHeatTransferModel(PetscReal time, PetscReal dt, PetscInt dim, const ablate::boundarySolver::BoundarySolver::BoundaryFVFaceGeom *fg,
231  const PetscFVCellGeom *boundaryCell, const PetscInt uOff[], PetscScalar *boundaryValues, const PetscScalar *stencilValues[],
232  const PetscInt aOff[], PetscScalar *auxValues, const PetscScalar *stencilAuxValues[], PetscInt stencilSize, const PetscInt stencil[],
233  const PetscScalar stencilWeights[], void *ctx);
234 };
235 
236 } // namespace ablate::boundarySolver::physics
237 
238 #endif // ABLATELIBRARY_SUBLIMATION_HPP
Definition: boundaryProcess.hpp:8
Definition: boundarySolver.hpp:13
Definition: sublimation.hpp:18
SerializerType Serialize() const override
Definition: sublimation.hpp:129
void Initialize(ablate::boundarySolver::BoundarySolver &bSolver) override
Definition: sublimation.cpp:107
static PetscErrorCode SublimationFunction(PetscInt dim, const ablate::boundarySolver::BoundarySolver::BoundaryFVFaceGeom *fg, const PetscFVCellGeom *boundaryCell, const PetscInt uOff[], const PetscScalar *boundaryValues, const PetscScalar *stencilValues[], const PetscInt aOff[], const PetscScalar *auxValues, const PetscScalar *stencilAuxValues[], PetscInt stencilSize, const PetscInt stencil[], const PetscScalar stencilWeights[], const PetscInt sOff[], PetscScalar source[], void *ctx)
Definition: sublimation.cpp:137
PetscErrorCode Restore(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override
Definition: sublimation.hpp:149
PetscErrorCode Save(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override
Definition: sublimation.hpp:137
const std::string & GetId() const override
Definition: sublimation.hpp:123
static PetscErrorCode UpdateBoundaryHeatTransferModel(PetscReal time, PetscReal dt, PetscInt dim, const ablate::boundarySolver::BoundarySolver::BoundaryFVFaceGeom *fg, const PetscFVCellGeom *boundaryCell, const PetscInt uOff[], PetscScalar *boundaryValues, const PetscScalar *stencilValues[], const PetscInt aOff[], PetscScalar *auxValues, const PetscScalar *stencilAuxValues[], PetscInt stencilSize, const PetscInt stencil[], const PetscScalar stencilWeights[], void *ctx)
Definition: sublimation.cpp:429
static PetscErrorCode SublimationOutputFunction(PetscInt dim, const ablate::boundarySolver::BoundarySolver::BoundaryFVFaceGeom *fg, const PetscFVCellGeom *boundaryCell, const PetscInt uOff[], const PetscScalar *boundaryValues, const PetscScalar *stencilValues[], const PetscInt aOff[], const PetscScalar *auxValues, const PetscScalar *stencilAuxValues[], PetscInt stencilSize, const PetscInt stencil[], const PetscScalar stencilWeights[], const PetscInt sOff[], PetscScalar source[], void *ctx)
Definition: sublimation.cpp:280
static PetscErrorCode SublimationPreRHS(BoundarySolver &, TS ts, PetscReal time, bool initialStage, Vec locX, void *ctx)
Definition: sublimation.cpp:359
void Setup(ablate::boundarySolver::BoundarySolver &bSolver) override
Definition: sublimation.cpp:31
Definition: serializable.hpp:13
SerializerType
Definition: serializable.hpp:18
Definition: solver.hpp:17