Fork me on GitHub

Full API Documentation

Quick Summary

These are the things you will use a lot.

util.describe(f[, verbose]) try to extract function arguements name
Minuit Minuit(fcn, throw_nan=False, pedantic=True, frontend=None, forced_parameters=None, print_level=1, errordef=None, **kwds)
Minuit.migrad(self, int ncall=10000[, resume]) Run migrad.
Minuit.minos(self[, var, sigma]) Run minos for parameter var
Minuit.values values: object
Minuit.args args: object
Minuit.errors errors: object
Minuit.get_merrors(self) Returns a dictionary of varname-> MinosError Struct
Minuit.fval fval: ‘double’
Minuit.fitarg fitarg: object
Minuit.mnprofile(self, vname[, bins, bound, ...]) calculate minos profile around specify range ie. Migrad results
Minuit.draw_mnprofile(self, vname[, bins, ...]) Draw minos profile around specify range ie.
Minuit.mncontour(self, x, y, int numpoints=20) Minos contour scan.
Minuit.mncontour_grid(self, x, y[, bins, ...]) compute gridded minos contour.
Minuit.draw_mncontour(self, x, y[, bins, ...]) draw minos contour

Minuit

class iminuit.Minuit

Minuit(fcn, throw_nan=False, pedantic=True, frontend=None, forced_parameters=None, print_level=1, errordef=None, **kwds)

Construct minuit object from given fcn

Arguments:

  • fcn: function to optimized. Minuit automagically find how to call the function and each parameters. More information about how Minuit detects function signature can be found in Function Signature Extraction Ordering

