ABLATE Source Documentation  0.12.33
twoPhaseEulerAdvection.hpp
1 #ifndef ABLATELIBRARY_TWOPHASEEULERADVECTION_HPP
2 #define ABLATELIBRARY_TWOPHASEEULERADVECTION_HPP
3 
4 #include <petsc.h>
5 #include "eos/perfectGas.hpp"
6 #include "eos/stiffenedGas.hpp"
7 #include "eos/twoPhase.hpp"
8 #include "finiteVolume/compressibleFlowFields.hpp"
9 #include "finiteVolume/fluxCalculator/fluxCalculator.hpp"
10 #include "process.hpp"
11 
12 namespace ablate::finiteVolume::processes {
13 
15  public:
16  inline const static std::string VOLUME_FRACTION_FIELD = eos::TwoPhase::VF;
17  inline const static std::string DENSITY_VF_FIELD = ablate::finiteVolume::CompressibleFlowFields::CONSERVED + VOLUME_FRACTION_FIELD;
18 
19  struct TimeStepData {
20  PetscReal cfl;
21  eos::ThermodynamicFunction computeSpeedOfSound;
22  };
23  TimeStepData timeStepData;
24 
25  private:
26  struct DecodeDataStructGas {
27  PetscReal etot;
28  PetscReal rhotot;
29  PetscReal Yg;
30  PetscReal Yl;
31  PetscReal gam1;
32  PetscReal gam2;
33  PetscReal cvg;
34  PetscReal cpl;
35  PetscReal p0l;
36  };
37  struct DecodeDataStructStiff {
38  PetscReal etot;
39  PetscReal rhotot;
40  PetscReal Yg;
41  PetscReal Yl;
42  PetscReal gam1;
43  PetscReal gam2;
44  PetscReal cpg;
45  PetscReal cpl;
46  PetscReal p0g;
47  PetscReal p0l;
48  };
49  static PetscErrorCode FormFunctionGas(SNES snes, Vec x, Vec F, void *ctx);
50  static PetscErrorCode FormJacobianGas(SNES snes, Vec x, Mat J, Mat P, void *ctx);
51  static PetscErrorCode FormFunctionStiff(SNES snes, Vec x, Vec F, void *ctx);
52  static PetscErrorCode FormJacobianStiff(SNES snes, Vec x, Mat J, Mat P, void *ctx);
53 
54  PetscErrorCode MultiphaseFlowPreStage(TS flowTs, ablate::solver::Solver &flow, PetscReal stagetime);
58  class TwoPhaseDecoder {
59  public:
60  virtual void DecodeTwoPhaseEulerState(PetscInt dim, const PetscInt *uOff, const PetscReal *conservedValues, const PetscReal *normal, PetscReal *density, PetscReal *densityG,
61  PetscReal *densityL, PetscReal *normalVelocity, PetscReal *velocity, PetscReal *internalEnergy, PetscReal *internalEnergyG, PetscReal *internalEnergyL,
62  PetscReal *aG, PetscReal *aL, PetscReal *MG, PetscReal *ML, PetscReal *p, PetscReal *T, PetscReal *alpha) = 0;
63  virtual ~TwoPhaseDecoder() = default;
64  };
65 
69  class PerfectGasPerfectGasDecoder : public TwoPhaseDecoder {
70  const std::shared_ptr<eos::PerfectGas> eosGas;
71  const std::shared_ptr<eos::PerfectGas> eosLiquid;
72 
76  std::vector<PetscReal> gasEulerFieldScratch;
77  std::vector<PetscReal> liquidEulerFieldScratch;
78 
82  eos::ThermodynamicFunction gasComputeTemperature;
83  eos::ThermodynamicTemperatureFunction gasComputeInternalEnergy;
84  eos::ThermodynamicTemperatureFunction gasComputeSpeedOfSound;
85  eos::ThermodynamicTemperatureFunction gasComputePressure;
86 
87  eos::ThermodynamicFunction liquidComputeTemperature;
88  eos::ThermodynamicTemperatureFunction liquidComputeInternalEnergy;
89  eos::ThermodynamicTemperatureFunction liquidComputeSpeedOfSound;
90  eos::ThermodynamicTemperatureFunction liquidComputePressure;
91 
92  public:
93  PerfectGasPerfectGasDecoder(PetscInt dim, const std::shared_ptr<eos::PerfectGas> &perfectGasEos1, const std::shared_ptr<eos::PerfectGas> &perfectGasEos2);
94  void DecodeTwoPhaseEulerState(PetscInt dim, const PetscInt *uOff, const PetscReal *conservedValues, const PetscReal *normal, PetscReal *density, PetscReal *densityG, PetscReal *densityL,
95  PetscReal *normalVelocity, PetscReal *velocity, PetscReal *internalEnergy, PetscReal *internalEnergyG, PetscReal *internalEnergyL, PetscReal *aG, PetscReal *aL,
96  PetscReal *MG, PetscReal *ML, PetscReal *p, PetscReal *T, PetscReal *alpha) override;
97  };
98 
102  class PerfectGasStiffenedGasDecoder : public TwoPhaseDecoder {
103  const std::shared_ptr<eos::PerfectGas> eosGas;
104  const std::shared_ptr<eos::StiffenedGas> eosLiquid;
105 
109  std::vector<PetscReal> gasEulerFieldScratch;
110  std::vector<PetscReal> liquidEulerFieldScratch;
111 
115  eos::ThermodynamicFunction gasComputeTemperature;
116  eos::ThermodynamicTemperatureFunction gasComputeInternalEnergy;
117  eos::ThermodynamicTemperatureFunction gasComputeSpeedOfSound;
118  eos::ThermodynamicTemperatureFunction gasComputePressure;
119 
120  eos::ThermodynamicFunction liquidComputeTemperature;
121  eos::ThermodynamicTemperatureFunction liquidComputeInternalEnergy;
122  eos::ThermodynamicTemperatureFunction liquidComputeSpeedOfSound;
123  eos::ThermodynamicTemperatureFunction liquidComputePressure;
124 
125  public:
126  PerfectGasStiffenedGasDecoder(PetscInt dim, const std::shared_ptr<eos::PerfectGas> &perfectGasEos1, const std::shared_ptr<eos::StiffenedGas> &perfectGasEos2);
127 
128  void DecodeTwoPhaseEulerState(PetscInt dim, const PetscInt *uOff, const PetscReal *conservedValues, const PetscReal *normal, PetscReal *density, PetscReal *densityG, PetscReal *densityL,
129  PetscReal *normalVelocity, PetscReal *velocity, PetscReal *internalEnergy, PetscReal *internalEnergyG, PetscReal *internalEnergyL, PetscReal *aG, PetscReal *aL,
130  PetscReal *MG, PetscReal *ML, PetscReal *p, PetscReal *T, PetscReal *alpha) override;
131  };
132 
136  class StiffenedGasStiffenedGasDecoder : public TwoPhaseDecoder {
137  const std::shared_ptr<eos::StiffenedGas> eosGas;
138  const std::shared_ptr<eos::StiffenedGas> eosLiquid;
139 
143  std::vector<PetscReal> gasEulerFieldScratch;
144  std::vector<PetscReal> liquidEulerFieldScratch;
145 
149  eos::ThermodynamicFunction gasComputeTemperature;
150  eos::ThermodynamicTemperatureFunction gasComputeInternalEnergy;
151  eos::ThermodynamicTemperatureFunction gasComputeSpeedOfSound;
152  eos::ThermodynamicTemperatureFunction gasComputePressure;
153 
154  eos::ThermodynamicFunction liquidComputeTemperature;
155  eos::ThermodynamicTemperatureFunction liquidComputeInternalEnergy;
156  eos::ThermodynamicTemperatureFunction liquidComputeSpeedOfSound;
157  eos::ThermodynamicTemperatureFunction liquidComputePressure;
158 
159  public:
160  StiffenedGasStiffenedGasDecoder(PetscInt dim, const std::shared_ptr<eos::StiffenedGas> &perfectGasEos1, const std::shared_ptr<eos::StiffenedGas> &perfectGasEos2);
161 
162  void DecodeTwoPhaseEulerState(PetscInt dim, const PetscInt *uOff, const PetscReal *conservedValues, const PetscReal *normal, PetscReal *density, PetscReal *densityG, PetscReal *densityL,
163  PetscReal *normalVelocity, PetscReal *velocity, PetscReal *internalEnergy, PetscReal *internalEnergyG, PetscReal *internalEnergyL, PetscReal *aG, PetscReal *aL,
164  PetscReal *MG, PetscReal *ML, PetscReal *p, PetscReal *T, PetscReal *alpha) override;
165  };
166 
167  const std::shared_ptr<eos::EOS> eosTwoPhase;
168  std::shared_ptr<eos::EOS> eosGas;
169  std::shared_ptr<eos::EOS> eosLiquid;
170  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorGasGas;
171  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorGasLiquid;
172  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorLiquidGas;
173  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorLiquidLiquid;
174 
178  std::shared_ptr<TwoPhaseDecoder> decoder;
179 
180  public:
181  static PetscErrorCode UpdateAuxTemperatureField2Gas(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[],
182  PetscScalar *auxField, void *ctx);
183 
184  static PetscErrorCode UpdateAuxPressureField2Gas(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[],
185  PetscScalar *auxField, void *ctx);
186 
187  static PetscErrorCode UpdateAuxVelocityField2Gas(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[],
188  PetscScalar *auxField, void *ctx);
189 
190  TwoPhaseEulerAdvection(std::shared_ptr<eos::EOS> eosTwoPhase, const std::shared_ptr<parameters::Parameters> &parameters, std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorGasGas,
191  std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorGasLiquid, std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorLiquidGas,
192  std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculatorLiquidLiquid);
193  void Setup(ablate::finiteVolume::FiniteVolumeSolver &flow) override;
194 
195  private:
196  // static function to compute time step for twoPhase euler advection
197  static double ComputeCflTimeStep(TS ts, ablate::finiteVolume::FiniteVolumeSolver &flow, void *ctx);
198 
199  static PetscErrorCode CompressibleFlowComputeEulerFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[],
200  const PetscInt aOff[], const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar *flux, void *ctx);
201  static PetscErrorCode CompressibleFlowComputeVFFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[], const PetscInt aOff[],
202  const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar *flux, void *ctx);
203 
204  public:
212  static std::shared_ptr<TwoPhaseDecoder> CreateTwoPhaseDecoder(PetscInt dim, const std::shared_ptr<eos::EOS> &eosGas, const std::shared_ptr<eos::EOS> &eosLiquid);
213 };
214 
215 } // namespace ablate::finiteVolume::processes
216 #endif // ABLATELIBRARY_TWOPHASEEULERADVECTION_HPP
static const std::string CONSERVED
The conserved prefix used for fields that have a conserved and non conserved form.
Definition: compressibleFlowFields.hpp:22
Definition: finiteVolumeSolver.hpp:28
Definition: twoPhaseEulerAdvection.hpp:14
void Setup(ablate::finiteVolume::FiniteVolumeSolver &flow) override
Definition: twoPhaseEulerAdvection.cpp:226
static std::shared_ptr< TwoPhaseDecoder > CreateTwoPhaseDecoder(PetscInt dim, const std::shared_ptr< eos::EOS > &eosGas, const std::shared_ptr< eos::EOS > &eosLiquid)
Definition: twoPhaseEulerAdvection.cpp:800
Definition: solver.hpp:17