The macro creates a simple signal model and two background models, which are added to a RooWorkspace. The macro creates a toy dataset, and then uses a RooStats ProfileLikleihoodCalculator to do a hypothesis test of the background-only and signal+background hypotheses. In this example, shape uncertainties are not taken into account, but normalization uncertainties are.
␛[1mRooFit v3.60 -- Developed by Wouter Verkerke and David Kirkby␛[0m
Copyright (C) 2000-2013 NIKHEF, University of California & Stanford University
All rights reserved, please read http://roofit.sourceforge.net/license.txt
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooAddPdf::model
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooGaussian::sigModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::invMass
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::mH
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::sigma1
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooProduct::fsig
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::mu
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::ratioSigEff
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::fsigExpected
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooGaussian::zjjModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::mZ
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::sigma1_z
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::fzjj
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooChebychev::qcdModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::a0
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::a1
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::a2
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing dataset modelData
[#1] INFO:ObjectHandling -- RooWorkSpace::import(myWS) changing name of dataset from modelData to data
[#1] INFO:Minimization -- createConstraintTerm: caching constraint set under name CACHE_CONSTR_OF_PDF_model_FOR_OBS_invMass with 0 entries
[#0] PROGRESS:Minimization -- ProfileLikelihoodCalcultor::DoGLobalFit - find MLE
[#0] PROGRESS:Minimization -- ProfileLikelihoodCalcultor::DoMinimizeNLL - using Minuit / Migrad with strategy 1
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sigModel,zjjModel,qcdModel)
[#1] INFO:Minimization --
RooFitResult: minimized FCN value: 717.039, estimated distance to minimum: 8.90226e-10
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
fzjj 3.1152e-01 +/- 5.03e-02
mu 1.0968e+00 +/- 3.03e-01
[#0] PROGRESS:Minimization -- ProfileLikelihoodCalcultor::GetHypoTest - do conditional fit
[#0] PROGRESS:Minimization -- ProfileLikelihoodCalcultor::DoMinimizeNLL - using Minuit / Migrad with strategy 1
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization --
RooFitResult: minimized FCN value: 723.97, estimated distance to minimum: 2.09862e-09
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
fzjj 2.6213e-01 +/- 5.18e-02
-------------------------------------------------
The p-value for the null is 9.83108e-05
Corresponding to a significance of 3.72332
-------------------------------------------------
[#1] INFO:Minimization -- createConstraintTerm picked up cached constraints from workspace with 0 entries
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sigModel,zjjModel,qcdModel)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (sigModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (zjjModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (qcdModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
[#1] INFO:Minimization -- createConstraintTerm picked up cached constraints from workspace with 0 entries
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (sigModel,zjjModel,qcdModel)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (zjjModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (qcdModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: ()
#include <string>
void rs102_hypotestwithshapes()
{
AddModel(wspace);
AddData(wspace);
DoHypothesisTest(wspace);
MakePlots(wspace);
delete wspace;
}
{
Double_t lowRange = 60, highRange = 200;
RooRealVar invMass(
"invMass",
"M_{inv}", lowRange, highRange,
"GeV");
RooRealVar sigma1(
"sigma1",
"Width of Gaussian", 12., 2, 100);
RooGaussian sigModel(
"sigModel",
"Signal Model", invMass, mH, sigma1);
mH.setConstant();
sigma1.setConstant();
RooRealVar sigma1_z(
"sigma1_z",
"Width of Gaussian", 10., 6, 100);
RooGaussian zjjModel(
"zjjModel",
"Z+jets Model", invMass, mZ, sigma1_z);
mZ.setConstant();
sigma1_z.setConstant();
a0.setConstant();
a1.setConstant();
a2.setConstant();
RooRealVar fzjj(
"fzjj",
"fraction of zjj background events", .4, 0., 1);
RooRealVar fsigExpected(
"fsigExpected",
"expected fraction of signal events", .2, 0., 1);
fsigExpected.setConstant();
RooRealVar mu(
"mu",
"signal strength in units of SM expectation", 1, 0., 2);
RooRealVar ratioSigEff(
"ratioSigEff",
"ratio of signal efficiency to nominal signal efficiency", 1., 0., 2);
ratioSigEff.setConstant(
kTRUE);
RooProduct fsig(
"fsig",
"fraction of signal events",
RooArgSet(mu, ratioSigEff, fsigExpected));
RooAddPdf model(
"model",
"sig+zjj+qcd background shapes",
RooArgList(sigModel, zjjModel, qcdModel),
}
{
}
{
cout << "-------------------------------------------------" << endl;
cout <<
"The p-value for the null is " << htr->
NullPValue() << endl;
cout <<
"Corresponding to a significance of " << htr->
Significance() << endl;
cout << "-------------------------------------------------\n\n" << endl;
}
{
frame->
SetTitle(
"An example fit to the signal + background model");
xframe2->
SetTitle(
"An example fit to the background-only model");
}
Bool_t setRealValue(const char *name, Double_t newVal=0, Bool_t verbose=kFALSE)
Set value of a RooAbsRealLValye stored in set with given name to newVal No error messages are printed...
RooAbsData is the common abstract base class for binned and unbinned datasets.
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
virtual RooFitResult * fitTo(RooAbsData &data, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none())
Fit PDF to given dataset.
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
virtual RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none(), const RooCmdArg &arg9=RooCmdArg::none(), const RooCmdArg &arg10=RooCmdArg::none()) const
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
RooPlot * frame(const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create a new RooPlot on the heap with a drawing frame initialized for this object,...
void setConstant(Bool_t value=kTRUE)
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Chebychev polynomial p.d.f.
RooDataSet is a container class to hold unbinned data.
A RooPlot is a plot frame and a container for graphics objects within that frame.
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
A RooProduct represents the product of a given set of RooAbsReal objects.
RooRealVar represents a variable that can be changed from the outside.
virtual void setVal(Double_t value)
Set value of variable to 'value'.
virtual void SetNullParameters(const RooArgSet &set)
set parameter values for the null if using a common PDF
virtual void SetModel(const ModelConfig &model)
set the model (in this case can set only the parameters for the null hypothesis)
virtual void SetData(RooAbsData &data)
Set the DataSet, add to the the workspace if not already there.
HypoTestResult is a base class for results from hypothesis tests.
virtual Double_t Significance() const
familiar name for the Null p-value in terms of 1-sided Gaussian significance
virtual Double_t NullPValue() const
Return p-value for null hypothesis.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
virtual void SetWorkspace(RooWorkspace &ws)
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
The ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface cla...
virtual HypoTestResult * GetHypoTest() const
Return the hypothesis test result obtained from the likelihood ratio of the maximum likelihood value ...
The RooWorkspace is a persistable container for RooFit projects.
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
Bool_t import(const RooAbsArg &arg, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg())
Import a RooAbsArg object, e.g.
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
RooCmdArg Rename(const char *suffix)
RooCmdArg Hesse(Bool_t flag=kTRUE)
RooCmdArg DataError(Int_t)
RooCmdArg Save(Bool_t flag=kTRUE)
RooCmdArg PrintLevel(Int_t code)
RooCmdArg Minos(Bool_t flag=kTRUE)
RooCmdArg Components(Args_t &&... argsOrArgSet)
RooCmdArg LineColor(Color_t color)
RooCmdArg LineStyle(Style_t style)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.