1 #ifndef ABLATELIBRARY_ZERORKEOS_HPP
2 #define ABLATELIBRARY_ZERORKEOS_HPP
7 #include "chemistryModel.hpp"
9 #include "eos/zerork/sourceCalculatorZeroRK.hpp"
10 #include "monitors/logs/log.hpp"
11 #include "parameters/parameters.hpp"
12 #include "utilities/intErrorChecker.hpp"
13 #include "utilities/vectorUtilities.hpp"
14 #include "zerork/mechanism.h"
15 #include "zerork/utilities.h"
16 #include "zerork_cfd_plugin.h"
18 namespace ablate::eos {
33 const char* cklogfilename =
"mech.cklog";
35 std::vector<double> stateVector;
44 explicit zerorkEOS(std::filesystem::path reactionFile, std::filesystem::path thermoFile,
const std::shared_ptr<ablate::parameters::Parameters>& options = {});
46 std::shared_ptr<zerork::mechanism> mech;
48 const std::filesystem::path reactionFile;
50 const std::filesystem::path thermoFile;
73 void View(std::ostream& stream)
const override;
86 [[nodiscard]]
virtual const std::vector<std::string>&
GetProgressVariables()
const override {
return ablate::utilities::VectorUtilities::Empty<std::string>; }
105 PetscInt eulerOffset;
106 PetscInt densityYiOffset;
108 std::shared_ptr<zerork::mechanism> mech;
114 inline const static double TREF = 298.15;
118 PetscErrorCode (*
function)(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx) =
nullptr;
158 [[nodiscard]] EOSFunction
GetFieldFunctionFunction(
const std::string& field, ThermodynamicProperty property1, ThermodynamicProperty property2,
159 std::vector<std::string> otherProperties)
const override;
178 [[nodiscard]] std::shared_ptr<FunctionContext>
BuildFunctionContext(ablate::eos::ThermodynamicProperty property,
const std::vector<domain::Field>& fields,
bool checkDensityYi =
true)
const;
188 static PetscErrorCode DensityFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
189 static PetscErrorCode TemperatureFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
190 static PetscErrorCode PressureFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
191 static PetscErrorCode InternalSensibleEnergyFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
192 static PetscErrorCode SensibleEnthalpyFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
193 static PetscErrorCode SpecificHeatConstantVolumeFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
194 static PetscErrorCode SpecificHeatConstantPressureFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
195 static PetscErrorCode SpeedOfSoundFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
196 static PetscErrorCode SpeciesSensibleEnthalpyFunction(
const PetscReal conserved[], PetscReal* property,
void* ctx);
208 static PetscErrorCode DensityTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
209 static PetscErrorCode TemperatureTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
210 static PetscErrorCode PressureTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
211 static PetscErrorCode InternalSensibleEnergyTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
212 static PetscErrorCode SensibleEnthalpyTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
213 static PetscErrorCode SpecificHeatConstantVolumeTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
214 static PetscErrorCode SpecificHeatConstantPressureTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
215 static PetscErrorCode SpeedOfSoundTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
216 static PetscErrorCode SpeciesSensibleEnthalpyTemperatureFunction(
const PetscReal conserved[], PetscReal T, PetscReal* property,
void* ctx);
227 static PetscErrorCode DensityMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
228 static PetscErrorCode TemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
229 static PetscErrorCode PressureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
230 static PetscErrorCode InternalSensibleEnergyMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
231 static PetscErrorCode SensibleEnthalpyMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
232 static PetscErrorCode SpecificHeatConstantVolumeMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
233 static PetscErrorCode SpecificHeatConstantPressureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
234 static PetscErrorCode SpeedOfSoundMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
235 static PetscErrorCode SpeciesSensibleEnthalpyMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
247 static PetscErrorCode DensityTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
248 static PetscErrorCode TemperatureTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
249 static PetscErrorCode PressureTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
250 static PetscErrorCode InternalSensibleEnergyTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
251 static PetscErrorCode SensibleEnthalpyTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
252 static PetscErrorCode SpecificHeatConstantVolumeTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
253 static PetscErrorCode SpecificHeatConstantPressureTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
254 static PetscErrorCode SpeedOfSoundTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
255 static PetscErrorCode SpeciesSensibleEnthalpyTemperatureMassFractionFunction(
const PetscReal conserved[],
const PetscReal yi[], PetscReal T, PetscReal* property,
void* ctx);
262 using ThermodynamicTemperatureStaticFunction = PetscErrorCode (*)(
const PetscReal conserved[], PetscReal temperature, PetscReal* property,
void* ctx);
263 std::map<ThermodynamicProperty, std::tuple<ThermodynamicStaticFunction, ThermodynamicTemperatureStaticFunction>> thermodynamicFunctions = {
264 {ThermodynamicProperty::Density, {DensityFunction, DensityTemperatureFunction}},
265 {ThermodynamicProperty::Pressure, {PressureFunction, PressureTemperatureFunction}},
266 {ThermodynamicProperty::Temperature, {TemperatureFunction, TemperatureTemperatureFunction}},
267 {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyFunction, InternalSensibleEnergyTemperatureFunction}},
268 {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyFunction, SensibleEnthalpyTemperatureFunction}},
269 {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeFunction, SpecificHeatConstantVolumeTemperatureFunction}},
270 {ThermodynamicProperty::SpecificHeatConstantPressure, {SpecificHeatConstantPressureFunction, SpecificHeatConstantPressureTemperatureFunction}},
271 {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundFunction, SpeedOfSoundTemperatureFunction}},
272 {ThermodynamicProperty::SpeciesSensibleEnthalpy,
273 {SpeciesSensibleEnthalpyFunction, SpeciesSensibleEnthalpyTemperatureFunction}}
276 using ThermodynamicStaticMassFractionFunction = PetscErrorCode (*)(
const PetscReal conserved[],
const PetscReal yi[], PetscReal* property,
void* ctx);
277 using ThermodynamicTemperatureStaticMassFractionFunction = PetscErrorCode (*)(
const PetscReal conserved[],
const PetscReal yi[], PetscReal temperature, PetscReal* property,
void* ctx);
278 std::map<ThermodynamicProperty, std::tuple<ThermodynamicStaticMassFractionFunction, ThermodynamicTemperatureStaticMassFractionFunction>> thermodynamicMassFractionFunctions = {
279 {ThermodynamicProperty::Density, {DensityMassFractionFunction, DensityTemperatureMassFractionFunction}},
280 {ThermodynamicProperty::Pressure, {PressureMassFractionFunction, PressureTemperatureMassFractionFunction}} ,
281 {ThermodynamicProperty::Temperature, {TemperatureMassFractionFunction, TemperatureTemperatureMassFractionFunction}},
282 {ThermodynamicProperty::InternalSensibleEnergy, {InternalSensibleEnergyMassFractionFunction, InternalSensibleEnergyTemperatureMassFractionFunction}},
283 {ThermodynamicProperty::SensibleEnthalpy, {SensibleEnthalpyMassFractionFunction, SensibleEnthalpyTemperatureMassFractionFunction}},
284 {ThermodynamicProperty::SpecificHeatConstantVolume, {SpecificHeatConstantVolumeMassFractionFunction, SpecificHeatConstantVolumeTemperatureMassFractionFunction}},
285 {ThermodynamicProperty::SpecificHeatConstantPressure,
286 {SpecificHeatConstantPressureMassFractionFunction, SpecificHeatConstantPressureTemperatureMassFractionFunction}} ,
287 {ThermodynamicProperty::SpeedOfSound, {SpeedOfSoundMassFractionFunction, SpeedOfSoundTemperatureMassFractionFunction}},
288 {ThermodynamicProperty::SpeciesSensibleEnthalpy,
289 {SpeciesSensibleEnthalpyMassFractionFunction, SpeciesSensibleEnthalpyTemperatureMassFractionFunction}}
Definition: chemistryModel.hpp:18
Definition: zerork.hpp:20
ThermodynamicTemperatureFunction GetThermodynamicTemperatureFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: zerork.cpp:68
PetscErrorCode(*)(const PetscReal conserved[], PetscReal *property, void *ctx) ThermodynamicStaticFunction
Definition: zerork.hpp:261
std::vector< std::string > species
hold a copy of the constrains that can be used for single or batch source calculation
Definition: zerork.hpp:24
ThermodynamicFunction GetThermodynamicFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: zerork.cpp:62
virtual const std::vector< std::string > & GetProgressVariables() const override
Definition: zerork.hpp:86
EOSFunction GetFieldFunctionFunction(const std::string &field, ThermodynamicProperty property1, ThermodynamicProperty property2, std::vector< std::string > otherProperties) const override
Definition: zerork.cpp:598
ThermodynamicMassFractionFunction GetThermodynamicMassFractionFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const
Definition: zerork.cpp:74
zerorkEOS(std::filesystem::path reactionFile, std::filesystem::path thermoFile, const std::shared_ptr< ablate::parameters::Parameters > &options={})
Definition: zerork.cpp:5
void View(std::ostream &stream) const override
Definition: zerork.cpp:33
static const std::array< std::string, 2 > validThermoFileExtensions
Definition: zerork.hpp:99
std::shared_ptr< FunctionContext > BuildFunctionContext(ablate::eos::ThermodynamicProperty property, const std::vector< domain::Field > &fields, bool checkDensityYi=true) const
Definition: zerork.cpp:38
static void FillreactorMassFracVectorFromMassFractions(int nSpc, const double *Yi, std::vector< double > &reactorYi)
Definition: zerork.cpp:579
std::shared_ptr< SourceCalculator > CreateSourceCalculator(const std::vector< domain::Field > &fields, const ablate::domain::Range &cellRange) override
Definition: zerork.cpp:751
std::map< std::string, std::map< std::string, int > > GetSpeciesElementalInformation() const override
Definition: zerork.cpp:735
std::map< std::string, double > GetSpeciesMolecularMass() const override
Definition: zerork.cpp:737
std::map< std::string, double > GetElementInformation() const override
Definition: zerork.cpp:733
ThermodynamicTemperatureMassFractionFunction GetThermodynamicTemperatureMassFractionFunction(ThermodynamicProperty property, const std::vector< domain::Field > &fields) const override
Definition: zerork.cpp:81
const std::set< ThermodynamicProperty > speciesSizedProperties
Definition: zerork.hpp:295
const std::vector< std::string > & GetSpeciesVariables() const override
Definition: zerork.hpp:80
static void FillreactorMassFracVectorFromDensityMassFractions(int nSpc, double density, const double *densityYi, std::vector< double > &reactorYi)
Definition: zerork.cpp:561
Definition: chemistryModel.hpp:66
Definition: zerork.hpp:102
Definition: zerork.hpp:116
PetscInt propertySize
the property size being set
Definition: zerork.hpp:122
std::shared_ptr< void > context
optional context to pass into the function
Definition: zerork.hpp:120
hold a struct that can be used for chemistry constraints
Definition: sourceCalculatorZeroRK.hpp:26