Strangeon EOS inference pipeline
This is an example notebook about how to use our tools to analysis a observation constraint on strangeon star equation of state.
Here in this notebook, we are using a strangeon matter EOS.
The following is the package that we need.
[1]:
import InferenceWorkflow.BayesianSampler as sampler
import InferenceWorkflow.Likelihood as likelihood
import InferenceWorkflow.prior as prior
import math
import numpy as np
import EOSgenerators.Strangeon_EOS as Strangeon
from TOVsolver.constant import oneoverfm_MeV, m_rho, m_w,G,c
import TOVsolver.main as main
/Users/a9/CompactOject/InferenceWorkflow/Likelihood.py:314: SyntaxWarning: invalid escape sequence '\('
"""
The following is the strangeon matter EOS function.
We note that because the parameter Nq is an integer in the strangeon matter EOS.
In the following, we use a possible value of this parameter, Nq=18, as an example to carry out the Bayesian analysis.
[2]:
Nq=18
Set up prior
Next step, we need to set up the prior, first use parameters array to specify the variable name, this process should consistent with what you need to call them.
Define a prior transform function to define prior. Cube are set of random number from 0 to 1. This prior setting is standard set-up of UltraNest package, since we are using UltraNest to do nest-sampling.
We provided two options call from prior:”normal_Prior” and “flat_prior”.
Then the Parameters prior should all set.
We note that since we are doing Equation of state Inference from mass-radius data of neutron star measurement. The center density of the star should be also sampled. Otherwise will be a partially-defined prior, did not span all parameters space, and proved to be different with full-scope inference.
This request as randomly generate a density from a EoS range, however, this process is not that trivial, since we need to determine the upper limit of the central density of compact star — different equation of state will predict different upper bound, so here we need to use the prior-setting EoS parameters computing the EOS by
strangeon.compute_EOS Compute out EOS, put into
main.OutputMR
find out Mass Radius of this equation of state, then find out the last stable point of this equation of state.(first mass points that give the direvative to be negative)
and find out that index by len() function, then reset this max_d to be upper limit of this density range.
[3]:
parameters = [ 'epsilon', 'ns', 'd1']
# for two or more MR measurements, define
# parameters = ['epsilon', 'ns', 'd1', 'd2']
def prior_transform(cube):
params = cube.copy()
params[0] = prior.flat_prior(10, 170,cube[0]) #epsilon=10-170MeV
params[1] = prior.flat_prior(0.17,0.36,cube[1]) #ns=0.17-0.36fm^-3
epsilon = params[0]
ns = params[1]
theta1 = np.array([Nq, epsilon, ns])
n_min = 3 * theta1[2] / theta1[0]
n_max = 0.16 * 8 * 3 / theta1[0]
n_values = np.linspace(n_min, n_max, 2000)
energy_density, pressure = Strangeon.compute_EOS(n_values, theta1)
eps_total=energy_density
pres_total=pressure
RFSU2R = []
MFSU2R = []
density = np.logspace(14.3, 15.6, 50)
if all(x<y for x,y in zip(eps_total[:], eps_total[1:])) and all(x<y for x, y in zip(pres_total[:], pres_total[1:])):
MR = main.OutputMR('',energy_density,pressure).T
else:
MR = []
if len(MR) == False:
params[2] = 0
#params[3] = 0
#this line for showing how to add one more observation
else:
for i in range(len(MR[1])):
RFSU2R.append(MR[0][i])
MFSU2R.append(MR[1][i])
if i > 20 and MR[1][i] - MR[1][i-1]< 0:
break
if len(MFSU2R)==False:
params[2] = 0
#params[3] = 0
#this line for showing how to add one more observation
else:
max_index = len(MFSU2R)
max_d = np.log10(density[max_index-1])
params[2] = 14.3 + (max_d - 14.3) * cube[2]
#params[3] = 14.3 + (max_d - 14.3) * cube[3]
#this line for showing how to add one more observation
return params
In the upper part, we define a flat (uniform) prior for the parameters in the strangeon matter equation of state, due to the lack of constraints from terrestrial experiments.
Note that the above code is an example of Bayesian analysis for a given mass and radius observation measurement. For example, if you use the NICER data for the measurements of J0030, then you should define another parameter, except the strangeon EOS parameters, e.g. “d1” for the centre density of this measurement, in the meantime add “params[2]” to this code.
If you further consider the adjoint analysis with J0030+J0740, then you should define the other two parameters, e.g. “d1” and “d2” for the centre density of these two measurements, in the meantime add “params[3]” to the above code.
Set up likehood
We need to set up a likelihood, Using standard definition way of UltraNest, that is below.
Here the likelihood is generated from a simulated mass radius measurement, which is 𝑀=1.4𝑀⊙ and 𝑅=13 km, With a 5% Mass radius measurement uncertainty,
so here
likelihood.MRlikihood_Gaussian
function will be use for our likelihood, please check likelihood.MRlikihood_Gaussian to see the original code, and more choice of likelihood. eg:
1.If we have some real mass-radius measurements, say PSR J0030 or PSR J0740, come from NICER, a KDE kernel could be trained to feed into
likelihood.MRlikihood_kernel(eps_total,pres_total,x,d1)
set the KDE kernel as a input for this function
2.If we gain measurement from radio-timing, say only measure the neutron star mass, then
likelihood.Masslikihood_Gaussian(eps_total,pres_total,x,d1)
Which will give the likelihood from single mass measurement, x is the parameters of that measurement, you should specify where this measurement mass is located and what is the sigma width of this mass measurement.
3.If we have nuclear measurements, and want to constrain this RMF model by nuclear properties like K(The Incompressibility of nuclear matter),J ( the symmetry energy at saturation density) and L( the slope of symmetry energy at saturation density). You can choose:
likelihood.Kliklihood(theta,K_low,K_up) likelihood.Jliklihood(theta,K_low,K_up) likelihood.Lliklihood(theta,K_low,K_up)
We are defaulting a hard-cut flat constrain, so if you don’t like this default hard cut, also could define the likelihood by youself with similiar style.
4.If we have a Tidal measurements from Gravitational wave detector, we can use it to do constraint:
likelihood.TidalLikihood_kernel(eps_total,pres_total,x,d1)
Where x is sampled distribution from real measurements, the standard is
kernel, chrip = x,
where the kernel is a whole set sampling from GW event, that is [chrip mass, M2/M1, tidal of M1, tidal of M2] four quantities. Chrip is the single smapling that comes only the chrip mass sampling.
[4]:
import scipy.stats as stats
def likelihood_transform(theta):
# This is a demonstration code for only introduce one constraint from one mass-radius observation.
# Could be very easy to implement more constraint from nuclear quantity, since that do not need to
# sample more central density of real neutron star. If user want to expand to two mass radius measurement
# the code could be:
epsilon, ns, d1 = theta
# comment this line if you need two measuremnts.
#epsilon, ns, d1, d2 = theta
####################################################################################################################
############ This is the block to compute out all the EoS you need based on your parameters#########################
theta1 = np.array([Nq, epsilon, ns])
n_min = 3 * theta1[2] / theta1[0]
n_max = 0.16 * 8 * 3 / theta1[0]
n_values = np.linspace(n_min, n_max, 2000)
energy_density, pressure = Strangeon.compute_EOS(n_values, theta1)
####################################################################################################################
probMRgaussian = likelihood.MRlikihood_Gaussian(energy_density,pressure,(1.4,13,0.07,0.65),d1)
prob = probMRgaussian
return prob
In the following, we will show how to modify the likehood_transform function when considering more observations.
def likelihood_transform(theta):
Nq, epsilon, ns, d1 = theta
theta = np.array([epsilon, ns])
n_min = 3 * theta[1] / Nq #n_min=3*ns/Nq
n_max = 0.16 * 8 * 3 / Nq #n_max=0.16*8*3/Nq
n_values = np.linspace(n_min, n_max, 100) # 100 points between n_min and n_max
energy_density_total, pressure_total = Strangeon_compute_EOS(n_values, theta)
#1. This line is to compute MR likelihood from a Simulated MR measurement:
probMRgaussian = MRlikihood_Gaussian(energy_density_total,pressure_total,(1.4,13,0.07,0.65),d1)
#probMRgaussian = likelihood.MRlikihood_Gaussian(energy_density_total,pressure_total,(2.08,13,0.08,0.65),d2)
#2. This is a block that constrain from given real MR measurement, say J0030:
#J0030 = numpy.loadtxt('data/PST_equal_sampled_MR.txt', delimiter=' ')
#J30R_list, J30M_list = zip(*J0030)
#J30R_list = numpy.array(J30R_list).T
#J30M_list = numpy.array(J30M_list).T
#Rmin = J30R_list.min()
#Rmax = J30R_list.max()
#Mmin = J30M_list.min()
#Mmax = J30M_list.max()
#X3, Y3 = numpy.mgrid[Rmin:Rmax:500j, Mmin:Mmax:100j]
#positions = numpy.vstack([X3.ravel(), Y3.ravel()])
#values = numpy.vstack([J30R_list, J30M_list])
#kernel3 = stats.gaussian_kde(values)
#probMRJ0030 = likelihood.MRlikelihhood_kernel(eps_total,pres_total,kernel3,d1)
#3. This is to compute the constraint from experiment of nuclearmatter
# 250<K<400, 25<J<38, 30<L<86:
# hint: since this K,J,L sampling don't need to sample central density so this
# theta should be redifined.
#probK = likelihood.Kliklihood(theta,250,400)
#probJ = likelihood.Jliklihood(theta,250,400)
#probL = likelihood.Lliklihood(theta,250,400)
#4. This is block to cosntrain from a GW event, say GW170817, here the file of this
# event is origanized by [chrip mass, M2/M1, tidal of M1, tidal of M2, sampling weight]:
#GW170817 = np.load('GW170817_McQL1L2weights.npy')
#chrip170817 = stats.gaussian_kde(GW170817[:,0],weights = GW170817[:,4])
#kernelGW = stats.gaussian_kde(GW170817.T[0:4],weights = GW170817[:,4])
#probGW = likelihood.TidalLikihood_kernel(eps_total,pres_total,(kernelGW,chrip170817),d1)
prob = probMRgaussian
#prob = probGW#+ probMRJ0030 + probK + probJ + probL + probGW
return prob
Set up sampler
Here next, we define sampler, there is two different sampler we provided for.
Considering where you need resume file:
sampler.UltranestSampler and sampler.UltranestSamplerResume
Here since it is our first run, so we only use first one. Some of the sampler parameters is requested, first is step number, our choice for UltraNest sampler is slicesampler, which could easily be sliced up your total computation load, and parallelize, speed up sampling. So step as suggested by documentation of UltraNest, we use 2*len(parameters).
live_point we set 2000, it will influence the sampling precision, We suggest for 7 dimension space, maybe 5000 is a better choice, however, since my computer only have limited resources, we set 2000.
max_calls set 10000, it is how many iteration after it will stop, we suggest to set this number significantly higher, otherwise maybe will broken before the inference converging to a definite value. That result will be un-phyiscal.
[ ]:
step = 2 * len(parameters)
live_point = 400
max_calls = 60000
samples = sampler.UltranestSampler(parameters,likelihood_transform,prior_transform,step,live_point,max_calls)
[ ]: