ABLATE Source Documentation  0.12.34
lodiBoundary.hpp
1 #ifndef ABLATELIBRARY_LODIBOUNDARY_HPP
2 #define ABLATELIBRARY_LODIBOUNDARY_HPP
3 
4 #include "boundarySolver/boundaryProcess.hpp"
5 #include "eos/eos.hpp"
6 #include "finiteVolume/compressibleFlowFields.hpp"
7 #include "finiteVolume/processes/flowProcess.hpp"
8 #include "finiteVolume/processes/navierStokesTransport.hpp"
9 #include "finiteVolume/processes/pressureGradientScaling.hpp"
10 
11 namespace ablate::boundarySolver::lodi {
12 class LODIBoundary : public BoundaryProcess {
13  protected:
14  typedef enum {
15  RHO = finiteVolume::CompressibleFlowFields::RHO,
16  RHOE = finiteVolume::CompressibleFlowFields::RHOE,
17  RHOVELN = finiteVolume::CompressibleFlowFields::RHOU,
18  RHOVELT1 = finiteVolume::CompressibleFlowFields::RHOV,
19  RHOVELT2 = finiteVolume::CompressibleFlowFields::RHOW
20  } BoundaryEulerComponents;
21 
22  // Global parameters
23  const std::shared_ptr<finiteVolume::processes::PressureGradientScaling> pressureGradientScaling;
24 
25  void GetVelAndCPrims(PetscReal velNorm, PetscReal speedOfSound, PetscReal Cp, PetscReal Cv, PetscReal& velNormPrim, PetscReal& speedOfSoundPrim);
26 
27  void GetEigenValues(PetscReal veln, PetscReal c, PetscReal velnprm, PetscReal cprm, PetscReal lamda[]) const;
28 
29  void GetmdFdn(const PetscInt sOff[], const PetscReal* velNormCord, PetscReal rho, PetscReal T, PetscReal Cp, PetscReal Cv, PetscReal C, PetscReal Enth, PetscReal velnprm, PetscReal Cprm,
30  const PetscReal* conserved, const PetscInt uOff[], const PetscReal* sL, const PetscReal transformationMatrix[3][3], PetscReal* mdFdn) const;
31 
32  // Compute known/shared values
33  PetscInt dims, nEqs, nSpecEqs, nEvEqs, eulerId, speciesId;
34 
36  std::vector<PetscInt> evIds;
37 
39  std::vector<PetscInt> nEvComps;
40 
41  // Keep track of the required fields
42  std::vector<std::string> fieldNames;
43 
44  // Store eos decode params
45  const std::shared_ptr<eos::EOS> eos;
46  eos::ThermodynamicFunction computeTemperature;
47  eos::ThermodynamicTemperatureFunction computePressureFromTemperature;
48  eos::ThermodynamicTemperatureFunction computeSpeedOfSound;
49  eos::ThermodynamicTemperatureFunction computeSpecificHeatConstantPressure;
50  eos::ThermodynamicTemperatureFunction computeSpecificHeatConstantVolume;
51  eos::ThermodynamicTemperatureFunction computeSensibleEnthalpyFunction;
52  eos::ThermodynamicFunction computePressure;
53 
54  public:
55  explicit LODIBoundary(std::shared_ptr<eos::EOS> eos, std::shared_ptr<finiteVolume::processes::PressureGradientScaling> pressureGradientScaling = {});
56 
57  void Setup(ablate::boundarySolver::BoundarySolver& bSolver) override;
58 
66  void Setup(PetscInt dims, PetscInt nEqs, PetscInt nSpecEqs = 0, std::vector<PetscInt> nEvEqs = {}, const std::vector<domain::Field>& fields = {});
67 
68  private:
69  eos::ThermodynamicTemperatureFunction computeTemperatureFunction;
70 };
71 
72 } // namespace ablate::boundarySolver::lodi
73 #endif // ABLATELIBRARY_LODIBOUNDARY_HPP
Definition: boundaryProcess.hpp:8
Definition: boundarySolver.hpp:13
Definition: lodiBoundary.hpp:12
std::vector< PetscInt > evIds
There maybe more than one ev field so keep a separate id for each.
Definition: lodiBoundary.hpp:36
std::vector< PetscInt > nEvComps
Track the number of components in each ev field.
Definition: lodiBoundary.hpp:39
void Setup(ablate::boundarySolver::BoundarySolver &bSolver) override
Definition: lodiBoundary.cpp:116