1 #ifndef ABLATELIBRARY_SPECIESTRANSPORT_HPP
2 #define ABLATELIBRARY_SPECIESTRANSPORT_HPP
4 #include "eos/transport/transportModel.hpp"
5 #include "finiteVolume/fluxCalculator/fluxCalculator.hpp"
6 #include "flowProcess.hpp"
8 namespace ablate::finiteVolume::processes {
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;
17 struct AdvectionData {
19 PetscInt numberSpecies;
28 ablate::finiteVolume::fluxCalculator::FluxCalculatorFunction fluxCalculatorFunction;
29 void* fluxCalculatorCtx;
31 AdvectionData advectionData;
33 struct DiffusionData {
38 PetscInt numberSpecies;
44 std::vector<PetscReal> speciesSpeciesSensibleEnthalpy;
46 std::vector<PetscReal> speciesDiffusionCoefficient;
48 DiffusionData diffusionData;
51 struct DiffusionTimeStepData {
53 PetscInt numberSpecies;
56 PetscReal stabilityFactor;
61 std::vector<PetscReal> speciesDiffusionCoefficient;
63 DiffusionTimeStepData diffusionTimeStepData;
66 PetscInt numberSpecies;
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 = {});
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);
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);
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);
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);
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);
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);
Definition: finiteVolumeSolver.hpp:28
Definition: flowProcess.hpp:8
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