#!/usr/bin/env python # coding: utf-8 # # rf604_constraints # Likelihood and minimization: fitting with constraints # # # # # **Author:** Clemens Lange, Wouter Verkerke (C++ version) # This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Wednesday, April 17, 2024 at 11:19 AM. # In[1]: from __future__ import print_function import ROOT # Create model and dataset # ---------------------------------------------- # Construct a Gaussian pdf # In[2]: x = ROOT.RooRealVar("x", "x", -10, 10) m = ROOT.RooRealVar("m", "m", 0, -10, 10) s = ROOT.RooRealVar("s", "s", 2, 0.1, 10) gauss = ROOT.RooGaussian("gauss", "gauss(x,m,s)", x, m, s) # Construct a flat pdf (polynomial of 0th order) # In[3]: poly = ROOT.RooPolynomial("poly", "poly(x)", x) # model = f*gauss + (1-f)*poly # In[4]: f = ROOT.RooRealVar("f", "f", 0.5, 0.0, 1.0) model = ROOT.RooAddPdf("model", "model", [gauss, poly], [f]) # Generate small dataset for use in fitting below # In[5]: d = model.generate({x}, 50) # Create constraint pdf # ----------------------------------------- # Construct Gaussian constraint pdf on parameter f at 0.8 with # resolution of 0.1 # In[6]: fconstraint = ROOT.RooGaussian("fconstraint", "fconstraint", f, 0.8, 0.1) # Method 1 - add internal constraint to model # ------------------------------------------------------------------------------------- # Multiply constraint term with regular pdf using ROOT.RooProdPdf Specify in # fitTo() that internal constraints on parameter f should be used # Multiply constraint with pdf # In[7]: modelc = ROOT.RooProdPdf("modelc", "model with constraint", [model, fconstraint]) # Fit model (without use of constraint term) # In[8]: r1 = model.fitTo(d, Save=True, PrintLevel=-1) # Fit modelc with constraint term on parameter f # In[9]: r2 = modelc.fitTo(d, Constrain={f}, Save=True, PrintLevel=-1) # Method 2 - specify external constraint when fitting # ------------------------------------------------------------------------------------------ # Construct another Gaussian constraint pdf on parameter f at 0.8 with # resolution of 0.1 # In[10]: fconstext = ROOT.RooGaussian("fconstext", "fconstext", f, 0.2, 0.1) # Fit with external constraint # In[11]: r3 = model.fitTo(d, ExternalConstraints={fconstext}, Save=True, PrintLevel=-1) # Print the fit results # In[12]: print("fit result without constraint (data generated at f=0.5)") r1.Print("v") print("fit result with internal constraint (data generated at f=0.5, is f=0.8+/-0.2)") r2.Print("v") print("fit result with (another) external constraint (data generated at f=0.5, is f=0.2+/-0.1)") r3.Print("v")