Logo ROOT  
Reference Guide
RooStatsUtils.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer 28/07/2008
3
4/*************************************************************************
5 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "Rtypes.h"
13
14#if !defined(R__ALPHA) && !defined(R__SOLARIS) && !defined(R__ACC) && !defined(R__FBSD)
16#endif
17
18#include "TTree.h"
19
20#include "RooUniform.h"
21#include "RooProdPdf.h"
22#include "RooExtendPdf.h"
23#include "RooSimultaneous.h"
26#include <typeinfo>
27
28using namespace std;
29
30namespace RooStats {
31
33 static RooStatsConfig theConfig;
34 return theConfig;
35 }
36
38 // Asimov significance
39 // formula [10] and [20] from https://www.pp.rhul.ac.uk/~cowan/stat/notes/medsigNote.pdf
40 // case we have a sigma_b
41 double sb2 = sigma_b*sigma_b;
42 // formula below has a large error when sigma_b becomes zero
43 // better to use the approximation for sigma_b=0 for very small values
44 double r = sb2/b;
45 if (r > 1.E-12) {
46 double bpsb2 = b + sb2;
47 double b2 = b*b;
48 double spb = s+b;
49 double za2 = 2.*( (spb)* std::log( ( spb)*(bpsb2)/(b2+ spb*sb2) ) -
50 (b2/sb2) * std::log(1. + ( sb2 * s)/(b * bpsb2) ) );
51 return sqrt(za2);
52
53 }
54 // case when the background (b) is known
55 double za2 = 2.*( (s+b) * std::log(1. + s/b) -s );
56 return std::sqrt(za2);
57 }
58
59 /// Use an offset in NLL calculations.
60 void UseNLLOffset(bool on) {
62 }
63
64 /// Test of RooStats should by default offset NLL calculations.
65 bool IsNLLOffset() {
67 }
68
69 void FactorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints) {
70 // utility function to factorize constraint terms from a pdf
71 // (from G. Petrucciani)
72 const std::type_info & id = typeid(pdf);
73 if (id == typeid(RooProdPdf)) {
74 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
75 RooArgList list(prod->pdfList());
76 for (int i = 0, n = list.getSize(); i < n; ++i) {
77 RooAbsPdf *pdfi = (RooAbsPdf *) list.at(i);
78 FactorizePdf(observables, *pdfi, obsTerms, constraints);
79 }
80 } else if (id == typeid(RooExtendPdf)) {
81 TIterator *iter = pdf.serverIterator();
82 // extract underlying pdf which is extended; first server is the pdf; second server is the number of events variable
83 RooAbsPdf *updf = dynamic_cast<RooAbsPdf *>(iter->Next());
84 assert(updf != 0);
85 delete iter;
86 FactorizePdf(observables, *updf, obsTerms, constraints);
87 } else if (id == typeid(RooSimultaneous)) { //|| id == typeid(RooSimultaneousOpt)) {
88 RooSimultaneous *sim = dynamic_cast<RooSimultaneous *>(&pdf);
89 assert(sim != 0);
91 for (int ic = 0, nc = cat->numBins((const char *)0); ic < nc; ++ic) {
92 cat->setBin(ic);
93 RooAbsPdf* catPdf = sim->getPdf(cat->getLabel());
94 // it is possible that a pdf is not defined for every category
95 if (catPdf != 0) FactorizePdf(observables, *catPdf, obsTerms, constraints);
96 }
97 delete cat;
98 } else if (pdf.dependsOn(observables)) {
99 if (!obsTerms.contains(pdf)) obsTerms.add(pdf);
100 } else {
101 if (!constraints.contains(pdf)) constraints.add(pdf);
102 }
103 }
104
105
106 void FactorizePdf(RooStats::ModelConfig &model, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints) {
107 // utility function to factorize constraint terms from a pdf
108 // (from G. Petrucciani)
109 if (!model.GetObservables() ) {
110 oocoutE((TObject*)0,InputArguments) << "RooStatsUtils::FactorizePdf - invalid input model: missing observables" << endl;
111 return;
112 }
113 return FactorizePdf(*model.GetObservables(), pdf, obsTerms, constraints);
114 }
115
116
117 RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name) {
118 // make a nuisance pdf by factorizing out all constraint terms in a common pdf
119 RooArgList obsTerms, constraints;
120 FactorizePdf(observables, pdf, obsTerms, constraints);
121 if(constraints.getSize() == 0) {
122 oocoutW((TObject *)0, Eval) << "RooStatsUtils::MakeNuisancePdf - no constraints found on nuisance parameters in the input model" << endl;
123 return 0;
124 } else if(constraints.getSize() == 1) {
125 return dynamic_cast<RooAbsPdf *>(constraints.first()->clone(name));
126 }
127 return new RooProdPdf(name,"", constraints);
128 }
129
131 // make a nuisance pdf by factorizing out all constraint terms in a common pdf
132 if (!model.GetPdf() || !model.GetObservables() ) {
133 oocoutE((TObject*)0, InputArguments) << "RooStatsUtils::MakeNuisancePdf - invalid input model: missing pdf and/or observables" << endl;
134 return 0;
135 }
136 return MakeNuisancePdf(*model.GetPdf(), *model.GetObservables(), name);
137 }
138
139 RooAbsPdf * StripConstraints(RooAbsPdf &pdf, const RooArgSet &observables) {
140 const std::type_info & id = typeid(pdf);
141
142 if (id == typeid(RooProdPdf)) {
143
144 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
145 RooArgList list(prod->pdfList()); RooArgList newList;
146
147 for (int i = 0, n = list.getSize(); i < n; ++i) {
148 RooAbsPdf *pdfi = (RooAbsPdf *) list.at(i);
149 RooAbsPdf *newPdfi = StripConstraints(*pdfi, observables);
150 if(newPdfi != NULL) newList.add(*newPdfi);
151 }
152
153 if(newList.getSize() == 0) return NULL; // only constraints in product
154 // return single component (no longer a product)
155 else if(newList.getSize() == 1) return dynamic_cast<RooAbsPdf *>(newList.at(0)->clone(TString::Format("%s_unconstrained",
156 newList.at(0)->GetName())));
157 else return new RooProdPdf(TString::Format("%s_unconstrained", prod->GetName()).Data(),
158 TString::Format("%s without constraints", prod->GetTitle()).Data(), newList);
159
160 } else if (id == typeid(RooExtendPdf)) {
161
162 TIterator *iter = pdf.serverIterator();
163 // extract underlying pdf which is extended; first server is the pdf; second server is the number of events variable
164 RooAbsPdf *uPdf = dynamic_cast<RooAbsPdf *>(iter->Next());
165 RooAbsReal *extended_term = dynamic_cast<RooAbsReal *>(iter->Next());
166 assert(uPdf != NULL); assert(extended_term != NULL); assert(iter->Next() == NULL);
167 delete iter;
168
169 RooAbsPdf *newUPdf = StripConstraints(*uPdf, observables);
170 if(newUPdf == NULL) return NULL; // only constraints in underlying pdf
171 else return new RooExtendPdf(TString::Format("%s_unconstrained", pdf.GetName()).Data(),
172 TString::Format("%s without constraints", pdf.GetTitle()).Data(), *newUPdf, *extended_term);
173
174 } else if (id == typeid(RooSimultaneous)) { //|| id == typeid(RooSimultaneousOpt)) {
175
176 RooSimultaneous *sim = dynamic_cast<RooSimultaneous *>(&pdf); assert(sim != NULL);
177 RooAbsCategoryLValue *cat = (RooAbsCategoryLValue *) sim->indexCat().Clone(); assert(cat != NULL);
178 RooArgList pdfList;
179
180 for (int ic = 0, nc = cat->numBins((const char *)NULL); ic < nc; ++ic) {
181 cat->setBin(ic);
182 RooAbsPdf* catPdf = sim->getPdf(cat->getLabel());
183 RooAbsPdf* newPdf = NULL;
184 // it is possible that a pdf is not defined for every category
185 if (catPdf != NULL) newPdf = StripConstraints(*catPdf, observables);
186 if (newPdf == NULL) { delete cat; return NULL; } // all channels must have observables
187 pdfList.add(*newPdf);
188 }
189
190 return new RooSimultaneous(TString::Format("%s_unconstrained", sim->GetName()).Data(),
191 TString::Format("%s without constraints", sim->GetTitle()).Data(), pdfList, *cat);
192
193 } else if (pdf.dependsOn(observables)) {
194 return (RooAbsPdf *) pdf.clone(TString::Format("%s_unconstrained", pdf.GetName()).Data());
195 }
196
197 return NULL; // just a constraint term
198 }
199
200 RooAbsPdf * MakeUnconstrainedPdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name) {
201 // make a clone pdf without all constraint terms in a common pdf
202 RooAbsPdf * unconstrainedPdf = StripConstraints(pdf, observables);
203 if(!unconstrainedPdf) {
204 oocoutE((TObject *)NULL, InputArguments) << "RooStats::MakeUnconstrainedPdf - invalid observable list passed (observables not found in original pdf) or invalid pdf passed (without observables)" << endl;
205 return NULL;
206 }
207 if(name != NULL) unconstrainedPdf->SetName(name);
208 return unconstrainedPdf;
209 }
210
212 // make a clone pdf without all constraint terms in a common pdf
213 if(!model.GetPdf() || !model.GetObservables()) {
214 oocoutE((TObject *)NULL, InputArguments) << "RooStatsUtils::MakeUnconstrainedPdf - invalid input model: missing pdf and/or observables" << endl;
215 return NULL;
216 }
217 return MakeUnconstrainedPdf(*model.GetPdf(), *model.GetObservables(), name);
218 }
219
220 // Helper class for GetAsTTree
221 class BranchStore {
222 public:
223 std::map<TString, Double_t> fVarVals;
224 double fInval;
225 TTree *fTree;
226
227 BranchStore(const vector <TString> &params = vector <TString>(), double _inval = -999.) : fTree(0) {
228 fInval = _inval;
229 for(unsigned int i = 0;i<params.size();i++)
230 fVarVals[params[i]] = _inval;
231 }
232
233 ~BranchStore() {
234 if (fTree) {
235 for(std::map<TString, Double_t>::iterator it = fVarVals.begin();it!=fVarVals.end();++it) {
236 TBranch *br = fTree->GetBranch( it->first );
237 if (br) br->ResetAddress();
238 }
239 }
240 }
241
242 void AssignToTTree(TTree &myTree) {
243 fTree = &myTree;
244 for(std::map<TString, Double_t>::iterator it = fVarVals.begin();it!=fVarVals.end();++it) {
245 const TString& name = it->first;
246 myTree.Branch( name, &fVarVals[name], TString::Format("%s/D", name.Data()));
247 }
248 }
249 void ResetValues() {
250 for(std::map<TString, Double_t>::iterator it = fVarVals.begin();it!=fVarVals.end();++it) {
251 const TString& name = it->first;
252 fVarVals[name] = fInval;
253 }
254 }
255 };
256
257 BranchStore* CreateBranchStore(const RooDataSet& data) {
258 if (data.numEntries() == 0) {
259 return new BranchStore;
260 }
261 vector <TString> V;
262 const RooArgSet* aset = data.get(0);
263 RooAbsArg *arg(0);
264 TIterator *it = aset->createIterator();
265 for(;(arg = dynamic_cast<RooAbsArg*>(it->Next()));) {
266 RooRealVar *rvar = dynamic_cast<RooRealVar*>(arg);
267 if (rvar == NULL)
268 continue;
269 V.push_back(rvar->GetName());
270 if (rvar->hasAsymError()) {
271 V.push_back(TString::Format("%s_errlo", rvar->GetName()));
272 V.push_back(TString::Format("%s_errhi", rvar->GetName()));
273 }
274 else if (rvar->hasError()) {
275 V.push_back(TString::Format("%s_err", rvar->GetName()));
276 }
277 }
278 delete it;
279 return new BranchStore(V);
280 }
281
282 void FillTree(TTree &myTree, const RooDataSet &data) {
283 BranchStore *bs = CreateBranchStore(data);
284 bs->AssignToTTree(myTree);
285
286 for(int entry = 0;entry<data.numEntries();entry++) {
287 bs->ResetValues();
288 const RooArgSet* aset = data.get(entry);
289 RooAbsArg *arg(0);
290 RooLinkedListIter it = aset->iterator();
291 for(;(arg = dynamic_cast<RooAbsArg*>(it.Next()));) {
292 RooRealVar *rvar = dynamic_cast<RooRealVar*>(arg);
293 if (rvar == NULL)
294 continue;
295 bs->fVarVals[rvar->GetName()] = rvar->getValV();
296 if (rvar->hasAsymError()) {
297 bs->fVarVals[TString::Format("%s_errlo", rvar->GetName())] = rvar->getAsymErrorLo();
298 bs->fVarVals[TString::Format("%s_errhi", rvar->GetName())] = rvar->getAsymErrorHi();
299 }
300 else if (rvar->hasError()) {
301 bs->fVarVals[TString::Format("%s_err", rvar->GetName())] = rvar->getError();
302 }
303 }
304 myTree.Fill();
305 }
306
307 delete bs;
308 }
309
311 TTree* myTree = new TTree(name, desc);
312 FillTree(*myTree, data);
313 return myTree;
314 }
315
316
317 // useful function to print in one line the content of a set with their values
318 void PrintListContent(const RooArgList & l, std::ostream & os ) {
319 bool first = true;
320 os << "( ";
321 for (int i = 0; i< l.getSize(); ++i) {
322 if (first) {
323 first=kFALSE ;
324 } else {
325 os << ", " ;
326 }
327 l[i].printName(os);
328 os << " = ";
329 l[i].printValue(os);
330 }
331 os << ")\n";
332 }
333}
ROOT::R::TRInterface & r
Definition: Object.C:4
#define b(i)
Definition: RSha256.hxx:100
#define oocoutW(o, a)
Definition: RooMsgService.h:48
#define oocoutE(o, a)
Definition: RooMsgService.h:49
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
#define NamespaceImp(name)
Definition: Rtypes.h:380
char name[80]
Definition: TGX11.cxx:109
double sqrt(double)
double log(double)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:71
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) const
Test whether we depend on (ie, are served by) any object in the specified collection.
Definition: RooAbsArg.cxx:734
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition: RooAbsArg.h:83
TIterator * serverIterator() const R__SUGGEST_ALTERNATIVE("Use servers() and begin()
virtual TObject * clone(const char *newname=0) const =0
void SetName(const char *name)
Set the name of the TNamed.
Definition: RooAbsArg.cxx:2224
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
virtual void setBin(Int_t ibin, const char *rangeName=0)
Set category to i-th fit bin, which is the i-th registered state.
virtual Int_t numBins(const char *rangeName) const
Returm the number of fit bins ( = number of types )
virtual const char * getLabel() const
Return label string of current state.
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
RooAbsArg * first() const
RooLinkedListIter iterator(Bool_t dir=kIterForward) const R__SUGGEST_ALTERNATIVE("begin()
TIterator-style iteration over contained elements.
TIterator * createIterator(Bool_t dir=kIterForward) const R__SUGGEST_ALTERNATIVE("begin()
TIterator-style iteration over contained elements.
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:307
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Definition: RooExtendPdf.h:22
A wrapper around TIterator derivatives.
TObject * Next() override
RooProdPdf is an efficient implementation of a product of PDFs of the form.
Definition: RooProdPdf.h:31
const RooArgList & pdfList() const
Definition: RooProdPdf.h:68
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
Double_t getAsymErrorLo() const
Definition: RooRealVar.h:60
Bool_t hasAsymError(Bool_t allowZero=kTRUE) const
Definition: RooRealVar.h:62
Double_t getAsymErrorHi() const
Definition: RooRealVar.h:61
virtual Double_t getValV(const RooArgSet *nset=0) const
Return value of variable.
Definition: RooRealVar.cxx:248
Bool_t hasError(Bool_t allowZero=kTRUE) const
Definition: RooRealVar.h:57
Double_t getError() const
Definition: RooRealVar.h:56
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.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
Definition: ModelConfig.h:249
RooAbsPdf * GetPdf() const
get model PDF (return NULL if pdf has not been specified or does not exist)
Definition: ModelConfig.h:234
A TTree is a list of TBranches.
Definition: TBranch.h:91
virtual void ResetAddress()
Reset the address of the branch.
Definition: TBranch.cxx:2488
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition: TString.cxx:2311
A TTree represents a columnar dataset.
Definition: TTree.h:72
virtual Int_t Fill()
Fill all branches.
Definition: TTree.cxx:4487
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:5170
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition: TTree.h:341
const Int_t n
Definition: legend1.C:16
@ InputArguments
Definition: RooGlobalFunc.h:68
Namespace for the RooStats classes.
Definition: Asimov.h:20
Double_t AsimovSignificance(Double_t s, Double_t b, Double_t sigma_b=0.0)
Compute the Asimov Median significance for a Poisson process with s = expected number of signal event...
TTree * GetAsTTree(TString name, TString desc, const RooDataSet &data)
RooAbsPdf * MakeUnconstrainedPdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name=NULL)
void FillTree(TTree &myTree, const RooDataSet &data)
BranchStore * CreateBranchStore(const RooDataSet &data)
RooAbsPdf * StripConstraints(RooAbsPdf &pdf, const RooArgSet &observables)
RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name)
void FactorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints)
RooStatsConfig & GetGlobalRooStatsConfig()
Retrieve the config object which can be used to set flags for things like offsetting the likelihood o...
void UseNLLOffset(bool on)
Use an offset in NLL calculations.
bool IsNLLOffset()
Test of RooStats should by default offset NLL calculations.
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
static constexpr double s
Definition: first.py:1
auto * l
Definition: textangle.C:4