jesterTOV.eos.spectral.SpectralDecomposition_EOS_model#
- class SpectralDecomposition_EOS_model(
- crust_name='SLy',
- n_points_high=500,
- reparametrized=False,
- sigma_scale=1.0,
Bases:
Interpolate_EOS_modelSpectral decomposition equation of state model.
This class implements the spectral decomposition EOS parametrization following Lindblom (2010) PRD 82, 103011 and Lackey & Wade (2018) PRD 98, 063004. The implementation matches LALSuite’s XLALSimNeutronStarEOS4ParameterSpectralDecomposition exactly.
The adiabatic index is parametrized as:
\[\log \Gamma(x) = \gamma_0 + \gamma_1 x + \gamma_2 x^2 + \gamma_3 x^3\]where \(x = \log(p/p_0)\) is the dimensionless log-pressure.
- Reference values (geometric units):
\(\varepsilon_0 = 9.54629006 \times 10^{-11} \, \mathrm{m}^{-2}\) \(p_0 = 4.43784199 \times 10^{-13} \, \mathrm{m}^{-2}\)
Thermodynamic relations from PRD 82, 103011 (2010):
\[\begin{split}\mu(x) &= \exp\left[-\int_0^x \frac{1}{\Gamma(x')} \, dx'\right] \\ \varepsilon(x) &= \varepsilon_0 \cdot \exp\left[\int_0^x \frac{\mu(x')}{1+\mu(x')} \, dx'\right] \\ p(x) &= p_0 \cdot \exp(x)\end{split}\]Reparametrization (
reparametrized=True):When this flag is set the model expects tilde parameters \(\tilde\gamma = [\tilde\gamma_0, \tilde\gamma_1, \tilde\gamma_2, \tilde\gamma_3]\) drawn from a 4-D unit normal, and maps them to physical spectral coefficients via:
\[\boldsymbol{\gamma} = \boldsymbol{\mu}_\text{radio} + L_\text{wide}\,\tilde{\boldsymbol{\gamma}}\]where \(\boldsymbol{\mu}_\text{radio}\) is the posterior mean from a radio-timing inference run and \(L_\text{wide} = \sigma_\text{scale}\,L\) is the Cholesky factor of the radio posterior covariance scaled by
sigma_scale(default 1.0). Increasingsigma_scalebroadens the reparametrized prior around the radio posterior. This concentrates the sampler around physically reasonable EOS.- __init__(
- crust_name='SLy',
- n_points_high=500,
- reparametrized=False,
- sigma_scale=1.0,
Initialize spectral decomposition EOS model.
- Parameters:
crust_name (
str) – Name of crust model to use (‘SLy’, ‘DH’, ‘BPS’). Default ‘SLy’ for LALSuite compatibility.n_points_high (
int) – Number of high-density points to generate. Default 500 to match LALSuite although this is a bit high.reparametrized (
bool) – If True, expect tilde parameters \(\tilde\gamma_i\) drawn from a 4-D unit normal and apply the affine map \(\gamma = \mu_\text{radio} + L_\text{wide}\,\tilde\gamma\) insideconstruct_eos(). The_tildesuffix in the parameter names signals to the inference system that the whitened parametrization is active. Default False (original parametrization).sigma_scale (
float) – Multiplicative factor applied directly to the base Cholesky factor \(L\) to form \(L_\text{wide} = \sigma_\text{scale}\,L\). This scales the standard deviation of the reparametrized distribution. Only used whenreparametrized=True. Values greater than 1 broaden the reparametrized prior around the radio posterior; the default 1.0 uses the posterior standard deviation exactly. Must be strictly positive.
Methods
__init__([crust_name, n_points_high, ...])Initialize spectral decomposition EOS model.
construct_eos(params)Construct full EOS from spectral parameters.
Return list of spectral parameters required for this EOS.
interpolate_eos(n, p, e)Convert physical EOS quantities to geometric units and compute auxiliary quantities.
Attributes
Lower-triangular Cholesky factor \(L\) of the radio posterior covariance (\(\Sigma_\text{radio} = L\,L^\top\)).
Posterior mean \(\boldsymbol{\mu}_\text{radio}\) of the spectral coefficients.
- REPARAM_L_BASE: Float[Array, '4 4'] = Array([[ 3.41174267e-01, 0.00000000e+00, 0.00000000e+00, 0.00000000e+00], [-2.12029883e-01, 1.15271551e-01, 0.00000000e+00, 0.00000000e+00], [ 3.18696840e-02, -3.73794902e-02, 8.27531880e-03, 0.00000000e+00], [-1.38603940e-03, 2.39274018e-03, -7.76593593e-04, 1.86425321e-04]], dtype=float64)#
Lower-triangular Cholesky factor \(L\) of the radio posterior covariance (\(\Sigma_\text{radio} = L\,L^\top\)). The effective transform uses \(L_\text{wide} = \sigma_\text{scale}\,L\) where
sigma_scaleis set at initialisation time.
- REPARAM_MEAN: Float[Array, '4'] = Array([ 0.92232738, 0.34177628, -0.08307007, 0.00413816], dtype=float64)#
Posterior mean \(\boldsymbol{\mu}_\text{radio}\) of the spectral coefficients.
- construct_eos(params)[source]#
Construct full EOS from spectral parameters.
This method: 1. Validates parameters and checks gamma bounds 2. Loads low-density crust (preprocessed via Crust class) 3. Generates high-density spectral region (500 points by default) 4. Stitches them together 5. Converts to geometric units and computes auxiliary quantities
- e0_geom = 9.54629006e-11#
- p0_geom = 4.43784199e-13#
- xmax = 12.3081#