ABLATE Source Documentation  0.12.33
speciesTransport.hpp
1 #ifndef ABLATELIBRARY_SPECIESTRANSPORT_HPP
2 #define ABLATELIBRARY_SPECIESTRANSPORT_HPP
3 
4 #include "eos/transport/transportModel.hpp"
5 #include "finiteVolume/fluxCalculator/fluxCalculator.hpp"
6 #include "flowProcess.hpp"
7 
8 namespace ablate::finiteVolume::processes {
9 
10 class SpeciesTransport : public FlowProcess {
11  private:
12  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculator;
13  const std::shared_ptr<eos::EOS> eos;
14  const std::shared_ptr<eos::transport::TransportModel> transportModel;
15 
16  // Store ctx needed for static function advection function passed to PETSc
17  struct AdvectionData {
18  /* number of gas species */
19  PetscInt numberSpecies;
20 
21  // EOS function calls
22  eos::ThermodynamicFunction computeTemperature;
23  eos::ThermodynamicTemperatureFunction computeInternalEnergy;
24  eos::ThermodynamicTemperatureFunction computeSpeedOfSound;
26 
27  /* store method used for flux calculator */
28  ablate::finiteVolume::fluxCalculator::FluxCalculatorFunction fluxCalculatorFunction;
29  void* fluxCalculatorCtx;
30  };
31  AdvectionData advectionData;
32 
33  struct DiffusionData {
34  /* diffusivity */
36 
37  /* number of gas species */
38  PetscInt numberSpecies;
39 
40  /* functions to compute species enthalpy */
41  eos::ThermodynamicTemperatureFunction computeSpeciesSensibleEnthalpyFunction;
42 
43  /* store a scratch space for speciesSpeciesSensibleEnthalpy */
44  std::vector<PetscReal> speciesSpeciesSensibleEnthalpy;
45  /* store an optional scratch space for individual species diffusion */
46  std::vector<PetscReal> speciesDiffusionCoefficient;
47  };
48  DiffusionData diffusionData;
49 
51  struct DiffusionTimeStepData {
52  /* number of gas species */
53  PetscInt numberSpecies;
54 
56  PetscReal stabilityFactor;
57 
58  /* diffusivity */
60  /* store an optional scratch space for individual species diffusion */
61  std::vector<PetscReal> speciesDiffusionCoefficient;
62  };
63  DiffusionTimeStepData diffusionTimeStepData;
64 
65  // Store ctx needed for static function diffusion function passed to PETSc
66  PetscInt numberSpecies;
67 
68  public:
69  explicit SpeciesTransport(std::shared_ptr<eos::EOS> eos, std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalcIn = {}, std::shared_ptr<eos::transport::TransportModel> transportModel = {},
70  const std::shared_ptr<parameters::Parameters>& parametersIn = {});
71 
77 
81  static PetscErrorCode UpdateAuxMassFractionField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt aOff[],
82  PetscScalar* auxField, void* ctx);
83 
88  static void NormalizeSpecies(TS ts, ablate::solver::Solver&);
89 
90  private:
99  static PetscErrorCode DiffusionEnergyFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[],
100  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void* ctx);
101 
110  static PetscErrorCode DiffusionEnergyFluxVariableDiffusionCoefficient(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[],
111  const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[],
112  const PetscScalar gradAux[], PetscScalar flux[], void* ctx);
113 
122  static PetscErrorCode DiffusionSpeciesFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[],
123  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void* ctx);
124 
133  static PetscErrorCode DiffusionSpeciesFluxVariableDiffusionCoefficient(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[],
134  const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[],
135  const PetscScalar gradAux[], PetscScalar flux[], void* ctx);
136 
143  static PetscErrorCode AdvectionFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[], const PetscInt aOff[],
144  const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar* flux, void* ctx);
145 
146  // static function to compute the conduction based time step
147  static double ComputeViscousDiffusionTimeStep(TS ts, ablate::finiteVolume::FiniteVolumeSolver& flow, void* ctx);
148 };
149 
150 } // namespace ablate::finiteVolume::processes
151 #endif // ABLATELIBRARY_SPECIESTRANSPORT_HPP
Definition: finiteVolumeSolver.hpp:28
Definition: speciesTransport.hpp:10
static PetscErrorCode UpdateAuxMassFractionField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[], PetscScalar *auxField, void *ctx)
Definition: speciesTransport.cpp:94
static void NormalizeSpecies(TS ts, ablate::solver::Solver &)
Definition: speciesTransport.cpp:339
void Setup(ablate::finiteVolume::FiniteVolumeSolver &flow) override
Definition: speciesTransport.cpp:35
Definition: solver.hpp:17