ABLATE Source Documentation  0.12.33
evTransport.hpp
1 #ifndef ABLATELIBRARY_EVTRANSPORT_HPP
2 #define ABLATELIBRARY_EVTRANSPORT_HPP
3 
4 #include <eos/transport/transportModel.hpp>
5 #include "eos/transport/transportModel.hpp"
6 #include "finiteVolume/fluxCalculator/fluxCalculator.hpp"
7 #include "flowProcess.hpp"
8 
9 namespace ablate::finiteVolume::processes {
16 class EVTransport : public FlowProcess {
17  private:
18  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculator;
19  const std::shared_ptr<eos::EOS> eos;
20  const std::shared_ptr<eos::transport::TransportModel> transportModel;
21 
22  // Store ctx needed for static function advection function passed to PETSc
23  struct AdvectionData {
24  /* number of extra species */
25  PetscInt numberEV;
26 
27  // EOS function calls
28  eos::ThermodynamicFunction computeTemperature;
29  eos::ThermodynamicTemperatureFunction computeInternalEnergy;
30  eos::ThermodynamicTemperatureFunction computeSpeedOfSound;
32 
33  /* store method used for flux calculator */
34  ablate::finiteVolume::fluxCalculator::FluxCalculatorFunction fluxCalculatorFunction;
35  void* fluxCalculatorCtx;
36  };
37 
38  struct DiffusionData {
39  /* number of extra species */
40  PetscInt numberEV;
41 
42  /* functions to compute diffusion */
43  eos::ThermodynamicFunction diffFunction;
44 
45  /* store a scratch space for evDiffusionCoefficient */
46  std::vector<PetscReal> evDiffusionCoefficient;
47  };
48 
49  // Store an AdvectionData, diffusionData, and numberEV for each ev field
50  std::vector<AdvectionData> advectionDatas;
51  std::vector<DiffusionData> diffusionDatas;
52  std::vector<PetscInt> numberEVs;
53 
54  public:
55  explicit EVTransport(std::shared_ptr<eos::EOS> eos, std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalcIn = {}, std::shared_ptr<eos::transport::TransportModel> transportModel = {});
56 
61  static void BoundExtraVariables(TS ts, ablate::solver::Solver& solver, const std::string& field);
62 
67  static void BoundExtraVariablesMinusOneToOne(TS ts, ablate::solver::Solver& solver, const std::string& field);
68 
73  static void PositiveExtraVariables(TS ts, ablate::solver::Solver& solver, const std::string& field);
74 
80 
84  static PetscErrorCode UpdateEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt* aOff,
85  PetscScalar* auxField, void* ctx);
86 
90  static PetscErrorCode UpdateBoundEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt* aOff,
91  PetscScalar* auxField, void* ctx);
92 
96  static PetscErrorCode UpdateMinusOneToOneBoundEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues,
97  const PetscInt* aOff, PetscScalar* auxField, void* ctx);
98 
102  static PetscErrorCode UpdatePositiveEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt* aOff,
103  PetscScalar* auxField, void* ctx);
104 
105  private:
114  static PetscErrorCode DiffusionEVFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[],
115  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void* ctx);
116 
125  static PetscErrorCode DiffusionEVFluxVariableDiffusionCoefficient(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[],
126  const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[],
127  PetscScalar flux[], void* ctx);
128 
135  static PetscErrorCode AdvectionFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[], const PetscInt aOff[],
136  const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar* flux, void* ctx);
137 };
138 
139 } // namespace ablate::finiteVolume::processes
140 #endif // ABLATELIBRARY_SPECIESTRANSPORT_HPP
Definition: finiteVolumeSolver.hpp:28
Definition: evTransport.hpp:16
static PetscErrorCode UpdatePositiveEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt *aOff, PetscScalar *auxField, void *ctx)
Definition: evTransport.cpp:125
static void PositiveExtraVariables(TS ts, ablate::solver::Solver &solver, const std::string &field)
Definition: evTransport.cpp:355
static PetscErrorCode UpdateEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt *aOff, PetscScalar *auxField, void *ctx)
Definition: evTransport.cpp:83
static void BoundExtraVariablesMinusOneToOne(TS ts, ablate::solver::Solver &solver, const std::string &field)
Definition: evTransport.cpp:319
static PetscErrorCode UpdateBoundEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt *aOff, PetscScalar *auxField, void *ctx)
Definition: evTransport.cpp:97
void Setup(ablate::finiteVolume::FiniteVolumeSolver &flow) override
Definition: evTransport.cpp:9
static void BoundExtraVariables(TS ts, ablate::solver::Solver &solver, const std::string &field)
Definition: evTransport.cpp:284
static PetscErrorCode UpdateMinusOneToOneBoundEVField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt *aOff, PetscScalar *auxField, void *ctx)
Definition: evTransport.cpp:111
Definition: solver.hpp:17