ABLATE Source Documentation  0.12.33
navierStokesTransport.hpp
1 #ifndef ABLATELIBRARY_NAVIERSTOKESTRANSPORT_HPP
2 #define ABLATELIBRARY_NAVIERSTOKESTRANSPORT_HPP
3 
4 #include <petsc.h>
5 #include "eos/transport/transportModel.hpp"
6 #include "finiteVolume/fluxCalculator/fluxCalculator.hpp"
7 #include "flowProcess.hpp"
8 #include "pressureGradientScaling.hpp"
9 
10 namespace ablate::finiteVolume::processes {
11 
13  public:
14  // Store ctx needed for static function advection function passed to PETSc
15  struct AdvectionData {
16  // flow CFL
17  PetscReal cfl;
18 
19  /* number of gas species */
20  PetscInt numberSpecies;
21 
22  // EOS function calls
23  eos::ThermodynamicFunction computeTemperature;
24  eos::ThermodynamicTemperatureFunction computeInternalEnergy;
25  eos::ThermodynamicTemperatureFunction computeSpeedOfSound;
27 
28  /* store method used for flux calculator */
29  ablate::finiteVolume::fluxCalculator::FluxCalculatorFunction fluxCalculatorFunction;
30  void* fluxCalculatorCtx;
31  };
32 
33  // Store ctx needed for static function diffusion function passed to PETSc
34  struct DiffusionData {
35  /* thermal conductivity*/
37  /* dynamic viscosity*/
39  };
40 
41  private:
42  const std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalculator;
43  const std::shared_ptr<eos::EOS> eos;
44  const std::shared_ptr<eos::transport::TransportModel> transportModel;
45  AdvectionData advectionData;
46 
47  eos::ThermodynamicTemperatureFunction computeTemperatureFunction;
48 
49  DiffusionData diffusionData;
50 
51  eos::ThermodynamicFunction computePressureFunction;
52 
53  // Store the required ctx for time stepping
54  struct CflTimeStepData {
55  /* thermal conductivity*/
56  AdvectionData* advectionData;
57 
61  std::shared_ptr<ablate::finiteVolume::processes::PressureGradientScaling> pgs;
62  };
63  CflTimeStepData timeStepData;
64 
66  struct DiffusionTimeStepData {
68  PetscReal conductionStabilityFactor;
69 
71  PetscReal viscousStabilityFactor;
72 
73  /* thermal conductivity*/
75  /* dynamic viscosity*/
77  /* specific heat*/
79  /* density */
81  };
82  DiffusionTimeStepData diffusionTimeStepData;
83 
84  // static function to compute time step for euler advection
85  static double ComputeCflTimeStep(TS ts, ablate::finiteVolume::FiniteVolumeSolver& flow, void* ctx);
86 
87  // static function to compute the conduction based time step
88  static double ComputeConductionTimeStep(TS ts, ablate::finiteVolume::FiniteVolumeSolver& flow, void* ctx);
89 
90  // static function to compute the conduction based time step
91  static double ComputeViscousDiffusionTimeStep(TS ts, ablate::finiteVolume::FiniteVolumeSolver& flow, void* ctx);
92 
93  public:
97  static PetscErrorCode UpdateAuxTemperatureField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt aOff[],
98  PetscScalar* auxField, void* ctx);
102  static PetscErrorCode UpdateAuxVelocityField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt aOff[],
103  PetscScalar* auxField, void* ctx);
104 
108  static PetscErrorCode UpdateAuxPressureField(PetscReal time, PetscInt dim, const PetscFVCellGeom* cellGeom, const PetscInt uOff[], const PetscScalar* conservedValues, const PetscInt aOff[],
109  PetscScalar* auxField, void* ctx);
110 
115  NavierStokesTransport(const std::shared_ptr<parameters::Parameters>& parameters, std::shared_ptr<eos::EOS> eos, std::shared_ptr<fluxCalculator::FluxCalculator> fluxCalcIn = {},
116  std::shared_ptr<eos::transport::TransportModel> transportModel = {}, std::shared_ptr<ablate::finiteVolume::processes::PressureGradientScaling> = {});
117 
122  void Setup(ablate::finiteVolume::FiniteVolumeSolver& flow) override;
123 
131  static PetscErrorCode AdvectionFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[], const PetscInt aOff[],
132  const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar* flux, void* ctx);
133 
141  static PetscErrorCode DiffusionFlux(PetscInt dim, const PetscFVFaceGeom* fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[],
142  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void* ctx);
143 
153  static PetscErrorCode CompressibleFlowComputeStressTensor(PetscInt dim, PetscReal mu, const PetscReal* gradVel, PetscReal* tau);
154 };
155 
156 } // namespace ablate::finiteVolume::processes
157 #endif // ABLATELIBRARY_NAVIERSTOKESTRANSPORT_HPP
Definition: finiteVolumeSolver.hpp:28
Definition: navierStokesTransport.hpp:12
NavierStokesTransport(const std::shared_ptr< parameters::Parameters > &parameters, std::shared_ptr< eos::EOS > eos, std::shared_ptr< fluxCalculator::FluxCalculator > fluxCalcIn={}, std::shared_ptr< eos::transport::TransportModel > transportModel={}, std::shared_ptr< ablate::finiteVolume::processes::PressureGradientScaling >={})
Definition: navierStokesTransport.cpp:10
static PetscErrorCode UpdateAuxVelocityField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[], PetscScalar *auxField, void *ctx)
Definition: navierStokesTransport.cpp:528
static PetscErrorCode UpdateAuxPressureField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[], PetscScalar *auxField, void *ctx)
Definition: navierStokesTransport.cpp:551
static PetscErrorCode UpdateAuxTemperatureField(PetscReal time, PetscInt dim, const PetscFVCellGeom *cellGeom, const PetscInt uOff[], const PetscScalar *conservedValues, const PetscInt aOff[], PetscScalar *auxField, void *ctx)
Definition: navierStokesTransport.cpp:541
void Setup(ablate::finiteVolume::FiniteVolumeSolver &flow) override
Definition: navierStokesTransport.cpp:35
static PetscErrorCode AdvectionFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscScalar fieldL[], const PetscScalar fieldR[], const PetscInt aOff[], const PetscScalar auxL[], const PetscScalar auxR[], PetscScalar *flux, void *ctx)
Definition: navierStokesTransport.cpp:97
static PetscErrorCode CompressibleFlowComputeStressTensor(PetscInt dim, PetscReal mu, const PetscReal *gradVel, PetscReal *tau)
Definition: navierStokesTransport.cpp:504
static PetscErrorCode DiffusionFlux(PetscInt dim, const PetscFVFaceGeom *fg, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar field[], const PetscScalar grad[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar aux[], const PetscScalar gradAux[], PetscScalar flux[], void *ctx)
Definition: navierStokesTransport.cpp:448