IMinuit.jl
Julia wrapper of the Python package iminuit
, which is the interface to the C++ MINUIT2, widely used in fitting in the high-energy physics community.
IMinuit.ArrayFit
IMinuit.Data
IMinuit.Fit
IMinuit.Minuit
IMinuit.chisq
IMinuit.contour_df
IMinuit.contour_df_given_parameter
IMinuit.contour_df_samples
IMinuit.func_argnames
IMinuit.get_contours
IMinuit.get_contours_all
IMinuit.get_contours_given_parameter
IMinuit.get_contours_samples
IMinuit.method_argnames
IMinuit.model_fit
IMinuit.@model_fit
IMinuit.@plt_best
IMinuit.@plt_data
IMinuit.ArrayFit
— TypeArrayFit <: AbstractFit
struct for fit with parameters collected in an Array
.
IMinuit.Data
— TypeData(x::T, y::T, err::T) where {T<:Vector{Real}}
Data(df::DataFrame)
Fields: x, y, err, ndata
This defines a type for data with three columns:x, y, err
; ndata
is the number of data rows. Different Data
sets can be concatenated as vcat(dat1, dat2, dat3)
.
Only symmetric errors (of y
) are supported.
IMinuit.Fit
— TypeFit <: AbstractFit
struct for fit with individual parameters.
IMinuit.Minuit
— MethodMinuit(fcn; kwds...)
Minuit(fcn, start; kwds...)
Minuit(fcn, m::AbatractFit; kwds...)
Wrapper of the iminuit
function Minuit
.
fcn
is the function to be optimized.start
: an array/tuple of the starting values of the parameters.kwds
is the list of keyword arguments ofMinuit
. For more information, refer to theiminuit
manual.m
: a fit that was previously defined; its parameters at the latest stage can be passed to the new fit.
Example:
fit = Minuit(fcn, [1, 0]; name = ["a", "b"], error = 0.1*ones(2),
fix_a = true, limit_b = (0, 50) )
migrad(fit)
where the parameters are collected in an array par
which is the argument of fcn(par)
. In this case, one can use external code (e.g., using ForwardDiff: gradient
) to compute the gradient as gradfun(par) = gradient(fcn, par)
, and include grad = gradfun
as a keyword argument.
If fcn
is defined as fcn(a, b)
, then the starting values need to be set as Minuit(fcn, a = 1, b = 0)
.
From iminuit
:
Minuit(fcn, throw_nan=False, pedantic=True, forced_parameters=None, print_level=0, errordef=None, grad=None, use_array_call=False, **kwds)
IMinuit.chisq
— Methodchisq(dist::Function, data, par; fitrange = ())
defines the $\chi^2$ function: fun
the function to be fitted to the data given by data
. The parameters are collected into par
, given as an array or a tuple.
data
can be either given as theData
type, or of the form(xdata, ydata [, err])
.
If no err
is given explicitly, the errors are assumed to be 1 for all data points.
fitrange
: default to the whole data set; may given as, e.g.,2:10
,
which means only fitting to the 2nd to the 10th data points.
IMinuit.contour_df
— Methodcontour_df(fit::AbstractFit, χsq; npts=20, limits=true, sigma = 1.0)
parameters in the form of a dataframe.
IMinuit.contour_df_given_parameter
— Methodcontour_df_given_parameter(fit::AbstractFit, χsq, para::T, range; limits = true) where {T <: Union{Symbol, String}}
return parameter sets in one sigma for a given parameter constrained in a range as a DataFrame
.
See also get_contours_given_parameter
.
IMinuit.contour_df_samples
— Methodcontour_df_samples(fit::AbstractFit, χsq, paras, ranges; nsamples = 100, MNbounds=true)
gives 1σ parameter sets as a DataFrame
for given parameters constrained in ranges
:
- if
paras
is a single parameter, then take equally spacednsamples
inranges
given in the form of(min, max)
; - if
paras
contain more parameters, thenparas
should be of the form(:para1, :para2)
,ranges
should be of the form((min1, max1), (min2, max2))
; paras
can be more than 2. Values for the parameters given inparas
are randomly sampled in the givenranges
.- is
MNbounds
is true, then constrain the parameters in the range provided byMINOS
no matter whether that is valid or not (to be improved by checking the validity) - if
igrad
is true, then useForwardDiff.gradient
to compute the gradient
IMinuit.func_argnames
— Methodfunc_argnames(f::Function)
Extracting the argument names of a function as an array.
IMinuit.get_contours
— Methodget_contours(fit::AbstractFit, χsq, parameters_combination::Vector{Int}; npts::Int=20, limits=true, sigma = 1.0)
For a given fit fit
and the $\chi^2$ function χsq
, gives an array of parameter arrays, with each array corresponding to a set of parameters obtained from calculating the MINOS
$1σ$ contour (try to find npts
points in the contour) for the two parameters in parameters_combination
.
parameters_combination
is an Int
array of the numbers of that two parameters, e.g. it is [1, 2]
for the first two parameters and [2, 3]
or the second and third parameters.
If limits
is true
, then fix one parameter to its bounds from MINOS
of the best fit and get the values for the other parameters; this runs over all parameters.
IMinuit.get_contours_all
— Methodget_contours_all(fit::AbstractFit, χsq; npts=20, limits=true, sigma = 1.0)
For a given fit fit
and the $\chi^2$ function χsq
, gives a list of parameters sets which are at the edge of $1σ$ MINOS
contours for all combinations of varying parameters. The case of limits
being true
runs only once.
IMinuit.get_contours_given_parameter
— Methodget_contours_given_parameter(fit::AbstractFit, χsq, para::T, range)
where {T <: Union{Symbol, String}}
gives parameter sets in one sigma for a given parameter constrained in a range. If fit
is an ArrayFit
and no user-defined names have been given to the parameters, then para
is "x0"
or :x0
for the 1st parameter, "x1"
or :x1
for the 2nd parameter, ...
IMinuit.get_contours_samples
— Methodget_contours_samples(fit::AbstractFit, χsq, paras, ranges; nsamples = 100, MNbounds = true)
return 1σ parameter sets as an Array
(the latter returns a DataFrame
) for given parameters constrained in ranges
:
- if
paras
is a single parameter, then take equally spacednsamples
inranges
given in the form of(min, max)
; - if
paras
contain more ($\geq 2$) parameters, thenparas
should be of the form(:para1, :para2)
,ranges
should be of the form((min1, max1), (min2, max2))
,
and values for the parameters given in paras
are randomly sampled in the given ranges
;
- if
MNbounds
is true, then constrain the parameters in the range provided byMINOS
no matter whether that is valid or not (to be improved by checking the validity) - if
igrad
is true, then useForwardDiff.gradient
to compute the gradient. - For using array parameters, if no user-defined names have been given to the parameters,
paras
should be given such that "x0"
or :x0
for the 1st parameter, "x1"
or :x1
for the 2nd parameter, ...
IMinuit.method_argnames
— Methodmethod_argnames(m::Method)
Extracting the argument names of a method as an array.
Modified from [`methodshow.jl`](https://github.com/JuliaLang/julia/blob/master/base/methodshow.jl) (`Vector{Any}`` changed to `Vector{Symbol}`)
IMinuit.model_fit
— Methodmodel_fit(model::Function, data::Data, start_values; kws...)
convenient wrapper for fitting a model
to data
; the returning stype is ArrayFit
, which can be passed to migrad
, minos
etc.
model
is the function to be fitted todata
; it should be of the formmodel(x, params)
withparams
given either as an array or a tuple.
IMinuit.@model_fit
— Macro@model_fit model data start_values kws...
convenient wrapper for fitting a model
to data
; see model_fit
.
IMinuit.@plt_best
— Macro@plt_best(dist, fit, data, kws...)
@plt_best!(dist, fit, data, kws...)
A convenient macro for comparing the best-fit result with the data; all combinations of keyword settings for plot
in Plots
can be used for the optional arguments kws...
The ordering of dist
, fit
, and data
does not matter.
IMinuit.@plt_data
— Macro@plt_data(data, kws...)
@plt_data!(data, kws...)
Convenient mascros to make an errorbar plot of the data
; all combinations of keyword settings for scatter
in Plots
can be used for the optional arguments kws...