Files
SERiF/src/eos/public/EOS.h

267 lines
13 KiB
C++

#pragma once
#include "EOSio.h"
#include "helm.h"
#include <string>
#include "composition.h"
#include <iomanip>
namespace serif::eos {
/**
* @brief Input parameters for an EOS calculation.
*
* This struct holds the necessary physical conditions (composition, density, temperature)
* required to query the Equation of State.
*/
struct EOSInput {
serif::composition::Composition composition; ///< The composition of the system.
double density; ///< The density of the system in cgs (g/cm^3).
double temperature; ///< The temperature of the system in cgs (K).
friend std::ostream& operator<<(std::ostream& os, const EOSInput& input) {
os << "<EOSInput: "
<< "Density: " << input.density << " g/cm^3, "
<< "Temperature: " << input.temperature << " K, "
<< "Composition: " << input.composition << ">";
return os;
}
};
/**
* @brief Represents a thermodynamic parameter and its derivatives.
*
* This struct stores a specific thermodynamic quantity (e.g., pressure, energy, entropy),
* its breakdown into gas and radiation components, and its partial derivatives
* with respect to density, temperature, mean atomic mass number, and mean atomic number.
* All values are in cgs units unless otherwise specified.
*/
struct EOSParameter {
explicit EOSParameter(std::string name_)
: total(), gas(), radiation(),
dDensity(), dTemperature(),
dMeanAtomicMassNumber(),
dMeanAtomicNumber(),
name(std::move(name_)) {}
double total; ///< Total value of the parameter (gas + radiation) (cgs).
double gas; ///< Gas contribution to the parameter (cgs).
double radiation; ///< Radiation contribution to the parameter (cgs).
double dDensity; ///< Derivative of the total parameter with respect to density (cgs units / (g/cm^3)).
double dTemperature; ///< Derivative of the total parameter with respect to temperature (cgs units / K).
double dMeanAtomicMassNumber; ///< Derivative of the total parameter with respect to mean atomic mass number (Abar) (cgs units / (g/mol)).
double dMeanAtomicNumber; ///< Derivative of the total parameter with respect to mean atomic number (Zbar) (cgs units / dimensionless).
std::string name; ///< Name of the parameter (e.g., "Pressure", "Energy", "Entropy").
friend std::ostream& operator<<(std::ostream& os, const EOSParameter& param) {
os << std::setprecision(5) << "<EOSParameter (" << param.name << "): " << param.total << " (gas: " << param.gas
<< ", radiation: " << param.radiation << ") "
<< "d/dRho: " << param.dDensity << ", d/dT: " << param.dTemperature
<< ", d/dAbar: " << param.dMeanAtomicMassNumber
<< ", d/dZbar: " << param.dMeanAtomicNumber << ">";
return os;
}
};
/**
* @brief Output from an EOS calculation.
*
* This struct contains various thermodynamic quantities and their derivatives
* calculated by the EOS for a given set of input conditions.
* It includes fundamental properties like electron fraction and chemical potential,
* as well as detailed breakdowns of pressure, energy, and entropy.
* Additionally, it provides methods to calculate derived quantities like
* susceptibilities, sound speed, adiabatic gradients, and specific heats.
*/
struct EOSOutput {
double electronFraction{}; ///< Electron fraction (ye), dimensionless.
double electronChemicalPotential{}; ///< Electron chemical potential (eta_e) in cgs (erg/g).
double neutronExcessFraction{}; ///< Neutron excess fraction (xnefer), dimensionless.
EOSParameter pressure {"pressure"}; ///< Pressure output data, including total, gas, radiation, and derivatives.
EOSParameter energy {"energy"}; ///< Internal energy output data, including total, gas, radiation, and derivatives.
EOSParameter entropy {"entropy"}; ///< Entropy output data, including total, gas, radiation, and derivatives.
/**
* @brief Calculates the temperature susceptibility (chi_T).
* @return Temperature susceptibility, dimensionless.
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, chi_T = (d ln P / d ln T)_rho.
*/
double chiTemperature();
/**
* @brief Calculates the density susceptibility (chi_rho).
* @return Density susceptibility, dimensionless.
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, chi_rho = (d ln P / d ln rho)_T.
*/
double chiRho();
/**
* @brief Calculates the adiabatic sound speed.
* @return Sound speed in cgs (cm/s).
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, c_s^2 = gamma1 * P / rho.
*/
double soundSpeed();
/**
* @brief Calculates the adiabatic temperature gradient (nabla_ad).
* @return Adiabatic gradient, dimensionless.
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, nabla_ad = (P * chi_T) / (rho * T * c_p * chi_rho).
*/
double adiabaticGradient();
/**
* @brief Calculates the first adiabatic index (Gamma1).
* @return First adiabatic index, dimensionless.
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, Gamma1 = (d ln P / d ln rho)_S.
*/
double gamma1();
/**
* @brief Calculates the second adiabatic index (Gamma2).
* @return Second adiabatic index, dimensionless.
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, Gamma2 / (Gamma2 - 1) = (d ln P / d ln T)_S.
*/
double gamma2();
/**
* @brief Calculates the third adiabatic index (Gamma3).
* @return Third adiabatic index, dimensionless.
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, Gamma3 - 1 = (d ln T / d ln rho)_S.
*/
double gamma3();
/**
* @brief Calculates the specific heat capacity at constant volume (c_v).
* @return Specific heat capacity at constant volume in cgs (erg/K/g).
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, c_v = (dE/dT)_rho.
*/
double specificHeatCapacityAtConstantVolume();
/**
* @brief Calculates the specific heat capacity at constant pressure (c_p).
* @return Specific heat capacity at constant pressure in cgs (erg/K/g).
* @note Placeholder: Actual calculation needs to be implemented based on available EOS derivatives.
* Typically, c_p = c_v + (T / rho^2) * ( (d P / d T)_rho^2 / (d P / d rho)_T ).
*/
double specificHeatCapacityAtConstantPressure();
/**
* @brief Returns the format of the EOS data used to generate this output.
* @return The EOSFormat enum value (currently only EOSFormat::HELM).
*/
EOSFormat EOSFormat() const;
friend std::ostream& operator<<(std::ostream& os, const EOSOutput& output) {
os << "EOSOutput:\n"
<< "\tElectron Fraction: " << output.electronFraction << "\n"
<< "\tElectron Chemical Potential: " << output.electronChemicalPotential << "\n"
<< "\tNeutron Excess Fraction: " << output.neutronExcessFraction << "\n\t"
<< output.pressure << "\n\t"
<< output.energy << "\n\t"
<< output.entropy;
return os;
}
};
/**
* @class EOS
* @brief Main class for accessing Equation of State data.
*
* This class provides an interface to an underlying EOS table (e.g., Helmholtz EOS).
* It handles loading the EOS data and provides a method to retrieve thermodynamic
* properties for given physical conditions.
*
* @example
* @code
* #include "EOS.h"
* #include "composition.h" // For serif::composition::Composition
* #include <iostream>
*
* int main() {
* try {
* // Initialize EOS from a Helmholtz table file
* serif::eos::EOS helmEOS("path/to/helm_table.dat", serif::eos::EOSFormat::HELM);
*
* // Define input conditions
* serif::eos::EOSInput input;
* input.density = 1.0e6; // g/cm^3
* input.temperature = 1.0e7; // K
* // Assuming a simple composition (e.g., pure Helium-4 for demonstration)
* // In a real scenario, initialize Composition properly.
* // For example, if Composition has a constructor like:
* // Composition(const std::map<std::pair<int, int>, double>& mass_fractions);
* // std::map<std::pair<int, int>, double> he4_mass_fraction = {{{2, 4}, 1.0}};
* // input.composition = serif::composition::Composition(he4_mass_fraction);
* // For now, let's assume Composition can be default constructed or set up simply:
* input.composition.addSpecies(2, 4, 1.0); // Z=2, A=4 (He-4), mass fraction 1.0
*
* // Get EOS output
* serif::eos::EOSOutput output = helmEOS.get(input);
*
* // Access results
* std::cout << "Pressure (total): " << output.pressure.total << " dyne/cm^2" << std::endl;
* std::cout << "Energy (total): " << output.energy.total << " erg/g" << std::endl;
* std::cout << "Entropy (total): " << output.entropy.total << " erg/K/g" << std::endl;
* std::cout << "Electron fraction: " << output.electronFraction << std::endl;
*
* // Example of accessing derivatives
* std::cout << "dP/dRho: " << output.pressure.dDensity << std::endl;
* std::cout << "dE/dT: " << output.energy.dTemperature << std::endl;
*
* } catch (const std::exception& e) {
* std::cerr << "An error occurred: " << e.what() << std::endl;
* return 1;
* }
* return 0;
* }
* @endcode
*/
class EOS {
public:
/**
* @brief Constructs an EOS object by loading data from a file.
* @param filename The path to the EOS data file.
* @param format The format of the EOS data file (e.g., EOSFormat::HELM).
* @throw std::runtime_error If the file cannot be opened or read, or if the format is unsupported.
*/
explicit EOS(const std::string& filename, EOSFormat format=EOSFormat::HELM);
/**
* @brief Constructs an EOS object from an existing EOSio reader.
* @param reader An EOSio object that has already loaded the EOS data.
*/
explicit EOS(const EOSio& reader);
/**
* @brief Default destructor.
*/
~EOS() = default;
/**
* @brief Retrieves thermodynamic properties for the given input conditions.
* @param in An EOSInput struct containing the density, temperature, and composition.
* @return An EOSOutput struct containing the calculated thermodynamic properties.
* @throw std::runtime_error If the underlying EOS calculation fails (e.g., out of table bounds for Helmholtz).
*
* This method queries the loaded EOS table (e.g., Helmholtz) using the provided
* density, temperature, and composition (mean atomic mass Abar, mean atomic number Zbar).
* It populates and returns an EOSOutput struct with various thermodynamic quantities
* such as pressure, energy, entropy, their derivatives, electron fraction, etc.
*/
[[nodiscard]] EOSOutput get(const EOSInput& in);
/**
* @brief Gets the format of the loaded EOS data.
* @return The EOSFormat enum value.
*/
[[nodiscard]] EOSFormat getFormat() const;
/**
* @brief Gets a constant reference to the internal EOSio reader.
* @return A const reference to the EOSio object.
*/
[[nodiscard]] const EOSio& getReader() const;
private:
EOSio m_reader; ///< The EOS I/O handler responsible for reading and storing EOS table data.
};
}