Public Documentation
Documentation for Minuit2.jl
public interface.
Index
Minuit2.BinnedNLL
Minuit2.Constant
Minuit2.CostSum
Minuit2.ExtendedBinnedNLL
Minuit2.ExtendedUnbinnedNLL
Minuit2.LeastSquares
Minuit2.Minuit
Minuit2.Minuit
Minuit2.UnbinnedNLL
Minuit2.RooFit.AbstractPdf
Minuit2.RooFit.AddPdf
Minuit2.RooFit.AddPdf
Minuit2.RooFit.ArgusPdf
Minuit2.RooFit.Chebyshev
Minuit2.RooFit.DataSet
Minuit2.RooFit.Exponential
Minuit2.RooFit.FitResult
Minuit2.RooFit.Gaussian
Minuit2.RooFit.RealVar
Minuit2.RooFit.RealVar
Minuit2.FCN
Minuit2.FCN
Minuit2.chi2
Minuit2.chi2_grad
Minuit2.contour
Minuit2.get_nargs
Minuit2.hesse!
Minuit2.matrix
Minuit2.migrad!
Minuit2.minos!
Minuit2.mncontour
Minuit2.mnprofile
Minuit2.multinomial_chi2
Minuit2.multinomial_chi2_grad
Minuit2.poisson_chi2
Minuit2.poisson_chi2_grad
Minuit2.profile
Minuit2.robust_low_level_fit
Minuit2.scan!
Minuit2.simplex!
Minuit2.soft_l1_cost
Minuit2.soft_l1_cost_grad
Minuit2.RooFit.ConstVar
Minuit2.RooFit.fitTo
Minuit2.RooFit.generate
Types
This is the list of all types and functions defined for Minuit2
Minuit2.BinnedNLL
— MethodBinnedNLL(bincounts::AbstractArray, binedges::Union{AbstractArray, Tuple}, cdf::Function; use_pdf=:none, verbose=0, grad=nothing, names=())
Binned negative log-likelihood.
Use this if only the shape of the fitted PDF is of interest and the data is binned. This cost function works with normal and weighted histograms. The histogram can be one- or multi-dimensional.
Arguments
bincounts::AbstractArray
: Histogram counts. If this is an array with dimension D, where D is the number of histogram axes.xe::Union{AbstractArray, Tuple}
: Bin edge locations, must be len(n) + 1, where n is the number of bins. If the histogram has more than one axis, xe must be a collection of the bin edge locations along each axis.cdf::Function
: Cumulative density function of the form f(xe, par0, par1, ..., parN), where xe is a bin edge and par0, ... are model parameters. The corresponding density must be normalized to unity over the space covered by the histogram. If the model is multivariate, xe must be an array-like with shape (D, N), where D is the dimension and N is the number of points where the model is evaluated.verbose::Int
: Verbosity level. 0: is no output.grad::Union{Function, Nothing} : Optionally pass the gradient of the
cdf``. Has the same calling signature like the cdf, but must return an array with the shape (K,), where K is the number of parameters. The gradient can be used by Minuit to improve or speed up convergence.use_pdf::Symbol
: Either:none
,:numerical
, or:approximate
. If the model cdf is not available, but the model pdf is, this option can be set to "numerical" or "approximate" to compute the integral of the pdf over the bin. The option "numerical" uses numerical integration, which is accurate but computationally expensive and only supported for 1D histograms. The option "approximate" uses the zero-order approximation of evaluating the pdf at the bin center, multiplied with the bin area. This is fast and works in higher dimensions, but can lead to biased results if the curvature of the pdf inside the bin is significant.names
: Optional names for each parameter of the model (in order). Must have the same length as there are model parameters.
Minuit2.Constant
— MethodConstant(value::Float64; verbose::Int=0)
Constant cost function with fixed value.
Minuit2.CostSum
— MethodCostSum(costs::CostFunction...; verbose::Int=0)
Minuit2.ExtendedBinnedNLL
— Method ExtendedBinnedNLL(bincounts::AbstractArray, binedges::Union{AbstractArray, Tuple}, cdf::Function; use_pdf=:none, verbose=0, grad=nothing, names=())
Binned extended negative log-likelihood.
Use this if shape and normalization of the fitted PDF are of interest and the data is binned. This cost function works with normal and weighted histograms. The histogram can be one- or multi-dimensional.
The cost function works for both weighted data. The cost function assumes that the weights are independent of the data.
The cost function has a minimum value that is asymptotically chi2-distributed. It is constructed from the log-likelihood assuming a poisson distribution and using the saturated model as a reference.
Arguments
bincounts::AbstractArray
: Histogram counts. If this is an array with dimension D, where D is the number of histogram axes.xe::Union{AbstractArray, Tuple}
: Bin edge locations, must be len(n) + 1, where n is the number of bins. If the histogram has more than one axis, xe must be a collection of the bin edge locations along each axis.cdf::Function
: Cumulative density function of the form f(xe, par0, par1, ..., parN), where xe is a bin edge and par0, ... are model parameters. The corresponding density must be normalized to unity over the space covered by the histogram. If the model is multivariate, xe must be an array-like with shape (D, N), where D is the dimension and N is the number of points where the model is evaluated.verbose::Int
: Verbosity level. 0: is no output.grad::Union{Function, Nothing} : Optionally pass the gradient of the
cdf``. Has the same calling signature like the cdf, but must return an array with the shape (K,), where K is the number of parameters. The gradient can be used by Minuit to improve or speed up convergence.use_pdf::Symbol
: Either:none
,:numerical
, or:approximate
. If the model cdf is not available, but the model pdf is, this option can be set to "numerical" or "approximate" to compute the integral of the pdf over the bin. The option "numerical" uses numerical integration, which is accurate but computationally expensive and only supported for 1D histograms. The option "approximate" uses the zero-order approximation of evaluating the pdf at the bin center, multiplied with the bin area. This is fast and works in higher dimensions, but can lead to biased results if the curvature of the pdf inside the bin is significant.names
: Optional names for each parameter of the model (in order). Must have the same length as there are model parameters.
Minuit2.ExtendedUnbinnedNLL
— MethodExtendedUnbinnedNLL(data::AbstractArray, scaled_pdf::Function; log=false, verbose=0, mask=nothing, grad=nothing, names=())
Unbinned extended negative log-likelihood.
Use this if shape and normalization of the fitted PDF are of interest and the original unbinned data is available. The data can be one- or multi-dimensional.
Arguments
data::AbstractArray
: Sample of observations. If the observations are multidimensional, data must have the shape (D, N), where D is the number of dimensions and N the number of data points.scaled_pdf::Function
: Probability density function of the form f(data, par0, [par1, ...]), where data is the sample and par0, ... are model parameters. Must return a tuple (<integral over f in data window>, <f evaluated at data points>). The first value is the density integrated over the data window, the interval that we consider for the fit. For example, if the data are exponentially distributed, but we fit only the interval (0, 5), then the first value is the density integrated from 0 to 5. If the data are multivariate, data passed to f has shape (D,N), where D is the number of dimensions and N the number of data points.verbose::Int
: Verbosity level. 0: is no output.log::Bool=false
: Distributions of the exponential family (normal, exponential, poisson, ...) allow one to compute the logarithm of the pdf directly, which is more accurate and efficient than numerically computing $log(pdf)$. Set this totrue
, if the model returns the logpdf instead of the pdf.mask::Union{Vector{Bool}, BitVector, Nothing}
: Optional mask to select a subset of the data. Must have the same length as data.grad::Union{Function, Nothing}
: Optionally pass the gradient of the pdf. Has the same calling signature like the pdf, but must return an array with the shape (K, N), where N is the number of data points and K is the number of parameters.names
: Optional names for each parameter of the model (in order). Must have the same length as there are model parameters.
Minuit2.LeastSquares
— MethodLeastSquares(x::AbstractArray, y::AbstractVector, yerror, model::Function;
loss=:linear, verbose=0, model_grad=nothing, names=(), mask=nothing)
Least-squares cost function (aka chisquare function).
Use this if you have data of the form (x, y +/- yerror), where x can be one-dimensional or multi-dimensional, but y is always one-dimensional.
Arguments
x::AbstractArray
: Locations where the model is evaluated. If the model is multivariate, x must have shape (D, N), where D is the number of dimensions and N the number of data points.y::AbstractVector
: Observed values. Must have the same length as x.yerror
: Estimated uncertainty of observed values. Must have same shape as y or be a scalar, which is then broadcasted to same shape as y.model::Function
: Function of the form f(x, par0, [par1, ...]) whose output is compared to observed values, where x is the location and par0, ... are model parameters. If the model is multivariate, x has shape (D,), where D is the N the number of data points.loss::Union{Symbol, Function}
: The loss function can be modified to make the fit robust against outliers. Only `:linear
and:soft_l1
are currently implemented, but users can pass any loss function as this argument. It should be a monotonic, twice differentiable function, which accepts the squared residual and returns a modified squared residual.verbose::Int
: Verbosity level. 0: is no output.model_grad::Union{Function, Nothing}
: Optionally pass the gradient of the model. Has the same calling signature like the model, but must return an array with the shape (K,), where K is the number of parameters. The gradient can be used by Minuit to improve or speed up convergence.names
: Optional names for each parameter of the model (in order). Must have the same length as there are model parameters.mask::Union{Vector{Bool}, BitVector, Nothing}
: Optional mask to select a subset of the data. Must have the same length as x.
Notes
Alternative loss functions make the fit more robust against outliers by weakening the pull of outliers. The mechanical analog of a least-squares fit is a system with attractive forces. The loss function can be modified to make the fit robust against outliers.
Minuit2.Minuit
— TypeMinuit structure
Direct or calculated fields
funcname::String
: Name of the functionx0::AbstractVector
: Initial parameters valuesmethod::Symbol
: The minimization algorithm to use. Possible values are:migrad
,:simplex
tolerance::Real
: Tolerance for the minimization. If set to 0, Minuit will use a default value.precision::Union{Real,Nothing}
: The precision for the minimizationstrategy::Int
: The strategy for the minimization (0,1(default),2). See the manual for details.values
: The values of the parameters at the minimumerrors
: The errors of the parameters at the minimumfixed
: The fixed status of the parameterslimits
: The limits of the parametersis_valid
: If the minimization was successfulfval
: The function value at the minimumedm
: The estimated distance to minimumnfcn
: The number of function callsngrad
: The number of gradient callsniter
: The number of iterationsnpars
: The number of parametersndof
: Number of degrees of freedomcovariance
: The covariance matrix of the parametersis_above_max_edm
: If the estimated distance to minimum is above the maximumhas_parameters_at_limit
: If any of the parameters are at the limitshas_accurate_covar
: If the covariance matrix is accuratehas_posdef_covar
: If the covariance matrix is positive definitehas_made_posdef_covar
: If the covariance matrix was made positive definitehesse_failed
: If the Hesse algorithm failedhas_covariance
: If the covariance matrix is availablecovariance
: The covariance matrix of the parametershas_accurate_covar
: If the covariance matrix is accuratehas_valid_parameters
: If the parameters are validhas_reached_call_limit
: If the maximum number of function calls was reachedminos
: The Minos errorsparameters
: The parameters values and errors
Minuit2.Minuit
— MethodMinuit(fcn, x0...; grad=nothing, error=(), errordef=1.0, names=(), limits=(), fixed=(), method=:migrad, maxfcn=0,
tolerance=0.1, precision=nothing, strategy=1, kwargs...)
Initialize a Minuit object.
This does not start the minimization or perform any other work yet. Algorithms are started by calling the corresponding methods.
Arguments
fcn::Union{Function,CostFunction}
: Function to minimize. See notes for details on what kind of functions are accepted.x0::AbstractArray
: Starting values for the minimization. See notes for details on how to set starting values.grad::Union{Function,Nothing, Bool}
: Ifgrad
is a function, it must be a function that calculates the gradient and returns an iterable object with one entry for each parameter, which is the derivative offcn
for that parameter. Ifnothing
(default), Minuit will use the gradient of the provided cost function. If it does not exists Minuit will compute the gradient numerically. Ifgrad
isfalse
, Minuit will not use the gradient.error::AbstractArray
: Starting values for the errors of the parameters. If not provided, Minuit will use 0.1 for all parameters.errordef::Real
: Error definition of the function. Minuit defines parameter errors as the change in parameter value required to change the function value byerrordef
. Normally, for chisquared fits it is 1, and for negative log likelihood, its value is 0.5. If the user wants instead the 2-sigma errors for chisquared fits, it becomes 4, asChi2(x+n*sigma) = Chi2(x) + n*n
.names::Vector{String}
: Names of the parameters. If not provided, Minuit will try to extract the names from the function signature.limits::AbstractArray
: Limits for the parameters. If not provided, Minuit will use no limits for all parameters.fixed::AbstractArray
: Fixed status of the parameters. If not provided, Minuit will usefalse
for all parameters. The limits are a tuple of two values, the lower and upper limit. If the parameter is fixed, the limits are set to(-Inf, Inf)
.method::Symbol
: The minimization algorithm to use. Possible values are:migrad
,:simplex
maxfcn::Int
: Maximum number of function calls. If set to 0, Minuit will use a default value.tolerance::Real
: Tolerance for the minimization. If set to 0, Minuit will use a default value.arraycall::Bool
: If the function takes a single argument which can be an array or takes multiple arguments. If not provided, Minuit will try to detect it via reflection. Seeget_nargs
.kwargs
: Additional keyword arguments. Starting values for the minimization as keyword arguments. See notes for details on how to set starting values.
Notes
Function to minimize
By default, Minuit assumes that the callable fcn
behaves like chi-square function, meaning that the function minimum in repeated identical random experiments is chi-square distributed up to an arbitrary additive constant. This is important for the correct error calculation. If fcn
returns a log-likelihood, one should multiply the result with -2 to adapt it. If the function returns the negated log-likelihood, one can alternatively set the attribute errordef
to make Minuit calculate errors properly.
Minuit reads the function signature of fcn
to detect the number and names of the function parameters. Two kinds of function signatures are understood.
a. Function with positional arguments.
The function has positional arguments, one for each fit parameter. Example:
fcn(a, b, c) = ...
The parameters a, b, c must accept a real number Minuit automatically detects the parameters names in this case.
b. Function with arguments passed as a single AbstractArray.
The function has a single argument which is an AbstractArray. Example:
function fcn_v(x) = ...
To use this form, starting values (x0
) needs to be passed to Minuit in form of a Tuple
or Vector
. In some cases, the detection may fail, and will be necessary to use the names
keyword to set the parameter names.
Parameter initialization
Initial values for the minimization can be set with positional arguments or via keywords. This is best explained through an example:
fcn(x, y) = (x - 2)^2 + (y - 3)^2
The following ways of passing starting values are equivalent:
Minuit(fcn, x=1, y=2)
Minuit(fcn, y=2, x=1) # order is irrelevant when keywords are used ...
Minuit(fcn, [1,2]) # ... but here the order matters
Positional arguments can also be used if the function has no signature:
fcn_no_sig(args...) = ...
Minuit(fcn_no_sig, [1,2])
If the arguments are explicitly named with the names
keyword described further below, keywords can be used for initialization:
Minuit(fcn_no_sig, x=1, y=2, names=("x", "y")) # this also works
If the function accepts a single AbstractVector, then the initial values must be passed as a single array-like object:
fcn_v(x) = return (x[1] - 2) ** 2 + (x[2] - 3) ** 2
Minuit(fcn_v, (1, 2))
Setting the values with keywords is not possible in this case. Minuit deduces the number of parameters from the length of the initialization sequence.
Minuit2.UnbinnedNLL
— MethodUnbinnedNLL(data::AbstractArray, pdf::Function; log=false, verbose=0, mask=nothing, grad=nothing, names=())
Unbinned negative log-likelihood cost function.
Arguments
data::AbstractArray
: Sample of observations. If the observations are multidimensional, data must have the shape (D, N), where D is the number of dimensions and N the number of data points.pdf::Function
: Probability density function of the form f(data, par0, [par1, ...]), where data is the data sample and par0, ... are model parameters. If the data are multivariate, data passed to f has shape (D,), where D is the number of dimensions and N the number of data points.verbose::Int
: Verbosity level. 0: is no output (default). 1: print current args and negative log-likelihood value.log::Bool=false
: Distributions of the exponential family (normal, exponential, poisson, ...) allow one to compute the logarithm of the pdf directly, which is more accurate and efficient than numerically computing $log(pdf)$. Set this totrue
, if the model returns the logpdf instead of the pdf.mask::Union{Vector{Bool}, BitVector, Nothing}
: Optional mask to select a subset of the data. Must have the same length as data.grad::Union{Function, Nothing}
: Optionally pass the gradient of the pdf. Has the same calling signature like the pdf, but must return an array with the shape (K, N), where N is the number of data points and K is the number of parameters. Iflog
is True, the function must return the gradient of the logpdf instead of the pdf. The gradient can be used by Minuit to improve or speed up convergence and to compute the sandwich estimator for the variance of the parameter estimates.names
: Optional names for each parameter of the model (in order). Must have the same length as there are model parameters.
Returns
- Cost function object.
Minuit2.RooFit.AbstractPdf
— TypeAbstractPdf
Abstract type for all RooFit distributions.
Minuit2.RooFit.AddPdf
— MethodAddPdf(name, pdfs, fractions)
Construct an AddPdf distribution with a name, a vector of AbstractPdf pdfs and a vector of RealVar fractions.
Arguments
name::Symbol
: Name of the distribution.pdfs::Vector{<:AbstractPdf}
: Vector of AbstractPdf pdfs.fractions::Vector{<:RealVar}
: Vector of RealVar fractions or coefficients.- If the number of pdfs is equal to the number of fractions, the fractions represent absolute populations.
- If the number of pdfs is equal to the number of fractions + 1, the fractions represent recursive fractions.
Returns
- An AddPdf distribution.
Minuit2.RooFit.AddPdf
— MethodAddPdf(name, pdf1::AbstractPdf, pdf2::AbstractPdf, fraction::RealVar)
Construct an AddPdf distribution with a name, two AbstractPdf pdfs and a RealVar fraction.
Minuit2.RooFit.ArgusPdf
— TypeArgusPdf(name, m, m₀, c, p)
ArgusPdf describes the ARGUS background shape.
$\mathrm{Argus}(m, m_0, c, p) = \mathcal{N} \cdot m \cdot \left[ 1 - \left( \frac{m}{m_0} \right)^2 \right]^p \cdot \exp\left[ c \cdot \left(1 - \left(\frac{m}{m_0}\right)^2 \right) \right]$
Minuit2.RooFit.Chebyshev
— MethodChebyshev(name, x, coeffs)
Construct a Chebyshev distribution with a name, a RealVar x and a vector of coefficients.
Minuit2.RooFit.DataSet
— TypeDataSet
DataSet is a container for a set of data points.
Minuit2.RooFit.Exponential
— MethodExponential(name, x, c)
Construct an Exponential distribution with a name
, a RealVar x
and a RealVar c
.
$\mathrm{Exponential}(x, c) = \mathcal{N} \cdot \exp(c\cdot x)$, where $\mathcal{N}$ is a normalization constant that depends on the range and values of the arguments
Minuit2.RooFit.FitResult
— TypeFitResult
FitResult is a container for the result of a fit.
Minuit2.RooFit.Gaussian
— MethodGaussian(name, x, μ, σ)
Construct a Gaussian distribution with a name, a RealVar x, a RealVar μ and a RealVar σ.
Minuit2.RooFit.RealVar
— TypeRealVar{T<:Real}
Entity to represent a real variable with a name, value, limits and number of bins.
Minuit2.RooFit.RealVar
— MethodRealVar(name, value=0.; limits=(-Inf,Inf), nbins=0)
Construct a RealVar with a name, value, limits and number of bins.
Functions
Minuit2.FCN
— FunctionFCN(cost::CostFunction, grad=true)
Create a JuliaFcn object from a CostFunction.
Arguments
fnc::CostFunction
: The CostFunction to minimize.grad::Bool=true
: Iftrue
, the gradient of the cost function is used. Iffalse
, the gradient is not used.
Returns
JuliaFcn
: A JuliaFcn object inheriting from the abstract C++ classMinuit::FCNBase
that can be used in Minuit.
Minuit2.FCN
— FunctionFCN(fnc, grad=nothing, arraycall=false, errordef=1.0)
Create a JuliaFcn object from a Julia function fnc
and its gradient grad
.
Arguments
fnc::Function
: The Julia function to minimize. It can either accept a set of discrete arguments or a single argument of typeAbstractVector
. This is decided in conjunction with the argumentarraycall
.grad::Function=nothing
; Gradient Julia function. The input arguments follow the same as forfcn
and it returns aVector
, of length the number of parameters, with the gradients.arraycall::Bool=false
: Iftrue
, the functionfcn
accepts a single argument of typeAbstractVector
. Iffalse
, the function accepts a set of discrete arguments.errordef::Real=1.0
: Error definition of the function. Minuit defines parameter errors as the change in parameter value required to change the function value byerrordef
. Normally, for chisquared fits it is 1, and for negative log likelihood, its value is 0.5.
Returns
JuliaFcn
: A JuliaFcn object inheriting from the abstract C++ classMinuit::FCNBase
that can be used in Minuit.
Usage
fcn(x, y) = (x - 2)^2 + (y - 3)^2
grad(x, y) = [2*(x - 2), 2*(y - 3)]
jf = FCN(fcn, grad)
jf(1.0, 1.0) # returns 5.0
jf.grad(1.0, 1.0) # returns [-2.0, -4.0]
jf.nfcn # returns the number of function calls
jf.ngrad # returns the number of gradient calls
jf.has_gradient # returns true
Minuit2.chi2
— Methodchi2(y, ye, ym)
Compute (potentially) chi2-distributed cost.
The value returned by this function is chi2-distributed, if the observed values are normally distributed around the expected values with the provided standard deviations.
Arguments
y
: Observed values.ye
: Uncertainties of values.ym
: Expected values.
Returns
- Value of cost function.
Minuit2.chi2_grad
— Methodchi2_grad(y, ye, ym, gym)
Compute gradient of function chi2
.
Arguments
y
: Observed values.ye
: Uncertainties of values.ym
: Expected values.gym
: Gradient of ym with respect to K model parameters.
Returns
- Gradient of cost function with respect to model parameters.
Minuit2.contour
— Methodcontour(m::Minuit, x, y; size=50, bound=2, grid=nothing, subtract_min=false)
Get a 2D contour of the function around the minimum.
It computes the contour via a function scan over two parameters, while keeping all other parameters fixed. The related :meth:mncontour
works differently: for each pair of parameter values in the scan, it minimises the function with the respect to all other parameters.
This method is useful to inspect the function near the minimum to detect issues (the contours should look smooth). It is not a confidence region unless the function only has two parameters. Use :meth:mncontour
to compute confidence regions.
Arguments
x
: First parameter for scan (name or index).y
`: Second parameter for scan (name or index).size=50
: Number of scanning points per parameter (Default: 50). It can be tuple(nx,ny)
as the number of scanning points per parameter. Ignored ifgrid
is set.bound=2
: Either ((v1min,v1max),(v2min,v2max)) or the number ofsigma
s to scan symmetrically from minimum.grid::Tuple{AbstractVector,AbstractVector} : Grid points to scan over. If
grid$is set, `size$ andbound
are ignored.subtract_min::Bool=false
: Subtract minimum from return values
Returns
- Tuple(
xv
,yv
,zv
) : Tuple of 1D arrays with the x and y values and a 2D array with the function values.
Minuit2.get_nargs
— Methodget_nargs(f)::Int
Returns the number of arguments of the function f
. The first argument is the function itself, so the number of arguments is length(m.sig.parameters)-1
.
Minuit2.hesse!
— Methodhesse!(m::Minuit, ncall::Int=0)
Run Hesse algorithm to compute asymptotic errors.
The Hesse method estimates the covariance matrix by inverting the matrix of second derivatives (Hesse matrix) at the minimum. To get parameters correlations, you need to use this. The Minos algorithm is another way to estimate parameter uncertainties, see function minos
.
Arguments
ncall::Int=0
: Approximate upper limit for the number of calls made by the Hesse algorithm. If set to 0, use the adaptive heuristic from the Minuit2 library.
Notes
The covariance matrix is asymptotically (in large samples) valid. By valid we mean that confidence intervals constructed from the errors contain the true value with a well-known coverage probability (68 % for each interval). In finite samples, this is likely to be true if your cost function looks like a hyperparabola around the minimum.
In practice, the errors very likely have correct coverage if the results from Minos and Hesse methods agree. It is possible to construct artificial functions where this rule is violated, but in practice it should always work.
Minuit2.matrix
— Methodmatrix(m::Minuit; correlation=false)
Get the covariance matrix of the parameters.
Minuit2.migrad!
— Functionmigrad!(m::Minuit, strategy=1; ncall=0)
Run Migrad minimization.
Migrad from the Minuit2 library is a robust minimisation algorithm which earned its reputation in 40+ years of almost exclusive usage in high-energy physics. How Migrad works is described in the Minuit paper. It uses first and approximate second derivatives to achieve quadratic convergence near the minimum.
Parameters
m::Minuit
: The Minuit object to minimize.strategy::Int
: The minimization strategy. The default value is 1, which is the recommended value for most cases. The value 0 is faster, but less reliable. The value 2 is slower, but more reliable. The value 3 or higher is slower, but even more reliable.
Keyword Parameters
ncall::Int=0
: Approximate upper limit for the number of calls. If set to 0, use the adaptive heuristic from the Minuit2 library.
The following two parameters controls the behavior of robust_low_level_fit
, which essentially tries to restart Migrad procedure when it fails to converge after the first iteration:
iterate::Int=5
: Number of iterations to run the minimization. Default is 5.use_simplex::Bool=true
: Iftrue
, use the simplex algorithm in retry to find the minimum. Iffalse
, use the migrad algorithm.
Minuit2.minos!
— Methodminos!(m::Minuit, name::String)
Run Minos algorithm to compute asymmetric errors for a single parameter.
The Minos algorithm uses the profile likelihood method to compute (generally asymmetric) confidence intervals. It scans the negative log-likelihood or (equivalently) the least-squares cost function around the minimum to construct a confidence interval.
Arguments
m::Minuit
: The Minuit object to minimize.parameters::AbstractVector{String}
: Names of the parameters to compute the Minos errors for.cl::Number
: Confidence level of the interval. If not set, a standard 68 % interval is computed (default). If 0 < cl < 1, the value is interpreted as the confidence level (a probability). For convenience, values cl >= 1 are interpreted as the probability content of a central symmetric interval covering that many standard deviations of a normal distribution. For example, cl=1 is interpreted as 68.3 %, and cl=2 is 84.3 %, and so on. Using values other than 0.68, 0.9, 0.95, 0.99, 1, 2, 3, 4, 5 require the scipy module.ncall::Int
: Limit the number of calls made by Minos. If 0, an adaptive internal heuristic of the Minuit2 library is used (Default: 0).
Notes
Asymptotically (large samples), the Minos interval has a coverage probability equal to the given confidence level. The coverage probability is the probability for the interval to contain the true value in repeated identical experiments.
The interval is invariant to transformations and thus not distorted by parameter limits, unless the limits intersect with the confidence interval. As a rule-of-thumb: when the confidence intervals computed with the Hesse and Minos algorithms differ strongly, the Minos intervals are preferred. Otherwise, Hesse intervals are preferred.
Running Minos is computationally expensive when there are many fit parameters. Effectively, it scans over one parameter in small steps and runs a full minimisation for all other parameters of the cost function for each scan point. This requires many more function evaluations than running the Hesse algorithm.
Minuit2.mncontour
— Methodnmcontour(x, y; cl=0.68, size=50, interpolated=0, ncall=0, iterate=5, use_simplex=true)
Get 2D Minos confidence region.
This scans over two parameters and minimises all other free parameters for each scan point. This scan produces a statistical confidence region according to the profile likelihood method with a confidence level cl
, which is asymptotically equal to the coverage probability of the confidence region according to Wilks' theorem. Note that 1D projections of the 2D confidence region are larger than 1D Minos intervals computed for the same confidence level. This is not an error, but a consequence of Wilks'theorem.
The calculation is expensive since a numerical minimisation has to be performed at various points.
Arguments
x
: Variable name of the first parameter.y
: Variable name of the second parameter.cl::Real=0.68
: Confidence level of the contour. If not set a standard 68 % contour is computed (default). If 0 < cl < 1, the value is interpreted as the confidence level (a probability). For convenience, values cl >= 1 are interpreted as the probability content of a central symmetric interval covering that many standard deviations of a normal distribution.size::Int=50
: Number of points on the contour to find. Increasing this makes the contour smoother, but requires more computation time.interpolated::Int=0
: Number of interpolated points on the contour. If you set this to a value larger than size, cubic spline interpolation is used to generate a smoother curve and the interpolated coordinates are returned. Values smaller than size are ignored. Good results can be obtained with size=20, interpolated=200.
Returns
contour::Vector(Tuple{Float64,Float64})
: Contour points of the form [(x1, y1)...(xn, yn)].
Minuit2.mnprofile
— Methodmnprofile(m::Minuit, var; size=30, bound=2, grid=nothing, subtract_min=false,
ncall=0, iterate=5, use_simplex=true)
Get Minos profile over a specified interval.
Scans over one parameter and minimises the function with respect to all other parameters for each scan point.
Arguments
var
: Parameter to scan over.size=30
: Number of scanning points. Ignored if grid is set.bound=2
: If bound is a tuple, (left, right) scanning bound, or the number of sigmas to scan symmetrically around the minimum. Ignored if grid is set.grid
: Parameter values on which to compute the profile. Ifgrid
is set,size
andbound
are ignored.subtract_min=false
: If true, subtract offset so that smallest value is zero.ncall=0
: Approximate maximum number of calls before minimization will be aborted. If set to 0, use the adaptive heuristic from the Minuit2 library. Note: The limit may be slightly violated, because the condition is checked only after a full iteration of the algorithm, which usually performs several function calls. iterate : int, optional Automatically call Migrad up to N times if convergence was not reached (Default: 5). This simple heuristic makes Migrad converge more often even if the numerical precision of the cost function is low. Setting this to 1 disables the feature. use_simplex: bool, optional If we have to iterate, set this to True to call the Simplex algorithm before each call to Migrad (Default: True). This may improve convergence in pathological cases (which we are in when we have to iterate).
Returns
- Tuple(
x
,y
,ok
) : Tuple of 1D arrays with the parameter values, function values and booleans whether the fit succeeded or not.
Minuit2.multinomial_chi2
— Methodmultinomial_chi2(n, mu)
Compute asymptotically chi2-distributed cost for multinomially-distributed data. See Baker & Cousins, NIM 221 (1984) 437-442.
Arguments
n
: Observed counts.mu
Expected counts. Must satisfy sum(mu) == sum(n).
Returns
- Cost function value.
Notes
The implementation makes the result asymptotically chi2-distributed, which helps to maximise the numerical accuracy for Minuit.
Minuit2.multinomial_chi2_grad
— Methodmultinomial_chi2_grad(n, mu, gmu)
Compute gradient of function multinomial_chi2
.
Minuit2.poisson_chi2
— Methodpoisson_chi2(n, mu)
Compute asymptotically chi2-distributed cost for Poisson-distributed data. See Baker & Cousins, NIM 221 (1984) 437-442.
Arguments
n
: Observed counts.mu
: Expected counts per bin.
Returns
- Cost function value.
Notes
The implementation makes the result asymptotically chi2-distributed, which helps to maximise the numerical accuracy for Minuit.
If sum(mu) == sum(n), the result is equal to multinomial_chi2
.
Minuit2.poisson_chi2_grad
— Methodpoisson_chi2_grad(n, mu, gmu)
Compute gradient of function poisson_chi2
.
Minuit2.profile
— Methodprofile(m::Minuit, var; size=100, bound=2, grid=nothing, subtract_min=false)
Calculate 1D cost function profile over a range.
A 1D scan of the cost function around the minimum, useful to inspect the minimum. For a fit with several free parameters this is not the same as the Minos profile computed by mncontour
.
Arguments
m::Minuit
: The Minuit object to minimize.var
: The parameter to scan over (name or index).size=100
: Number of scanning points. Ignored ifgrid
is set.bound=2
: Number ofsigma
s to scan symmetrically around the minimum. Ignored ifgrid
is set.grid::AbstractVector
: Grid points to scan over. Ifgrid
is set,size
andbound
are ignored.subtract_min::Bool=false
: Subtract minimum from return values.
Returns
- Tuple(
x
,y
) : Tuple of 1D arrays with the parameter values and the function values.
Minuit2.robust_low_level_fit
— Methodrobust_low_level_fit(fcn, state, ncall, strategy, tolerance, precision, iterate, use_simplex)
A meta algorithm that:
- runs Migrad with user-specified configs (such as strategy and tolerance)
- checks if we converged
- if not, sets strategy=2, and start running N more iterations
- in each additional iteration, if
use_simplex
, runs SIMPLEX (again, think NelderMead) first, followed by a Migrad
Minuit2.scan!
— Functionscan!(m::Minuit, maxfcn = 0, strategy=1)
Run Scan algorithm to find the minimum. Thus is a brute force algorithm that scans the function in a hypercube around the initial values. The values must be within the limits.
Arguments
m::Minuit
: The Minuit object to minimize.maxfcn::Int=0
: Maximum number of function calls. If set to 0, Minuit will use a default value.strategy::Int=1
: The minimization strategy. The default value is 1, which is the recommended value for most cases. The value 0 is faster, but less reliable. The value 2 is slower, but more reliable. The value 3 or higher is slower, but even more reliable.
Minuit2.simplex!
— Functionsimplex!(m::Minuit, strategy=1)
Run Simplex minimization.
Parameters
m::Minuit
: The Minuit object to minimize.strategy::Int
: The minimization strategy. The default value is 1, which is the recommended value for most cases. The value 0 is faster, but less reliable. The value 2 is slower, but more reliable. The value 3 or higher is slower, but even more reliable.
Minuit2.soft_l1_cost
— Methodsoft_l1_cost(y, ye, ym)
Minuit2.soft_l1_cost_grad
— Methodsoft_l1_cost_grad(y, ye, ym, gym)
Minuit2.RooFit.ConstVar
— MethodConstVar(name, value=0.)
Construct a ConstVar with a name ans value.
Minuit2.RooFit.fitTo
— MethodfitTo(d::AbstractPdf, data)
Fit a distribution d
to a data set data
.
Minuit2.RooFit.generate
— Methodgenerate(d::AbstractPdf, n::Integer=1000; nbins=0)
Generate n
random numbers from a distribution d
.