models
Module contents
- xga.models.convert_to_odr_compatible(model_func, new_par_name='β', new_data_name='x_values')[source]
This is a bit of a weird one; its meant to convert model functions from the standard XGA setup (i.e. pass x values, then parameters as individual variables), into the form expected by Scipy’s ODR. I’d recommend running a check to compare results from the original and converted functions where-ever this function is called - I don’t completely trust it.
- Parameters:
model_func (FunctionType) – The original model function to be converted.
new_par_name (str) – The name we want to use for the new list/array of fit parameters.
new_data_name (str) – The new name we want to use for the x_data.
- Returns:
A successfully converted model function (hopefully) which can be used with ODR.
- Return type:
FunctionType
- xga.models.derivative(func, x0, dx=1.0, n=1, args=(), order=3)[source]
Find the nth derivative of a function at a point.
Given a function, use a central difference formula with spacing dx to compute the nth derivative at x0.
This is intended as a drop-in replacement for Scipy’s misc.derivative function, which was deprecated in Scipy v1.10.0 and removed after Scipy v1.14.1. It has been directly copied/reconstructed from Scipy code.
- Parameters:
func (FunctionType) – Input function
x0 – The point at which the nth derivative is found.
dx – Spacing.
n – Order of the derivative. Default is 1.
args – Arguments
order – Number of points to use, must be odd.
models.base module
- class xga.models.base.BaseModel1D(x_unit, y_unit, start_pars, par_priors, model_name, model_pub_name, par_pub_names, describes, info, x_lims=None)[source]
Bases:
objectThe superclass of XGA’s 1D models, with base functionality implemented, including the numerical methods for calculating derivatives and abel transforms which can be overwritten by subclasses if analytical solutions are available. The BaseModel class shouldn’t be instantiated by itself, as it won’t do anything.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values.
par_priors (List[Dict[str, Quantity/str]]) – The priors on the model parameters, for any fitting function that uses them.
model_name (str) – The simple name of the particular model, e.g. ‘beta’.
model_pub_name (str) – A smart name for the model that might be used in a publication plot, e.g. ‘Beta Profile’
par_pub_names (List[str]) – The names of the parameters of the model, as they should be used in plots for publication. As matplotlib supports LaTeX formatting for labels these should use $$ syntax.
describes (str) – An identifier for the type of data this model describes, e.g. ‘Surface Brightness’.
info (dict) – A dictionary with information about the model, used by the info() method. Can hold a general description, reference, authors etc.
x_lims (Quantity) – Upper and lower limits outside of which the model may not be valid, to be passed as a single non-scalar astropy quantity, with the first entry being the lower limit and the second entry being the upper limit. Default is None.
- get_realisations(x)[source]
This method uses the parameter distributions added to this model by a fitting process to generate random realisations of this model at a given x-position (or positions).
- Parameters:
x (Quantity) – The x-position(s) at which realisations of the model should be generated from the associated parameter distributions. If multiple, non-distribution, radii are to be used, make sure to pass them as an (M,), single dimension, astropy quantity, where M is the number of separate radii to generate realisations for. To marginalise over a radius distribution when generating realisations, pass a multi-dimensional astropy quantity; i.e. for a single set of realisations pass a (1, N) quantity, where N is the number of samples in the parameter posteriors, for realisations for M different radii pass a (M, N) quantity.
- Returns:
The model realisations, in a Quantity with shape (len(x), num_samples) if x has multiple radii in it (num_samples is the number of samples in the parameter distributions), and (num_samples,) if only a single x value is passed.
- Return type:
Quantity
- abstract static model(x, pars)[source]
This is where the model function is actually defined, this MUST be overridden by every subclass model, hence why I’ve used the abstract method decorator.
- Parameters:
x (Quantity) – The x-position at which the model should be evaluated.
pars (List[Quantity]) – The parameters of model to be evaluated.
- Returns:
The y-value of the model at x.
- derivative(x, dx, use_par_dist=False)[source]
Calculates a numerical derivative of the model at the specified x value, using the specified dx value. This method will be overridden in models that have an analytical solution to their first derivative, in which case the dx value will become irrelevant.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured.
dx (Quantity) – The dx value to use during the calculation.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- nth_derivative(x, dx, order, use_par_dist=False)[source]
A method to calculate the nth order derivative of the model using a numerical method.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured.
dx (Quantity) – The dx value to use during the calculation.
order (int) – The order of the desired derivative.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The value(s) of the nth order derivative of the model at x, either calculated from the current best fit parameters, or a distribution.
- Return type:
Quantity
- inverse_abel(x, use_par_dist=False, method='direct')[source]
This method uses numerical methods to generate the inverse abel transform of the model. It may be overridden by models that have analytical solutions to the inverse abel transform. All numerical inverse abel transform methods are from the pyabel module, and please be aware that in my (limited) experience the numerical solutions tend to diverge from analytical solutions at large radii.
- Parameters:
x (Quantity) – The x location(s) at which to calculate the value of the inverse abel transform.
use_par_dist (bool) – Should the parameter distributions be used to calculate a inverse abel transform distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
method (str) – The method that should be used to calculate the values of this inverse abel transform. You may pass ‘direct’, ‘basex’, ‘hansenlaw’, ‘onion_bordas’, ‘onion_peeling’, ‘two_point’, or ‘three_point’.
- Returns:
The inverse abel transform result.
- Return type:
Quantity
- volume_integral(outer_radius, inner_radius=None, use_par_dist=False)[source]
Calculates a numerical value for the volume integral of the function over a sphere of radius outer_radius. The scipy quad function is used. This method can either return a single value calculated using the current model parameters, or a distribution of values using the parameter distributions (assuming that this model has had a fit run on it).
This method may be overridden if there is an analytical solution to a particular model’s volume integration over a sphere.
The results of calculations with single values of outer and inner radius are stored in the model object to reduce processing time if they are needed again, but if a distribution of radii are passed then the results will not be stored and will be re-calculated each time.
- Parameters:
outer_radius (Quantity) – The radius to integrate out to. Either a single value or, if you want to marginalise over a radius distribution when ‘use_par_dist=True’, a non-scalar quantity of the same length as the number of samples in the parameter posteriors.
inner_radius (Quantity) – The inner bound of the radius integration. Default is None, which results in an inner radius of 0 in the units of outer_radius being used.
use_par_dist (bool) – Should the parameter distributions be used to calculate a volume integral distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value
- Returns:
The result of the integration, either a single value or a distribution.
- Return type:
Quantity
- allowed_prior_types(table_format='fancy_grid')[source]
Simple method to display the allowed prior types and their expected formats. :param str table_format: The desired format of the allowed models table. This is passed to the
tabulate module (allowed formats can be found here - https://pypi.org/project/tabulate/), and alters the way the printed table looks.
- static compare_units(check_pars, good_pars)[source]
Simple method that will be used in the inits of subclasses to make sure that any custom start values passed in by the user match the expected units of the default start parameters for that model.
- Parameters:
check_pars (List[Quantity]) – The first list of parameters, these are being checked.
good_pars (List[Quantity]) – The second list of parameters, these are taken as having ‘correct’ units.
- Returns:
Only if the check pars pass the tests. We return the check pars list but with all elements converted to EXACTLY the same units as good_pars, not just equivalent.
- Return type:
List[Quantity]
- info(table_format='fancy_grid')[source]
A method that gives some information about this particular model. :param str table_format: The desired format of the allowed models table. This is passed to the
tabulate module (allowed formats can be found here - https://pypi.org/project/tabulate/), and alters the way the printed table looks.
- predicted_dist_view(radius, bins='auto', colour='lightslategrey', figsize=(6, 5))[source]
A simple view method, to visualise the predicted value distribution at a particular radius. Only usable if this model has had parameter distributions assigned to it.
- Parameters:
radius (Quantity) – The radius at which you wish to evaluate this model and view the predicted distribution.
bins (Union[str, int]) – Equivelant to the plt.hist bins argument, set either the number of bins or the algorithm to decide on the number of bins.
colour (str) – Set the colour of the histogram.
figsize (tuple) – The desired dimensions of the figure.
- par_dist_view(bins='auto', colour='lightseagreen')[source]
Very simple method that allows you to view the parameter distributions that have been added to this model. The model parameter and uncertainties are indicated with red lines, highlighting the value and enclosing the 1sigma confidence region.
- Parameters:
bins (Union[str, int]) – Equivalent to the plt.hist bins argument, set either the number of bins or the algorithm to decide on the number of bins.
colour (str) – Set the colour of the histogram.
- view(radii=None, xscale='log', yscale='log', figsize=(8, 8), colour='black')[source]
Very simple view method to visualise XGA models with the current parameters.
- Parameters:
radii (Quantity) – Radii at which to calculate points to plot, doesn’t need to be set if the model has x limits defined.
xscale (str) – The scale to apply to the x-axis, default is log.
yscale (str) – The scale to apply to the y-axis, default is log.
figsize (tuple) – The size of figure to be set up.
colour (str) – The colour that the line in the plot should be.
- property model_pars
Property that returns the current parameters of the model, by default they are the same as the parameter start values.
- Returns:
A list of astropy quantities representing the values of the parameters of this model.
- Return type:
List[Quantity]
- property model_par_errs
Property that returns the uncertainties on the current parameters of the model, by default these will be zero as the default model_pars are the same as the start_pars.
- Returns:
A list of astropy quantities representing the uncertainties on the parameters of this model.
- Return type:
List[Quantity]
- property start_pars
Property that returns the current start parameters of the model, by which I mean the values that certain types of fitting function will use to start their fit.
- Returns:
A list of astropy quantities representing the values of the start parameters of this model.
- Return type:
List[Quantity]
- property unitless_start_pars
Returns sanitised start parameters which are floats rather than astropy quantities, sometimes necessary for fitting methods.
- Returns:
Array of floats representing model start parameters.
- Return type:
np.ndarray
- property par_priors
Property that returns the current priors on parameters of the model, these will be used by any fitting function that sets priors on parameters. Each entry in this list will be a dictionary with two keys ‘prior’ and ‘type’. The value for prior will be an astropy quantity, and the value for type will be a prior type (so uniform, gaussian, etc.)
- Returns:
A list of astropy quantities representing the values of the start parameters of this model.
- Return type:
List[Quantity]
- property x_unit
Property to access the expected x-unit of this model.
- Returns:
Astropy unit of the x values of the model.
- Return type:
Unit
- property y_unit
Property to access the expected y-unit of this model.
- Returns:
Astropy unit of the y values of the model.
- Return type:
Unit
- property x_lims
Property to access the x limits within which the model is considered valid, the default is None if no x limits were set for the model on instantiation.
- Returns:
A non-scalar astropy quantity with two entries, the first is a lower limit, and the second an upper limit. The default is None if no x limits were set.
- Return type:
Quantity
- property par_units
A list of units for the parameters of this model.
- Returns:
A list of astropy units.
- Return type:
List[Unit]
- property name
Property getter for the simple name of the model, which the user would enter when requesting a particular model to be fit to a profile, for instance.
- Returns:
String representation of the simple name of the model.
- Return type:
str
- property publication_name
Property getter for the publication name of the model, which is what would be added in a plot meant for publication, for instance.
- Returns:
String representation of the publication (i.e. pretty) name of the model.
- Return type:
str
- property par_publication_names
Property getter for the publication names of the model parameters. These would be used in a plot for instance, and so can make use of Matplotlib’s ability to render LaTeX math code.
- Returns:
List of string representation of the publication (i.e. pretty) names of the model parameters.
- Return type:
List[str]
- property describes
A one or two word description of the type of data this model describes.
- Returns:
A string description.
- Return type:
str
- property num_pars
Property getter for the number of parameters associated with this model.
- Returns:
Number of parameters.
- Return type:
int
- property par_dists
A property that returns the currently stored distributions for the model parameters, by default these will be empty quantities as no fitting will have occurred. Once a fit has been performed involving the model however, the distributions can be set externally.
- Returns:
A list of astropy quantities containing parameter distributions for all model parameters.
- Return type:
List[Quantity]
- property fit_warning
Returns any warnings generated by a fitting function that acted upon this model.
- Returns:
A string containing warnings.
- Return type:
str
- property acceptance_fraction
Property getter for the acceptance fraction of an MCMC fit (if one has been associated with this model).
- Returns:
The acceptance fraction.
- Return type:
int
- property emcee_sampler
Property getter for the emcee sampler used to fit this model, if applicable. By default this will be None, as the it has to be set externally, as the model won’t necessarily be fit by emcee
- Returns:
The emcee sampler used to fit this model.
- Return type:
em.EnsembleSampler
- property cut_off
- Property getter for the number of steps that an MCMC fitting method decided should be removed for burn-in.
By default this will be None, as the it has to be set externally, as the model won’t necessarily be fit by emcee
- Returns:
The number of steps to be removed for burn-in.
- Return type:
int
- property fit_method
Property getter for the method used to fit the model instance, this will be None if no fit has been run using this model.
- Returns:
The fit method.
- Return type:
str
- property par_names
The names of the parameters as they appear in the signature of the model python function.
- Returns:
A list of parameter names.
- Return type:
List[str]
- property success
If an fit has been run using this model then this property will tell you whether the fit method considered it to be ‘successful’ or not. If no fit has been run using this model then the value is None.
- Returns:
Was the fit successful?
- Return type:
bool
- property profile
The profile that this model has been fit to.
- Returns:
The profile object that this has been fit to, if no fit has been performed then this property will return None.
- Return type:
models.density module
- class xga.models.density.KingProfile1D(x_unit='kpc', y_unit=Unit('solMass / Mpc3'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of the King profile, describing an isothermal sphere. This describes a radial density profile and assumes spherical symmetry.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, beta, r_core, norm)[source]
The model function for the king profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
beta (Quantity) – The beta slope parameter of the model.
r_core (Quantity) – The core radius.
norm (Quantity) – The normalisation of the model.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the king profile at a given point, overriding the numerical method implemented in the BaseModel1D class, as this simple model has an easily derivable first derivative.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured. If multiple, non-distribution, radii are to be used, make sure to pass them as an (M,), single dimension, astropy quantity, where M is the number of separate radii to generate realisations for. To marginalise over a radius distribution when generating realisations, pass a multi-dimensional astropy quantity; i.e. for a single set of realisations pass a (1, N) quantity, where N is the number of samples in the parameter posteriors, for realisations for M different radii pass a (M, N) quantity.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- class xga.models.density.DoubleKingProfile1D(x_unit='kpc', y_unit=Unit('solMass / Mpc3'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of the double King profile, simply the sum of two King profiles. This describes a radial density profile and assumes spherical symmetry.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, beta_one, r_core_one, norm_one, beta_two, r_core_two, norm_two)[source]
The model function for the double King profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
beta_one (Quantity) – The beta slope parameter of the first King model.
r_core_one (Quantity) – The core radius of the first King model.
norm_one (Quantity) – The normalisation of the first King model.
beta_two (Quantity) – The beta slope parameter of the second King model.
r_core_two (Quantity) – The core radius of the second King model.
norm_two (Quantity) – The normalisation of the second King model.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the double King profile at a given point, overriding the numerical method implemented in the BaseModel1D class, as this simple model has an easily derivable first derivative.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured. If multiple, non-distribution, radii are to be used, make sure to pass them as an (M,), single dimension, astropy quantity, where M is the number of separate radii to generate realisations for. To marginalise over a radius distribution when generating realisations, pass a multi-dimensional astropy quantity; i.e. for a single set of realisations pass a (1, N) quantity, where N is the number of samples in the parameter posteriors, for realisations for M different radii pass a (M, N) quantity.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- class xga.models.density.SimpleVikhlininDensity1D(x_unit='kpc', y_unit=Unit('solMass / Mpc3'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of a simplified version of Vikhlinin’s full density model. Used relatively recently in https://doi.org/10.1051/0004-6361/201833325 by Ghirardini et al., a simplified form of Vikhlinin’s full density model, which can be found in https://doi.org/10.1086/500288.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, beta, r_core, alpha, r_s, epsilon, norm)[source]
The model function for the simplified Vikhlinin density profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
beta (Quantity) – The beta parameter of the model.
r_core (Quantity) – The core radius of the model.
alpha (Quantity) – The alpha parameter of the model.
r_s (Quantity) – The radius near where a change of slope by epsilon occurs.
epsilon (Quantity) – The epsilon parameter of the model.
norm (Quantity) – The overall normalisation of the model.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the simple Vikhlinin density profile at a given point, overriding the numerical method implemented in the BaseModel1D class.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured. If multiple, non-distribution, radii are to be used, make sure to pass them as an (M,), single dimension, astropy quantity, where M is the number of separate radii to generate realisations for. To marginalise over a radius distribution when generating realisations, pass a multi-dimensional astropy quantity; i.e. for a single set of realisations pass a (1, N) quantity, where N is the number of samples in the parameter posteriors, for realisations for M different radii pass a (M, N) quantity.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- class xga.models.density.VikhlininDensity1D(x_unit='kpc', y_unit=Unit('solMass / Mpc3'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of Vikhlinin’s full density model for galaxy cluster intra-cluster medium, which can be found in https://doi.org/10.1086/500288. It is a radial profile, so an assumption of spherical symmetry is baked in.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, beta_one, r_core_one, alpha, r_s, epsilon, gamma, norm_one, beta_two, r_core_two, norm_two)[source]
The model function for the full Vikhlinin density profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
beta_one (Quantity) – The beta parameter of the model.
r_core_one (Quantity) – The core radius of the model.
alpha (Quantity) – The alpha parameter of the model.
r_s (Quantity) – The radius near where a change of slope by epsilon occurs.
epsilon (Quantity) – The epsilon parameter of the model.
gamma (Quantity) – Width of slope change transition region.
norm_one (Quantity) – The normalisation of the model first part of the model.
beta_two (Quantity) – The beta parameter slope of the small core part of the model.
:param Quantity r_core_two:The core radius of the small core part of the model. :param Quantity norm_two: The normalisation of the additive, small core part of the model.
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the full Vikhlinin density profile at a given point, overriding the numerical method implemented in the BaseModel1D class.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured. If multiple, non-distribution, radii are to be used, make sure to pass them as an (M,), single dimension, astropy quantity, where M is the number of separate radii to generate realisations for. To marginalise over a radius distribution when generating realisations, pass a multi-dimensional astropy quantity; i.e. for a single set of realisations pass a (1, N) quantity, where N is the number of samples in the parameter posteriors, for realisations for M different radii pass a (M, N) quantity.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
models.entropy module
- class xga.models.entropy.CoreConstPowerEntropy(x_unit='kpc', y_unit=Unit('cm2 keV'), cust_start_pars=None)[source]
Bases:
BaseModel1DA model to fit galaxy cluster radial entropy profiles, made up of a constant core combined with a power law that dominates in the outskirts. The functional form is simple, but was found to be the best model for fitting many different entropy profiles by Cavagnolo et al. (https://doi.org/10.1088/0067-0049/182/1/12).
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, k_zero, k_100kpc, alpha)[source]
The model function for the constant-core and power-law entropy model.
- Parameters:
x (Quantity) – The radii to calculate y values for.
k_zero (Quantity) – Parameter quantifying the typical excess of core entropy above the best fitting power-law found at larger radii.
k_100kpc (Quantity) – A normalization for entropy at 100kpc.
alpha (Quantity) – The power law index.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- class xga.models.entropy.PowerEntropy(x_unit='kpc', y_unit=Unit('cm2 keV'), cust_start_pars=None)[source]
Bases:
BaseModel1DA model to fit galaxy cluster radial entropy profiles, made up of a power law. This is another model that was fit to a large sample of clusters by Cavagnolo et al. (https://doi.org/10.1088/0067-0049/182/1/12). This is an even simpler version of the CoreConstPowerEntropy model.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, k_100kpc, alpha)[source]
The model function for the power-law entropy model.
- Parameters:
x (Quantity) – The radii to calculate y values for.
k_zero (Quantity) – Parameter quantifying the typical excess of core entropy above the best fitting power-law found at larger radii.
k_100kpc (Quantity) – A normalization for entropy at 100kpc.
alpha (Quantity) – The power law index.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
models.fitting module
- xga.models.fitting.log_likelihood(theta, r, y, y_err, m_func)[source]
Uses a simple Gaussian likelihood function, returns the logged value.
- Parameters:
theta (np.ndarray) – The knowledge we have (think theta in Bayesian parlance) - gets fed into the model we’ve chosen.
r (np.ndarray) – The radii at which we have measured profile values.
y (np.ndarray) – The values we have measured for the profile.
y_err (np.ndarray) – The uncertainties on the measured profile values.
m_func – The model function that is being fit to.
- Returns:
The log-likelihood value.
- Return type:
np.ndarray
- xga.models.fitting.log_uniform_prior(theta, pr)[source]
This function acts as a uniform prior. Using the limits for the parameters in the chosen model (either user defined or default), the function checks whether the passed theta values sit within those limits. If they do then of course probability is 1, so we return the natural log (as this is a log prior), otherwise the probability is 0, so return -infinity.
- Parameters:
theta (np.ndarray) – The knowledge we have (think theta in Bayesian parlance) - gets fed into the model we’ve chosen.
pr (List) – A list of upper and lower limits for the parameters in theta, the limits of the uniform, uninformative priors.
- Returns:
The log prior value.
- Return type:
float
- xga.models.fitting.log_prob(theta, r, y, y_err, m_func, pr)[source]
The combination of the log prior and log likelihood.
- Parameters:
theta (np.ndarray) – The knowledge we have (think theta in Bayesian parlance) - gets fed into the model we’ve chosen.
r (np.ndarray) – The radii at which we have measured profile values.
y (np.ndarray) – The values we have measured for the profile.
y_err (np.ndarray) – The uncertainties on the measured profile values.
m_func – The model function that is being fit to.
pr (List) – A list of upper and lower limits for the parameters in theta, the limits of the uniform, uninformative priors.
- Returns:
The log probability value.
- Return type:
np.ndarray
models.mass module
- class xga.models.mass.NFW(x_unit='kpc', y_unit=Unit('solMass'), cust_start_pars=None)[source]
Bases:
BaseModel1DA simple model to fit galaxy cluster mass profiles (https://ui.adsabs.harvard.edu/abs/1997ApJ…490..493N/abstract) - a cumulative mass version of the Navarro-Frenk-White profile. Typically, the NFW is formulated in terms of mass density, but one can derive a mass profile from it (https://ui.adsabs.harvard.edu/abs/2006MNRAS.368..518V/abstract).
The NFW is extremely widely used, though generally for dark matter mass profiles, but will act as a handy functional form to fit to data-driven mass profiles derived from X-ray observations of clusters.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, Msun for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, r_scale, rho_zero)[source]
The model function for the constant-core and power-law entropy model.
- Parameters:
x (Quantity) – The radii to calculate y values for.
r_scale (Quantity) – The scale radius parameter.
rho_zero (Quantity) – A density normalization parameter.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
models.misc module
- xga.models.misc.straight_line(x_values, gradient, intercept)[source]
As simple a model as you can get, a straight line. Possible uses include fitting very simple scaling relations.
- Parameters:
x_values (np.ndarray/float) – The x_values to retrieve corresponding y values for.
gradient (float) – The gradient of the straight line.
intercept (float) – The intercept of the straight line.
- Returns:
The y values corresponding to the input x values.
- Return type:
Union[np.ndarray, float]
- xga.models.misc.power_law(x_values, slope, norm)[source]
A simple power law model, with slope and normalisation parameters. This is the standard model for fitting cluster scaling relations in XGA.
- Parameters:
x_values (np.ndarray/float) – The x_values to retrieve corresponding y values for.
slope (float) – The slope parameter of the power law.
norm (float) – The normalisation parameter of the power law.
- Returns:
The y values corresponding to the input x values.
- Return type:
Union[np.ndarray, float]
models.pressure module
- class xga.models.pressure.SimpleGNFWThermalPressure(x_unit='kpc', y_unit=Unit('keV / cm3'), cust_start_pars=None)[source]
Bases:
BaseModel1DA model to fit galaxy cluster radial thermal pressure profiles, based on the generalized NFW profile equation proposed by https://ui.adsabs.harvard.edu/abs/2007ApJ…668….1N/abstract and used in https://ui.adsabs.harvard.edu/abs/2010A%26A…517A..92A/abstract - this one has had several parameters removed or frozen.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV/cm^3 for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, p_zero, r_scale, alpha, beta)[source]
The model function for the simplified generalized NFW pressure profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
p_zero (Quantity) – The pressure normalization for the model.
r_scale (Quantity) – The scale radius for the radial pressure profile.
alpha (Quantity) – The alpha slope parameter, for the intermediate slope.
beta (Quantity) – The beta slope parameter, for the outer slope.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- class xga.models.pressure.GNFWThermalPressure(x_unit='kpc', y_unit=Unit('keV / cm3'), cust_start_pars=None)[source]
Bases:
BaseModel1DThe full generalized NFW profile equation used to fit galaxy cluster radial thermal pressure profiles, based on https://ui.adsabs.harvard.edu/abs/2007ApJ…668….1N/abstract, but first used for pressure in https://ui.adsabs.harvard.edu/abs/2010A%26A…517A..92A/abstract.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV/cm^3 for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, p_zero, r_scale, alpha, beta, gamma)[source]
The model function for the simplified generalized NFW pressure profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
p_zero (Quantity) – The pressure normalization for the model.
r_scale (Quantity) – The scale radius for the radial pressure profile.
alpha (Quantity) – The alpha slope parameter, for the intermediate slope.
beta (Quantity) – The beta slope parameter, for the outer slope.
gamma (Quantity) – The beta slope parameter, for the inner slope.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
models.sb module
- class xga.models.sb.BetaProfile1D(x_unit='kpc', y_unit=Unit('ct / (s arcmin2)'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of the beta profile, essentially a projected isothermal king profile, it can be used to describe a simple galaxy cluster radial surface brightness profile.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, beta, r_core, norm)[source]
The model function for the beta profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
beta (Quantity) – The beta slope parameter of the model.
r_core (Quantity) – The core radius.
norm (Quantity) – The normalisation of the model.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the beta profile at a given point, overriding the numerical method implemented in the BaseModel1D class, as this simple model has an easily derivable first derivative.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- inverse_abel(x, use_par_dist=False, method='analytical')[source]
This overrides the inverse abel method of the model superclass, as there is an analytical solution to the inverse abel transform of the single beta model. The form of the inverse abel transform is that of the king profile, but with an extra transformation applied to the normalising parameter. This method can either return a single value calculated using the current model parameters, or a distribution of values using the parameter distributions (assuming that this model has had a fit run on it).
- Parameters:
x (Quantity) – The x location(s) at which to calculate the value of the inverse abel transform.
use_par_dist (bool) – Should the parameter distributions be used to calculate a inverse abel transform distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
method (str) – The method that should be used to calculate the values of this inverse abel transform. Default for this overriding method is ‘analytical’, in which case the analytical solution is used. You may pass ‘direct’, ‘basex’, ‘hansenlaw’, ‘onion_bordas’, ‘onion_peeling’, ‘two_point’, or ‘three_point’ to calculate the transform numerically.
- Returns:
The inverse abel transform result.
- Return type:
Quantity
- class xga.models.sb.DoubleBetaProfile1D(x_unit='kpc', y_unit=Unit('ct / (s arcmin2)'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of the double beta profile, a summation of two single beta models. Often thought to deal better with peaky cluster cores that you might get from a cool-core cluster, this model can be used to describe a galaxy cluster radial surface brightness profile.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, beta_one, r_core_one, norm_one, beta_two, r_core_two, norm_two)[source]
The model function for the double beta profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
norm_one (Quantity) – The normalisation of the first beta profile.
beta_one (Quantity) – The beta slope parameter of the first component beta profile.
r_core_one (Quantity) – The core radius of the first component beta profile.
norm_two (Quantity) – The normalisation of the second beta profile.
beta_two (Quantity) – The beta slope parameter of the second component beta profile.
r_core_two (Quantity) – The core radius of the second component beta profile.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the double beta profile at a given point, overriding the numerical method implemented in the BaseModel1D class, as this simple model has an easily derivable first derivative.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- inverse_abel(x, use_par_dist=False, method='analytical')[source]
This overrides the inverse abel method of the model superclass, as there is an analytical solution to the inverse abel transform of the double beta model. The form of the inverse abel transform is that of two summed king profiles, but with extra transformations applied to the normalising parameters. This method can either return a single value calculated using the current model parameters, or a distribution of values using the parameter distributions (assuming that this model has had a fit run on it).
- Parameters:
x (Quantity) – The x location(s) at which to calculate the value of the inverse abel transform.
use_par_dist (bool) – Should the parameter distributions be used to calculate a inverse abel transform distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
method (str) – The method that should be used to calculate the values of this inverse abel transform. Default for this overriding method is ‘analytical’, in which case the analytical solution is used. You may pass ‘direct’, ‘basex’, ‘hansenlaw’, ‘onion_bordas’, ‘onion_peeling’, ‘two_point’, or ‘three_point’ to calculate the transform numerically.
- Returns:
The inverse abel transform result.
- Return type:
Quantity
models.temperature module
- class xga.models.temperature.SimpleVikhlininTemperature1D(x_unit='kpc', y_unit=Unit('keV'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of the simplified version of Vikhlinin’s temperature model. This is for the description of 3D temperature profiles (also can be fit to projected temperature) of galaxy clusters.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, r_cool, a_cool, t_min, t_zero, r_tran, c_power)[source]
The model function for the simplified Vikhlinin temperature profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
r_cool (Quantity) – Parameter describing the radius of the cooler core region.
a_cool (Quantity) – Power law parameter for the cooler core region.
t_min (Quantity) – A minimum temperature parameter for the model.
t_zero (Quantity) – A normalising temperature parameter for the model.
r_tran (Quantity) – The radius of the transition region of this broken power law model.
c_power (Quantity) – The power law index for the part of the model which describes the outer region of the cluster.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity
- derivative(x, dx=<Quantity 0.>, use_par_dist=False)[source]
Calculates the gradient of the simple Vikhlinin temperature profile at a given point, overriding the numerical method implemented in the BaseModel1D class.
- Parameters:
x (Quantity) – The point(s) at which the slope of the model should be measured.
dx (Quantity) – This makes no difference here, as this is an analytical derivative. It has been left in so that the inputs for this method don’t vary between models.
use_par_dist (bool) – Should the parameter distributions be used to calculate a derivative distribution; this can only be used if a fit has been performed using the model instance. Default is False, in which case the current parameters will be used to calculate a single value.
- Returns:
The calculated slope of the model at the supplied x position(s).
- Return type:
Quantity
- class xga.models.temperature.VikhlininTemperature1D(x_unit='kpc', y_unit=Unit('keV'), cust_start_pars=None)[source]
Bases:
BaseModel1DAn XGA model implementation of the full version of Vikhlinin’s temperature model. This is for the description of 3D temperature profiles (also can be fit to projected temperature) of galaxy clusters.
- Parameters:
x_unit (Unit/str) – The unit of the x-axis of this model, kpc for instance. May be passed as a string representation or an astropy unit object.
y_unit (Unit/str) – The unit of the output of this model, keV for instance. May be passed as a string representation or an astropy unit object.
cust_start_pars (List[Quantity]) – The start values of the model parameters for any fitting function that used start values. The units are checked against default start values.
- static model(x, r_cool, a_cool, t_min, t_zero, r_tran, a_power, b_power, c_power)[source]
The model function for the full Vikhlinin temperature profile.
- Parameters:
x (Quantity) – The radii to calculate y values for.
r_cool (float) – Parameter describing the radius of the cooling region (I THINK - NOT CERTAIN YET).
a_cool (float) – Power law parameter for the cooling region (I THINK - NOT CERTAIN YET).
t_min (float) – A minimum temperature parameter for the model (I THINK - NOT CERTAIN YET).
t_zero (float) – A normalising temperature parameter for the model (I THINK - NOT CERTAIN YET).
r_tran (float) – The radius of the transition region of this broken power law model.
a_power (float) – The first power law index.
b_power (float) – The second power law index.
c_power (float) – the third power law index.
- Returns:
The y values corresponding to the input x values.
- Return type:
Quantity