Builtin Keyword Arguments:

  • throw_nan: fcn can be set to raise RuntimeError when it encounters nan. (Default False)

  • pedantic: warns about parameters that do not have initial value or initial error/stepsize set.

  • frontend: Minuit frontend. There are two builtin frontends.

    1. ConsoleFrontend which is design to print out to terminal.
    2. HtmlFrontend which is designed to give a nice output in IPython notebook session.

    By Default, Minuit switch to HtmlFrontend automatically if it is called in IPython session. It uses ConsoleFrontend otherwise.

  • forced_parameters: tell Minuit not to do function signature detection and use this argument instead. (Default None (automagically detect signature)

  • print_level: set the print_level for this Minuit. 0 is quiet. 1 print out at the end of migrad/hesse/minos. The reason it has this cAmEl case is to keep it compatible with PyMinuit.

  • errordef: Optionals. Amount of increase in fcn to be defined as 1 \(\sigma\). If None is given, it will look at fcn.default_errordef(). If fcn.default_errordef() is not defined or not callable iminuit will give a warning and set errordef to 1. Default None(which means errordef=1 with a warning).

Parameter Keyword Arguments:

Similar to PyMinuit. iminuit allows user to set initial value, initial stepsize/error, limits of parameters and whether parameter should be fixed or not by passing keyword arguments to Minuit. This is best explained through an example:

def f(x,y):
    return (x-2)**2 + (y-3)**2
  • Initial value(varname):

    #initial value for x and y
    m = Minuit(f, x=1, y=2)
    
  • Initial step size/error(fix_varname):

    #initial step size for x and y
    m = Minuit(f, error_x=0.5, error_y=0.5)
    
  • Limits (limit_varname=tuple):

    #limits x and y
    m = Minuit(f, limit_x=(-10,10), limit_y=(-20,20))
    
  • Fixing parameters:

    #fix x but vary y
    m = Minuit(f, fix_x=True)
    

Note

Tips: you can use python dictionary expansion to programatically change fitting argument.

kwdarg = dict(x=1., error_x=0.5)
m = Minuit(f, **kwdarg)

You can obtain also obtain fit arguments from Minuit object to reuse it later too. Note that fitarg will be automatically updated to minimum value and corresponding error when you ran migrad/hesse:

m = Minuit(f, x=1, error_x=0.5)
my_fitarg = m.fitarg
another_fit = Minuit(f, **my_fitarg)
migrad(self, int ncall=10000, resume=True, int nsplit=1)

Run migrad.

Migrad is an age-tested(over 40 years old, no kidding), super robust and stable minimization algorithm. It even has wiki page. You can read how it does the magic at here.

Arguments:

  • ncall: integer (approximate) maximum number of call before migrad stop trying. Default 10000
  • resume: boolean indicating whether migrad should resume from previous minimizer attempt(True) or should start from the beginning(False). Default True.
  • split: split migrad in to split runs. Max fcn call for each run is ncall/nsplit. Migrad stops when it found the function minimum to be valid or ncall is reached. This is useful for getting progress. However, you need to make sure that ncall/nsplit is large enough. Otherwise, migrad will think that the minimum is invalid due to exceeding max call (ncall/nsplit). Default 1(no split).

Return:

hesse(self)

Run HESSE.

HESSE estimates error matrix by the second derivative at the minimim. This error matrix is good if your \(\chi^2\) or likelihood profile is parabolic at the minimum. From my experience, most of simple fits are.

minos() makes no parabolic assumption and scan the likelihood and give the correct error asymmetric error in all cases(Unless your likelihood profile is utterly discontinuous near the minimum). But it is much more computationally expensive.

Returns
list of Minuit Parameter Struct
minos(self, var=None, sigma=1.0, unsigned int maxcall=1000)

Run minos for parameter var

If var is None it runs minos for all parameters

Arguments:

  • var: optional variable name. Default None.(run minos for every variable)
  • sigma: number of \(\sigma\) error. Default 1.0.

Returns

Dictionary of varname to Minos Error Struct if minos is requested for all parameters.
args

args: object

current values of parameters in form of tuple(var1,var2, ...)

values

values: object

curret values of parameters in form of dictionary(name -> value)

errors

errors: object

curret values of parameters in form of dictionary(name -> error)

fitarg

fitarg: object

current value state in form of dictionary(name->value, error_name->error, fix_name->fix , limit_name->(lower, upper)). This is very useful when you wan to save the fit parameters and use it later. For example,:

m = Minuit(f,x=1)
m.migrad()
fitarg = m.fitarg

m2 = Minuit(f,**fitarg)
merrors

merrors: object

last minos error in form of dictionary((name,1.0)->upper, (name,-1.0) ->lower). It is here for PyMinuit compatibility.

This is not recommend use get_merrors() which return dictionary of name-> Minos Error Struct instead.

fval

fval: ‘double’

current function value

See also

get_fmin()

edm

edm: ‘double’

estimated distance to minimum

See also

get_fmin()

covariance

covariance: object

covariance matrix in form of dictionary(v1,v2)->covariance

See also

matrix()

gcc

gcc: object

global correlation coefficiencts in form of dict(varname->gcc)

errordef

errordef: ‘double’

Amount of change in fcn that defines 1 \(sigma\) error. Default value is 1.0. errordef should be 1.0 for \(\chi^2\) cost funciton and 0.5 for negative log likelihood function.

tol

tol: ‘double’

Tolerance. One fo migrad convergence criteria is edm < maxedm. Maxedm is calculated by 0.0001*tol*UP.

contour(self, x, y, bins=20, bound=2, args=None, subtract_min=False)

2D countour scan.

The contour returns is obtained by fixing all others parameters and varying x and y.

Arguments

  • x variable name for X axis of scan

  • y variable name for Y axis of scan

  • bound

    If bound is 2x2 array [[v1min,v1max],[v2min,v2max]]. If bound is a number, it specifies how many \(\sigma\) symmetrically from minimum (minimum+- bound*:math:sigma). Default 2

  • subtract_min subtract_minimum off from return value. This

    makes it easy to label confidence interval. Default False.

Returns

x_bins, y_bins, values

values[y, x] <– this choice is so that you can pass it to through matplotlib contour()

See also

mncontour()

Note

If subtract_min=True, the return value has the minimum subtracted off. The value on the contour can be interpreted loosely as \(i^2 \times \textrm{up}\) where i is number of standard deviation away from the fitted value WITHOUT taking into account correlation with other parameters that’s fixed.

draw_contour(self, x, y, bins=20, bound=2, args=None, show_sigma=False)

Convenient wrapper for drawing contour. The argument is the same as contour(). If show_sigma=True`(Default), the label on the contour lines will show how many :math:sigma` away from the optimal value instead of raw value.

Note

Like contour(). The error shown on the plot is not strictly 1 \(\sigma\) contour since the other parameters are fixed.

draw_mncontour(self, x, y, bins=100, nsigma=2, numpoints=20, sigma_res=4)

draw minos contour

Arguments

  • x, y parameter name
  • bins number of bin in contour grid.
  • nsigma number of sigma contour to draw
  • numpoints number of points to calculate for each contour
  • sigma_res number of sigma level to calculate MnContours. Default 4.

returns

x, y, gridvalue, contour

gridvalue is interorlated nsigma

draw_mnprofile(self, vname, bins=30, bound=2, subtract_min=False, band=True, text=True)

Draw minos profile around specify range ie. Migrad results with vname fixed at various places within bound

Arguments

  • vname variable name to scan

  • bins number of scanning bin. Default 30.

  • bound

    If bound is tuple, (left, right) scanning bound. If bound is a number, it specifies how many \(\sigma\) symmetrically from minimum (minimum+- bound* \(\sigma\)). Default 2

  • subtract_min subtract_minimum off from return value. This

    makes it easy to label confidence interval. Default False.

  • band show green band to indicate the increase of fcn by errordef. Default True.

  • text show text for the location where the fcn is increased by errordef. This is less accurate than minos(). Default True.

Returns

bins(center point), value, migrad results
from iminuit import Minuit


def f(x, y, z):
    return (x-1)**2 + (y-x)**2 + (z-2)**2

m = Minuit(f, print_level=0, pedantic=False)
m.migrad()
m.draw_mnprofile('y')

(Source code, png, hires.png, pdf)

_images/draw_mnprofile.png
draw_profile(self, vname, bins=100, bound=2, args=None, subtract_min=False, band=True, text=True)

A convenient wrapper for drawing profile using matplotlib

Note

This is not a real minos profile. It’s just a simple 1D scan. The number shown on the plot is taken from the green band. They are not minos error. To get a real minos profile call mnprofile() or draw_mnprofile()

Arguments

In addition to argument listed on profile(). draw_profile take these addition argument:

  • band show green band to indicate the increase of fcn by errordef. Note again that this is NOT minos error in general. Default True.
  • text show text for the location where the fcn is increased by errordef. This is less accurate than minos() Note again that this is NOT minos error in general. Default True.
get_fmin(self)

return current FunctionMinimum Struct

get_initial_param_state(self)

get initiail setting inform of MinuitParameter Struct

get_merrors(self)

Returns a dictionary of varname-> MinosError Struct

get_num_call_fcn(self)

return number of total call to fcn(not just the last operation)

get_param_states(self)

Return a list of current MinuitParameter Struct for all parameters

is_clean_state(self)

check if minuit is at clean state ie. no migrad call

is_fixed(self, vname)

check if variable vname is (initialy) fixed

latex_initial_param(self)

build iminuit.latex.LatexTable for initial parameter

latex_matrix(self)

Build LatexFactory object that contains correlation matrix

latex_param(self)

build iminuit.latex.LatexTable for current parameter

list_of_fixed_param(self)

return list of (initially) fixed parameters

list_of_vary_param(self)

return list of (initially) float vary parameters

matrix(self, correlation=False, skip_fixed=True)

return error/correlation matrix in tuple or tuple format.

matrix_accurate(self)

check if covariance(of the last migrad) is accurate.

merrors

merrors: object

merrors_struct

merrors_struct: object

migrad_ok(self)

check if minimum is valid

mncontour(self, x, y, int numpoints=20, sigma=1.0)

Minos contour scan. A proper n sigma contour scan. This is line where the minimum of fcn with x,y is fixed at points on the line and letting the rest of variable varied is change by sigma * errordef^2 . The calculation is very very expensive since it has to run migrad at various points.

Arguments

  • x string variable name of the first parameter
  • y string variable name of the second parameter
  • numpoints number of points on the line to find. Default 20.
  • sigma number of sigma for the contour line. Default 1.0.

Returns

x minos error struct, y minos error struct, contour line

contour line is a list of the form [[x1,y1]...[xn,yn]]

mncontour_grid(self, x, y, bins=100, nsigma=2, numpoints=20, int sigma_res=4, edges=False)

compute gridded minos contour.

Arguments

  • x,**y** parameter name
  • bins number of bins in the grid. The boundary of the grid is selected automatically by the minos error computed. Default 100.
  • nsigma number of sigma to draw. Default 2
  • numpoints number of points to calculate mncontour for each sigma points(there are sigma_res*nsigma total)
  • sigma_res number of sigma level to calculate MnContours
  • edges Return bin edges instead of mid value(pass True if you want to draw it using pcolormesh)

Returns

xgrid, ygrid, sigma, rawdata

rawdata is tuple of (x,y,sigma_level)

See also

draw_mncontour()

from iminuit import Minuit


def f(x, y, z):
    return (x-1)**2 + (y-x)**2 + (z-2)**2

m = Minuit(f, print_level=0, pedantic=False)
m.migrad()
m.draw_mncontour('x', 'y')

(Source code, png, hires.png, pdf)

_images/draw_mncontour.png
mnprofile(self, vname, bins=30, bound=2, subtract_min=False)

calculate minos profile around specify range ie. Migrad results with vname fixed at various places within bound

Arguments

  • vname variable name to scan

  • bins number of scanning bin. Default 30.

  • bound

    If bound is tuple, (left, right) scanning bound. If bound is a number, it specifies how many \(\sigma\) symmetrically from minimum (minimum+- bound* \(\sigma\)). Default 2

  • subtract_min subtract_minimum off from return value. This

    makes it easy to label confidence interval. Default False.

Returns

bins(center point), value, migrad results
narg

narg: object

ncalls

ncalls: ‘int’

np_matrix(self, correlation=False)

return error/correlation matrix in numpy array format.

parameters

parameters: object

pos2var

pos2var: object

print_all_minos(self)

print all minos errors(and its states)

print_fmin(self)

print current function minimum state

print_initial_param(self, **kwds)

Print initial parameters

print_level

print_level: object

print_matrix(self, **kwds)

show error_matrix

print_param(self, **kwds)

print current parameter state. Extra keyword arguments will be passed to frontend.print_param

profile(self, vname, bins=100, bound=2, args=None, subtract_min=False)

calculate cost function profile around specify range.

Arguments

  • vname variable name to scan

  • bins number of scanning bin. Default 100.

  • bound If bound is tuple, (left, right) scanning bound. If bound is a number, it specifies how many \(\sigma\) symmetrically from minimum (minimum+- bound* \(\sigma\)). Default 2

  • subtract_min subtract_minimum off from return value. This

    makes it easy to label confidence interval. Default False.

Returns

bins(center point), value

See also

mnprofile()

set_errordef(self, double errordef)

set error parameter 1 for \(\chi^2\) and 0.5 for log likelihood

set_print_level(self, lvl)

set printlevel 0 quiet, 1 normal, 2 paranoid, 3 really paranoid

set_strategy(self, stra)

set strategy 0=fast , 1=default, 2=slow but accurate

set_up(self, double errordef)

alias for set_errordef()

strategy

strategy: ‘unsigned int’

var2pos

var2pos: object

Utility Functions

:mod:util module provides describe function and various function to manipulate fitarguments.

class iminuit.util.Struct(**kwds)
iminuit.util.arguments_from_docstring(doc)

Parse first line of docstring for argument name

Docstring should be of the form ‘min(iterable[, key=func])
‘.
It can also parse cython docstring of the form Minuit.migrad(self[, int ncall_me =10000, resume=True, int nsplit=1])
iminuit.util.describe(f, verbose=False)

try to extract function arguements name

iminuit.util.fitarg_rename(fitarg, ren)

rename variable names in fitarg with rename function

#simple renaming
fitarg_rename({'x':1, 'limit_x':1, 'fix_x':1, 'error_x':1},
    lambda pname: 'y' if pname=='x' else pname)
#{'y':1, 'limit_y':1, 'fix_y':1, 'error_y':1},

#prefixing
figarg_rename({'x':1, 'limit_x':1, 'fix_x':1, 'error_x':1},
    lambda pname: 'prefix_'+pname)
#{'prefix_x':1, 'limit_prefix_x':1, 'fix_prefix_x':1, 'error_prefix_x':1}
iminuit.util.true_param(p)

check if p is parameter name not a limit/error/fix attributes

iminuit.util.param_name(p)

extract parameter name from attributes eg

fix_x -> x error_x -> x limit_x -> x

iminuit.util.extract_iv(b)

extract initial value from fitargs dictionary

iminuit.util.extract_limit(b)

extract limit from fitargs dictionary

iminuit.util.extract_error(b)

extract error from fitargs dictionary

iminuit.util.extract_fix(b)

extract fix attribute from fitargs dictionary

iminuit.util.remove_var(b, exclude)

exclude variable in exclude list from b

Return Value Struct:

iminuit uses various structs as return value. This section lists the struct and all it’s field

Function Minimum Struct

They are usually return value from Minuit.get_fmin() and Minuit.migrad() Function Mimum Struct has the following attributes:

  • fval: FCN minimum value

  • edm: Estimated Distance to Minimum.

  • nfcn: Number of function call in last mimizier call

  • up: UP parameter. This determine how minimizer define 1 \(\sigma\) error

  • is_valid: Validity of function minimum. This is defined as
    • has_valid_parameters
    • and not has_reached_call_limit
    • and not is_above_max_edm
  • has_valid_parameters: Validity of parameters. This means:
    1. The parameters must have valid error(if it’s not fixed). Valid error is not necessarily accurate.
    2. The parameters value must be valid
  • has_accurate_covariance: Boolean indicating whether covariance matrix is accurate.

  • has_pos_def_covar: Positive definiteness of covariance

  • has_made_posdef_covar: Whether minimizer has to force covariance matrix to be positive definite by adding diagonal matrix.

  • hesse_failed: Successfulness of the hesse call after minimizer.

  • has_covaraince: Has Covariance.

  • is_above_max_edm: Is EDM above 0.0001*tolerance*up? The convergence of migrad is defined by EDM being below this number.

  • has_reached_call_limit: Whether the last minimizer exceeds number of FCN calls it is allowd.

Minos Error Struct

Minos Error Struct is used in return value from Minuit.minos(). You can also call Minuit.get_merrors() to get accumulated dictionary all minos errors that has been calculated. It contains various minos status:

  • lower: lower error value
  • upper: upper error value
  • is_valid: Validity of minos error value. This means lower_valid and upper_valid
  • lower_valid: Validity of lower error
  • upper_valid: Validity of upper error
  • at_lower_limit: minos calculation hits the lower limit on parameters
  • at_upper_limit: minos calculation hits the upper limit on parameters
  • lower_new_min: found a new minimum while scanning cost function for lower error value
  • upper_new_min: found a new minimum while scanning cost function for upper error value
  • nfn: number of call to FCN in the last minos scan
  • min: the value of the parameter at the minimum

Minuit Parameter Struct

Minuit Parameter Struct is return value from Minuit.hesse() You can, however, access the latest parameter by calling Minuit.get_param_states(). Minuit Parameter Struct has the following attrubutes:

  • number: parameter number
  • name: parameter name
  • value: parameter value
  • error: parameter parabolic error(like those from hesse)
  • is_fixed: is the parameter fixed
  • is_const: is the parameter a constant(We do not support const but you can alway use fixing parameter instead)
  • has_limits: parameter has limits set
  • has_lower_limit: parameter has lower limit set. We do not support one sided limit though.
  • has_upper_limit: parameter has upper limit set.
  • lower_limit: value of lower limit for this parameter
  • upper_limit: value of upper limit for this parameter

Function Signature Extraction Ordering

  1. Using f.func_code.co_varnames, f.func_code.co_argcount All functions that is defined like:

    def f(x,y):
        return (x-2)**2+(y-3)**2
    

    or:

    f = lambda x,y: (x-2)**2+(y-3)**2
    

    Has these two attributes.

  2. Using f.__call__.func_code.co_varnames, f.__call__.co_argcount. Minuit knows how to dock off self parameter. This allow you to do things like encapsulate your data with in fitting algorithm:

    class MyChi2:
        def __init__(self, x, y):
            self.x, self.y = (x,y)
        def f(self, x, m, c):
            return m*x + c
        def __call__(self,m,c):
            return sum([(self.f(x,m,c)-y)**2
                       for x,y in zip(self.x ,self.y)])
    
  3. If all fail, Minuit will call inspect.getargspec for function signature. Builtin C functions will fall into this category since they have no signature information. inspect.getargspec will parse docstring to get function signature.

This order is very similar to PyMinuit signature detection. Actually, it is a superset of PyMinuit signature detection. The difference is that it allows you to fake function signature by having func_code attribute in the object. This allows you to make a generic functor of your custom cost function. This is how probfit was written:

f = lambda x,m,c: m*x+c
#the beauty here is that all you need to build
#a Chi^2 is just your function and data
class GenericChi2:
    def __init__(self, f, x, y):
        self.f = f
        args = describe(f)#extract function signature
        self.func_code = Struct(
                co_varnames = args[1:],#dock off independent param
                co_argcount = len(args)-1
            )
    def __call__(self, *arg):
        return sum((self.f(x,*arg)-y)**2 for x,y in zip(self.x, self.y))

m = Minuit(GenericChi2(f,x,y))

Note

If you are unsure what minuit will parse your function signature as , you can use describe() which returns tuple of arument names minuit will use as call signature.