pdf#
PDF and PDF ratio implementation.
This module provides the PDF modeling and evaluation framework for csky. Here is an overview of the key concepts:
- PDF ratio models
In csky, it is computationally convenient and efficient to work in terms of log likelihood ratios rather than directly in terms of the likelihood itself. This enables us to represent all PDF calculations in terms of PDF ratios. In the standard case, it looks like:
\[\Lambda(\vec\mu) = -2 \ln \left[ \frac{L(n_s=0)}{L(\vec\mu)} \right] = 2\sum_\text{events} \ln \left[ \left(\frac{S(\text{event}|\vec\mu)}{B(\text{event}|\vec\mu)} - 1 \right) \cdot \frac{n_s}{N} + 1 \right],\]where \(S\) and \(B\) are the signal and background PDFs, and the fit parameters \(\vec\mu\) always start with \(n_s\) and typically also include a spectral index \(\gamma\). Because we have only the ratio \(S/B\), we can define PDF ratios \(W(\text{event}|\vec\mu)=S(\text{event}|\vec\mu)/B(\text{event}|\vec\mu)\).
Because the \(\text{events}\) vary with each trial, and \(\vec\mu\) changes with each log likelihood ratio evaluation, we want a piecemeal prescription for evaluating \(W\). The first step is the PDF ratio model, which you can think of as \(W(?|?)\). In csky, classes providing such a definition are derived from
PDFRatioModel
. These are the classes that the user is most likely to interact with directly, as they must be configured prior to running any trials.- PDF ratio evaluators
We have discussed trials so far without providing a definition. In csky, a trial is simply a list of one or more event ensembles (structured as
utils.Events
) along with a count of events that are present in principle but masked out of the detailed calculations in practice.Once we have a trial, we want to be able to evaluate \(\Lambda(\vec\mu)\). In general, many evaluations will be required to find the parameters \(\hat{\vec\mu}\) which maximize \(\Lambda\). Usually, significant optimizations can be achieved by caching certain values in between evaluations. Therefore, we do not directly ask PDF ratio models to evaluate PDF ratios. Instead we ask them for PDF ratio evaluators, which you can think of as \(W(\text{events}|?)\). In csky, classes providing this functionality are derived from
PDFRatioEvaluator
.- Acceptance Models
One of the crucial observations during the design of this framework is that we need to care about variable signal acceptance at many steps in the analysis process:
For source stacking, \(S(\text{event}|\vec\mu)\) becomes a weighted average over the values it would take for each individual source. The weighting in this average includes a factor \(w_\text{acc}\), the relative signal acceptance for each source. This typically depends on declination, but in the case of time-dependent analysis it also depends on the relative flux expectations for ecah source at the event arrival times.
For multi-dataset analysis, the overall likelihood is simply a product over the per-dataset likelihoods. However, we do not fit \((n_s)_d\) independently for each dataset \(d\). Instead, these values are constrained such that the ratios \((n_s)_1/(n_s)_2\), etc. are equal to the ratios of the absolute signal acceptances for each dataset, \((w_\text{acc})_1/(w_\text{acc})_2\) and so on.
For counts <–> flux/fluence conversion, we need a factor \(n_s/w_\text{acc}\). If we have a signal injector for the specific spectrum under consideration, this \(w_\text{acc}\) can be obtained directly from the signal MC. However, in general \(w_\text{acc}=w_\text{acc}(\vec\mu)\), and we want to be able to calculate its value for arbitrary \(\hat{\vec\mu}\).
Because variable signal acceptance is relevant in so many places, in csky we treat it as roughly on par with PDF ratio models and evaluators. Specifically, we only evaluate likelihoods in terms of acceptance weighted PDF ratio model/evaluator pairs:
AccWeightedPDFRatioModel
andAccWeightedPDFRatioEvaluator
. Each class implementingAccWeightedPDFRatioModel
must provide an acceptance parameterization, which by convention is an instance of a nested class (e.g.PointSourceSpacePDFRatioModel.AccModel
) that implements theAccModel
interface.- Parameterizations
Often, we can pre-configure some of the functions needed by PDF ratio evaluators, so that they can be re-used for many different analyses. The most common examples are the background space PDF (e.g.
BgSinDecParameterization
) and the declination- and spectrum-dependent signal acceptance (e.g.SinDecAccParameterization
). These parameterizations do not constitute PDF ratio models or evaluators themselves, but they are defined in this module because they are tools used by those classes. Note that for most applications we can also pre-compute anEnergyPDFRatioModel
, which is indeed a properPDFRatioModel
.
TODO:
discuss signal subtraction
document
TransientTimePDFRatioModel
andTransientTimePDFRatioEvaluator
Classes:
|
Traditional parameterization of the background space PDF. |
|
KDE-like parameterization of the background space PDF. |
|
Background space PDF that accounts for azimuthal uncertainty. |
Uniform background space PDF. |
|
|
Signal acceptance parameterization vs sin(dec) over a range of power law spectra. |
|
Signal acceptance parameterization vs sin(dec) for a specific spectrum. |
|
Parameterization of the dec and ra angular resolution in reco sin(dec) and log(energy) bins. |
|
Acceptance model base class. |
PDF ratio model base class. |
|
|
Model of energy PDF ratios. |
|
Model of energy PDF ratios for a custom flux. |
Uniform energy PDF ratio model. |
|
|
Generic PDF ratio model. |
Base class for acceptance-weighted PDF ratio models. |
|
|
|
|
Model of a product of separable PDF ratio models. |
|
Space PDF ratio model for one or more point or point-like sources. |
|
Space PDF ratio model that fits the sources. |
|
Template space PDF ratio model. |
|
Generic mixture model of any PDFRatio, also MultiPDFRatioModel. |
|
Space PDF ratio model for one or more point or point-like sources whose |
Base class for time PDF ratio models. |
|
|
Untriggered flare time PDF ratio model. |
|
Binned lightcurve time PDF ratio model. |
|
|
|
Space * Time acceptance weighted PDF ratio model. |
PDF ratio evaluator base class. |
|
|
|
|
|
|
|
|
|
|
Evaluator for SegmentedPDFRatioModel |
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
- class csky.pdf.BgSinDecParameterization(ev, dOmega_corr=[], hkw={}, skw={}, bg_mc_weight='', _other=None)[source]#
Bases:
object
Traditional parameterization of the background space PDF.
This class implements the traditional background space PDF as a function of the sin of the declination. The PDF satisfies the condition \(\int d\Omega f = \int d\alpha \int d(\sin(\delta)) f = 1\).
Methods:
__init__
(ev[, dOmega_corr, hkw, skw, ...])Construct a BgSinDecParameterization.
__call__
([ev, sindec])Evaluate for given events or sindecs.
- __init__(ev, dOmega_corr=[], hkw={}, skw={}, bg_mc_weight='', _other=None)[source]#
Construct a BgSinDecParameterization.
- Parameters:
ev (
utils.Events
) – the background-like eventshkw (mapping)
skw (mapping) –
histlite.Hist.spline_fit()
kwargs_other (BgSinDecParameterization) – another similarly-constructed BgSinDecKDEParameterization to be added to the first one.
- class csky.pdf.BgSinDecKDEParameterization(ev, hkw={}, skw={}, _other=None)[source]#
Bases:
BgSinDecParameterization
KDE-like parameterization of the background space PDF.
This class implements a background space PDF using a pseudo-KDE, where the kernel widths are determined by the event angular uncertainties.
Methods:
__init__
(ev[, hkw, skw, _other])Construct a BgSinDecParameterization.
__call__
([ev, sindec])Evaluate for given events or sindecs.
- __init__(ev, hkw={}, skw={}, _other=None)[source]#
Construct a BgSinDecParameterization.
- Parameters:
ev (
utils.Events
) – the background-like eventshkw (mapping)
skw (mapping) –
histlite.Hist.spline_fit()
kwargs_other (BgSinDecParameterization) – another similarly-constructed BgSinDecKDEParameterization to be added to the first one.
- __call__(ev=[], sindec=None)#
Evaluate for given events or sindecs.
- Parameters:
ev (
utils.Events
) – the eventssindec (array-like) – the sin(dec) array)
- Returns:
ndarray
- class csky.pdf.BgAzimuthSinDecParameterization(ev, hkw2={}, skw2={}, hkw={}, skw={}, smooth=None, _other=None)[source]#
Bases:
object
Background space PDF that accounts for azimuthal uncertainty.
This class implements a zenith-dependentn azimuthal asymmetry that modulates the expectation measured by a BgSinDecParameterization.
Methods:
- class csky.pdf.NullBgSpaceParameterization[source]#
Bases:
object
Uniform background space PDF.
This class is intended to implement a uniform background space PDF. It has not been tested thoroughly… if at all.
Methods:
- __init__()#
- class csky.pdf.SinDecAccParameterization(sig_ev, hkw={}, skw={}, gammas=array([1., 1.125, 1.25, 1.375, 1.5, 1.625, 1.75, 1.875, 2., 2.125, 2.25, 2.375, 2.5, 2.625, 2.75, 2.875, 3., 3.125, 3.25, 3.375, 3.5, 3.625, 3.75, 3.875, 4.]))[source]#
Bases:
object
Signal acceptance parameterization vs sin(dec) over a range of power law spectra.
This class implements a spline lookup of the signal acceptance \(A(\sin(\delta),\gamma)\). This is the time-independent acceptance – it must be integrated over time to obtain, e.g., the relative acceptance ratios across multiple data taking configurations.
Methods:
__init__
(sig_ev[, hkw, skw, gammas])Construct a SinDecAccParameterization.
__call__
(a, **params)Compute acceptance.
- __init__(sig_ev, hkw={}, skw={}, gammas=array([1., 1.125, 1.25, 1.375, 1.5, 1.625, 1.75, 1.875, 2., 2.125, 2.25, 2.375, 2.5, 2.625, 2.75, 2.875, 3., 3.125, 3.25, 3.375, 3.5, 3.625, 3.75, 3.875, 4.]))[source]#
Construct a SinDecAccParameterization.
- Parameters:
sig_ev (utils.Events) – the signal-like events
hkw (mapping) –
histlite.hist()
kwargsskw (mapping) –
histlite.Hist.spline_fit()
kwargsgammas (array-like) – set of gammas to evaluate
- class csky.pdf.SinDecCustomFluxAccParameterization(sig_ev, flux, hkw={}, skw={})[source]#
Bases:
object
Signal acceptance parameterization vs sin(dec) for a specific spectrum.
- Parameters:
sig_ev (utils.Events) – the signal-like events
flux (hyp.Flux) – the spectrum of interest
hkw (mapping) –
histlite.hist()
kwargsskw (mapping) –
histlite.Hist.spline_fit()
kwargs
Methods:
- class csky.pdf.Gauss2DAngResParameterization(sig, flux=PowerLawFlux(gamma=2.5), hkw={}, smooth=None, smooth_bins=None, fit=True)[source]#
Bases:
object
Parameterization of the dec and ra angular resolution in reco sin(dec) and log(energy) bins.
This is mainly here for reproducing the MESC-2yr analysis. You probably don’t want or need to use it.
Methods:
- class csky.pdf.AccModel[source]#
Bases:
object
Acceptance model base class.
This is a base class for acceptance models. Any acceptance-weighted PDFRatioModel must implement an acceptance model that derives from this class.
Methods:
get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
__init__
()- abstract get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)[source]#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- __init__()#
- class csky.pdf.PDFRatioModel[source]#
Bases:
object
PDF ratio model base class.
This is the base class for PDF ratio models. Models based on e.g. space, energy, and/or time information must derive from this class. These models essentially define the signal/background discrimination that will be used in the likelihood, and they should be aware of everything that is needed to compute the PDF ratios, once given a set of events.
PDFRatioEvaluator
instances, created automatically by the trial running machinery intrial
, are responsible for the actual evaluation based on a given set of events.Methods:
__call__
(ev, i)Apply the model to a set of events to obtain a PDFRatioEvaluator.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
__init__
()- abstract __call__(ev, i)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- set_ra(ra)[source]#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- __init__()#
- class csky.pdf.EnergyPDFRatioModel(bg_ev, sig_ev, features=['sindec', 'log10energy'], normalize_axes=[1], hkw={}, skw={}, gammas=array([1., 1.125, 1.25, 1.375, 1.5, 1.625, 1.75, 1.875, 2., 2.125, 2.25, 2.375, 2.5, 2.625, 2.75, 2.875, 3., 3.125, 3.25, 3.375, 3.5, 3.625, 3.75, 3.875, 4.]), bg_mc_weight='')[source]#
Bases:
PDFRatioModel
Model of energy PDF ratios.
This class implements the classic energy PDF ratio. Specifically, it computes 1D or 2D histograms for an ensemble of background-like and signal-like events, and then fits a spline (by default, in terms of log-value) tto the resulting histogram ratio. This is done on a grid of spectral indices \(\gamma`(see also :class:`CustomFluxEnergyPDFRatioModel\) for a fixed, non-power-law spectrum). The dimensions in practice are typically \(\sin(\delta,\log_{10)(E)})\). The axes along which to normalize can be specified; typically this is done along the \(\log_{10}(E)\) axis.
Formally, we can say that \(S=S(\log_{10}(E)|\vec{\mu})\), and similarly for \(B\).
In principle, this class can be used for arbitrary 1D or 2D spline-of-histogram PDF ratios, and it could be readily generalized to a signal weighting parameter other than \(\gamma\). If these cases arise in practice, this class should be refactored accordingly.
The 2D limitation is driven by difficulties fitting higher dimensional splines with scipy. A motivated developer could in principle migrate these fits to photospline, in which case arbitrary dimensionality should be possible.
Methods:
__init__
(bg_ev, sig_ev[, features, ...])Construct an EnergyPDFRatioModel.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __init__(bg_ev, sig_ev, features=['sindec', 'log10energy'], normalize_axes=[1], hkw={}, skw={}, gammas=array([1., 1.125, 1.25, 1.375, 1.5, 1.625, 1.75, 1.875, 2., 2.125, 2.25, 2.375, 2.5, 2.625, 2.75, 2.875, 3., 3.125, 3.25, 3.375, 3.5, 3.625, 3.75, 3.875, 4.]), bg_mc_weight='')[source]#
Construct an EnergyPDFRatioModel.
- Parameters:
bg_ev (utils.Events) – background-like events
sig_ev (utils.Events) – signal-like events
features (list of str) – the event features to histogram
normalize_axes (list of int) – the axes along which to normalize
hkw (mapping) –
histlite.hist()
kwargsskw (mapping) –
histlite.Hist.spline_fit()
kwargsgammas (array of float) – spectral indices to consider
bg_mc_weights (str) – if given, background histogram is weighted by bg_ev[bg_mc_weight]. otherwise, the default is equal weights per event.
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.CustomFluxEnergyPDFRatioModel(bg_ev, sig_ev, flux, features=['sindec', 'log10energy'], normalize_axes=[1], hkw={}, skw={}, bg_mc_weight='', keep_pdfs=False)[source]#
Bases:
PDFRatioModel
Model of energy PDF ratios for a custom flux.
This class implements the energy PDF ratio similar to
EnergyPDFRatioModel
, except for a specific flux which need not be a simple unbroken power law.Methods:
__init__
(bg_ev, sig_ev, flux[, features, ...])Construct a CustomFluxEnergyPDFRatioModel.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __init__(bg_ev, sig_ev, flux, features=['sindec', 'log10energy'], normalize_axes=[1], hkw={}, skw={}, bg_mc_weight='', keep_pdfs=False)[source]#
Construct a CustomFluxEnergyPDFRatioModel.
- Parameters:
flux (hyp.Flux) – the spectral hypothesis.
parameters (Other) – see
EnergyPDFRatioModel.__init__()
.
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.NullEnergyPDFRatioModel[source]#
Bases:
PDFRatioModel
Uniform energy PDF ratio model.
This class is intended to implement a uniform energy PDF ratio, i.e. simply \(S/B=1\). It has not been tested thoroughly… if at all. It also is likely useless; if you don’t want an energy PDF to be used, you should build your likelihood without one.
Methods:
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
__init__
()get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- __init__()#
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.GenericPDFRatioModel(func, features, fits, ana=None, src=None)[source]#
Bases:
PDFRatioModel
Generic PDF ratio model.
This class allows an arbitrary callable to be used as a PDFRatioModel.
Methods:
__init__
(func, features, fits[, ana, src])Construct a GenericPDFRatioModel.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the GenericPDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __init__(func, features, fits, ana=None, src=None)[source]#
Construct a GenericPDFRatioModel.
- Parameters:
func (callable) – python object that is callable and returns S/B. If update_bg setting is set to True, then the object must also implement a get_updated method.
features (str->str mapping) – mapping of func argname => Events column name
fits (str->number-or-sequence mapping) – mapping with str keys and number (for fixed values) or sequence (for fittable values) pairs. sequence form should be like (min_value, seed1, seed2, …, max_value)
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_updated(evs)[source]#
Update the GenericPDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.AccWeightedPDFRatioModel[source]#
Bases:
PDFRatioModel
Base class for acceptance-weighted PDF ratio models.
This is the base class for PDF ratio models which track acceptance weighting. Examples include:
PointSourceSpacePDFRatioModel
Tracks declination-dependent signal acceptance.
SpaceTimePDFRatioModel
Integrates declination-dependent signal acceptance over per-source time PDFs.
MultiPDFRatioModel
Combines exactly one
AccWeightedPDFRatioModel
with one or more additional non-acceptance-weightedPDFRatioModels`(such as :class:`EnergyPDFRatioModel
).
Attributes:
Acceptance model.
Methods:
__call__
(ev, i)Apply the model to a set of events to obtain a PDFRatioEvaluator.
__init__
()get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- abstract __call__(ev, i)#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- __init__()#
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.AccWeightedGenericPDFRatioModel(func, features, fits, ana=None, src=None, acc_param=None, bg_param=None, kent_min=0.12217304763960307, cut_n_sigma=5, sigsub=False)[source]#
Bases:
AccWeightedPDFRatioModel
Methods:
__init__
(func, features, fits[, ana, src, ...])Construct an AccGenericPDFRatioModel.
__call__
(ev[, i])Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the GenericPDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
Attributes:
The acceptance model (just a formality in this case).
Classes:
AccModel
(src, acc_param, livetime)Acceptance model for PointSourceSpacePDFRatioModel.
- __init__(func, features, fits, ana=None, src=None, acc_param=None, bg_param=None, kent_min=0.12217304763960307, cut_n_sigma=5, sigsub=False)[source]#
Construct an AccGenericPDFRatioModel. Intended for use with funcs that already take acceptance weighting into account.
- Parameters:
func (callable) – python object that is callable and returns S/B. If update_bg setting is set to True, then the object must also implement a get_updated method.
features (str->str mapping) – mapping of func argname => Events column name
fits (str->number-or-sequence mapping) – mapping with str keys and number (for fixed values) or sequence (for fittable values) pairs. sequence form should be like (min_value, seed1, seed2, …, max_value)
- __call__(ev, i=(None, None))[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_updated(evs)[source]#
Update the GenericPDFRatioModel in light of a set of events.
- Parameters:
evs (utils.Events) – the event ensembles
- abstract property acc_model#
The acceptance model (just a formality in this case).
- class AccModel(src, acc_param, livetime)[source]#
Bases:
AccModel
Acceptance model for PointSourceSpacePDFRatioModel. This acceptance model evaluates the time-independent acceptance for each source, weights these by the per-source intrinsic weights, and finally multiplies by the livetime.
Methods:
__init__
(src, acc_param, livetime)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.MultiPDFRatioModel(acc_weighted_model, *other_models, **kwargs)[source]#
Bases:
AccWeightedPDFRatioModel
Model of a product of separable PDF ratio models.
Methods:
__init__
(acc_weighted_model, *other_models, ...)Construct a MultiPDFRatioModel.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
Attributes:
Acceptance model.
- __init__(acc_weighted_model, *other_models, **kwargs)[source]#
Construct a MultiPDFRatioModel.
- Parameters:
acc_weighted_model (
AccWeightedPDFRatioModel
) – the acceptance-weighted model*other_models (
PDFRatioModel
) – one or more additional models
- set_ra(ra)[source]#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- class csky.pdf.PointSourceSpacePDFRatioModel(ana, src, bg_param=None, acc_param=None, cut_n_sigma=5, sigsub=False, mc_pdf=False, kent_min=0.12217304763960307)[source]#
Bases:
AccWeightedPDFRatioModel
Space PDF ratio model for one or more point or point-like sources.
This class implements the space PDF ratio; it is one of the most important workhorses of csky. The source list should give the coordinates \((\alpha,\delta)\) for one or more sources; multiple sources will result in a stacking analysis. Per-source extensions and intrinsic weights may be given as well.
A declination band cut is applied at a given \(n\sigma\) from each event, using the per-event angular error estimates. If signal subtraction is not enabled, an additional right ascension cut is applied(this is the so-called “box” cut).
Methods:
__init__
(ana, src[, bg_param, acc_param, ...])Construct a PointSourceSpacePDFRatioModel.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
__call__
(ev[, i])Apply the model to a set of events to obtain a PDFRatioEvaluator.
Attributes:
Acceptance model.
Classes:
AccModel
(src, acc_param, livetime)Acceptance model for PointSourceSpacePDFRatioModel.
- __init__(ana, src, bg_param=None, acc_param=None, cut_n_sigma=5, sigsub=False, mc_pdf=False, kent_min=0.12217304763960307)[source]#
Construct a PointSourceSpacePDFRatioModel.
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysissrc (
utils.Sources
) – the source list. required keys:ra
anddec
. optional keys:weight
,extension
bg_param – the background space PDF parameterization
acc_param – the signal acceptance parameterization
cut_n_sigma (float) – number of sigmas away from a source an event must be within to be included
sigsub (bool) – whether signal subtraction will be enabled
- set_ra(ra)[source]#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- class AccModel(src, acc_param, livetime)[source]#
Bases:
AccModel
Acceptance model for PointSourceSpacePDFRatioModel.
This acceptance model evaluates the time-independent acceptance for each source, weights these by the per-source intrinsic weights, and finally multiplies by the livetime.
Methods:
__init__
(src, acc_param, livetime)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- __call__(ev, i=(None, None))[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- class csky.pdf.FitPointSourceSpacePDFRatioModel(ana, src, bg_param=None, acc_param=None, cut_n_sigma=5, cut_extension=0.05235987755982989, sigma=None, sigsub=False)[source]#
Bases:
AccWeightedPDFRatioModel
Space PDF ratio model that fits the sources.
This model, which is implemented internally in terms of PointSourceSpacePDFRatioModel, adds per source coordinates \((\alpha,\delta)\) to the fitter parameters. This is mainly useful for finding the true hottest spot in a small, already-identified region, in which case the user can ask for such refinement without mentioning this class explicitly by name.
Do not expect good performance for multiple sources. This class may be broken if a time PDF ratio model is used.
Methods:
__init__
(ana, src[, bg_param, acc_param, ...])Construct a FitPointSourceSpacePDFRatioModel
__call__
(ev[, i])Apply the model to a set of events to obtain a PDFRatioEvaluator.
params_to_src
(src, **params)Extract source coordinate fit parameters into a source list.
src_to_params
(src)Convert a source list into a fit parameters mapping.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
Attributes:
Acceptance model.
Classes:
AccModel
(model)Acceptance model for FitPointSourceSpacePDFRatioModel.
- __init__(ana, src, bg_param=None, acc_param=None, cut_n_sigma=5, cut_extension=0.05235987755982989, sigma=None, sigsub=False)[source]#
Construct a FitPointSourceSpacePDFRatioModel
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysissrc (
utils.Sources
) – the source list. required keys:ra
anddec
. optional keys:weight
,extension
bg_param – the background space PDF parameterization
acc_param – the signal acceptance parameterization
cut_n_sigma (float) – number of sigmas away from a source an event must be within to be included
cut_extension (float) – additional Gaussian smearing width(in radians) to apply before pre-computing band and/or box cuts
sigma (float) – offset(in radians) of additional seed grid points in (ra,dec)
sigsub (bool) – whether signal subtraction will be enabled
- class AccModel(model)[source]#
Bases:
AccModel
Acceptance model for FitPointSourceSpacePDFRatioModel.
This acceptance model functions similarly to
PointSourceSpacePDFRatioModel.AccModel
, but it requires the per-source coordinates to be specified in order to perform the evaluation.Methods:
__init__
(model)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- __call__(ev, i=(None, None))[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- static params_to_src(src, **params)[source]#
Extract source coordinate fit parameters into a source list.
- Parameters:
src (
utils.Sources
) – base source list(may include weights and extensions)- Returns:
the source list
- Return type:
utils.Sources
- static src_to_params(src)[source]#
Convert a source list into a fit parameters mapping.
- Parameters:
src (
utils.Sources
) – source list- Returns:
the fit parameters mapping
- Return type:
dict
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.TemplateSpacePDFRatioModel(ana, template, bg_param=None, acc_param=None, flux=None, bins_energy=None, sigmas=None, cut_n_sigma=5, sigsub=False, fast_weight=False, sindec_bandwidth=0.017453292519943295, dir='', quiet=False, floor=1e-12, keep_pdfs=False)[source]#
Bases:
AccWeightedPDFRatioModel
Template space PDF ratio model.
This class implements a template space PDF ratio using healpy [1]. The template should be specified either as a 1D array of shape
(npix,)
or as a 2D array of shape(npix,nEbins)
. In the former case, a spectrum must be given, wheres in the latter case the energy bin edges must also be given. Because the Gaussian smoothing calculations can be computationally intense(despite the relatively efficient implementation in healpy), it is possible to specify a directory in which to keep a cache.If per-pixel binned spectra are given, the multiplication of the map by the signal acceptance is computed in declination bands using the signal MC directly unless requested otherwise(in which case the pixel-averaged spectrum is used).
[1] https://healpy.readthedocs.io/en/latest/
Methods:
__init__
(ana, template[, bg_param, ...])Construct a TemplateSpacePDFRatioModel.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
Attributes:
Acceptance model.
Classes:
AccModel
(src, acc_param, livetime, sindec_range)- __init__(ana, template, bg_param=None, acc_param=None, flux=None, bins_energy=None, sigmas=None, cut_n_sigma=5, sigsub=False, fast_weight=False, sindec_bandwidth=0.017453292519943295, dir='', quiet=False, floor=1e-12, keep_pdfs=False)[source]#
Construct a TemplateSpacePDFRatioModel.
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysistemplate (ndarray) – (npix,) array of per-pixel weights, or (npix,nEbins) array of per-pixel and per-energy-bin dN/dE values
bg_param – the background space PDF parameterization
acc_param – the signal acceptance parameterization
flux (
hyp.Flux
) – the signal spectrumbins_energy (kind) – The energy bin edges. shape: (nEbins+1,) Alternatively (for backwards-compatibility) only the left edges of the energy bins. The right edge of the last bin is logarithmically averaged from the previous bins in that case.
sigmas (ndarray) – if given, the Gaussian smoothing grid to consider. if not given, reasonable defaults are chosen depending on whether
ana
includes tracks or cascadescut_n_sigma (float) – number of sigmas away from a source an event must be within to be included
sigsub (bool) – whether signal subtraction will be enabled
fast_weight (bool) – whether to compute the acceptance-weighted maps very fast using the
acc_param
; otherwise, by default small declination bands of signal MC are used for each row of constant declination pixelssindec_bandwidth (float) – width of declination bands for computing the acceptance-weighted map.
dir (str) – path to a directory in which to cache this template’s state(filename will include
ana.key
)quiet (bool) – if true, suppress status printouts during construction
floor (float) – minimum allowed value for the normalized template and smoothed signal space PDFs
keep_pdfs (bool) – when have more pdfs have option of keep bkg pdfs for custom fluxes
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- class AccModel(src, acc_param, livetime, sindec_range)[source]#
Bases:
AccModel
Methods:
__init__
(src, acc_param, livetime, sindec_range)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.SegmentedPDFRatioModel(ana, models, logits=[-10, -1, 0, 1, 10], weights=[0, 10, 100, 1000], logit_fit=True, dir='', global_params={}, sigsub=True)[source]#
Bases:
AccWeightedPDFRatioModel
Generic mixture model of any PDFRatio, also MultiPDFRatioModel.
This class implements a generic additive mixture model of PDFRatioModel objects, fitting a weight w_i for each model. The result is a signal PDF of the form sum_i w_i PDF_i(). The relative weights of each model are fit parameters (which can also be kept constant). The weights will be normalized. Internally this is handled with a softmax transformation. The fitted parameters are the logits $v_i$, and they have to be converted back to the physical weights via the get_weights method. The fit parameters will be propagated through so that each PDF region
will retain the original parameters.
The physical interpretation is, that $w_i$s are the fraction of the total amount of signal events ns in that segment. There is a second (less efficient) way of using this class, if setting logit_fit to False (default is True): No softmax transformation is used than and directly the w_i are fit. simultaneously, the overall ns is kept at 1, so that (if not manually overwritten otherwise) the w_i can directly be interpreted as the number of signal events in that template. This is needed, if we want to fit negative ns in each segment individually.
Methods:
__init__
(ana, models[, logits, weights, ...])Construct a SegmentedPDFRatioModel
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
get_model_name
(idx[, name])Returns the correctly formatted name of either the model (default) with index idx or if name is given the name of e.g. a parameter which belongs to a model with index idx.
compute_weights_from_logits
(n_models, **params)Compute the weights w_i from the logits v_i.
get_weights
(**params)Compute the weights w_i from the logits v_i, if the pdf is in logit/softmax mode, or return the individual ns right away.
split_params
(n_models, params[, remove_weights])Turns a dict with param names p_i to a list with n_model dicts with param names p.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
Attributes:
Acceptance model.
Classes:
AccModel
(acc_models)- __init__(ana, models, logits=[-10, -1, 0, 1, 10], weights=[0, 10, 100, 1000], logit_fit=True, dir='', global_params={}, sigsub=True)[source]#
Construct a SegmentedPDFRatioModel
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysismodels (list of models) – A list of individual models which will be combined. The models have to be instantiated. The models need to be acceptance corrected.
logits (array of float) – The template model logits seeds and range for the fitter. This is ignored if logit_fit is false. The first and last values are taken as bounds and the central values as seeds.
weights (array of float) – The different weight range and seed for the fitter. This is only used if logit_fit is false. The first and last values are taken as bounds and the central values as seeds.
logit_fit (bool) – If True (default), combine the models with weights, which are calculated via the softmax transform so they are normalized. If False, the models are combined via some generic factor, which can then be interpreted as the ns in that segment. This introduces a redundancy if csky fits the overall ns, which is why ns must be fixed to 1 in that case.
dir (str) – path to a directory in which to cache this template’s state (filename will include
ana.key
)global_params (dict) – A dictionary of global parameters which will all be passed to each segment model in each call. Usefull to e.g. pass gamma=2.7 to each model.
sigsub (bool) – If sigsub is enabled. This options checks consistency between all provided models. All of them must be consistent in sigsub.
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- static get_model_name(idx, name='model')[source]#
Returns the correctly formatted name of either the model (default) with index idx or if name is given the name of e.g. a parameter which belongs to a model with index idx.
Using this method ensures consistency with SegmentedPDFRatioModel.name_pattern.
- static compute_weights_from_logits(n_models, **params)[source]#
Compute the weights w_i from the logits v_i.
- get_weights(**params)[source]#
Compute the weights w_i from the logits v_i, if the pdf is in logit/softmax mode, or return the individual ns right away.
- static split_params(n_models, params, remove_weights=True)[source]#
Turns a dict with param names p_i to a list with n_model dicts with param names p. Useful to pass to original models.
- class AccModel(acc_models)[source]#
Bases:
AccModel
Methods:
__init__
(acc_models)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.PriorSpacePDFRatioModel(ana, src, bg_param=None, acc_param=None, cut_n_sigma=5, sigsub=False, mc_pdf=False, kent_min=0.12217304763960307)[source]#
Bases:
PointSourceSpacePDFRatioModel
- Space PDF ratio model for one or more point or point-like sources whose
position is described by a spatial prior
The source list should contain the healpix prior for each of the sources.
Methods:
__init__
(ana, src[, bg_param, acc_param, ...])Construct a PointSourceSpacePDFRatioModel.
__call__
(ev[, i])Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
Classes:
AccModel
(src, acc_param, livetime)Acceptance model for PriorSpacePDFRatioModel.
Attributes:
Acceptance model.
- __init__(ana, src, bg_param=None, acc_param=None, cut_n_sigma=5, sigsub=False, mc_pdf=False, kent_min=0.12217304763960307)[source]#
Construct a PointSourceSpacePDFRatioModel.
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysissrc (
utils.Sources
) – the source list. required keys:ra
anddec
. optional keys:weight
,extension
bg_param – the background space PDF parameterization
acc_param – the signal acceptance parameterization
cut_n_sigma (float) – number of sigmas away from a source an event must be within to be included
sigsub (bool) – whether signal subtraction will be enabled
- __call__(ev, i=(None, None))[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- class AccModel(src, acc_param, livetime)[source]#
Bases:
AccModel
Acceptance model for PriorSpacePDFRatioModel.
This acceptance model evaluates the time-independent acceptance for each source, weights these by the per-source weights given the fitted positions, and finally multiplies by the livetime.
Methods:
__init__
(src, acc_param, livetime)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.TimePDFRatioModel[source]#
Bases:
PDFRatioModel
Base class for time PDF ratio models.
Time PDF ratio models must provide an additional method,
TimePDFRatioModel.get_frac_during_livetime()
, which gives the fraction of the total integrated signal lightcurve that is covered by active livetime for a givenanalysis.SubAnalysis
.Methods:
get_frac_during_livetime
(**params)Get the fraction of the total lightcurve covered by a given dataset.
__call__
(ev, i)Apply the model to a set of events to obtain a PDFRatioEvaluator.
__init__
()get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- abstract get_frac_during_livetime(**params)[source]#
Get the fraction of the total lightcurve covered by a given dataset.
- Returns:
- fraction of the total integrated signal lightcurve that is covered by
active livetime for a given
analysis.SubAnalysis
- Return type:
float
- abstract __call__(ev, i)#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- __init__()#
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.UntriggeredTimePDFRatioModel(ana, src, t0_min=None, t0_max=None, dt_min=1.1574074074074074e-11, dt_max=200, livetime=None, n_seeds_t0=6, n_seeds_dt=5, box=False, box_mode='center', use_grl=True, cut_n_sigma=5)[source]#
Bases:
TimePDFRatioModel
Untriggered flare time PDF ratio model.
This class implements a Gaussian or box-time-window flare signal hypothesis for which the peak time and flare duration are free to float. Gaussian is the default, but
box=True
can easily be selected. The anchor time is the center of the flare by default, butbox_mode='pre'
andbox_mode='post'
allow precursor and afterglow fits, respectively. The user can specify the allowed range for each of these values.In principle it is possible to fit for flares from multiple source candidates simultaneously in a stacking analysis, though this is not well tested if at all.
Methods:
__init__
(ana, src[, t0_min, t0_max, dt_min, ...])Construct a UntriggeredTimePDFRatioModel.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_t0_dt
([self, n_src, box_mode])Extract timing parameters as aligned arrays
(t0,dt)
get_frac_during_livetime
([t0_array, dt_array])Get the fraction of the total lightcurve covered by a given dataset.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __init__(ana, src, t0_min=None, t0_max=None, dt_min=1.1574074074074074e-11, dt_max=200, livetime=None, n_seeds_t0=6, n_seeds_dt=5, box=False, box_mode='center', use_grl=True, cut_n_sigma=5)[source]#
Construct a UntriggeredTimePDFRatioModel.
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysissrc (
utils.Sources
) – the source list(this class just needs to know it’s length)t0_min (float) – earliest allowed t0 (default: time of first event in the sub subanalysis
t0_max (float) – latest allowed t0 (default: time of first event in the sub subanalysis
dt_min (float) – smallest allowed dt
dt_max (float) – largest allowed dt
n_seed_t0 (int) – number of t0 seeds to try
n_seeds_dt (int) – number of dt seeds to try
cut_n_sigma (float) – number of sigmas away from a source an event must be within to be included
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- static get_t0_dt(self=None, n_src=None, box_mode='center', **params)[source]#
Extract timing parameters as aligned arrays
(t0,dt)
- Returns:
the t0 and dt arrays
- Return type:
tuple of (ndarray,ndarray)
- get_frac_during_livetime(t0_array=None, dt_array=None, **params)[source]#
Get the fraction of the total lightcurve covered by a given dataset.
- Returns:
- fraction of the total integrated signal lightcurve that is covered by
active livetime for a given
analysis.SubAnalysis
- Return type:
float
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.BinnedTimePDFRatioModel(ana, src, lcs, n_seeds_thresh=10, n_seeds_lag=10, range_lag=(-1, 1), thresh_seeding='quantiles', lag_seeding='uniform', use_grl=True, cut_n_sigma=None, thresh_seeds=None, lag_seeds=None)[source]#
Bases:
TimePDFRatioModel
Binned lightcurve time PDF ratio model.
This class implements a binned lightcurve signal time PDF. The PDFs are specified as histlite histograms. By default the flux threshold is a fit parameter. In principle, stacking of multiple sources should be supported, though this has not been well tested, if at all. Note that in the stacking case, this implementation fits the threshold for each source independently; some refactoring would be required to constrain the thresholds to match.
Todo
- probably want to (optionally?) correct
BinnedTimePDFRatioModel.get_frac_during_livetime()
to integrate only over true livetime.
Methods:
__init__
(ana, src, lcs[, n_seeds_thresh, ...])Construct a BinnedTimePDFRatioModel
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_threshs
(**params)Extract per-source thresholds.
get_lags
(**params)Extract per-source lagolds.
get_lc_pdfs
([threshs_array, lags_array])Generate normalized per-source time PDFs from lightcurves, given the per-source thresholds.
get_lc_pdfs_cropped
(**params)Get cropped signal time PDFs spanning only the livetime of the originally specified sub analysis.
get_frac_during_livetime
(**params)Get the fraction of the total lightcurve covered by a given dataset.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __init__(ana, src, lcs, n_seeds_thresh=10, n_seeds_lag=10, range_lag=(-1, 1), thresh_seeding='quantiles', lag_seeding='uniform', use_grl=True, cut_n_sigma=None, thresh_seeds=None, lag_seeds=None)[source]#
Construct a BinnedTimePDFRatioModel
- Parameters:
ana (
analysis.SubAnalysis
) – the sub analysissrc (
utils.Sources
) – the source list(this class just needs to know it’s length)lcs (array of
histlite.Hist
) – the lightcurves. values are flux; bin edges are in units of MJDn_seeds_thresh (int) – number of flux threshold seeds to try
n_seeds_lag (int) – number of lag seeds to try
range_lag (tuple of (float, float)) – range of lag times to fit
thresh_seeding (str) – one of ‘quantiles’ (recommended), ‘log’, ‘linear’,’custom’; determines how threshold seeds are selected
lag_seeding (str) – one of ‘uniform’,’custom’; determines how lag seeds are selected
cut_n_sigma (float) – ignored (but present for consistency with other time models)
thresh_seeds (list) – list of threshold seeds if thresh_seeding is ‘custom’
lag_seeds (list) – list of lag seeds if lag_seeding is ‘custom’
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_threshs(**params)[source]#
Extract per-source thresholds.
- Returns:
the flux thresholds for each source
- Return type:
ndarray of float
- get_lags(**params)[source]#
Extract per-source lagolds.
- Returns:
the flux lags for each source
- Return type:
ndarray of float
- get_lc_pdfs(threshs_array=None, lags_array=None, **params)[source]#
Generate normalized per-source time PDFs from lightcurves, given the per-source thresholds.
- Parameters:
threshs_array (ndarray of float) – the thresholds
lags_array (ndarray of float) – the lags in MJD
params – if
threshs_array
and/or lags_array are not given, they are extracted from the fit parameters which are given as kwargs
- Returns:
the signal time PDFs
- Return type:
ndarray of
histlite.Hist
- get_lc_pdfs_cropped(**params)[source]#
Get cropped signal time PDFs spanning only the livetime of the originally specified sub analysis.
Crop to sub analysis, or edge of lc pdf
- Returns:
the cropped signal time PDFs
- Return type:
ndarray of
histlite.Hist
The PDFs are not renormalized after being cropped.
- get_frac_during_livetime(**params)[source]#
Get the fraction of the total lightcurve covered by a given dataset.
- Returns:
- fraction of the total integrated signal lightcurve that is covered by
active livetime for a given
analysis.SubAnalysis
- Return type:
float
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.TransientTimePDFRatioModel(ana, src, use_grl=True, cut_n_sigma=4)[source]#
Bases:
TimePDFRatioModel
Methods:
__init__
(ana, src[, use_grl, cut_n_sigma])__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_frac_during_livetime
(**params)Get the fraction of the total lightcurve covered by a given dataset.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_frac_during_livetime(**params)[source]#
Get the fraction of the total lightcurve covered by a given dataset.
- Returns:
- fraction of the total integrated signal lightcurve that is covered by
active livetime for a given
analysis.SubAnalysis
- Return type:
float
- get_updated(evs)#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- set_ra(ra)#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class csky.pdf.SpaceTimePDFRatioModel(space_model, time_model, transient=False)[source]#
Bases:
AccWeightedPDFRatioModel
Space * Time acceptance weighted PDF ratio model.
This class implements a space * time PDF ratio model, where the acceptance weight is computed by integrating the declination-dependent signal acceptance over the per-source time PDFs.
Methods:
__init__
(space_model, time_model[, transient])Consstruct a SpaceTimePDFRatioModel.
set_ra
(ra)Set the right ascension of the(presumably one-and-only) source.
__call__
(ev)Apply the model to a set of events to obtain a PDFRatioEvaluator.
get_updated
(evs)Update the PDFRatioModel in light of a set of events.
Attributes:
Acceptance model.
Classes:
AccModel
(spacetime_model)Acceptance model for SpaceTimePDFRatioModel.
- __init__(space_model, time_model, transient=False)[source]#
Consstruct a SpaceTimePDFRatioModel.
- Parameters:
space_model (
AccWeightedPDFRatioModel
) – the space model(should be a point-like source model, not a template model)time_model (
TimePDFRatioModel
) – the time modeltransient (bool) – whether the time model should be understood as a transient source model. TODO: why is this here? probably a never-finished idea for some optimization?
Todo
why is transient here? probably a never-finished idea for some optimization?
- set_ra(ra)[source]#
Set the right ascension of the(presumably one-and-only) source.
- Parameters:
ra (float) – the right ascension
Users do not need to call this method explicitly.
This method sets the right ascension of a single source, enabling a powerful optimization when computing skymaps. Specifically, space PDFs can be recomputed at a new coordinate without repeating “band” and/or “box” masking calculations. The default implementation here does nothing, but any space PDFRatioModel’s should update when this method is called.
- class AccModel(spacetime_model)[source]#
Bases:
AccModel
Acceptance model for SpaceTimePDFRatioModel.
This acceptance model corrects the acceptance given by the space model’s
acc_model
, accounting for the fraction of the total livetime covered by the signal time PDF.Methods:
__init__
(spacetime_model)get_acc_per_source
(**params)Get the absolute acceptance on a per-source basis.
get_acc_total
(**params)Get the absolute acceptance, summing over sources.
- get_acc_per_source(**params)[source]#
Get the absolute acceptance on a per-source basis.
- Parameters:
params – fit arguments
- get_acc_total(**params)#
Get the absolute acceptance, summing over sources.
- Parameters:
params – fit arguments
- __call__(ev)[source]#
Apply the model to a set of events to obtain a PDFRatioEvaluator.
- Parameters:
ev (utils.Events) – the events
i (tuple of ([array of int], [array of int])) – i_ev, i_src: aligned arrays indicating pairs of events and sources to which the evaluator should apply
- Returns:
PDFRatioEvaluator
- get_updated(evs)[source]#
Update the PDFRatioModel in light of a set of events.
- Parameters:
evs (sequence of utils.Events) – the event ensembles
Users do not need to call this method explicitly.
This method updates the PDFRatioModel, given a sequence of event ensembles in the form of utils.Events. This allows an optional mode in which the background PDFs should be recomputed to account not only for data events(
evs[0]
) but also for injected signal events(evs[1:]
).
- class csky.pdf.PDFRatioEvaluator[source]#
Bases:
object
PDF ratio evaluator base class.
This is the base class for PDF ratio evaluators. Evaluators based on e.g. space, energy, and/or time information must derive from this class. Subclasses apply corresponding PDFRatioModels to a given set of events such that, once further given a set of fit parameter kwargs, they can provide the S/B PDF ratios, as well as the S-tilde/B scrambled signal PDF ratios used in the signal subtracting likelihood.
Users should not need to make direct, explicit use of any PDFRatioEvaluators; they are created automatically by the trial running machinery in
trial
.Methods:
- __init__()#
- class csky.pdf.AccWeightedPDFRatioEvaluator[source]#
Bases:
PDFRatioEvaluator
Methods:
- abstract __call__(mask=None, **params)#
Call self as a function.
- __init__()#
- class csky.pdf.MultiPDFRatioEvaluator(ev, model)[source]#
Bases:
AccWeightedPDFRatioEvaluator
Methods:
- class csky.pdf.PointSourceSpacePDFRatioEvaluator(ev, model, i=(None, None))[source]#
Bases:
AccWeightedPDFRatioEvaluator
Methods:
- class csky.pdf.MCPointSourceSpacePDFRatioEvaluator(ev, model, i=(None, None))[source]#
Bases:
AccWeightedPDFRatioEvaluator
Methods:
- class csky.pdf.Gauss2DPointSourceSpacePDFRatioEvaluator(ev, model, i=(None, None))[source]#
Bases:
PointSourceSpacePDFRatioEvaluator
Methods:
- __call__(_mask=None, **params)#
Call self as a function.
- class csky.pdf.TemplateSpacePDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.SegmentedPDFRatioEvaluator(ev, model)[source]#
Bases:
AccWeightedPDFRatioEvaluator
Evaluator for SegmentedPDFRatioModel
Methods:
- class csky.pdf.PriorSpacePDFRatioEvaluator(ev, model, i=(None, None))[source]#
Bases:
PointSourceSpacePDFRatioEvaluator
Methods:
- class csky.pdf.EnergyPDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.CustomFluxEnergyPDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.NullEnergyPDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.GenericPDFRatioEvaluator(ev, model, i=(None, None))[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.FitPointSourceSpacePDFRatioEvaluator(ev, model, i=(None, None))[source]#
Bases:
AccWeightedPDFRatioEvaluator
Methods:
- class csky.pdf.UntriggeredTimePDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.BinnedTimePDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods:
- class csky.pdf.TransientTimePDFRatioEvaluator(ev, model)[source]#
Bases:
PDFRatioEvaluator
Methods: