API
Types
EFTfitter.Correlation
EFTfitter.EFTfitterModel
EFTfitter.HighestDensityRegion
EFTfitter.Measurement
EFTfitter.MeasurementDistribution
EFTfitter.MeasurementRanks
EFTfitter.NuisanceCorrelation
EFTfitter.Observable
EFTfitter.SmallestInterval
EFTfitter.SumOfSmallestIntervals
EFTfitter.UncertaintyRanks
Functions
EFTfitter.BLUE
EFTfitter.apply_criterion
EFTfitter.cov_to_cor
EFTfitter.criterion_value
EFTfitter.get_correlations
EFTfitter.get_measurement_distributions
EFTfitter.get_measurements
EFTfitter.get_nuisance_correlations
EFTfitter.get_observables
EFTfitter.get_parameters
EFTfitter.get_smallest_interval_edges
EFTfitter.measurement_models
EFTfitter.rank_measurements
EFTfitter.rank_uncertainties
EFTfitter.to_correlation_matrix
EFTfitter.uncertainty_models
Documentation
EFTfitter.Correlation
— Typestruct Correlation
Fields:
matrix::Array{Float64, 2}
: Observables that are measured.active::Bool
: Use this uncertainty category in fit. Defaults totrue
.
Constructors:
Correlation(matrix::Array{<:Real, 2}; active::Bool = true)
EFTfitter.EFTfitterModel
— Typestruct EFTfitterModel
This is the central type for using EFTfitter. It comprises all information necessary for performing an analysis. Only active Measurement
and Correlation
objects will be considered.
Fields:
parameters::BAT.NamedTupleDist
measurements::NamedTuple{<:Any, <:Tuple{Vararg{Measurement}}}
measurementdistributions::NamedTuple{<:Any, <:Tuple{Vararg{MeasurementDistribution}}}
correlations::NamedTuple{<:Any, <:Tuple{Vararg{Correlation}}}
nuisances::Union{NamedTuple{<:Any, <:Tuple{Vararg{NuisanceCorrelation}}}, Nothing}
Constructors:
EFTfitterModel(
parameters::BAT.NamedTupleDist,
measurements::NamedTuple{<:Any, <:Tuple{Vararg{AbstractMeasurement}}},
correlations::NamedTuple{<:Any, <:Tuple{Vararg{AbstractCorrelation}}},
nuisances::Union{NamedTuple{<:Any, <:Tuple{Vararg{NuisanceCorrelation}}}, Nothing} = nothing
)
Examples:
model = EFTfitterModel(parameters, measurements, correlations) # no nuisance correlations
model = EFTfitterModel(parameters, measurements, correlations, nuisances) # with nuisance correlations
)
EFTfitter.HighestDensityRegion
— Typestruct HighestDensityRegion <: AbstractRankingCriterion
Type for specifying the ranking criterion to be the area of the two-dimensional marginal posterior region containing a probability mass p
for two certain parameters.
Fields:
keys::NTuple{2, Symbol}
: Names of the parameters that are used for the marginalized regions.p::Float64 = 0.9
: Probability mass to be enclosed in the smallest interval.bins::T = 200
: Number of bins for the histograms to calculate interval widths.
Constructors:
HighestDensityRegion(keys::NTuple{2, Symbol}, p=0.9, bins=200)
EFTfitter.Measurement
— Typestruct Measurement
Fields:
observable::Observable
: Observable that is measured.value::Float64;
: Measured value.uncertainties::NamedTuple{<:Any, <:Tuple{Vararg{Real}}}
: Uncertainties of the measurement as NamedTuple.active::Bool
: Use or exclude measurement in fit. Defaults totrue
.
Constructors:
Measurement(
observable::Observable,
value::Float64;
uncertainties::NamedTuple{<:Any, <:Tuple{Vararg{Real}}},
active::Bool = true
)
Measurement(
observable::Function,
value::Float64;
uncertainties::NamedTuple{<:Any, <:Tuple{Vararg{Real}}},
active::Bool = true
)
EFTfitter.MeasurementDistribution
— Typestruct MeasurementDistribution
Fields: * observable::Array{Observable, 1}
: Observables that are measured. * value::Array{Float64, 1}
: Measured values. * uncertainties::NamedTuple{<:Any, <:Tuple{Vararg{Array{Float64, 1}}}}
: Uncertainties of the measurement as NamedTuple. * active::Array{Bool, 1}
: Use or exclude bins in fit. Defaults to true
for all bins. * bin_names::Array{Symbol, 1}
: Suffixes that will be appended to the name of the measurement distribution for the individual bins. Defaults to [_bin1, _bin2, ...].
Constructors:
MeasurementDistribution(
observable::Array{Observable, 1},
vals::Array{<:Real, 1};
uncertainties::NamedTuple{<:Any, <:Tuple{Vararg{Union{Vector{Float64}, Vector{Int64}}}}},
active::Union{Bool, Array{Bool, 1}} = [true for i in 1:length(vals)],
bin_names::Array{Symbol, 1} = [Symbol("bin$i") for i in 1:length(vals)]
)
MeasurementDistribution(
observable::Array{Function, 1},
vals::Array{<:Real, 1};
uncertainties::NamedTuple{<:Any, <:Tuple{Vararg{Union{Vector{Float64}, Vector{Int64}}}}},
active::Union{Bool, Array{Bool, 1}} = [true for i in 1:length(vals)],
bin_names::Array{Symbol, 1}
)
EFTfitter.NuisanceCorrelation
— Typestruct NuisanceCorrelation
Fields:
unc_key::Symbol
: Name of uncertainty category.meas1::Symbol
: Name of first measurement.meas2::Symbol
: Name of second measurement.prior::Distribution
: Prior distribution. Accepts the typeDistribution
and all other types accepted by BAT.NamedTupleDist, e.g.Interval
orReal
.
Constructors:
NuisanceCorrelation(unc_key::Symbol, meas1::Symbol, meas2::Symbol, prior::Any)
EFTfitter.Observable
— Typestruct Observable
Fields:
func::Function
: Function returning the predicted value of the observable as a function of the parametersmin::Float64
: Minimum boundary for values of the observable. Defaults to-Inf
.max::Float64
: Maximum boundary for values of the observable. Defaults toInf
.
Constructors:
Observable(
func::Function;
min::Float64 = -Inf
max::Float64 = Inf
)
EFTfitter.SmallestInterval
— Typestruct SmallestInterval <: AbstractRankingCriterion
Type for specifying the ranking criterion to be the total width of the one-dimensional marginalized smallest interval containing a probability mass p
for a certain parameter.
Fields:
key::Symbol
: Name of the parameter that is used for the marginalized intervals.p::Float64 = 0.9
: Probability mass to be enclosed in the smallest interval.
Constructors:
SmallestInterval(key::Symbol, p=0.9, bins=200)
EFTfitter.SumOfSmallestIntervals
— Typestruct SumOfSmallestIntervals <: AbstractRankingCriterion
Type for specifying the ranking criterion to be the summed width of all one-dimensional marginalized smallest intervals containing a probability mass p
.
Fields:
p::Float64 = 0.9
: Probability mass to be enclosed in the smallest intervals.bins::T = 200
: Number of bins for the histograms to calculate interval widths.
Constructors:
SumOfSmallestIntervals(p=0.9, bins=200)
EFTfitter.BLUE
— MethodBLUE(model::EFTfitterModel)
Calculates the best linear unbiased estimator (BLUE) for multiple measurements of the same observable, according to https://www.sciencedirect.com/science/article/pii/0168900288900186.
Note: Works only for an EFTfitterModel
where all measurements have the same observable. If this is not the case, an error is thrown.
Returns a NamedTuple
with the fields:
:value
: BLUE value:unc
: BLUE uncertainty:weights
: Array with the weights for each measurement
Example:
blue = BLUE(model)
println(blue.value, blue.unc, blue.weights)
EFTfitter.cov_to_cor
— Methodcov_to_cor(cov::Array{<:Real, 2})
Convert a covariance matrix cov
to a correlation matrix and a vector of uncertainty values.
Returns a matrix and a vector. Throws a warning when the covariance matrix is not positive definite.
Example:
cor, unc = cov_to_cor(cov)
EFTfitter.get_correlations
— Methodget_correlations(m::EFTfitterModel)
Returns a NamedTuple
with the Correlations
s in the EFTfitterModel
.
EFTfitter.get_measurement_distributions
— Methodget_measurement_distributions(m::EFTfitterModel)
Returns a NamedTuple
with the MeasurementDistribution
s in the EFTfitterModel
.
EFTfitter.get_measurements
— Methodget_measurements(m::EFTfitterModel)
Returns a NamedTuple
with the Measurement
s in the EFTfitterModel
.
EFTfitter.get_nuisance_correlations
— Methodget_nuisance_correlations(m::EFTfitterModel)
Returns a NamedTuple
with the NuisanceCorrelation
s in the EFTfitterModel
.
EFTfitter.get_observables
— Methodget_observables(m::EFTfitterModel)
Returns a NamedTuple
with the Observable
s in the EFTfitterModel
. Note: The upper and lower limits are ignored and for each unique Function
s only one Observable
is returned.
EFTfitter.get_parameters
— Methodget_parameters(m::EFTfitterModel)
Returns model parameters.
EFTfitter.get_smallest_interval_edges
— Methodget_smallest_interval_edges(
samples::DensitySampleVector,
key::Union{Symbol, Real},
p::Real;
bins=200,
atol=0.0)
Calculates the edges of the smallest intervals containing the fraction p
of the probability of the marginal distribution of parameter key
.
Returns a NamedTuple
with the keys lower
and upper
which both contain Vector
s with the corresponding bin edges. Keywords:
bins=200
: The number of bins used fpr calculating the intervals and edges.atol=0.0
: Intervals are joined together when they are seperated by less than this value.
EFTfitter.to_correlation_matrix
— Methodto_correlation_matrix(
measurements::NamedTuple{<:Any, <:Tuple{Vararg{AbstractMeasurement}}},
correlations::Tuple{Symbol,Symbol, Union{<:Real, Array{<:Real, 2}}}...
)
Returns a `Matrix{Float64}` as the correlation matrix for the measurements with all diagonal-elements being unity.
The correlations are specified by passing tuples of two `Symbols` (the keys of the `Measurements`) with a value or matrix for the correlation coefficients.
If the matrix is not positive-definite, a warning is shown but the matrix is still returned.
Example:
```julia
dist_corr = [1.0 0.5 0.0;
0.5 1.0 0.0;
0.0 0.0 1.0]
another_corr_matrix = to_correlation_matrix(
measurements,
(:Meas1, :Meas2, 0.4),
(:Meas1, :MeasDist, 0.1),
(:MeasDist, :MeasDist, dist_corr),
(:MeasDist_bin2, :MeasDist_bin3, 0.3),
)
```
EFTfitter.MeasurementRanks
— Typestruct MeasurementRanks
Type for the results of the measurement ranking.
Fields:
names::Vector{Symbol}
: Names of the measurements.values::Vector{Float64}
: Values of the ranking, calculated according to the criterion.criterion<:AbstractRankingCriterion
: Criterion that was used for ranking.
Constructors:
function MeasurementRanks(
names::Vector{Symbol},
values::Vector{Float64},
criterion::AbstractRankingCriterion;
order::Symbol = :values,
rev::Bool = false
)
Keyword arguments:
order::Symbol = :values
: Specifies how thenames
andvalues
are sorted. Further options are
:names
for alphabetical sorting based on the measurement names or :none
for keeping the initial order of measurements.
rev::Bool=false
: Switch to invert the order.
EFTfitter.UncertaintyRanks
— Typestruct UncertaintyRanks
Type for the results of the unceertainty ranking.
Fields:
names::Vector{Symbol}
: Names of the uncertainty categories.values::Vector{Float64}
: Values of the ranking, calculated according to the criterion.criterion<:AbstractRankingCriterion
: Criterion that was used for ranking.
Constructors:
function UncertaintyRanks(
names::Vector{Symbol},
values::Vector{Float64},
criterion::AbstractRankingCriterion;
order::Symbol = :values,
rev::Bool = false
)
Keyword arguments:
order::Symbol = :values
: Specifies how thenames
andvalues
are sorted. Further options are
:names
for alphabetical sorting based on the uncertainty names or :none
for keeping the initial order of uncertainty categories.
rev::Bool=false
: Switch to invert the order.
EFTfitter.apply_criterion
— Methodapply_criterion(criterion<:AbstractRankingCriterion, samples::DensitySampleVector)
Applies the criterion
to the samples
and returns the corresponding value.
Available ranking criteria:
SumOfSmallestIntervals(p=0.9)
: summed width of all one-dimensional marginalized smallest intervals containing p=90% of the posterior probabilitySmallestInterval(key=:C1, p=0.9)
: width of the one-dimensional marginalized smallest interval of parameter:C1
containing p=90% of the posterior probability.- `HighestDensityRegion(keys::NTuple{2, Symbol}, p=0.9, bins=200): volume of the highest density regione
EFTfitter.criterion_value
— Methodfunction criterion_value(
criterion::AbstractRankingCriterion,
model::EFTfitterModel;
sampling_algorithm::BAT.AbstractSamplingAlgorithm=MCMCSampling()
)
Computes the value of the criterion
for the model
using the sampling_algorithm
.
Available ranking criteria:
SmallestIntervalsSum(p=0.9)
: summed width of all one-dimensional marginalized smallest intervals containing p=90% of the posterior probabilitySmallestInterval(key=:C1, p=0.9)
: width of the one-dimensional marginalized smallest interval of parameterkey=:C1
containing p=90% of the posterior probability
EFTfitter.measurement_models
— Methodfunction measurement_models(model::EFTfitterModel)
Creates a Vector
of EFTfitterModel
where always one of the initially active measurements is deactivated at a time. Returns a Vector
of EFTfitterModel
and a Vector
of Symbol
with the names of the currently deactivated measurement.
EFTfitter.rank_measurements
— Methodrank_measurements(
model::EFTfitterModel;
sampling_algorithm::BAT.AbstractSamplingAlgorithm = MCMCSampling(),
criterion = SmallestIntervalsSum(p=0.9),
order = :values,
rev = true
)
Computes a ranking of the individual measurement in the EFTfitterModel
based on a certain ranking criterion by performing the sampling according to the specified sampling_algorithm
. Returns a MeasurementRank
object.
By default, the summed width of all one-dimensional marginalized smallest intervals containing 90% of the posterior probability is used as the ranking criterion. The measurements are ranked by the relative increase of this value.
Available ranking criteria:
SmallestIntervalsSum(p=0.9)
: summed width of all one-dimensional marginalized smallest intervals containing p=90% of the posterior probabilitySmallestInterval(key=:C1, p=0.9)
: width of the one-dimensional marginalized smallest interval of parameter:C1
containing p=90% of the posterior probability.
EFTfitter.rank_uncertainties
— Methodrank_uncertainties(
model::EFTfitterModel;
sampling_algorithm::BAT.AbstractSamplingAlgorithm = MCMCSampling(),
criterion = SmallestIntervalsSum(p=0.9),
order = :values,
rev = true
)
Computes a ranking of the types of uncertainty in the EFTfitterModel
based on a certain ranking criterion by performing the sampling according to the specified sampling_algorithm
. Returns a MeasurementRank
object.
By default, the summed width of all one-dimensional marginalized smallest intervals containing 90% of the posterior probability is used as the ranking criterion. The uncertainty types are ranked by the relative decrease of this value.
Available ranking criteria:
SmallestIntervalsSum(p=0.9)
: summed width of all one-dimensional marginalized smallest intervals containing p=90% of the posterior probabilitySmallestInterval(key=:C1, p=0.9)
: width of the one-dimensional marginalized smallest interval of parameterkey=:C1
containing p=90% of the posterior probability
EFTfitter.uncertainty_models
— Methodfunction uncertainty_models(model::EFTfitterModel)
Creates a Vector
of EFTfitterModel
where always one of the initially active uncertainty types is deactivated at a time. Returns a Vector
of EFTfitterModel
and a Vector
of Symbol
with the names of the currently deactivated uncertainty types.