ABLATE Source Documentation  0.12.33
arrheniusSublimation.hpp
1 #ifndef ABLATELIBRARY_ARRHENIUSSUBLIMATION_HPP
2 #define ABLATELIBRARY_ARRHENIUSSUBLIMATION_HPP
3 
4 #include <map>
5 #include <memory>
6 #include "oneDimensionHeatTransfer.hpp"
7 #include "solver/cellSolver.hpp"
8 #include "solver/timeStepper.hpp"
9 #include "sublimationModel.hpp"
10 
11 namespace ablate::boundarySolver::physics::subModels {
12 
14  private:
16  std::map<PetscInt, std::shared_ptr<OneDimensionHeatTransfer>> oneDimensionHeatTransfer;
17 
19  const std::shared_ptr<ablate::parameters::Parameters> properties;
20 
22  const std::shared_ptr<ablate::mathFunctions::MathFunction> initialization;
23 
25  const std::shared_ptr<ablate::parameters::Parameters> options;
26 
28  const PetscReal latentHeatOfFusion;
29 
31  const PetscReal solidDensity;
32 
34  const PetscReal preExponentialFactor;
35  const PetscReal activationEnergy;
36  const PetscReal parameterB;
37 
39  constexpr inline static PetscReal ugc = 8.31446261815324; // J/K-mol
40 
46  [[nodiscard]] PetscReal ComputeMassFluxRate(PetscReal temperature) const;
47 
48  public:
49  ArrheniusSublimation(const std::shared_ptr<ablate::parameters::Parameters> &properties, const std::shared_ptr<ablate::mathFunctions::MathFunction> &initialization,
50  const std::shared_ptr<ablate::parameters::Parameters> &options = {});
51 
56  void Initialize(ablate::boundarySolver::BoundarySolver &bSolver) override;
57 
63  bool RequiresUpdate() override { return true; };
64 
70  PetscErrorCode Update(PetscInt faceId, PetscReal dt, PetscReal heatFluxToSurface, PetscReal &temperature) override;
71 
76  PetscErrorCode Compute(PetscInt faceId, PetscReal heatFluxToSurface, SurfaceState &) override;
77 
82  [[nodiscard]] SerializerType Serialize() const override { return SerializerType::serial; }
83 
90  PetscErrorCode Save(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override;
91 
98  PetscErrorCode Restore(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override;
99 };
100 
101 } // namespace ablate::boundarySolver::physics::subModels
102 #endif // ABLATELIBRARY_SUBLIMATIONTEMPERATURE_HPP
Definition: boundarySolver.hpp:13
PetscErrorCode Restore(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override
Definition: arrheniusSublimation.cpp:92
bool RequiresUpdate() override
Definition: arrheniusSublimation.hpp:63
SerializerType Serialize() const override
Definition: arrheniusSublimation.hpp:82
void Initialize(ablate::boundarySolver::BoundarySolver &bSolver) override
Definition: arrheniusSublimation.cpp:17
PetscErrorCode Save(PetscViewer viewer, PetscInt sequenceNumber, PetscReal time) override
Definition: arrheniusSublimation.cpp:84
PetscErrorCode Compute(PetscInt faceId, PetscReal heatFluxToSurface, SurfaceState &) override
Definition: arrheniusSublimation.cpp:72
PetscErrorCode Update(PetscInt faceId, PetscReal dt, PetscReal heatFluxToSurface, PetscReal &temperature) override
Definition: arrheniusSublimation.cpp:55
SerializerType
Definition: serializable.hpp:18