Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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#include "TBranch.h"
20
21#include "RooArgSet.h"
22#include "RooWorkspace.h"
23#include "RooAbsPdf.h"
24#include "RooUniform.h"
25#include "RooProdPdf.h"
26#include "RooExtendPdf.h"
27#include "RooSimultaneous.h"
30
31using std::endl, std::vector;
32
33
34namespace {
35 template<class listT, class stringT> void getParameterNames(const listT* l,std::vector<stringT>& names){
36 // extract the parameter names from a list
37 if(!l) return;
38 for (auto const *obj : *l) {
39 names.push_back(obj->GetName());
40 }
41 }
42 void getArgs(RooWorkspace* ws, const std::vector<TString> names, RooArgSet& args){
43 for(const auto& p:names){
44 RooRealVar* v = ws->var(p.Data());
45 if(v){
46 args.add(*v);
47 }
48 }
49 }
50 }
51
52namespace RooStats {
53
55 static RooStatsConfig theConfig;
56 return theConfig;
57 }
58
59 double AsimovSignificance(double s, double b, double sigma_b ) {
60 // Asimov significance
61 // formula [10] and [20] from https://www.pp.rhul.ac.uk/~cowan/stat/notes/medsigNote.pdf
62 // case we have a sigma_b
63 double sb2 = sigma_b*sigma_b;
64 // formula below has a large error when sigma_b becomes zero
65 // better to use the approximation for sigma_b=0 for very small values
66 double r = sb2/b;
67 if (r > 1.E-12) {
68 double bpsb2 = b + sb2;
69 double b2 = b*b;
70 double spb = s+b;
71 double za2 = 2.*( (spb)* std::log( ( spb)*(bpsb2)/(b2+ spb*sb2) ) -
72 (b2/sb2) * std::log(1. + ( sb2 * s)/(b * bpsb2) ) );
73 return sqrt(za2);
74
75 }
76 // case when the background (b) is known
77 double za2 = 2.*( (s+b) * std::log(1. + s/b) -s );
78 return std::sqrt(za2);
79 }
80
81 /// Use an offset in NLL calculations.
82 void UseNLLOffset(bool on) {
84 }
85
86 /// Test of RooStats should by default offset NLL calculations.
87 bool IsNLLOffset() {
89 }
90
91 void FactorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints) {
92 // utility function to factorize constraint terms from a pdf
93 // (from G. Petrucciani)
94 if (auto prod = dynamic_cast<RooProdPdf *>(&pdf)) {
95 RooArgList list(prod->pdfList());
96 for (int i = 0, n = list.size(); i < n; ++i) {
97 RooAbsPdf *pdfi = static_cast<RooAbsPdf *>(list.at(i));
98 FactorizePdf(observables, *pdfi, obsTerms, constraints);
99 }
100 } else if (dynamic_cast<RooExtendPdf *>(&pdf)) {
101 // extract underlying pdf which is extended; first server is the pdf; second server is the number of events variable
102 auto iter = pdf.servers().begin();
103 assert(iter != pdf.servers().end());
104 assert(dynamic_cast<RooAbsPdf*>(*iter));
105 FactorizePdf(observables, static_cast<RooAbsPdf&>(**iter), obsTerms, constraints);
106 } else if (auto sim = dynamic_cast<RooSimultaneous *>(&pdf)) { //|| dynamic_cast<RooSimultaneousOpt>(&pdf)) {
107 assert(sim != nullptr);
108 RooAbsCategoryLValue *cat = static_cast<RooAbsCategoryLValue *>(sim->indexCat().clone(sim->indexCat().GetName()));
109 for (int ic = 0, nc = cat->numBins((const char *)nullptr); ic < nc; ++ic) {
110 cat->setBin(ic);
111 RooAbsPdf* catPdf = sim->getPdf(cat->getCurrentLabel());
112 // it is possible that a pdf is not defined for every category
113 if (catPdf != nullptr) FactorizePdf(observables, *catPdf, obsTerms, constraints);
114 }
115 delete cat;
116 } else if (pdf.dependsOn(observables)) {
117 if (!obsTerms.contains(pdf)) obsTerms.add(pdf);
118 } else {
119 if (!constraints.contains(pdf)) constraints.add(pdf);
120 }
121 }
122
123
124 void FactorizePdf(RooStats::ModelConfig &model, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints) {
125 // utility function to factorize constraint terms from a pdf
126 // (from G. Petrucciani)
127 if (!model.GetObservables() ) {
128 oocoutE(nullptr,InputArguments) << "RooStatsUtils::FactorizePdf - invalid input model: missing observables" << endl;
129 return;
130 }
131 return FactorizePdf(*model.GetObservables(), pdf, obsTerms, constraints);
132 }
133
134
135 RooAbsPdf * MakeNuisancePdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name) {
136 // make a nuisance pdf by factorizing out all constraint terms in a common pdf
137 RooArgList obsTerms;
138 RooArgList constraints;
139 FactorizePdf(observables, pdf, obsTerms, constraints);
140 if(constraints.empty()) {
141 oocoutW(nullptr, Eval) << "RooStatsUtils::MakeNuisancePdf - no constraints found on nuisance parameters in the input model" << endl;
142 return nullptr;
143 }
144 return new RooProdPdf(name,"", constraints);
145 }
146
148 // make a nuisance pdf by factorizing out all constraint terms in a common pdf
149 if (!model.GetPdf() || !model.GetObservables() ) {
150 oocoutE(nullptr, InputArguments) << "RooStatsUtils::MakeNuisancePdf - invalid input model: missing pdf and/or observables" << endl;
151 return nullptr;
152 }
153 return MakeNuisancePdf(*model.GetPdf(), *model.GetObservables(), name);
154 }
155
156 RooAbsPdf * StripConstraints(RooAbsPdf &pdf, const RooArgSet &observables) {
157
158 if (auto prod = dynamic_cast<RooProdPdf *>(&pdf)) {
159
160 RooArgList list(prod->pdfList()); RooArgList newList;
161
162 for (int i = 0, n = list.size(); i < n; ++i) {
163 RooAbsPdf *pdfi = static_cast<RooAbsPdf *>(list.at(i));
164 RooAbsPdf *newPdfi = StripConstraints(*pdfi, observables);
165 if(newPdfi != nullptr) newList.add(*newPdfi);
166 }
167
168 if (newList.empty()) {
169 return nullptr; // only constraints in product
170 // return single component (no longer a product)
171 } else if (newList.size() == 1) {
172 return dynamic_cast<RooAbsPdf *>(
173 newList.at(0)->clone(TString::Format("%s_unconstrained", newList.at(0)->GetName())));
174 } else {
175 return new RooProdPdf(TString::Format("%s_unconstrained", prod->GetName()).Data(),
176 TString::Format("%s without constraints", prod->GetTitle()).Data(), newList);
177 }
178
179 } else if (dynamic_cast<RooExtendPdf*>(&pdf)) {
180
181 auto iter = pdf.servers().begin();
182 // extract underlying pdf which is extended; first server is the pdf; second server is the number of events variable
183 auto uPdf = dynamic_cast<RooAbsPdf *>(*(iter++));
184 auto extended_term = dynamic_cast<RooAbsReal *>(*(iter++));
185 assert(uPdf != nullptr);
186 assert(extended_term != nullptr);
187 assert(iter == pdf.servers().end());
188
189 RooAbsPdf *newUPdf = StripConstraints(*uPdf, observables);
190 if(newUPdf == nullptr) return nullptr; // only constraints in underlying pdf
191 else return new RooExtendPdf(TString::Format("%s_unconstrained", pdf.GetName()).Data(),
192 TString::Format("%s without constraints", pdf.GetTitle()).Data(), *newUPdf, *extended_term);
193
194 } else if (auto sim = dynamic_cast<RooSimultaneous *>(&pdf)) { //|| dynamic_cast<RooSimultaneousOpt *>(&pdf)) {
195
196 assert(sim != nullptr);
197 RooAbsCategoryLValue *cat = static_cast<RooAbsCategoryLValue *>(sim->indexCat().Clone()); assert(cat != nullptr);
198 RooArgList pdfList;
199
200 for (int ic = 0, nc = cat->numBins((const char *)nullptr); ic < nc; ++ic) {
201 cat->setBin(ic);
202 RooAbsPdf* catPdf = sim->getPdf(cat->getCurrentLabel());
203 RooAbsPdf* newPdf = nullptr;
204 // it is possible that a pdf is not defined for every category
205 if (catPdf != nullptr) newPdf = StripConstraints(*catPdf, observables);
206 if (newPdf == nullptr) { delete cat; return nullptr; } // all channels must have observables
207 pdfList.add(*newPdf);
208 }
209
210 return new RooSimultaneous(TString::Format("%s_unconstrained", sim->GetName()).Data(),
211 TString::Format("%s without constraints", sim->GetTitle()).Data(), pdfList, *cat);
212
213 } else if (pdf.dependsOn(observables)) {
214 return static_cast<RooAbsPdf *>(pdf.clone(TString::Format("%s_unconstrained", pdf.GetName()).Data()));
215 }
216
217 return nullptr; // just a constraint term
218 }
219
220 RooAbsPdf * MakeUnconstrainedPdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name) {
221 // make a clone pdf without all constraint terms in a common pdf
222 RooAbsPdf * unconstrainedPdf = StripConstraints(pdf, observables);
223 if(!unconstrainedPdf) {
224 oocoutE(nullptr, InputArguments) << "RooStats::MakeUnconstrainedPdf - invalid observable list passed (observables not found in original pdf) or invalid pdf passed (without observables)" << endl;
225 return nullptr;
226 }
227 if(name != nullptr) unconstrainedPdf->SetName(name);
228 return unconstrainedPdf;
229 }
230
232 // make a clone pdf without all constraint terms in a common pdf
233 if(!model.GetPdf() || !model.GetObservables()) {
234 oocoutE(nullptr, InputArguments) << "RooStatsUtils::MakeUnconstrainedPdf - invalid input model: missing pdf and/or observables" << endl;
235 return nullptr;
236 }
237 return MakeUnconstrainedPdf(*model.GetPdf(), *model.GetObservables(), name);
238 }
239
240 // Helper class for GetAsTTree
242 public:
243 std::map<TString, double> fVarVals;
244 double fInval;
245 TTree *fTree = nullptr;
246
247 BranchStore(const vector<TString> &params = vector<TString>(), double _inval = -999.) : fInval(_inval)
248 {
249 for (unsigned int i = 0; i < params.size(); i++)
250 fVarVals[params[i]] = _inval;
251 }
252
254 if (fTree) {
255 for(std::map<TString, double>::iterator it = fVarVals.begin();it!=fVarVals.end();++it) {
256 TBranch *br = fTree->GetBranch( it->first );
257 if (br) br->ResetAddress();
258 }
259 }
260 }
261
262 void AssignToTTree(TTree &myTree) {
263 fTree = &myTree;
264 for(std::map<TString, double>::iterator it = fVarVals.begin();it!=fVarVals.end();++it) {
265 const TString& name = it->first;
266 myTree.Branch( name, &fVarVals[name], TString::Format("%s/D", name.Data()));
267 }
268 }
269 void ResetValues() {
270 for(std::map<TString, double>::iterator it = fVarVals.begin();it!=fVarVals.end();++it) {
271 const TString& name = it->first;
273 }
274 }
275 };
276
278 if (data.numEntries() == 0) {
279 return new BranchStore;
280 }
281 vector <TString> V;
282 for (auto *rvar : dynamic_range_cast<RooRealVar *>(* data.get(0))) {
283 if (rvar == nullptr)
284 continue;
285 V.push_back(rvar->GetName());
286 if (rvar->hasAsymError()) {
287 V.push_back(TString::Format("%s_errlo", rvar->GetName()));
288 V.push_back(TString::Format("%s_errhi", rvar->GetName()));
289 }
290 else if (rvar->hasError()) {
291 V.push_back(TString::Format("%s_err", rvar->GetName()));
292 }
293 }
294 return new BranchStore(V);
295 }
296
297 void FillTree(TTree &myTree, const RooDataSet &data) {
299 bs->AssignToTTree(myTree);
300
301 for(int entry = 0;entry<data.numEntries();entry++) {
302 bs->ResetValues();
303 for (auto const *rvar : dynamic_range_cast<RooRealVar *>(*data.get(entry))) {
304 if (rvar == nullptr)
305 continue;
306 bs->fVarVals[rvar->GetName()] = rvar->getValV();
307 if (rvar->hasAsymError()) {
308 bs->fVarVals[TString::Format("%s_errlo", rvar->GetName())] = rvar->getAsymErrorLo();
309 bs->fVarVals[TString::Format("%s_errhi", rvar->GetName())] = rvar->getAsymErrorHi();
310 }
311 else if (rvar->hasError()) {
312 bs->fVarVals[TString::Format("%s_err", rvar->GetName())] = rvar->getError();
313 }
314 }
315 myTree.Fill();
316 }
317
318 delete bs;
319 }
320
322 TTree* myTree = new TTree(name, desc);
323 FillTree(*myTree, data);
324 return myTree;
325 }
326
327
328 // useful function to print in one line the content of a set with their values
329 void PrintListContent(const RooArgList & l, std::ostream & os ) {
330 bool first = true;
331 os << "( ";
332 for (std::size_t i = 0; i< l.size(); ++i) {
333 if (first) {
334 first=false ;
335 } else {
336 os << ", " ;
337 }
338 l[i].printName(os);
339 os << " = ";
340 l[i].printValue(os);
341 }
342 os << ")\n";
343 }
344
345 // clone a workspace, copying all needed components and discarding all others
346 // start off with the old workspace
347 RooWorkspace* MakeReducedWorkspace(RooWorkspace *oldWS, const char *newName,
348 bool copySnapshots, const char *mcname,
349 const char *newmcname, bool copyData) {
350 auto objects = oldWS->allGenericObjects();
351 RooStats::ModelConfig *oldMC =
352 dynamic_cast<RooStats::ModelConfig *>(oldWS->obj(mcname));
353 for (auto it : objects) {
354 if (!oldMC) {
355 oldMC = dynamic_cast<RooStats::ModelConfig *>(it);
356 }
357 }
358 if (!oldMC)
359 throw std::runtime_error("unable to retrieve ModelConfig");
360
361 RooAbsPdf *origPdf = oldMC->GetPdf();
362
363 // start off with the old modelconfig
364 std::vector<TString> poilist;
365 std::vector<TString> nplist;
366 std::vector<TString> obslist;
367 std::vector<TString> globobslist;
368 RooAbsPdf *pdf = nullptr;
369 if (oldMC) {
370 pdf = oldMC->GetPdf();
371 ::getParameterNames(oldMC->GetParametersOfInterest(), poilist);
372 ::getParameterNames(oldMC->GetNuisanceParameters(), nplist);
373 ::getParameterNames(oldMC->GetObservables(), obslist);
374 ::getParameterNames(oldMC->GetGlobalObservables(), globobslist);
375 }
376 if (!pdf) {
377 if (origPdf)
378 pdf = origPdf;
379 }
380 if (!pdf) {
381 return nullptr;
382 }
383
384 // create them anew
385 RooWorkspace *newWS = new RooWorkspace(newName ? newName : oldWS->GetName());
386 newWS->autoImportClassCode(true);
387 RooStats::ModelConfig *newMC = new RooStats::ModelConfig(newmcname, newWS);
388
389 // Copy snapshots
390 if (copySnapshots) {
391 for (auto* snap : static_range_cast<RooArgSet const*>(oldWS->getSnapshots())) {
392 newWS->saveSnapshot(snap->GetName(), *snap, /*importValuesdf*/true);
393 }
394 }
395
396 newWS->import(*pdf, RooFit::RecycleConflictNodes());
397 RooAbsPdf *newPdf = newWS->pdf(pdf->GetName());
398 newMC->SetPdf(*newPdf);
399
400 if (copyData) {
401 for (auto d : oldWS->allData()) {
402 newWS->import(*d);
403 }
404 }
405
406 RooArgSet poiset;
407 ::getArgs(newWS, poilist, poiset);
408 RooArgSet npset;
409 ::getArgs(newWS, nplist, npset);
410 RooArgSet obsset;
411 ::getArgs(newWS, obslist, obsset);
412 RooArgSet globobsset;
413 ::getArgs(newWS, globobslist, globobsset);
414
415 newMC->SetParametersOfInterest(poiset);
416 newMC->SetNuisanceParameters(npset);
417 newMC->SetObservables(obsset);
418 newMC->SetGlobalObservables(globobsset);
419 newWS->import(*newMC);
420
421 return newWS;
422 }
423
424}
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define oocoutW(o, a)
#define oocoutE(o, a)
#define NamespaceImp(name)
Definition Rtypes.h:398
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void on
char name[80]
Definition TGX11.cxx:110
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
virtual TObject * clone(const char *newname=nullptr) const =0
void SetName(const char *name) override
Set the name of the TNamed.
const RefCountList_t & servers() const
List of all servers of this object.
Definition RooAbsArg.h:180
Abstract base class for objects that represent a discrete value that can be set from the outside,...
Int_t numBins(const char *rangeName=nullptr) const override
Return the number of fit bins ( = number of types )
void setBin(Int_t ibin, const char *rangeName=nullptr) override
Set category to i-th fit bin, which is the i-th registered state.
virtual const char * getCurrentLabel() const
Return label string of current state.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Container class to hold unbinned data.
Definition RooDataSet.h:33
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
Efficient implementation of a product of PDFs of the form.
Definition RooProdPdf.h:33
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
std::map< TString, double > fVarVals
void AssignToTTree(TTree &myTree)
BranchStore(const vector< TString > &params=vector< TString >(), double _inval=-999.)
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
virtual void SetObservables(const RooArgSet &set)
Specify the observables.
const RooArgSet * GetGlobalObservables() const
get RooArgSet for global observables (return nullptr if not existing)
virtual void SetParametersOfInterest(const RooArgSet &set)
Specify parameters of interest.
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
virtual void SetNuisanceParameters(const RooArgSet &set)
Specify the nuisance parameters (parameters that are not POI).
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the workspace if not already there.
Definition ModelConfig.h:88
virtual void SetGlobalObservables(const RooArgSet &set)
Specify the global observables.
Persistable container for RooFit projects.
TObject * obj(RooStringView name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
RooAbsPdf * pdf(RooStringView name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
bool saveSnapshot(RooStringView, const char *paramNames)
Save snapshot of values and attributes (including "Constant") of given parameters.
std::list< RooAbsData * > allData() const
Return list of all dataset in the workspace.
RooLinkedList const & getSnapshots() const
std::list< TObject * > allGenericObjects() const
Return list of all generic objects in the workspace.
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.
RooRealVar * var(RooStringView name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
static void autoImportClassCode(bool flag)
If flag is true, source code of classes not the ROOT distribution is automatically imported if on obj...
A TTree is a list of TBranches.
Definition TBranch.h:93
virtual void ResetAddress()
Reset the address of the branch.
Definition TBranch.cxx:2651
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
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:2378
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Int_t Fill()
Fill all branches.
Definition TTree.cxx:4603
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition TTree.cxx:5294
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:353
RooCmdArg RecycleConflictNodes(bool flag=true)
const Int_t n
Definition legend1.C:16
Namespace for the RooStats classes.
Definition Asimov.h:19
TTree * GetAsTTree(TString name, TString desc, const RooDataSet &data)
Create a TTree with the given name and description. All RooRealVars in the RooDataSet are represented...
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)
extract constraint terms from pdf
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...
double AsimovSignificance(double s, double b, double sigma_b=0.0)
Compute the Asimov Median significance for a Poisson process with s = expected number of signal event...
void UseNLLOffset(bool on)
function to set a global flag in RooStats to use NLL offset when performing nll computations Note tha...
RooAbsPdf * MakeUnconstrainedPdf(RooAbsPdf &pdf, const RooArgSet &observables, const char *name=nullptr)
remove constraints from pdf and return the unconstrained pdf
bool IsNLLOffset()
function returning if the flag to check if the flag to use NLLOffset is set
RooWorkspace * MakeReducedWorkspace(RooWorkspace *oldWS, const char *newName, bool copySnapshots, const char *mcname, const char *newmcname, bool copyData=true)
function that clones a workspace, copying all needed components and discarding all others
void PrintListContent(const RooArgList &l, std::ostream &os=std::cout)
useful function to print in one line the content of a set with their values
bool useLikelihoodOffset
Offset the likelihood by passing RooFit::Offset to fitTo().
TLine l
Definition textangle.C:4