With default parameters the macro will attempt to run the standard hist2workspace example and read the ROOT file that it produces.
The macro will scan through all the categories in a simPdf find the corresponding observable. For each category, it will loop through each of the nuisance parameters and plot
You can specify how many sigma to vary by changing nSigmaToVary. You can also change the signal rate by changing muVal.
␛[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) 0x564a6639ee40 RooRealVar:: alpha_syst2 = 0 L(-5 - 5) "alpha_syst2"
2) 0x564a663a5c60 RooRealVar:: alpha_syst3 = 0 L(-5 - 5) "alpha_syst3"
3) 0x564a65b051c0 RooRealVar:: gamma_stat_channel1_bin_0 = 1 L(0 - 1.25) "gamma_stat_channel1_bin_0"
4) 0x564a651c32a0 RooRealVar:: gamma_stat_channel1_bin_1 = 1 L(0 - 1.5) "gamma_stat_channel1_bin_1"
check expectedData by category
Is a simultaneous PDF
on type channel1
channelCat==channelCat::channel1
channel1 channel1
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 234 events out of 234 total events
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
channelCat==channelCat::channel1
channel1 channel1
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 234 events out of 234 total events
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
channelCat==channelCat::channel1
channel1 channel1
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 234 events out of 234 total events
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
channelCat==channelCat::channel1
channel1 channel1
[#1] INFO:Plotting -- RooTreeData::plotOn: plotting 234 events out of 234 total events
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(channel1_model_Int[obs_x_channel1]) using numeric integrator RooBinIntegrator to calculate Int(obs_x_channel1)
using namespace std;
void StandardHistFactoryPlotsWithCategories(const char *infile = "", const char *workspaceName = "combined",
const char *modelConfigName = "ModelConfig",
const char *dataName = "obsData")
{
double nSigmaToVary = 5.;
double muVal = 0;
bool doFit = false;
const char *filename = "";
if (!strcmp(infile, "")) {
filename = "results/example_combined_GaussExample_model.root";
if (!fileExist) {
#ifdef _WIN32
cout << "HistFactory file cannot be generated on Windows - exit" << endl;
return;
#endif
cout << "will run standard hist2workspace example" << endl;
gROOT->ProcessLine(
".! prepareHistFactory .");
gROOT->ProcessLine(
".! hist2workspace config/example.xml");
cout << "\n\n---------------------" << endl;
cout << "Done creating example input" << endl;
cout << "---------------------\n\n" << endl;
}
} else
filename = infile;
cout << "StandardRooStatsDemoMacro: Input file " << filename << " is not found" << endl;
return;
}
if (!w) {
cout << "workspace not found" << endl;
return;
}
ModelConfig *mc = (ModelConfig *)w->
obj(modelConfigName);
if (!data || !mc) {
cout << "data or ModelConfig was not found" << endl;
return;
}
if (doFit) {
mc->GetPdf()->fitTo(*data);
}
mc->GetNuisanceParameters()->Print("v");
int nPlotsMax = 1000;
cout << " check expectedData by category" << endl;
if (strcmp(mc->GetPdf()->ClassName(), "RooSimultaneous") == 0) {
cout << "Is a simultaneous PDF" << endl;
} else {
cout << "Is not a simultaneous PDF" << endl;
}
if (doFit) {
cout <<
tt->GetName() <<
" " << channelCat->
getLabel() << endl;
cout <<
"expected events = " << mc->GetPdf()->expectedEvents(*data->
get()) << endl;
return;
}
int nPlots = 0;
if (!simPdf) {
TIterator *it = mc->GetNuisanceParameters()->createIterator();
}
} else {
cout <<
"on type " <<
tt->GetName() <<
" " << endl;
TIterator *it = mc->GetNuisanceParameters()->createIterator();
cout <<
tt->GetName() <<
" " << channelCat->
getLabel() << endl;
if (strcmp(var->
GetName(),
"Lumi") == 0) {
cout << "working on lumi" << endl;
} else {
}
if (strcmp(var->
GetName(),
"Lumi") == 0) {
cout << "working on lumi" << endl;
} else {
}
if (strcmp(var->
GetName(),
"Lumi") == 0) {
cout << "working on lumi" << endl;
} else {
}
if (strcmp(var->
GetName(),
"Lumi") == 0) {
cout << "working on lumi" << endl;
} else {
}
++nPlots;
}
}
}
} else
for (
int i = 0; i < list->
GetSize(); ++i) {
}
}
char * Form(const char *fmt,...)
R__EXTERN TSystem * gSystem
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
virtual void Print(Option_t *options=0) const
Print TNamed name and title.
TIterator * typeIterator() const
Return iterator over all defined states.
RooAbsArg * first() const
RooAbsData is the common abstract base class for binned and unbinned datasets.
virtual const RooArgSet * get() const
virtual Double_t sumEntries() const =0
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
Calls RooPlot* plotOn(RooPlot* frame, const RooLinkedList& cmdList) const ;.
virtual Double_t expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
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,...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooCatType is an auxilary class for RooAbsCategory and defines a a single category state.
RooCategory represents a fundamental (non-derived) discrete value object.
virtual const char * getLabel() const
Return label string of current state.
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 SetName(const char *name)
Set the name of the RooPlot to 'name'.
void SetYTitle(const char *title)
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
RooRealVar represents a fundamental (non-derived) real valued object.
virtual void setVal(Double_t value)
Set value of variable to 'value'.
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const RooAbsCategoryLValue & indexCat() const
RooAbsPdf * getPdf(const char *catName) const
Return the p.d.f associated with the given index category name.
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.
void Print(Option_t *opts=0) const
Print contents of the workspace.
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseGeneralPurpose, Int_t netopt=0)
Create / open a file.
Iterator abstract base class.
virtual TObject * Next()=0
virtual void Add(TObject *obj)
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
virtual const char * GetName() const
Returns name of object.
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Template specialisation used in RooAbsArg:
RooCmdArg LineWidth(Width_t width)
RooCmdArg MarkerSize(Size_t size)
RooCmdArg Cut(const char *cutSpec)
RooCmdArg DataError(Int_t)
RooCmdArg LineColor(Color_t color)
RooCmdArg Normalization(Double_t scaleFactor)
RooCmdArg LineStyle(Style_t style)
Namespace for the RooStats classes.
Int_t CeilNint(Double_t x)