This tutorial shows an example of using SPlot to unfold two distributions. The physics context for the example is that we want to know the isolation distribution for real electrons from Z events and fake electrons from QCD. Isolation is our 'control' variable. To unfold them, we need a model for an uncorrelated variable that discriminates between Z and QCD. To do this, we use the invariant mass of two electrons. We model the Z with a Gaussian and the QCD with a falling exponential.
Note, since we don't have real data in this tutorial, we need to generate toy data. To do that we need a model for the isolation variable for both Z and QCD. This is only used to generate the toy data, and would not be needed if we had real data.
#include <iomanip>
{
AddModel(wspace);
AddData(wspace);
DoSPlot(wspace);
MakePlots(wspace);
}
{
Double_t lowRange = 0., highRange = 200.;
RooRealVar invMass(
"invMass",
"M_{inv}", lowRange, highRange,
"GeV");
RooRealVar isolation(
"isolation",
"isolation", 0., 20.,
"GeV");
std::cout << "make z model" << std::endl;
RooRealVar mZ(
"mZ",
"Z Mass", 91.2, lowRange, highRange);
RooRealVar sigmaZ(
"sigmaZ",
"Width of Gaussian", 2, 0, 10,
"GeV");
RooGaussian mZModel(
"mZModel",
"Z+jets Model", invMass, mZ, sigmaZ);
mZ.setConstant();
RooConstVar zIsolDecayConst(
"zIsolDecayConst",
"z isolation decay constant", -1);
RooExponential zIsolationModel(
"zIsolationModel",
"z isolation model", isolation, zIsolDecayConst);
std::cout << "make qcd model" << std::endl;
RooRealVar qcdMassDecayConst(
"qcdMassDecayConst",
"Decay const for QCD mass spectrum", -0.01, -100, 100,
"1/GeV");
RooExponential qcdMassModel(
"qcdMassModel",
"qcd Mass Model", invMass, qcdMassDecayConst);
RooConstVar qcdIsolDecayConst(
"qcdIsolDecayConst",
"Et resolution constant", -.1);
RooExponential qcdIsolationModel(
"qcdIsolationModel",
"QCD isolation model", isolation, qcdIsolDecayConst);
RooProdPdf qcdModel(
"qcdModel",
"2-d model for QCD",
RooArgSet(qcdMassModel, qcdIsolationModel));
RooRealVar zYield(
"zYield",
"fitted yield for Z", 500, 0., 5000);
RooRealVar qcdYield(
"qcdYield",
"fitted yield for QCD", 1000, 0., 10000);
std::cout << "make full model" << std::endl;
RooAddPdf model(
"model",
"z+qcd background models", {zModel, qcdModel}, {zYield, qcdYield});
RooAddPdf massModel(
"massModel",
"z+qcd invariant mass model", {mZModel, qcdMassModel}, {zYield, qcdYield});
std::cout << "import model" << std::endl;
}
{
std::cout << "make data set and import to workspace" << std::endl;
std::unique_ptr<RooDataSet>
data{model->
generate({*invMass, *isolation})};
}
{
std::cout << "Calculate sWeights" << std::endl;
std::cout << "\n\n------------------------------------------\nThe dataset before creating sWeights:\n";
std::cout << "\n\nThe dataset after creating sWeights:\n";
std::cout << "\n\n------------------------------------------\n\nCheck SWeights:" << std::endl;
std::cout << std::endl
<<
"Yield of Z is\t" << zYield->
getVal() <<
". From sWeights it is "
<< sData.GetYieldFromSWeight("zYield") << std::endl;
std::cout <<
"Yield of QCD is\t" << qcdYield->
getVal() <<
". From sWeights it is "
<< sData.GetYieldFromSWeight("qcdYield") << std::endl
<< std::endl;
for (
Int_t i = 0; i < 10; i++) {
std::cout << "z Weight for event " << i << std::right << std::setw(12) << sData.GetSWeight(i, "zYield") << " qcd Weight"
<< std::setw(12) << sData.GetSWeight(i, "qcdYield") << " Total Weight" << std::setw(12) << sData.GetSumOfEventSWeight(i)
<< std::endl;
}
std::cout << std::endl;
std::cout << "import new dataset with sWeights" << std::endl;
}
{
std::cout << "make plots" << std::endl;
RooPlot *frame2 = isolation->
frame(
Title(
"Isolation distribution with s weights to project out Z"));
RooPlot *frame3 = isolation->
frame(
Title(
"Isolation distribution with s weights to project out QCD"));
}
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
void graphVizTree(const char *fileName, const char *delimiter="\n", bool useTitle=false, bool useLatex=false)
Create a GraphViz .dot file visualizing the expression tree headed by this RooAbsArg object.
Abstract interface for all probability density functions.
RooPlot * plotOn(RooPlot *frame, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={}, const RooCmdArg &arg10={}) const override
Helper calling plotOn(RooPlot*, RooLinkedList&) const.
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
RooPlot * frame(const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create a new RooPlot on the heap with a drawing frame initialized for this object,...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
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.
Represents a constant real-valued object.
Container class to hold unbinned data.
static RooMsgService & instance()
Return reference to singleton instance.
void setGlobalKillBelow(RooFit::MsgLevel level)
void setSilentMode(bool flag)
Plot frame and a container for graphics objects within that frame.
TObject * findObject(const char *name, const TClass *tClass=nullptr) const
Find the named object in our list of items and return a pointer to it.
void Draw(Option_t *options=nullptr) override
Draw this plot and all of the elements it contains.
Efficient implementation of a product of PDFs of the form.
Variable that can be changed from the outside.
A class to calculate "sWeights" used to create an "sPlot".
Persistable container for RooFit projects.
RooAbsPdf * pdf(RooStringView name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
RooRealVar * var(RooStringView name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
RooAbsData * data(RooStringView name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
bool import(const RooAbsArg &arg, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={})
Import a RooAbsArg object, e.g.
TVirtualPad * cd(Int_t subpadnumber=0) override
Set current canvas & pad.
This class displays a legend box (TPaveText) containing several legend entries.
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
RooCmdArg RecycleConflictNodes(bool flag=true)
RooCmdArg Rename(const char *suffix)
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg DataError(Int_t)
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.
make z model
[#0] WARNING:InputArguments -- The parameter 'sigmaZ' with range [0, 10] of the RooGaussian 'mZModel' exceeds the safe range of (0, inf). Advise to limit its range.
make qcd model
make full model
import model
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooAddPdf::model
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooProdPdf::zModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooGaussian::mZModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::invMass
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::mZ
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::sigmaZ
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooExponential::zIsolationModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::isolation
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooConstVar::zIsolDecayConst
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::zYield
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooProdPdf::qcdModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooExponential::qcdMassModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::qcdMassDecayConst
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooExponential::qcdIsolationModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooConstVar::qcdIsolDecayConst
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooRealVar::qcdYield
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooGaussian::mZModel for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooRealVar::invMass for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooRealVar::mZ for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooRealVar::sigmaZ for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooRealVar::zYield for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooExponential::qcdMassModel for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooRealVar::qcdMassDecayConst for import of RooAddPdf::massModel
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) using existing copy of RooRealVar::qcdYield for import of RooAddPdf::massModel
make data set and import to workspace
[#1] INFO:ObjectHandling -- RooWorkspace::import(myWS) importing dataset modelData
[#1] INFO:ObjectHandling -- RooWorkSpace::import(myWS) changing name of dataset from modelData to data
Calculate sWeights
------------------------------------------
The dataset before creating sWeights:
RooDataSet::data[invMass,isolation] = 1500 entries
The dataset after creating sWeights:
RooDataSet::data[invMass,isolation,zYield_sw,L_zYield,qcdYield_sw,L_qcdYield] = 1500 entries
------------------------------------------
Check SWeights:
Yield of Z is 543.477. From sWeights it is 543.477
Yield of QCD is 956.917. From sWeights it is 956.917
z Weight for event 0 -0.0490677 qcd Weight 1.04942 Total Weight 1.00036
z Weight for event 1 -0.0490677 qcd Weight 1.04942 Total Weight 1.00036
z Weight for event 2 1.0038 qcd Weight -0.00368157 Total Weight 1.00012
z Weight for event 3 0.950485 qcd Weight 0.0496438 Total Weight 1.00013
z Weight for event 4 -0.0490677 qcd Weight 1.04942 Total Weight 1.00036
z Weight for event 5 0.994479 qcd Weight 0.00564054 Total Weight 1.00012
z Weight for event 6 -0.0490677 qcd Weight 1.04942 Total Weight 1.00036
z Weight for event 7 -0.0490677 qcd Weight 1.04942 Total Weight 1.00036
z Weight for event 8 1.04017 qcd Weight -0.0400605 Total Weight 1.00011
z Weight for event 9 1.04125 qcd Weight -0.0411384 Total Weight 1.00011
import new dataset with sWeights
make plots
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on invMass integrates over variables (isolation)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (zModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: (mZModel,zIsolationModel)
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on invMass integrates over variables (isolation)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) directly selected PDF components: (qcdModel)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(model) indirectly selected PDF components: (qcdMassModel,qcdIsolationModel)
[#1] INFO:Plotting -- RooAbsReal::plotOn(model) plot on invMass integrates over variables (isolation)
[#1] INFO:Plotting -- RooAbsReal::plotOn(zModel) plot on isolation integrates over variables (invMass)
[#1] INFO:Plotting -- RooAbsReal::plotOn(qcdModel) plot on isolation integrates over variables (invMass)