Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooAbsData.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooAbsData.cxx
19\class RooAbsData
20\ingroup Roofitcore
21
22RooAbsData is the common abstract base class for binned and unbinned
23datasets. The abstract interface defines plotting and tabulating entry
24points for its contents and provides an iterator over its elements
25(bins for binned data sets, data points for unbinned datasets).
26
27### Storing global observables in RooFit datasets
28
29RooFit groups model variables into *observables* and *parameters*, depending on
30if their values are stored in the dataset. For fits with parameter
31constraints, there is a third kind of variables, called *global observables*.
32These represent the results of auxiliary measurements that constrain the
33nuisance parameters. In the RooFit implementation, a likelihood is generally
34the sum of two terms:
35- the likelihood of the data given the parameters, where the normalization set
36 is the set of observables (implemented by RooNLLVar)
37- the constraint term, where the normalization set is the set of *global
38observables* (implemented by RooConstraintSum)
39
40Before this release, the global observable values were always taken from the
41model/pdf. With this release, a mechanism is added to store a snapshot of
42global observables in any RooDataSet or RooDataHist. For toy studies where the
43global observables assume a different values for each toy, the bookkeeping of
44the set of global observables and in particular their values is much easier
45with this change.
46
47Usage example for a model with global observables `g1` and `g2`:
48```
49auto data = model.generate(x, 1000); // data has only the single observables x
50data->setGlobalObservables(g1, g2); // now, data also stores a snapshot of g1 and g2
51
52// If you fit the model to the data, the global observables and their values
53// are taken from the dataset:
54model.fitTo(*data);
55
56// You can still define the set of global observables yourself, but the values
57// will be takes from the dataset if available:
58model.fitTo(*data, GlobalObservables(g1, g2));
59
60// To force `fitTo` to take the global observable values from the model even
61// though they are in the dataset, you can use the new `GlobalObservablesSource`
62// command argument:
63model.fitTo(*data, GlobalObservables(g1, g2), GlobalObservablesSource("model"));
64// The only other allowed value for `GlobalObservablesSource` is "data", which
65// corresponds to the new default behavior explained above.
66```
67
68In case you create a RooFit dataset directly by calling its constructor, you
69can also pass the global observables in a command argument instead of calling
70RooAbsData::setGlobalObservables() later:
71```
72RooDataSet data{"dataset", "dataset", x, RooFit::GlobalObservables(g1, g2)};
73```
74
75To access the set of global observables stored in a RooAbsData, call
76RooAbsData::getGlobalObservables(). It returns a `nullptr` if no global
77observable snapshots are stored in the dataset.
78**/
79
80#include "RooAbsData.h"
81
82#include "TBuffer.h"
83#include "TClass.h"
84#include "TMath.h"
85#include "TTree.h"
86
87#include "RooFormulaVar.h"
88#include "RooCmdConfig.h"
89#include "RooAbsRealLValue.h"
90#include "RooMsgService.h"
91#include "RooMultiCategory.h"
92#include "Roo1DTable.h"
93#include "RooAbsDataStore.h"
94#include "RooVectorDataStore.h"
95#include "RooTreeDataStore.h"
96#include "RooDataHist.h"
98#include "RooCategory.h"
99#include "RooTrace.h"
100#include "RooUniformBinning.h"
101#include "RooSimultaneous.h"
102
103#include "RooRealVar.h"
104#include "RooGlobalFunc.h"
105#include "RooPlot.h"
106#include "RooCurve.h"
107#include "RooHist.h"
108#include "RooHelpers.h"
109
110#include "ROOT/StringUtils.hxx"
111#include "TMatrixDSym.h"
112#include "TPaveText.h"
113#include "TH1.h"
114#include "TH2.h"
115#include "TH3.h"
116#include "Math/Util.h"
117
118#include <iostream>
119#include <memory>
120
121
122using namespace std;
123
125;
126
127static std::map<RooAbsData*,int> _dcc ;
128
130
131////////////////////////////////////////////////////////////////////////////////
132
134{
135 if (RooAbsData::Composite == s) {
136 cout << "Composite storage is not a valid *default* storage type." << endl;
137 } else {
139 }
140}
141
142////////////////////////////////////////////////////////////////////////////////
143
145{
146 return defaultStorageType;
147}
148
149////////////////////////////////////////////////////////////////////////////////
150
152{
153 _dcc[data]++ ;
154 //cout << "RooAbsData(" << data << ") claim incremented to " << _dcc[data] << endl ;
155}
156
157////////////////////////////////////////////////////////////////////////////////
158/// If return value is true variables can be deleted
159
161{
162 if (_dcc[data]>0) {
163 _dcc[data]-- ;
164 }
165
166 //cout << "RooAbsData(" << data << ") claim decremented to " << _dcc[data] << endl ;
167 return (_dcc[data]==0) ;
168}
169
170////////////////////////////////////////////////////////////////////////////////
171/// Default constructor
172
174{
175 claimVars(this) ;
176 _dstore = 0 ;
178
179 RooTrace::create(this) ;
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Constructor from a set of variables. Only fundamental elements of vars
184/// (RooRealVar,RooCategory etc) are stored as part of the dataset
185
186RooAbsData::RooAbsData(std::string_view name, std::string_view title, const RooArgSet& vars, RooAbsDataStore* dstore) :
187 TNamed(TString{name},TString{title}),
188 _vars("Dataset Variables"),
189 _cachedVars("Cached Variables"),
190 _dstore(dstore),
191 _namePtr(nullptr)
192{
193 if (dynamic_cast<RooTreeDataStore *>(dstore)) {
195 } else if (dynamic_cast<RooVectorDataStore *>(dstore)) {
197 } else {
199 }
200 // cout << "created dataset " << this << endl ;
201 claimVars(this);
202
203 // clone the fundamentals of the given data set into internal buffer
204 for (const auto var : vars) {
205 if (!var->isFundamental()) {
206 coutE(InputArguments) << "RooAbsDataStore::initialize(" << GetName()
207 << "): Data set cannot contain non-fundamental types, ignoring " << var->GetName()
208 << endl;
209 throw std::invalid_argument(std::string("Only fundamental variables can be placed into datasets. This is violated for ") + var->GetName());
210 } else {
211 _vars.addClone(*var);
212 }
213 }
214
215 // reconnect any parameterized ranges to internal dataset observables
216 for (auto var : _vars) {
217 var->attachArgs(_vars);
218 }
219
221
222 RooTrace::create(this);
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// Copy constructor
227
228RooAbsData::RooAbsData(const RooAbsData& other, const char* newname) :
229 TNamed(newname ? newname : other.GetName(),other.GetTitle()),
230 RooPrintable(other), _vars(),
231 _cachedVars("Cached Variables"),
232 _namePtr(newname ? RooNameReg::instance().constPtr(newname) : other._namePtr)
233{
234 //cout << "created dataset " << this << endl ;
235 claimVars(this) ;
236 _vars.addClone(other._vars) ;
238 // reconnect any parameterized ranges to internal dataset observables
239 for (auto var : _vars) {
240 var->attachArgs(_vars);
241 }
242
244 if (other._ownedComponents.size()>0) {
245
246 // copy owned components here
247
248 map<string,RooAbsDataStore*> smap ;
249 for (auto& itero : other._ownedComponents) {
250 RooAbsData* dclone = (RooAbsData*) itero.second->Clone();
251 _ownedComponents[itero.first] = dclone;
252 smap[itero.first] = dclone->store();
253 }
254
255 RooCategory* idx = (RooCategory*) _vars.find(*((RooCompositeDataStore*)other.store())->index()) ;
256 _dstore = new RooCompositeDataStore(newname?newname:other.GetName(),other.GetTitle(),_vars,*idx,smap) ;
258
259 } else {
260
261 // Convert to vector store if default is vector
262 _dstore = other._dstore->clone(_vars,newname?newname:other.GetName()) ;
263 storageType = other.storageType;
264 }
265
267
268 RooTrace::create(this) ;
269}
270
272 TNamed::operator=(other);
273 RooPrintable::operator=(other);
274
275 claimVars(this);
276 _vars.Clear();
277 _vars.addClone(other._vars);
278 _namePtr = other._namePtr;
279
280 // reconnect any parameterized ranges to internal dataset observables
281 for (const auto var : _vars) {
282 var->attachDataSet(*this) ;
283 }
284
285
286 if (other._ownedComponents.size()>0) {
287
288 // copy owned components here
289
290 map<string,RooAbsDataStore*> smap ;
291 for (auto& itero : other._ownedComponents) {
292 RooAbsData* dclone = (RooAbsData*) itero.second->Clone();
293 _ownedComponents[itero.first] = dclone;
294 smap[itero.first] = dclone->store();
295 }
296
297 RooCategory* idx = (RooCategory*) _vars.find(*((RooCompositeDataStore*)other.store())->index()) ;
298 _dstore = new RooCompositeDataStore(GetName(), GetTitle(), _vars, *idx, smap);
300
301 } else {
302
303 // Convert to vector store if default is vector
304 _dstore = other._dstore->clone(_vars);
305 storageType = other.storageType;
306 }
307
309
310 return *this;
311}
312
313
315 if (other._globalObservables) {
316 if(_globalObservables == nullptr) _globalObservables = std::make_unique<RooArgSet>();
317 else _globalObservables->clear();
318 other._globalObservables->snapshot(*_globalObservables);
319 } else {
320 _globalObservables.reset(nullptr);
321 }
322}
323
324
325////////////////////////////////////////////////////////////////////////////////
326/// Destructor
327
329{
330 if (releaseVars(this)) {
331 // will cause content to be deleted subsequently in dtor
332 } else {
334 }
335
336 // delete owned contents.
337 delete _dstore ;
338
339 // Delete owned dataset components
340 for(map<std::string,RooAbsData*>::iterator iter = _ownedComponents.begin() ; iter!= _ownedComponents.end() ; ++iter) {
341 delete iter->second ;
342 }
343
344 RooTrace::destroy(this) ;
345}
347////////////////////////////////////////////////////////////////////////////////
348/// Convert tree-based storage to vector-based storage
349
351{
352 if (auto treeStore = dynamic_cast<RooTreeDataStore*>(_dstore)) {
353 _dstore = new RooVectorDataStore(*treeStore, _vars, GetName());
354 delete treeStore;
356 }
357}
358
359////////////////////////////////////////////////////////////////////////////////
360
361Bool_t RooAbsData::changeObservableName(const char* from, const char* to)
362{
363 Bool_t ret = _dstore->changeObservableName(from,to) ;
364
365 RooAbsArg* tmp = _vars.find(from) ;
366 if (tmp) {
367 tmp->SetName(to) ;
368 }
369 return ret ;
370}
371
372////////////////////////////////////////////////////////////////////////////////
373
375{
376 _dstore->fill() ;
377}
378
379////////////////////////////////////////////////////////////////////////////////
380
382{
383 return nullptr != _dstore ? _dstore->numEntries() : 0;
384}
385
386////////////////////////////////////////////////////////////////////////////////
387
389{
390 _dstore->reset() ;
391}
392
393////////////////////////////////////////////////////////////////////////////////
394
395const RooArgSet* RooAbsData::get(Int_t index) const
396{
397 checkInit() ;
398 return _dstore->get(index) ;
399}
400
401////////////////////////////////////////////////////////////////////////////////
402/// Internal method -- Cache given set of functions with data
403
404void RooAbsData::cacheArgs(const RooAbsArg* cacheOwner, RooArgSet& varSet, const RooArgSet* nset, Bool_t skipZeroWeights)
405{
406 _dstore->cacheArgs(cacheOwner,varSet,nset,skipZeroWeights) ;
407}
408
409////////////////////////////////////////////////////////////////////////////////
410/// Internal method -- Remove cached function values
411
413{
416}
417
418////////////////////////////////////////////////////////////////////////////////
419/// Internal method -- Attach dataset copied with cache contents to copied instances of functions
420
421void RooAbsData::attachCache(const RooAbsArg* newOwner, const RooArgSet& cachedVars)
422{
423 _dstore->attachCache(newOwner, cachedVars) ;
424}
425
426////////////////////////////////////////////////////////////////////////////////
427
429{
430 _dstore->setArgStatus(set,active) ;
431}
432
433////////////////////////////////////////////////////////////////////////////////
434/// Control propagation of dirty flags from observables in dataset
435
437{
438 _dstore->setDirtyProp(flag) ;
439}
440
441////////////////////////////////////////////////////////////////////////////////
442/// Create a reduced copy of this dataset. The caller takes ownership of the returned dataset
443///
444/// The following optional named arguments are accepted
445/// <table>
446/// <tr><td> `SelectVars(const RooArgSet& vars)` <td> Only retain the listed observables in the output dataset
447/// <tr><td> `Cut(const char* expression)` <td> Only retain event surviving the given cut expression
448/// <tr><td> `Cut(const RooFormulaVar& expr)` <td> Only retain event surviving the given cut formula
449/// <tr><td> `CutRange(const char* name)` <td> Only retain events inside range with given name. Multiple CutRange
450/// arguments may be given to select multiple ranges
451/// <tr><td> `EventRange(int lo, int hi)` <td> Only retain events with given sequential event numbers
452/// <tr><td> `Name(const char* name)` <td> Give specified name to output dataset
453/// <tr><td> `Title(const char* name)` <td> Give specified title to output dataset
454/// </table>
455
456RooAbsData* RooAbsData::reduce(const RooCmdArg& arg1,const RooCmdArg& arg2,const RooCmdArg& arg3,const RooCmdArg& arg4,
457 const RooCmdArg& arg5,const RooCmdArg& arg6,const RooCmdArg& arg7,const RooCmdArg& arg8)
458{
459 // Define configuration for this method
460 RooCmdConfig pc(Form("RooAbsData::reduce(%s)",GetName())) ;
461 pc.defineString("name","Name",0,"") ;
462 pc.defineString("title","Title",0,"") ;
463 pc.defineString("cutRange","CutRange",0,"") ;
464 pc.defineString("cutSpec","CutSpec",0,"") ;
465 pc.defineObject("cutVar","CutVar",0,0) ;
466 pc.defineInt("evtStart","EventRange",0,0) ;
467 pc.defineInt("evtStop","EventRange",1,std::numeric_limits<int>::max()) ;
468 pc.defineObject("varSel","SelectVars",0,0) ;
469 pc.defineMutex("CutVar","CutSpec") ;
470
471 // Process & check varargs
472 pc.process(arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
473 if (!pc.ok(true)) {
474 return nullptr;
475 }
476
477 // Extract values from named arguments
478 const char* cutRange = pc.getString("cutRange",0,kTRUE) ;
479 const char* cutSpec = pc.getString("cutSpec",0,kTRUE) ;
480 RooFormulaVar* cutVar = static_cast<RooFormulaVar*>(pc.getObject("cutVar",0)) ;
481 Int_t nStart = pc.getInt("evtStart",0) ;
482 Int_t nStop = pc.getInt("evtStop",std::numeric_limits<int>::max()) ;
483 RooArgSet* varSet = static_cast<RooArgSet*>(pc.getObject("varSel")) ;
484 const char* name = pc.getString("name",0,kTRUE) ;
485 const char* title = pc.getString("title",0,kTRUE) ;
486
487 // Make sure varSubset doesn't contain any variable not in this dataset
488 RooArgSet varSubset ;
489 if (varSet) {
490 varSubset.add(*varSet) ;
491 for (const auto arg : varSubset) {
492 if (!_vars.find(arg->GetName())) {
493 coutW(InputArguments) << "RooAbsData::reduce(" << GetName() << ") WARNING: variable "
494 << arg->GetName() << " not in dataset, ignored" << endl ;
495 varSubset.remove(*arg) ;
496 }
497 }
498 } else {
499 varSubset.add(*get()) ;
500 }
501
502 RooAbsData* ret = 0 ;
503 if (cutSpec) {
504
505 RooFormulaVar cutVarTmp(cutSpec,cutSpec,*get()) ;
506 ret = reduceEng(varSubset,&cutVarTmp,cutRange,nStart,nStop,kFALSE) ;
507
508 } else if (cutVar) {
509
510 ret = reduceEng(varSubset,cutVar,cutRange,nStart,nStop,kFALSE) ;
511
512 } else {
513
514 ret = reduceEng(varSubset,0,cutRange,nStart,nStop,kFALSE) ;
515
516 }
517
518 if (!ret) return nullptr;
519
520 if (name) ret->SetName(name) ;
521 if (title) ret->SetTitle(title) ;
522
523 ret->copyGlobalObservables(*this);
524 return ret ;
525}
526
527////////////////////////////////////////////////////////////////////////////////
528/// Create a subset of the data set by applying the given cut on the data points.
529/// The cut expression can refer to any variable in the data set. For cuts involving
530/// other variables, such as intermediate formula objects, use the equivalent
531/// reduce method specifying the as a RooFormulVar reference.
532
534{
535 RooFormulaVar cutVar(cut,cut,*get()) ;
536 RooAbsData* ret = reduceEng(*get(),&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
537 ret->copyGlobalObservables(*this);
538 return ret;
539}
540
541////////////////////////////////////////////////////////////////////////////////
542/// Create a subset of the data set by applying the given cut on the data points.
543/// The 'cutVar' formula variable is used to select the subset of data points to be
544/// retained in the reduced data collection.
545
547{
548 RooAbsData* ret = reduceEng(*get(),&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
549 ret->copyGlobalObservables(*this);
550 return ret;
551}
552
553////////////////////////////////////////////////////////////////////////////////
554/// Create a subset of the data set by applying the given cut on the data points
555/// and reducing the dimensions to the specified set.
556///
557/// The cut expression can refer to any variable in the data set. For cuts involving
558/// other variables, such as intermediate formula objects, use the equivalent
559/// reduce method specifying the as a RooFormulVar reference.
560
561RooAbsData* RooAbsData::reduce(const RooArgSet& varSubset, const char* cut)
562{
563 // Make sure varSubset doesn't contain any variable not in this dataset
564 RooArgSet varSubset2(varSubset) ;
565 for (const auto arg : varSubset) {
566 if (!_vars.find(arg->GetName())) {
567 coutW(InputArguments) << "RooAbsData::reduce(" << GetName() << ") WARNING: variable "
568 << arg->GetName() << " not in dataset, ignored" << endl ;
569 varSubset2.remove(*arg) ;
570 }
571 }
572
573 RooAbsData* ret = nullptr;
574 if (cut && strlen(cut)>0) {
575 RooFormulaVar cutVar(cut, cut, *get(), false);
576 ret = reduceEng(varSubset2,&cutVar,0,0,std::numeric_limits<std::size_t>::max(),false);
577 } else {
578 ret = reduceEng(varSubset2,0,0,0,std::numeric_limits<std::size_t>::max(),false);
579 }
580 ret->copyGlobalObservables(*this);
581 return ret;
582}
583
584////////////////////////////////////////////////////////////////////////////////
585/// Create a subset of the data set by applying the given cut on the data points
586/// and reducing the dimensions to the specified set.
587///
588/// The 'cutVar' formula variable is used to select the subset of data points to be
589/// retained in the reduced data collection.
590
591RooAbsData* RooAbsData::reduce(const RooArgSet& varSubset, const RooFormulaVar& cutVar)
592{
593 // Make sure varSubset doesn't contain any variable not in this dataset
594 RooArgSet varSubset2(varSubset) ;
595 for(RooAbsArg * arg : varSubset) {
596 if (!_vars.find(arg->GetName())) {
597 coutW(InputArguments) << "RooAbsData::reduce(" << GetName() << ") WARNING: variable "
598 << arg->GetName() << " not in dataset, ignored" << endl ;
599 varSubset2.remove(*arg) ;
600 }
601 }
602
603 RooAbsData* ret = reduceEng(varSubset2,&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
604 ret->copyGlobalObservables(*this);
605 return ret;
606}
607
608
609RooPlot* RooAbsData::plotOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2,
610 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
611 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
612{
614 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
615 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
616 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
617 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
618 return plotOn(frame,l) ;
619}
620
621////////////////////////////////////////////////////////////////////////////////
622/// Create and fill a ROOT histogram TH1,TH2 or TH3 with the values of this dataset for the variables with given names
623/// The range of each observable that is histogrammed is always automatically calculated from the distribution in
624/// the dataset. The number of bins can be controlled using the [xyz]bins parameters. For a greater degree of control
625/// use the createHistogram() method below with named arguments
626///
627/// The caller takes ownership of the returned histogram
628
629TH1 *RooAbsData::createHistogram(const char* varNameList, Int_t xbins, Int_t ybins, Int_t zbins) const
630{
631 // Parse list of variable names
632 const auto varNames = ROOT::Split(varNameList, ",:");
633 RooLinkedList argList;
634 RooRealVar* vars[3] = {nullptr, nullptr, nullptr};
635
636 for (unsigned int i = 0; i < varNames.size(); ++i) {
637 if (i >= 3) {
638 coutW(InputArguments) << "RooAbsData::createHistogram(" << GetName() << "): Can only create 3-dimensional histograms. Variable "
639 << i << " " << varNames[i] << " unused." << std::endl;
640 continue;
641 }
642
643 vars[i] = static_cast<RooRealVar*>( get()->find(varNames[i].data()) );
644 if (!vars[i]) {
645 coutE(InputArguments) << "RooAbsData::createHistogram(" << GetName() << ") ERROR: dataset does not contain an observable named " << varNames[i] << std::endl;
646 return nullptr;
647 }
648 }
649
650 if (!vars[0]) {
651 coutE(InputArguments) << "RooAbsData::createHistogram(" << GetName() << "): No variable to be histogrammed in list '" << varNameList << "'" << std::endl;
652 return nullptr;
653 }
654
655 if (xbins<=0 || !vars[0]->hasMax() || !vars[0]->hasMin() ) {
656 argList.Add(RooFit::AutoBinning(xbins==0?vars[0]->numBins():abs(xbins)).Clone()) ;
657 } else {
658 argList.Add(RooFit::Binning(xbins).Clone()) ;
659 }
660
661 if (vars[1]) {
662 if (ybins<=0 || !vars[1]->hasMax() || !vars[1]->hasMin() ) {
663 argList.Add(RooFit::YVar(*vars[1],RooFit::AutoBinning(ybins==0?vars[1]->numBins():abs(ybins))).Clone()) ;
664 } else {
665 argList.Add(RooFit::YVar(*vars[1],RooFit::Binning(ybins)).Clone()) ;
666 }
667 }
668 if (vars[2]) {
669 if (zbins<=0 || !vars[2]->hasMax() || !vars[2]->hasMin() ) {
670 argList.Add(RooFit::ZVar(*vars[2],RooFit::AutoBinning(zbins==0?vars[2]->numBins():abs(zbins))).Clone()) ;
671 } else {
672 argList.Add(RooFit::ZVar(*vars[2],RooFit::Binning(zbins)).Clone()) ;
673 }
674 }
675
676
677 // Call implementation function
678 TH1* result = createHistogram(GetName(), *vars[0], argList);
679
680 // Delete temporary list of RooCmdArgs
681 argList.Delete() ;
682
683 return result ;
684}
685
686
688 const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
689 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
690{
692 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
693 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
694 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
695 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
696
697 return createHistogram(name,xvar,l) ;
698}
699
700////////////////////////////////////////////////////////////////////////////////
701///
702/// This function accepts the following arguments
703///
704/// \param[in] name Name of the ROOT histogram
705/// \param[in] xvar Observable to be mapped on x axis of ROOT histogram
706/// \return Histogram now owned by user.
707///
708/// <table>
709/// <tr><td> `AutoBinning(Int_t nbins, Double_y margin)` <td> Automatically calculate range with given added fractional margin, set binning to nbins
710/// <tr><td> `AutoSymBinning(Int_t nbins, Double_y margin)` <td> Automatically calculate range with given added fractional margin,
711/// with additional constraint that mean of data is in center of range, set binning to nbins
712/// <tr><td> `Binning(const char* name)` <td> Apply binning with given name to x axis of histogram
713/// <tr><td> `Binning(RooAbsBinning& binning)` <td> Apply specified binning to x axis of histogram
714/// <tr><td> `Binning(int nbins, double lo, double hi)` <td> Apply specified binning to x axis of histogram
715///
716/// <tr><td> `YVar(const RooAbsRealLValue& var,...)` <td> Observable to be mapped on y axis of ROOT histogram
717/// <tr><td> `ZVar(const RooAbsRealLValue& var,...)` <td> Observable to be mapped on z axis of ROOT histogram
718/// </table>
719///
720/// The YVar() and ZVar() arguments can be supplied with optional Binning() Auto(Sym)Range() arguments to control the binning of the Y and Z axes, e.g.
721/// ```
722/// createHistogram("histo",x,Binning(-1,1,20), YVar(y,Binning(-1,1,30)), ZVar(z,Binning("zbinning")))
723/// ```
724///
725/// The caller takes ownership of the returned histogram
726
727TH1 *RooAbsData::createHistogram(const char *name, const RooAbsRealLValue& xvar, const RooLinkedList& argListIn) const
728{
729 RooLinkedList argList(argListIn) ;
730
731 // Define configuration for this method
732 RooCmdConfig pc(Form("RooAbsData::createHistogram(%s)",GetName())) ;
733 pc.defineString("cutRange","CutRange",0,"",kTRUE) ;
734 pc.defineString("cutString","CutSpec",0,"") ;
735 pc.defineObject("yvar","YVar",0,0) ;
736 pc.defineObject("zvar","ZVar",0,0) ;
737 pc.allowUndefined() ;
738
739 // Process & check varargs
740 pc.process(argList) ;
741 if (!pc.ok(kTRUE)) {
742 return 0 ;
743 }
744
745 const char* cutSpec = pc.getString("cutString",0,kTRUE) ;
746 const char* cutRange = pc.getString("cutRange",0,kTRUE) ;
747
748 RooArgList vars(xvar) ;
749 RooAbsArg* yvar = static_cast<RooAbsArg*>(pc.getObject("yvar")) ;
750 if (yvar) {
751 vars.add(*yvar) ;
752 }
753 RooAbsArg* zvar = static_cast<RooAbsArg*>(pc.getObject("zvar")) ;
754 if (zvar) {
755 vars.add(*zvar) ;
756 }
757
758 pc.stripCmdList(argList,"CutRange,CutSpec") ;
759
760 // Swap Auto(Sym)RangeData with a Binning command
761 RooLinkedList ownedCmds ;
762 RooCmdArg* autoRD = (RooCmdArg*) argList.find("AutoRangeData") ;
763 if (autoRD) {
765 if (!getRange((RooRealVar&)xvar,xmin,xmax,autoRD->getDouble(0),autoRD->getInt(0))) {
766 RooCmdArg* bincmd = (RooCmdArg*) RooFit::Binning(autoRD->getInt(1),xmin,xmax).Clone() ;
767 ownedCmds.Add(bincmd) ;
768 argList.Replace(autoRD,bincmd) ;
769 }
770 }
771
772 if (yvar) {
773 RooCmdArg* autoRDY = (RooCmdArg*) ((RooCmdArg*)argList.find("YVar"))->subArgs().find("AutoRangeData") ;
774 if (autoRDY) {
776 if (!getRange((RooRealVar&)(*yvar),ymin,ymax,autoRDY->getDouble(0),autoRDY->getInt(0))) {
777 RooCmdArg* bincmd = (RooCmdArg*) RooFit::Binning(autoRDY->getInt(1),ymin,ymax).Clone() ;
778 //ownedCmds.Add(bincmd) ;
779 ((RooCmdArg*)argList.find("YVar"))->subArgs().Replace(autoRDY,bincmd) ;
780 }
781 delete autoRDY ;
782 }
783 }
784
785 if (zvar) {
786 RooCmdArg* autoRDZ = (RooCmdArg*) ((RooCmdArg*)argList.find("ZVar"))->subArgs().find("AutoRangeData") ;
787 if (autoRDZ) {
788 Double_t zmin,zmax ;
789 if (!getRange((RooRealVar&)(*zvar),zmin,zmax,autoRDZ->getDouble(0),autoRDZ->getInt(0))) {
790 RooCmdArg* bincmd = (RooCmdArg*) RooFit::Binning(autoRDZ->getInt(1),zmin,zmax).Clone() ;
791 //ownedCmds.Add(bincmd) ;
792 ((RooCmdArg*)argList.find("ZVar"))->subArgs().Replace(autoRDZ,bincmd) ;
793 }
794 delete autoRDZ ;
795 }
796 }
797
798
799 TH1* histo = xvar.createHistogram(name,argList) ;
800 fillHistogram(histo,vars,cutSpec,cutRange) ;
801
802 ownedCmds.Delete() ;
803
804 return histo ;
805}
806
807////////////////////////////////////////////////////////////////////////////////
808/// Construct table for product of categories in catSet
809
810Roo1DTable* RooAbsData::table(const RooArgSet& catSet, const char* cuts, const char* opts) const
811{
812 RooArgSet catSet2 ;
813
814 string prodName("(") ;
815 TIterator* iter = catSet.createIterator() ;
816 RooAbsArg* arg ;
817 while((arg=(RooAbsArg*)iter->Next())) {
818 if (dynamic_cast<RooAbsCategory*>(arg)) {
819 RooAbsCategory* varsArg = dynamic_cast<RooAbsCategory*>(_vars.find(arg->GetName())) ;
820 if (varsArg != 0) catSet2.add(*varsArg) ;
821 else catSet2.add(*arg) ;
822 if (prodName.length()>1) {
823 prodName += " x " ;
824 }
825 prodName += arg->GetName() ;
826 } else {
827 coutW(InputArguments) << "RooAbsData::table(" << GetName() << ") non-RooAbsCategory input argument " << arg->GetName() << " ignored" << endl ;
828 }
829 }
830 prodName += ")" ;
831 delete iter ;
832
833 RooMultiCategory tmp(prodName.c_str(),prodName.c_str(),catSet2) ;
834 return table(tmp,cuts,opts) ;
835}
836
837////////////////////////////////////////////////////////////////////////////////
838/// Print name of dataset
839
840void RooAbsData::printName(ostream& os) const
841{
842 os << GetName() ;
843}
844
845////////////////////////////////////////////////////////////////////////////////
846/// Print title of dataset
847
848void RooAbsData::printTitle(ostream& os) const
849{
850 os << GetTitle() ;
851}
852
853////////////////////////////////////////////////////////////////////////////////
854/// Print class name of dataset
855
856void RooAbsData::printClassName(ostream& os) const
857{
858 os << IsA()->GetName() ;
859}
860
861////////////////////////////////////////////////////////////////////////////////
862
863void RooAbsData::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
864{
865 _dstore->printMultiline(os,contents,verbose,indent) ;
866}
867
868////////////////////////////////////////////////////////////////////////////////
869/// Define default print options, for a given print style
870
872{
874}
875
876////////////////////////////////////////////////////////////////////////////////
877/// Calculate standardized moment.
878///
879/// \param[in] var Variable to be used for calculating the moment.
880/// \param[in] order Order of the moment.
881/// \param[in] cutSpec If specified, the moment is calculated on the subset of the data which pass the C++ cut specification expression 'cutSpec'
882/// \param[in] cutRange If specified, calculate inside the range named 'cutRange' (also applies cut spec)
883/// \return \f$ \frac{\left< \left( X - \left< X \right> \right)^n \right>}{\sigma^n} \f$, where n = order.
884
885Double_t RooAbsData::standMoment(const RooRealVar &var, Double_t order, const char* cutSpec, const char* cutRange) const
886{
887 // Hardwire invariant answer for first and second moment
888 if (order==1) return 0 ;
889 if (order==2) return 1 ;
890
891 return moment(var,order,cutSpec,cutRange) / TMath::Power(sigma(var,cutSpec,cutRange),order) ;
892}
893
894////////////////////////////////////////////////////////////////////////////////
895/// Calculate moment of requested order.
896///
897/// \param[in] var Variable to be used for calculating the moment.
898/// \param[in] order Order of the moment.
899/// \param[in] cutSpec If specified, the moment is calculated on the subset of the data which pass the C++ cut specification expression 'cutSpec'
900/// \param[in] cutRange If specified, calculate inside the range named 'cutRange' (also applies cut spec)
901/// \return \f$ \left< \left( X - \left< X \right> \right)^n \right> \f$ of order \f$n\f$.
902///
903
904Double_t RooAbsData::moment(const RooRealVar& var, Double_t order, const char* cutSpec, const char* cutRange) const
905{
906 Double_t offset = order>1 ? moment(var,1,cutSpec,cutRange) : 0 ;
907 return moment(var,order,offset,cutSpec,cutRange) ;
908
909}
910
911////////////////////////////////////////////////////////////////////////////////
912/// Return the 'order'-ed moment of observable 'var' in this dataset. If offset is non-zero it is subtracted
913/// from the values of 'var' prior to the moment calculation. If cutSpec and/or cutRange are specified
914/// the moment is calculated on the subset of the data which pass the C++ cut specification expression 'cutSpec'
915/// and/or are inside the range named 'cutRange'
916
917Double_t RooAbsData::moment(const RooRealVar& var, Double_t order, Double_t offset, const char* cutSpec, const char* cutRange) const
918{
919 // Lookup variable in dataset
920 auto arg = _vars.find(var.GetName());
921 if (!arg) {
922 coutE(InputArguments) << "RooDataSet::moment(" << GetName() << ") ERROR: unknown variable: " << var.GetName() << endl ;
923 return 0;
924 }
925
926 auto varPtr = dynamic_cast<const RooRealVar*>(arg);
927 // Check if found variable is of type RooRealVar
928 if (!varPtr) {
929 coutE(InputArguments) << "RooDataSet::moment(" << GetName() << ") ERROR: variable " << var.GetName() << " is not of type RooRealVar" << endl ;
930 return 0;
931 }
932
933 // Check if dataset is not empty
934 if(sumEntries(cutSpec, cutRange) == 0.) {
935 coutE(InputArguments) << "RooDataSet::moment(" << GetName() << ") WARNING: empty dataset" << endl ;
936 return 0;
937 }
938
939 // Setup RooFormulaVar for cutSpec if it is present
940 std::unique_ptr<RooFormula> select;
941 if (cutSpec) {
942 select.reset(new RooFormula("select",cutSpec,*get()));
943 }
944
945
946 // Calculate requested moment
948 for(Int_t index= 0; index < numEntries(); index++) {
949 const RooArgSet* vars = get(index) ;
950 if (select && select->eval()==0) continue ;
951 if (cutRange && vars->allInRange(cutRange)) continue ;
952
953 sum += weight() * TMath::Power(varPtr->getVal() - offset,order);
954 }
955
956 return sum/sumEntries(cutSpec, cutRange);
957}
958
959////////////////////////////////////////////////////////////////////////////////
960/// Internal method to check if given RooRealVar maps to a RooRealVar in this dataset
961
962RooRealVar* RooAbsData::dataRealVar(const char* methodname, const RooRealVar& extVar) const
963{
964 // Lookup variable in dataset
965 RooRealVar *xdata = (RooRealVar*) _vars.find(extVar.GetName());
966 if(!xdata) {
967 coutE(InputArguments) << "RooDataSet::" << methodname << "(" << GetName() << ") ERROR: variable : " << extVar.GetName() << " is not in data" << endl ;
968 return 0;
969 }
970 // Check if found variable is of type RooRealVar
971 if (!dynamic_cast<RooRealVar*>(xdata)) {
972 coutE(InputArguments) << "RooDataSet::" << methodname << "(" << GetName() << ") ERROR: variable : " << extVar.GetName() << " is not of type RooRealVar in data" << endl ;
973 return 0;
974 }
975 return xdata;
976}
977
978////////////////////////////////////////////////////////////////////////////////
979/// Internal method to calculate single correlation and covariance elements
980
981Double_t RooAbsData::corrcov(const RooRealVar &x, const RooRealVar &y, const char* cutSpec, const char* cutRange, Bool_t corr) const
982{
983 // Lookup variable in dataset
984 RooRealVar *xdata = dataRealVar(corr?"correlation":"covariance",x) ;
985 RooRealVar *ydata = dataRealVar(corr?"correlation":"covariance",y) ;
986 if (!xdata||!ydata) return 0 ;
987
988 // Check if dataset is not empty
989 if(sumEntries(cutSpec, cutRange) == 0.) {
990 coutW(InputArguments) << "RooDataSet::" << (corr?"correlation":"covariance") << "(" << GetName() << ") WARNING: empty dataset, returning zero" << endl ;
991 return 0;
992 }
993
994 // Setup RooFormulaVar for cutSpec if it is present
995 RooFormula* select = cutSpec ? new RooFormula("select",cutSpec,*get()) : 0 ;
996
997 // Calculate requested moment
998 Double_t xysum(0),xsum(0),ysum(0),x2sum(0),y2sum(0);
999 const RooArgSet* vars ;
1000 for(Int_t index= 0; index < numEntries(); index++) {
1001 vars = get(index) ;
1002 if (select && select->eval()==0) continue ;
1003 if (cutRange && vars->allInRange(cutRange)) continue ;
1004
1005 xysum += weight()*xdata->getVal()*ydata->getVal() ;
1006 xsum += weight()*xdata->getVal() ;
1007 ysum += weight()*ydata->getVal() ;
1008 if (corr) {
1009 x2sum += weight()*xdata->getVal()*xdata->getVal() ;
1010 y2sum += weight()*ydata->getVal()*ydata->getVal() ;
1011 }
1012 }
1013
1014 // Normalize entries
1015 xysum/=sumEntries(cutSpec, cutRange) ;
1016 xsum/=sumEntries(cutSpec, cutRange) ;
1017 ysum/=sumEntries(cutSpec, cutRange) ;
1018 if (corr) {
1019 x2sum/=sumEntries(cutSpec, cutRange) ;
1020 y2sum/=sumEntries(cutSpec, cutRange) ;
1021 }
1022
1023 // Cleanup
1024 if (select) delete select ;
1025
1026 // Return covariance or correlation as requested
1027 if (corr) {
1028 return (xysum-xsum*ysum)/(sqrt(x2sum-(xsum*xsum))*sqrt(y2sum-(ysum*ysum))) ;
1029 } else {
1030 return (xysum-xsum*ysum);
1031 }
1032}
1033
1034////////////////////////////////////////////////////////////////////////////////
1035/// Return covariance matrix from data for given list of observables
1036
1037TMatrixDSym* RooAbsData::corrcovMatrix(const RooArgList& vars, const char* cutSpec, const char* cutRange, Bool_t corr) const
1038{
1039 RooArgList varList ;
1040 TIterator* iter = vars.createIterator() ;
1041 RooRealVar* var ;
1042 while((var=(RooRealVar*)iter->Next())) {
1043 RooRealVar* datavar = dataRealVar("covarianceMatrix",*var) ;
1044 if (!datavar) {
1045 delete iter ;
1046 return 0 ;
1047 }
1048 varList.add(*datavar) ;
1049 }
1050 delete iter ;
1051
1052
1053 // Check if dataset is not empty
1054 if(sumEntries(cutSpec, cutRange) == 0.) {
1055 coutW(InputArguments) << "RooDataSet::covariance(" << GetName() << ") WARNING: empty dataset, returning zero" << endl ;
1056 return 0;
1057 }
1058
1059 // Setup RooFormulaVar for cutSpec if it is present
1060 RooFormula* select = cutSpec ? new RooFormula("select",cutSpec,*get()) : 0 ;
1061
1062 iter = varList.createIterator() ;
1063 TIterator* iter2 = varList.createIterator() ;
1064
1065 TMatrixDSym xysum(varList.getSize()) ;
1066 vector<double> xsum(varList.getSize()) ;
1067 vector<double> x2sum(varList.getSize()) ;
1068
1069 // Calculate <x_i> and <x_i y_j>
1070 for(Int_t index= 0; index < numEntries(); index++) {
1071 const RooArgSet* dvars = get(index) ;
1072 if (select && select->eval()==0) continue ;
1073 if (cutRange && dvars->allInRange(cutRange)) continue ;
1074
1075 RooRealVar* varx, *vary ;
1076 iter->Reset() ;
1077 Int_t ix=0,iy=0 ;
1078 while((varx=(RooRealVar*)iter->Next())) {
1079 xsum[ix] += weight()*varx->getVal() ;
1080 if (corr) {
1081 x2sum[ix] += weight()*varx->getVal()*varx->getVal() ;
1082 }
1083
1084 *iter2=*iter ; iy=ix ;
1085 vary=varx ;
1086 while(vary) {
1087 xysum(ix,iy) += weight()*varx->getVal()*vary->getVal() ;
1088 xysum(iy,ix) = xysum(ix,iy) ;
1089 iy++ ;
1090 vary=(RooRealVar*)iter2->Next() ;
1091 }
1092 ix++ ;
1093 }
1094
1095 }
1096
1097 // Normalize sums
1098 for (Int_t ix=0 ; ix<varList.getSize() ; ix++) {
1099 xsum[ix] /= sumEntries(cutSpec, cutRange) ;
1100 if (corr) {
1101 x2sum[ix] /= sumEntries(cutSpec, cutRange) ;
1102 }
1103 for (Int_t iy=0 ; iy<varList.getSize() ; iy++) {
1104 xysum(ix,iy) /= sumEntries(cutSpec, cutRange) ;
1105 }
1106 }
1107
1108 // Calculate covariance matrix
1109 TMatrixDSym* C = new TMatrixDSym(varList.getSize()) ;
1110 for (Int_t ix=0 ; ix<varList.getSize() ; ix++) {
1111 for (Int_t iy=0 ; iy<varList.getSize() ; iy++) {
1112 (*C)(ix,iy) = xysum(ix,iy)-xsum[ix]*xsum[iy] ;
1113 if (corr) {
1114 (*C)(ix,iy) /= sqrt((x2sum[ix]-(xsum[ix]*xsum[ix]))*(x2sum[iy]-(xsum[iy]*xsum[iy]))) ;
1115 }
1116 }
1117 }
1118
1119 if (select) delete select ;
1120 delete iter ;
1121 delete iter2 ;
1122
1123 return C ;
1124}
1125
1126////////////////////////////////////////////////////////////////////////////////
1127/// Create a RooRealVar containing the mean of observable 'var' in
1128/// this dataset. If cutSpec and/or cutRange are specified the
1129/// moment is calculated on the subset of the data which pass the C++
1130/// cut specification expression 'cutSpec' and/or are inside the
1131/// range named 'cutRange'
1132
1133RooRealVar* RooAbsData::meanVar(const RooRealVar &var, const char* cutSpec, const char* cutRange) const
1134{
1135 // Create a new variable with appropriate strings. The error is calculated as
1136 // RMS/Sqrt(N) which is generally valid.
1137
1138 // Create holder variable for mean
1139 TString name(var.GetName()),title("Mean of ") ;
1140 name.Append("Mean");
1141 title.Append(var.GetTitle());
1142 RooRealVar *meanv= new RooRealVar(name,title,0) ;
1143 meanv->setConstant(kFALSE) ;
1144
1145 // Adjust plot label
1146 TString label("<") ;
1147 label.Append(var.getPlotLabel());
1148 label.Append(">");
1149 meanv->setPlotLabel(label.Data());
1150
1151 // fill in this variable's value and error
1152 Double_t meanVal=moment(var,1,0,cutSpec,cutRange) ;
1153 Double_t N(sumEntries(cutSpec,cutRange)) ;
1154
1155 Double_t rmsVal= sqrt(moment(var,2,meanVal,cutSpec,cutRange)*N/(N-1));
1156 meanv->setVal(meanVal) ;
1157 meanv->setError(N > 0 ? rmsVal/sqrt(N) : 0);
1158
1159 return meanv;
1160}
1161
1162////////////////////////////////////////////////////////////////////////////////
1163/// Create a RooRealVar containing the RMS of observable 'var' in
1164/// this dataset. If cutSpec and/or cutRange are specified the
1165/// moment is calculated on the subset of the data which pass the C++
1166/// cut specification expression 'cutSpec' and/or are inside the
1167/// range named 'cutRange'
1168
1169RooRealVar* RooAbsData::rmsVar(const RooRealVar &var, const char* cutSpec, const char* cutRange) const
1170{
1171 // Create a new variable with appropriate strings. The error is calculated as
1172 // RMS/(2*Sqrt(N)) which is only valid if the variable has a Gaussian distribution.
1173
1174 // Create RMS value holder
1175 TString name(var.GetName()),title("RMS of ") ;
1176 name.Append("RMS");
1177 title.Append(var.GetTitle());
1178 RooRealVar *rms= new RooRealVar(name,title,0) ;
1179 rms->setConstant(kFALSE) ;
1180
1181 // Adjust plot label
1182 TString label(var.getPlotLabel());
1183 label.Append("_{RMS}");
1184 rms->setPlotLabel(label);
1185
1186 // Fill in this variable's value and error
1187 Double_t meanVal(moment(var,1,0,cutSpec,cutRange)) ;
1188 Double_t N(sumEntries(cutSpec, cutRange));
1189 Double_t rmsVal= sqrt(moment(var,2,meanVal,cutSpec,cutRange)*N/(N-1));
1190 rms->setVal(rmsVal) ;
1191 rms->setError(rmsVal/sqrt(2*N));
1192
1193 return rms;
1194}
1195
1196////////////////////////////////////////////////////////////////////////////////
1197/// Add a box with statistics information to the specified frame. By default a box with the
1198/// event count, mean and rms of the plotted variable is added.
1199///
1200/// The following optional named arguments are accepted
1201/// <table>
1202/// <tr><td> `What(const char* whatstr)` <td> Controls what is printed: "N" = count, "M" is mean, "R" is RMS.
1203/// <tr><td> `Format(const char* optStr)` <td> \deprecated Classing parameter formatting options, provided for backward compatibility
1204///
1205/// <tr><td> `Format(const char* what,...)` <td> Parameter formatting options.
1206/// <table>
1207/// <tr><td> const char* what <td> Controls what is shown:
1208/// - "N" adds name
1209/// - "E" adds error
1210/// - "A" shows asymmetric error
1211/// - "U" shows unit
1212/// - "H" hides the value
1213/// <tr><td> `FixedPrecision(int n)` <td> Controls precision, set fixed number of digits
1214/// <tr><td> `AutoPrecision(int n)` <td> Controls precision. Number of shown digits is calculated from error + n specified additional digits (1 is sensible default)
1215/// <tr><td> `VerbatimName(Bool_t flag)` <td> Put variable name in a \\verb+ + clause.
1216/// </table>
1217/// <tr><td> `Label(const chat* label)` <td> Add header label to parameter box
1218/// <tr><td> `Layout(Double_t xmin, Double_t xmax, Double_t ymax)` <td> Specify relative position of left,right side of box and top of box. Position of
1219/// bottom of box is calculated automatically from number lines in box
1220/// <tr><td> `Cut(const char* expression)` <td> Apply given cut expression to data when calculating statistics
1221/// <tr><td> `CutRange(const char* rangeName)` <td> Only consider events within given range when calculating statistics. Multiple
1222/// CutRange() argument may be specified to combine ranges.
1223///
1224/// </table>
1225
1226RooPlot* RooAbsData::statOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2,
1227 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
1228 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
1229{
1230 // Stuff all arguments in a list
1231 RooLinkedList cmdList;
1232 cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
1233 cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
1234 cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
1235 cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
1236
1237 // Select the pdf-specific commands
1238 RooCmdConfig pc(Form("RooTreeData::statOn(%s)",GetName())) ;
1239 pc.defineString("what","What",0,"MNR") ;
1240 pc.defineString("label","Label",0,"") ;
1241 pc.defineDouble("xmin","Layout",0,0.65) ;
1242 pc.defineDouble("xmax","Layout",1,0.99) ;
1243 pc.defineInt("ymaxi","Layout",0,Int_t(0.95*10000)) ;
1244 pc.defineString("formatStr","Format",0,"NELU") ;
1245 pc.defineInt("sigDigit","Format",0,2) ;
1246 pc.defineInt("dummy","FormatArgs",0,0) ;
1247 pc.defineString("cutRange","CutRange",0,"",kTRUE) ;
1248 pc.defineString("cutString","CutSpec",0,"") ;
1249 pc.defineMutex("Format","FormatArgs") ;
1250
1251 // Process and check varargs
1252 pc.process(cmdList) ;
1253 if (!pc.ok(kTRUE)) {
1254 return frame ;
1255 }
1256
1257 const char* label = pc.getString("label") ;
1258 Double_t xmin = pc.getDouble("xmin") ;
1259 Double_t xmax = pc.getDouble("xmax") ;
1260 Double_t ymax = pc.getInt("ymaxi") / 10000. ;
1261 const char* formatStr = pc.getString("formatStr") ;
1262 Int_t sigDigit = pc.getInt("sigDigit") ;
1263 const char* what = pc.getString("what") ;
1264
1265 const char* cutSpec = pc.getString("cutString",0,kTRUE) ;
1266 const char* cutRange = pc.getString("cutRange",0,kTRUE) ;
1267
1268 if (pc.hasProcessed("FormatArgs")) {
1269 RooCmdArg* formatCmd = static_cast<RooCmdArg*>(cmdList.FindObject("FormatArgs")) ;
1270 return statOn(frame,what,label,0,0,xmin,xmax,ymax,cutSpec,cutRange,formatCmd) ;
1271 } else {
1272 return statOn(frame,what,label,sigDigit,formatStr,xmin,xmax,ymax,cutSpec,cutRange) ;
1273 }
1274}
1275
1276////////////////////////////////////////////////////////////////////////////////
1277/// Implementation back-end of statOn() method with named arguments
1278
1279RooPlot* RooAbsData::statOn(RooPlot* frame, const char* what, const char *label, Int_t sigDigits,
1281 const char* cutSpec, const char* cutRange, const RooCmdArg* formatCmd)
1282{
1283 Bool_t showLabel= (label != 0 && strlen(label) > 0);
1284
1285 TString whatStr(what) ;
1286 whatStr.ToUpper() ;
1287 Bool_t showN = whatStr.Contains("N") ;
1288 Bool_t showR = whatStr.Contains("R") ;
1289 Bool_t showM = whatStr.Contains("M") ;
1290 Int_t nPar= 0;
1291 if (showN) nPar++ ;
1292 if (showR) nPar++ ;
1293 if (showM) nPar++ ;
1294
1295 // calculate the box's size
1296 Double_t dy(0.06), ymin(ymax-nPar*dy);
1297 if(showLabel) ymin-= dy;
1298
1299 // create the box and set its options
1300 TPaveText *box= new TPaveText(xmin,ymax,xmax,ymin,"BRNDC");
1301 if(!box) return 0;
1302 box->SetName(Form("%s_statBox",GetName())) ;
1303 box->SetFillColor(0);
1304 box->SetBorderSize(1);
1305 box->SetTextAlign(12);
1306 box->SetTextSize(0.04F);
1307 box->SetFillStyle(1001);
1308
1309 // add formatted text for each statistic
1310 RooRealVar N("N","Number of Events",sumEntries(cutSpec,cutRange));
1311 N.setPlotLabel("Entries") ;
1312 RooRealVar *meanv= meanVar(*(RooRealVar*)frame->getPlotVar(),cutSpec,cutRange);
1313 meanv->setPlotLabel("Mean") ;
1314 RooRealVar *rms= rmsVar(*(RooRealVar*)frame->getPlotVar(),cutSpec,cutRange);
1315 rms->setPlotLabel("RMS") ;
1316 TString *rmsText, *meanText, *NText ;
1317 if (options) {
1318 rmsText= rms->format(sigDigits,options);
1319 meanText= meanv->format(sigDigits,options);
1320 NText= N.format(sigDigits,options);
1321 } else {
1322 rmsText= rms->format(*formatCmd);
1323 meanText= meanv->format(*formatCmd);
1324 NText= N.format(*formatCmd);
1325 }
1326 if (showR) box->AddText(rmsText->Data());
1327 if (showM) box->AddText(meanText->Data());
1328 if (showN) box->AddText(NText->Data());
1329
1330 // cleanup heap memory
1331 delete NText;
1332 delete meanText;
1333 delete rmsText;
1334 delete meanv;
1335 delete rms;
1336
1337 // add the optional label if specified
1338 if(showLabel) box->AddText(label);
1339
1340 frame->addObject(box) ;
1341 return frame ;
1342}
1343
1344////////////////////////////////////////////////////////////////////////////////
1345/// Loop over columns of our tree data and fill the input histogram. Returns a pointer to the
1346/// input histogram, or zero in case of an error. The input histogram can be any TH1 subclass, and
1347/// therefore of arbitrary dimension. Variables are matched with the (x,y,...) dimensions of the input
1348/// histogram according to the order in which they appear in the input plotVars list.
1349
1350TH1 *RooAbsData::fillHistogram(TH1 *hist, const RooArgList &plotVars, const char *cuts, const char* cutRange) const
1351{
1352 // Do we have a valid histogram to use?
1353 if(0 == hist) {
1354 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: no valid histogram to fill" << endl;
1355 return 0;
1356 }
1357
1358 // Check that the number of plotVars matches the input histogram's dimension
1359 Int_t hdim= hist->GetDimension();
1360 if(hdim != plotVars.getSize()) {
1361 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: plotVars has the wrong dimension" << endl;
1362 return 0;
1363 }
1364
1365 // Check that the plot variables are all actually RooAbsReal's and print a warning if we do not
1366 // explicitly depend on one of them. Clone any variables that we do not contain directly and
1367 // redirect them to use our event data.
1368 RooArgSet plotClones,localVars;
1369 for(Int_t index= 0; index < plotVars.getSize(); index++) {
1370 const RooAbsArg *var= plotVars.at(index);
1371 const RooAbsReal *realVar= dynamic_cast<const RooAbsReal*>(var);
1372 if(0 == realVar) {
1373 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: cannot plot variable \"" << var->GetName()
1374 << "\" of type " << var->ClassName() << endl;
1375 return 0;
1376 }
1377 RooAbsArg *found= _vars.find(realVar->GetName());
1378 if(!found) {
1379 RooAbsArg *clone= plotClones.addClone(*realVar,kTRUE); // do not complain about duplicates
1380 assert(0 != clone);
1381 if(!clone->dependsOn(_vars)) {
1382 coutE(InputArguments) << ClassName() << "::" << GetName()
1383 << ":fillHistogram: Data does not contain the variable '" << realVar->GetName() << "'." << endl;
1384 return nullptr;
1385 }
1386 else {
1388 }
1389 localVars.add(*clone);
1390 }
1391 else {
1392 localVars.add(*found);
1393 }
1394 }
1395
1396 // Create selection formula if selection cuts are specified
1397 std::unique_ptr<RooFormula> select;
1398 if (cuts != nullptr && strlen(cuts) > 0) {
1399 select.reset(new RooFormula(cuts, cuts, _vars, false));
1400 if (!select || !select->ok()) {
1401 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: invalid cuts \"" << cuts << "\"" << endl;
1402 return 0 ;
1403 }
1404 }
1405
1406 // Lookup each of the variables we are binning in our tree variables
1407 const RooAbsReal *xvar = 0;
1408 const RooAbsReal *yvar = 0;
1409 const RooAbsReal *zvar = 0;
1410 switch(hdim) {
1411 case 3:
1412 zvar= dynamic_cast<RooAbsReal*>(localVars.find(plotVars.at(2)->GetName()));
1413 assert(0 != zvar);
1414 // fall through to next case...
1415 case 2:
1416 yvar= dynamic_cast<RooAbsReal*>(localVars.find(plotVars.at(1)->GetName()));
1417 assert(0 != yvar);
1418 // fall through to next case...
1419 case 1:
1420 xvar= dynamic_cast<RooAbsReal*>(localVars.find(plotVars.at(0)->GetName()));
1421 assert(0 != xvar);
1422 break;
1423 default:
1424 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: cannot fill histogram with "
1425 << hdim << " dimensions" << endl;
1426 break;
1427 }
1428
1429 // Parse cutRange specification
1430 const auto cutVec = ROOT::Split(cutRange ? cutRange : "", ",");
1431
1432 // Loop over events and fill the histogram
1433 if (hist->GetSumw2()->fN==0) {
1434 hist->Sumw2() ;
1435 }
1436 Int_t nevent= numEntries() ; //(Int_t)_tree->GetEntries();
1437 for(Int_t i=0; i < nevent; ++i) {
1438
1439 //Int_t entryNumber= _tree->GetEntryNumber(i);
1440 //if (entryNumber<0) break;
1441 get(i);
1442
1443 // Apply expression based selection criteria
1444 if (select && select->eval()==0) {
1445 continue ;
1446 }
1447
1448
1449 // Apply range based selection criteria
1450 Bool_t selectByRange = kTRUE ;
1451 if (cutRange) {
1452 for (const auto arg : _vars) {
1453 Bool_t selectThisArg = kFALSE ;
1454 for (auto const& cut : cutVec) {
1455 if (!cut.empty() && arg->inRange(cut.c_str())) {
1456 selectThisArg = kTRUE ;
1457 break ;
1458 }
1459 }
1460 if (!selectThisArg) {
1461 selectByRange = kFALSE ;
1462 break ;
1463 }
1464 }
1465 }
1466
1467 if (!selectByRange) {
1468 // Go to next event in loop over events
1469 continue ;
1470 }
1471
1472 Int_t bin(0);
1473 switch(hdim) {
1474 case 1:
1475 bin= hist->FindBin(xvar->getVal());
1476 hist->Fill(xvar->getVal(),weight()) ;
1477 break;
1478 case 2:
1479 bin= hist->FindBin(xvar->getVal(),yvar->getVal());
1480 static_cast<TH2*>(hist)->Fill(xvar->getVal(),yvar->getVal(),weight()) ;
1481 break;
1482 case 3:
1483 bin= hist->FindBin(xvar->getVal(),yvar->getVal(),zvar->getVal());
1484 static_cast<TH3*>(hist)->Fill(xvar->getVal(),yvar->getVal(),zvar->getVal(),weight()) ;
1485 break;
1486 default:
1487 assert(hdim < 3);
1488 break;
1489 }
1490
1491
1492 Double_t error2 = TMath::Power(hist->GetBinError(bin),2)-TMath::Power(weight(),2) ;
1494 if (we==0) we = weight() ;
1495 error2 += TMath::Power(we,2) ;
1496
1497
1498// Double_t we = weightError(RooAbsData::SumW2) ;
1499// Double_t error2(0) ;
1500// if (we==0) {
1501// we = weight() ; //sqrt(weight()) ;
1502// error2 = TMath::Power(hist->GetBinError(bin),2)-TMath::Power(weight(),2) + TMath::Power(we,2) ;
1503// } else {
1504// error2 = TMath::Power(hist->GetBinError(bin),2)-TMath::Power(weight(),2) + TMath::Power(we,2) ;
1505// }
1506 //hist->AddBinContent(bin,weight());
1507 hist->SetBinError(bin,sqrt(error2)) ;
1508
1509 //cout << "RooTreeData::fillHistogram() bin = " << bin << " weight() = " << weight() << " we = " << we << endl ;
1510
1511 }
1512
1513 return hist;
1514}
1515
1516////////////////////////////////////////////////////////////////////////////////
1517/// Split dataset into subsets based on states of given splitCat in this dataset.
1518/// A TList of RooDataSets is returned in which each RooDataSet is named
1519/// after the state name of splitCat of which it contains the dataset subset.
1520/// The observables splitCat itself is no longer present in the sub datasets.
1521/// If createEmptyDataSets is kFALSE (default) this method only creates datasets for states
1522/// which have at least one entry The caller takes ownership of the returned list and its contents
1523
1524TList* RooAbsData::split(const RooAbsCategory& splitCat, Bool_t createEmptyDataSets) const
1525{
1526 // Sanity check
1527 if (!splitCat.dependsOn(*get())) {
1528 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") ERROR category " << splitCat.GetName()
1529 << " doesn't depend on any variable in this dataset" << endl ;
1530 return 0 ;
1531 }
1532
1533 // Clone splitting category and attach to self
1534 RooAbsCategory* cloneCat =0;
1535 RooArgSet* cloneSet = 0;
1536 if (splitCat.isDerived()) {
1537 cloneSet = (RooArgSet*) RooArgSet(splitCat).snapshot(kTRUE) ;
1538 if (!cloneSet) {
1539 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") Couldn't deep-clone splitting category, abort." << endl ;
1540 return 0 ;
1541 }
1542 cloneCat = (RooAbsCategory*) cloneSet->find(splitCat.GetName()) ;
1543 cloneCat->attachDataSet(*this) ;
1544 } else {
1545 cloneCat = dynamic_cast<RooAbsCategory*>(get()->find(splitCat.GetName())) ;
1546 if (!cloneCat) {
1547 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") ERROR category " << splitCat.GetName()
1548 << " is fundamental and does not appear in this dataset" << endl ;
1549 return 0 ;
1550 }
1551 }
1552
1553 // Split a dataset in a series of subsets, each corresponding
1554 // to a state of splitCat
1555 TList* dsetList = new TList ;
1556
1557 // Construct set of variables to be included in split sets = full set - split category
1558 RooArgSet subsetVars(*get()) ;
1559 if (splitCat.isDerived()) {
1560 RooArgSet* vars = splitCat.getVariables() ;
1561 subsetVars.remove(*vars,kTRUE,kTRUE) ;
1562 delete vars ;
1563 } else {
1564 subsetVars.remove(splitCat,kTRUE,kTRUE) ;
1565 }
1566
1567 // Add weight variable explicitly if dataset has weights, but no top-level weight
1568 // variable exists (can happen with composite datastores)
1569 Bool_t addWV(kFALSE) ;
1570 RooRealVar newweight("weight","weight",-1e9,1e9) ;
1571 if (isWeighted() && !IsA()->InheritsFrom(RooDataHist::Class())) {
1572 subsetVars.add(newweight) ;
1573 addWV = kTRUE ;
1574 }
1575
1576 // If createEmptyDataSets is true, prepopulate with empty sets corresponding to all states
1577 if (createEmptyDataSets) {
1578 for (const auto& nameIdx : *cloneCat) {
1579 RooAbsData* subset = emptyClone(nameIdx.first.c_str(), nameIdx.first.c_str(), &subsetVars,(addWV?"weight":0)) ;
1580 dsetList->Add((RooAbsArg*)subset) ;
1581 }
1582 }
1583
1584
1585 // Loop over dataset and copy event to matching subset
1586 const bool propWeightSquared = isWeighted();
1587 for (Int_t i = 0; i < numEntries(); ++i) {
1588 const RooArgSet* row = get(i);
1589 RooAbsData* subset = (RooAbsData*) dsetList->FindObject(cloneCat->getCurrentLabel());
1590 if (!subset) {
1591 subset = emptyClone(cloneCat->getCurrentLabel(),cloneCat->getCurrentLabel(),&subsetVars,(addWV?"weight":0));
1592 dsetList->Add((RooAbsArg*)subset);
1593 }
1594 if (!propWeightSquared) {
1595 subset->add(*row, weight());
1596 } else {
1597 subset->add(*row, weight(), weightSquared());
1598 }
1599 }
1600
1601 delete cloneSet;
1602 return dsetList;
1603}
1604
1605////////////////////////////////////////////////////////////////////////////////
1606/// Split dataset into subsets based on the categorisation of the RooSimultaneous
1607/// A TList of RooDataSets is returned in which each RooDataSet is named
1608/// after the state name of splitCat of which it contains the dataset subset.
1609/// The observables splitCat itself is no longer present in the sub datasets, as well as the
1610/// observables of the other categories.
1611/// If createEmptyDataSets is kFALSE (default) this method only creates datasets for states
1612/// which have at least one entry The caller takes ownership of the returned list and its contents
1613
1614TList* RooAbsData::split(const RooSimultaneous& simpdf, Bool_t createEmptyDataSets) const
1615{
1616 auto& splitCat = const_cast<RooAbsCategoryLValue&>(simpdf.indexCat());
1617
1618 // Sanity check
1619 if (!splitCat.dependsOn(*get())) {
1620 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") ERROR category " << splitCat.GetName()
1621 << " doesn't depend on any variable in this dataset" << endl ;
1622 return nullptr;
1623 }
1624
1625 // Clone splitting category and attach to self
1626 RooAbsCategory* cloneCat =0;
1627 std::unique_ptr<RooArgSet> cloneSet;
1628 if (splitCat.isDerived()) {
1629 cloneSet.reset(static_cast<RooArgSet*>(RooArgSet(splitCat).snapshot(true)));
1630 if (!cloneSet) {
1631 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") Couldn't deep-clone splitting category, abort." << endl ;
1632 return nullptr;
1633 }
1634 cloneCat = (RooAbsCategory*) cloneSet->find(splitCat.GetName()) ;
1635 cloneCat->attachDataSet(*this) ;
1636 } else {
1637 cloneCat = dynamic_cast<RooAbsCategory*>(get()->find(splitCat.GetName())) ;
1638 if (!cloneCat) {
1639 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") ERROR category " << splitCat.GetName()
1640 << " is fundamental and does not appear in this dataset" << endl ;
1641 return nullptr;
1642 }
1643 }
1644
1645 // Split a dataset in a series of subsets, each corresponding
1646 // to a state of splitCat
1647 TList* dsetList = new TList ;
1648
1649 // Construct set of variables to be included in split sets = full set - split category
1650 RooArgSet subsetVars(*get()) ;
1651 if (splitCat.isDerived()) {
1652 std::unique_ptr<RooArgSet> vars{splitCat.getVariables()};
1653 subsetVars.remove(*vars,kTRUE,kTRUE) ;
1654 } else {
1655 subsetVars.remove(splitCat,kTRUE,kTRUE) ;
1656 }
1657
1658 // Add weight variable explicitly if dataset has weights, but no top-level weight
1659 // variable exists (can happen with composite datastores)
1660 Bool_t addWV(kFALSE) ;
1661 RooRealVar newweight("weight","weight",-1e9,1e9) ;
1662 if (isWeighted() && !IsA()->InheritsFrom(RooDataHist::Class())) {
1663 subsetVars.add(newweight) ;
1664 addWV = kTRUE ;
1665 }
1666
1667 // Get the observables for a given pdf in the RooSimultaneous, or an empty
1668 // RooArgSet if no pdf is set
1669 auto getPdfObservables = [this, &simpdf](const char * label) {
1670 RooArgSet obsSet;
1671 if(RooAbsPdf* catPdf = simpdf.getPdf(label)) {
1672 catPdf->getObservables(this->get(), obsSet);
1673 }
1674 return obsSet;
1675 };
1676
1677 // By default, remove all category observables from the subdatasets
1678 RooArgSet allObservables;
1679 for( const auto& catPair : splitCat) {
1680 allObservables.add(getPdfObservables(catPair.first.c_str()));
1681 }
1682 subsetVars.remove(allObservables, kTRUE, kTRUE);
1683
1684
1685 // If createEmptyDataSets is true, prepopulate with empty sets corresponding to all states
1686 if (createEmptyDataSets) {
1687 for (const auto& nameIdx : *cloneCat) {
1688 // Add in the subset only the observables corresponding to this category
1689 RooArgSet subsetVarsCat(subsetVars);
1690 subsetVarsCat.add(getPdfObservables(nameIdx.first.c_str()));
1691 RooAbsData* subset = emptyClone(nameIdx.first.c_str(), nameIdx.first.c_str(), &subsetVarsCat,(addWV?"weight":0)) ;
1692 dsetList->Add(subset) ;
1693 }
1694 }
1695
1696
1697 // Loop over dataset and copy event to matching subset
1698 const bool propWeightSquared = isWeighted();
1699 for (Int_t i = 0; i < numEntries(); ++i) {
1700 const RooArgSet* row = get(i);
1701 RooAbsData* subset = (RooAbsData*) dsetList->FindObject(cloneCat->getCurrentLabel());
1702 if (!subset) {
1703 // Add in the subset only the observables corresponding to this category
1704 RooArgSet subsetVarsCat(subsetVars);
1705 subsetVarsCat.add(getPdfObservables(cloneCat->getCurrentLabel()));
1706 subset = emptyClone(cloneCat->getCurrentLabel(),cloneCat->getCurrentLabel(),&subsetVarsCat,(addWV?"weight":0));
1707 dsetList->Add(subset);
1708 }
1709 if (!propWeightSquared) {
1710 subset->add(*row, weight());
1711 } else {
1712 subset->add(*row, weight(), weightSquared());
1713 }
1714 }
1715
1716 return dsetList;
1717}
1718
1719////////////////////////////////////////////////////////////////////////////////
1720/// Plot dataset on specified frame.
1721///
1722/// By default:
1723/// - An unbinned dataset will use the default binning of the target frame.
1724/// - A binned dataset will retain its intrinsic binning.
1725///
1726/// The following optional named arguments can be used to modify the behaviour:
1727/// \note Please follow the function links in the left column to learn about PyROOT specifics for a given option.
1728///
1729/// <table>
1730///
1731/// <tr><th> <th> Data representation options
1732/// <tr><td> RooFit::Asymmetry(const RooCategory& c)
1733/// <td> Show the asymmetry of the data in given two-state category [F(+)-F(-)] / [F(+)+F(-)].
1734/// Category must have two states with indices -1 and +1 or three states with indices -1,0 and +1.
1735/// <tr><td> RooFit::Efficiency(const RooCategory& c)
1736/// <td> Show the efficiency F(acc)/[F(acc)+F(rej)]. Category must have two states with indices 0 and 1
1737/// <tr><td> RooFit::DataError(Int_t)
1738/// <td> Select the type of error drawn:
1739/// - `Auto(default)` results in Poisson for unweighted data and SumW2 for weighted data
1740/// - `Poisson` draws asymmetric Poisson confidence intervals.
1741/// - `SumW2` draws symmetric sum-of-weights error ( \f$ \left( \sum w \right)^2 / \sum\left(w^2\right) \f$ )
1742/// - `None` draws no error bars
1743/// <tr><td> RooFit::Binning(int nbins, double xlo, double xhi)
1744/// <td> Use specified binning to draw dataset
1745/// <tr><td> RooFit::Binning(const RooAbsBinning&)
1746/// <td> Use specified binning to draw dataset
1747/// <tr><td> RooFit::Binning(const char* name)
1748/// <td> Use binning with specified name to draw dataset
1749/// <tr><td> RooFit::RefreshNorm()
1750/// <td> Force refreshing for PDF normalization information in frame.
1751/// If set, any subsequent PDF will normalize to this dataset, even if it is
1752/// not the first one added to the frame. By default only the 1st dataset
1753/// added to a frame will update the normalization information
1754/// <tr><td> RooFit::Rescale(Double_t f)
1755/// <td> Rescale drawn histogram by given factor.
1756/// <tr><td> RooFit::Cut(const char*)
1757/// <td> Only plot entries that pass the given cut.
1758/// Apart from cutting in continuous variables `Cut("x>5")`, this can also be used to plot a specific
1759/// category state. Use something like `Cut("myCategory == myCategory::stateA")`, where
1760/// `myCategory` resolves to the state number for a given entry and
1761/// `myCategory::stateA` resolves to the state number of the state named "stateA".
1762///
1763/// <tr><td> RooFit::CutRange(const char*)
1764/// <td> Only plot data from given range. Separate multiple ranges with ",".
1765/// \note This often requires passing the normalisation when plotting the PDF because RooFit does not save
1766/// how many events were being plotted (it will only work for cutting slices out of uniformly distributed
1767/// variables).
1768/// ```
1769/// data->plotOn(frame01, CutRange("SB1"));
1770/// const double nData = data->sumEntries("", "SB1");
1771/// // Make clear that the target normalisation is nData. The enumerator NumEvent
1772/// // is needed to switch between relative and absolute scaling.
1773/// model.plotOn(frame01, Normalization(nData, RooAbsReal::NumEvent),
1774/// ProjectionRange("SB1"));
1775/// ```
1776///
1777/// <tr><th> <th> Histogram drawing options
1778/// <tr><td> RooFit::DrawOption(const char* opt)
1779/// <td> Select ROOT draw option for resulting TGraph object
1780/// <tr><td> RooFit::LineStyle(Style_t style)
1781/// <td> Select line style by ROOT line style code, default is solid
1782/// <tr><td> RooFit::LineColor(Color_t color)
1783/// <td> Select line color by ROOT color code, default is black
1784/// <tr><td> RooFit::LineWidth(Width_t width)
1785/// <td> Select line with in pixels, default is 3
1786/// <tr><td> RooFit::MarkerStyle(Style_t style)
1787/// <td> Select the ROOT marker style, default is 21
1788/// <tr><td> RooFit::MarkerColor(Color_t color)
1789/// <td> Select the ROOT marker color, default is black
1790/// <tr><td> RooFit::MarkerSize(Size_t size)
1791/// <td> Select the ROOT marker size
1792/// <tr><td> RooFit::FillStyle(Style_t style)
1793/// <td> Select fill style, default is filled.
1794/// <tr><td> RooFit::FillColor(Color_t color)
1795/// <td> Select fill color by ROOT color code
1796/// <tr><td> RooFit::XErrorSize(Double_t frac)
1797/// <td> Select size of X error bar as fraction of the bin width, default is 1
1798///
1799/// <tr><th> <th> Misc. other options
1800/// <tr><td> RooFit::Name(const char* name)
1801/// <td> Give curve specified name in frame. Useful if curve is to be referenced later
1802/// <tr><td> RooFit::Invisible()
1803/// <td> Add curve to frame, but do not display. Useful in combination AddTo()
1804/// <tr><td> RooFit::AddTo(const char* name, double wgtSel, double wgtOther)
1805/// <td> Add constructed histogram to already existing histogram with given name and relative weight factors
1806///
1807/// </table>
1808
1809RooPlot* RooAbsData::plotOn(RooPlot* frame, const RooLinkedList& argList) const
1810{
1811 // New experimental plotOn() with varargs...
1812
1813 // Define configuration for this method
1814 RooCmdConfig pc(Form("RooAbsData::plotOn(%s)",GetName())) ;
1815 pc.defineString("drawOption","DrawOption",0,"P") ;
1816 pc.defineString("cutRange","CutRange",0,"",kTRUE) ;
1817 pc.defineString("cutString","CutSpec",0,"") ;
1818 pc.defineString("histName","Name",0,"") ;
1819 pc.defineObject("cutVar","CutVar",0) ;
1820 pc.defineObject("binning","Binning",0) ;
1821 pc.defineString("binningName","BinningName",0,"") ;
1822 pc.defineInt("nbins","BinningSpec",0,100) ;
1823 pc.defineDouble("xlo","BinningSpec",0,0) ;
1824 pc.defineDouble("xhi","BinningSpec",1,1) ;
1825 pc.defineObject("asymCat","Asymmetry",0) ;
1826 pc.defineObject("effCat","Efficiency",0) ;
1827 pc.defineInt("lineColor","LineColor",0,-999) ;
1828 pc.defineInt("lineStyle","LineStyle",0,-999) ;
1829 pc.defineInt("lineWidth","LineWidth",0,-999) ;
1830 pc.defineInt("markerColor","MarkerColor",0,-999) ;
1831 pc.defineInt("markerStyle","MarkerStyle",0,-999) ;
1832 pc.defineDouble("markerSize","MarkerSize",0,-999) ;
1833 pc.defineInt("fillColor","FillColor",0,-999) ;
1834 pc.defineInt("fillStyle","FillStyle",0,-999) ;
1835 pc.defineInt("errorType","DataError",0,(Int_t)RooAbsData::Auto) ;
1836 pc.defineInt("histInvisible","Invisible",0,0) ;
1837 pc.defineInt("refreshFrameNorm","RefreshNorm",0,1) ;
1838 pc.defineString("addToHistName","AddTo",0,"") ;
1839 pc.defineDouble("addToWgtSelf","AddTo",0,1.) ;
1840 pc.defineDouble("addToWgtOther","AddTo",1,1.) ;
1841 pc.defineDouble("xErrorSize","XErrorSize",0,1.) ;
1842 pc.defineDouble("scaleFactor","Rescale",0,1.) ;
1843 pc.defineMutex("DataError","Asymmetry","Efficiency") ;
1844 pc.defineMutex("Binning","BinningName","BinningSpec") ;
1845
1846 // Process & check varargs
1847 pc.process(argList) ;
1848 if (!pc.ok(kTRUE)) {
1849 return frame ;
1850 }
1851
1852 PlotOpt o ;
1853
1854 // Extract values from named arguments
1855 o.drawOptions = pc.getString("drawOption") ;
1856 o.cuts = pc.getString("cutString") ;
1857 if (pc.hasProcessed("Binning")) {
1858 o.bins = (RooAbsBinning*) pc.getObject("binning") ;
1859 } else if (pc.hasProcessed("BinningName")) {
1860 o.bins = &frame->getPlotVar()->getBinning(pc.getString("binningName")) ;
1861 } else if (pc.hasProcessed("BinningSpec")) {
1862 Double_t xlo = pc.getDouble("xlo") ;
1863 Double_t xhi = pc.getDouble("xhi") ;
1864 o.bins = new RooUniformBinning((xlo==xhi)?frame->getPlotVar()->getMin():xlo,
1865 (xlo==xhi)?frame->getPlotVar()->getMax():xhi,pc.getInt("nbins")) ;
1866 }
1867 const RooAbsCategoryLValue* asymCat = (const RooAbsCategoryLValue*) pc.getObject("asymCat") ;
1868 const RooAbsCategoryLValue* effCat = (const RooAbsCategoryLValue*) pc.getObject("effCat") ;
1869 o.etype = (RooAbsData::ErrorType) pc.getInt("errorType") ;
1870 o.histInvisible = pc.getInt("histInvisible") ;
1871 o.xErrorSize = pc.getDouble("xErrorSize") ;
1872 o.cutRange = pc.getString("cutRange",0,kTRUE) ;
1873 o.histName = pc.getString("histName",0,kTRUE) ;
1874 o.addToHistName = pc.getString("addToHistName",0,kTRUE) ;
1875 o.addToWgtSelf = pc.getDouble("addToWgtSelf") ;
1876 o.addToWgtOther = pc.getDouble("addToWgtOther") ;
1877 o.refreshFrameNorm = pc.getInt("refreshFrameNorm") ;
1878 o.scaleFactor = pc.getDouble("scaleFactor") ;
1879
1880 // Map auto error type to actual type
1881 if (o.etype == Auto) {
1883 if (o.etype == SumW2) {
1884 coutI(InputArguments) << "RooAbsData::plotOn(" << GetName()
1885 << ") INFO: dataset has non-integer weights, auto-selecting SumW2 errors instead of Poisson errors" << endl ;
1886 }
1887 }
1888
1889 if (o.addToHistName && !frame->findObject(o.addToHistName,RooHist::Class())) {
1890 coutE(InputArguments) << "RooAbsData::plotOn(" << GetName() << ") cannot find existing histogram " << o.addToHistName
1891 << " to add to in RooPlot" << endl ;
1892 return frame ;
1893 }
1894
1895 RooPlot* ret ;
1896 if (!asymCat && !effCat) {
1897 ret = plotOn(frame,o) ;
1898 } else if (asymCat) {
1899 ret = plotAsymOn(frame,*asymCat,o) ;
1900 } else {
1901 ret = plotEffOn(frame,*effCat,o) ;
1902 }
1903
1904 Int_t lineColor = pc.getInt("lineColor") ;
1905 Int_t lineStyle = pc.getInt("lineStyle") ;
1906 Int_t lineWidth = pc.getInt("lineWidth") ;
1907 Int_t markerColor = pc.getInt("markerColor") ;
1908 Int_t markerStyle = pc.getInt("markerStyle") ;
1909 Size_t markerSize = pc.getDouble("markerSize") ;
1910 Int_t fillColor = pc.getInt("fillColor") ;
1911 Int_t fillStyle = pc.getInt("fillStyle") ;
1912 if (lineColor!=-999) ret->getAttLine()->SetLineColor(lineColor) ;
1913 if (lineStyle!=-999) ret->getAttLine()->SetLineStyle(lineStyle) ;
1914 if (lineWidth!=-999) ret->getAttLine()->SetLineWidth(lineWidth) ;
1915 if (markerColor!=-999) ret->getAttMarker()->SetMarkerColor(markerColor) ;
1916 if (markerStyle!=-999) ret->getAttMarker()->SetMarkerStyle(markerStyle) ;
1917 if (markerSize!=-999) ret->getAttMarker()->SetMarkerSize(markerSize) ;
1918 if (fillColor!=-999) ret->getAttFill()->SetFillColor(fillColor) ;
1919 if (fillStyle!=-999) ret->getAttFill()->SetFillStyle(fillStyle) ;
1920
1921 if (pc.hasProcessed("BinningSpec")) {
1922 delete o.bins ;
1923 }
1924
1925 return ret ;
1926}
1927
1928////////////////////////////////////////////////////////////////////////////////
1929/// Create and fill a histogram of the frame's variable and append it to the frame.
1930/// The frame variable must be one of the data sets dimensions.
1931///
1932/// The plot range and the number of plot bins is determined by the parameters
1933/// of the plot variable of the frame (RooAbsReal::setPlotRange(), RooAbsReal::setPlotBins()).
1934///
1935/// The optional cut string expression can be used to select the events to be plotted.
1936/// The cut specification may refer to any variable contained in the data set.
1937///
1938/// The drawOptions are passed to the TH1::Draw() method.
1939/// \see RooAbsData::plotOn(RooPlot*,const RooLinkedList&) const
1941{
1942 if(0 == frame) {
1943 coutE(Plotting) << ClassName() << "::" << GetName() << ":plotOn: frame is null" << endl;
1944 return 0;
1945 }
1947 if(0 == var) {
1948 coutE(Plotting) << ClassName() << "::" << GetName()
1949 << ":plotOn: frame does not specify a plot variable" << endl;
1950 return 0;
1951 }
1952
1953 // create and fill a temporary histogram of this variable
1954 TString histName(GetName());
1955 histName.Append("_plot");
1956 std::unique_ptr<TH1> hist;
1957 if (o.bins) {
1958 hist.reset( var->createHistogram(histName.Data(), RooFit::AxisLabel("Events"), RooFit::Binning(*o.bins)) );
1959 } else if (!frame->getPlotVar()->getBinning().isUniform()) {
1960 hist.reset( var->createHistogram(histName.Data(), RooFit::AxisLabel("Events"),
1961 RooFit::Binning(frame->getPlotVar()->getBinning())) );
1962 } else {
1963 hist.reset( var->createHistogram(histName.Data(), "Events",
1964 frame->GetXaxis()->GetXmin(), frame->GetXaxis()->GetXmax(), frame->GetNbinsX()) );
1965 }
1966
1967 // Keep track of sum-of-weights error
1968 hist->Sumw2() ;
1969
1970 if(0 == fillHistogram(hist.get(), RooArgList(*var),o.cuts,o.cutRange)) {
1971 coutE(Plotting) << ClassName() << "::" << GetName()
1972 << ":plotOn: fillHistogram() failed" << endl;
1973 return 0;
1974 }
1975
1976 // If frame has no predefined bin width (event density) it will be adjusted to
1977 // our histograms bin width so we should force that bin width here
1978 Double_t nomBinWidth ;
1979 if (frame->getFitRangeNEvt()==0 && o.bins) {
1980 nomBinWidth = o.bins->averageBinWidth() ;
1981 } else {
1982 nomBinWidth = o.bins ? frame->getFitRangeBinW() : 0 ;
1983 }
1984
1985 // convert this histogram to a RooHist object on the heap
1986 RooHist *graph= new RooHist(*hist,nomBinWidth,1,o.etype,o.xErrorSize,o.correctForBinWidth,o.scaleFactor);
1987 if(0 == graph) {
1988 coutE(Plotting) << ClassName() << "::" << GetName()
1989 << ":plotOn: unable to create a RooHist object" << endl;
1990 return 0;
1991 }
1992
1993 // If the dataset variable has a wide range than the plot variable,
1994 // calculate the number of entries in the dataset in the plot variable fit range
1995 RooAbsRealLValue* dataVar = (RooAbsRealLValue*) _vars.find(var->GetName()) ;
1996 Double_t nEnt(sumEntries()) ;
1997 if (dataVar->getMin()<var->getMin() || dataVar->getMax()>var->getMax()) {
1998 RooAbsData* tmp = ((RooAbsData*)this)->reduce(*var) ;
1999 nEnt = tmp->sumEntries() ;
2000 delete tmp ;
2001 }
2002
2003 // Store the number of entries before the cut, if any was made
2004 if ((o.cuts && strlen(o.cuts)) || o.cutRange) {
2005 coutI(Plotting) << "RooTreeData::plotOn: plotting " << hist->GetSumOfWeights() << " events out of " << nEnt << " total events" << endl ;
2006 graph->setRawEntries(nEnt) ;
2007 }
2008
2009 // Add self to other hist if requested
2010 if (o.addToHistName) {
2011 RooHist* otherGraph = static_cast<RooHist*>(frame->findObject(o.addToHistName,RooHist::Class())) ;
2012
2013 if (!graph->hasIdenticalBinning(*otherGraph)) {
2014 coutE(Plotting) << "RooTreeData::plotOn: ERROR Histogram to be added to, '" << o.addToHistName << "',has different binning" << endl ;
2015 delete graph ;
2016 return frame ;
2017 }
2018
2019 RooHist* sumGraph = new RooHist(*graph,*otherGraph,o.addToWgtSelf,o.addToWgtOther,o.etype) ;
2020 delete graph ;
2021 graph = sumGraph ;
2022 }
2023
2024 // Rename graph if requested
2025 if (o.histName) {
2026 graph->SetName(o.histName) ;
2027 } else {
2028 TString hname(Form("h_%s",GetName())) ;
2029 if (o.cutRange && strlen(o.cutRange)>0) {
2030 hname.Append(Form("_CutRange[%s]",o.cutRange)) ;
2031 }
2032 if (o.cuts && strlen(o.cuts)>0) {
2033 hname.Append(Form("_Cut[%s]",o.cuts)) ;
2034 }
2035 graph->SetName(hname.Data()) ;
2036 }
2037
2038 // initialize the frame's normalization setup, if necessary
2039 frame->updateNormVars(_vars);
2040
2041
2042 // add the RooHist to the specified plot
2044
2045 return frame;
2046}
2047
2048////////////////////////////////////////////////////////////////////////////////
2049/// Create and fill a histogram with the asymmetry N[+] - N[-] / ( N[+] + N[-] ),
2050/// where N(+/-) is the number of data points with asymCat=+1 and asymCat=-1
2051/// as function of the frames variable. The asymmetry category 'asymCat' must
2052/// have exactly 2 (or 3) states defined with index values +1,-1 (and 0)
2053///
2054/// The plot range and the number of plot bins is determined by the parameters
2055/// of the plot variable of the frame (RooAbsReal::setPlotRange(), RooAbsReal::setPlotBins())
2056///
2057/// The optional cut string expression can be used to select the events to be plotted.
2058/// The cut specification may refer to any variable contained in the data set
2059///
2060/// The drawOptions are passed to the TH1::Draw() method
2061
2063{
2064 if(0 == frame) {
2065 coutE(Plotting) << ClassName() << "::" << GetName() << ":plotAsymOn: frame is null" << endl;
2066 return 0;
2067 }
2069 if(0 == var) {
2070 coutE(Plotting) << ClassName() << "::" << GetName()
2071 << ":plotAsymOn: frame does not specify a plot variable" << endl;
2072 return 0;
2073 }
2074
2075 // create and fill temporary histograms of this variable for each state
2076 TString hist1Name(GetName()),hist2Name(GetName());
2077 hist1Name.Append("_plot1");
2078 std::unique_ptr<TH1> hist1, hist2;
2079 hist2Name.Append("_plot2");
2080
2081 if (o.bins) {
2082 hist1.reset( var->createHistogram(hist1Name.Data(), "Events", *o.bins) );
2083 hist2.reset( var->createHistogram(hist2Name.Data(), "Events", *o.bins) );
2084 } else {
2085 hist1.reset( var->createHistogram(hist1Name.Data(), "Events",
2086 frame->GetXaxis()->GetXmin(), frame->GetXaxis()->GetXmax(),
2087 frame->GetNbinsX()) );
2088 hist2.reset( var->createHistogram(hist2Name.Data(), "Events",
2089 frame->GetXaxis()->GetXmin(), frame->GetXaxis()->GetXmax(),
2090 frame->GetNbinsX()) );
2091 }
2092
2093 assert(hist1 && hist2);
2094
2095 TString cuts1,cuts2 ;
2096 if (o.cuts && strlen(o.cuts)) {
2097 cuts1 = Form("(%s)&&(%s>0)",o.cuts,asymCat.GetName());
2098 cuts2 = Form("(%s)&&(%s<0)",o.cuts,asymCat.GetName());
2099 } else {
2100 cuts1 = Form("(%s>0)",asymCat.GetName());
2101 cuts2 = Form("(%s<0)",asymCat.GetName());
2102 }
2103
2104 if(! fillHistogram(hist1.get(), RooArgList(*var),cuts1.Data(),o.cutRange) ||
2105 ! fillHistogram(hist2.get(), RooArgList(*var),cuts2.Data(),o.cutRange)) {
2106 coutE(Plotting) << ClassName() << "::" << GetName()
2107 << ":plotAsymOn: createHistogram() failed" << endl;
2108 return 0;
2109 }
2110
2111 // convert this histogram to a RooHist object on the heap
2112 RooHist *graph= new RooHist(*hist1,*hist2,0,1,o.etype,o.xErrorSize,kFALSE,o.scaleFactor);
2113 graph->setYAxisLabel(Form("Asymmetry in %s",asymCat.GetName())) ;
2114
2115 // initialize the frame's normalization setup, if necessary
2116 frame->updateNormVars(_vars);
2117
2118 // Rename graph if requested
2119 if (o.histName) {
2120 graph->SetName(o.histName) ;
2121 } else {
2122 TString hname(Form("h_%s_Asym[%s]",GetName(),asymCat.GetName())) ;
2123 if (o.cutRange && strlen(o.cutRange)>0) {
2124 hname.Append(Form("_CutRange[%s]",o.cutRange)) ;
2125 }
2126 if (o.cuts && strlen(o.cuts)>0) {
2127 hname.Append(Form("_Cut[%s]",o.cuts)) ;
2128 }
2129 graph->SetName(hname.Data()) ;
2130 }
2131
2132 // add the RooHist to the specified plot
2134
2135 return frame;
2136}
2137
2138////////////////////////////////////////////////////////////////////////////////
2139/// Create and fill a histogram with the efficiency N[1] / ( N[1] + N[0] ),
2140/// where N(1/0) is the number of data points with effCat=1 and effCat=0
2141/// as function of the frames variable. The efficiency category 'effCat' must
2142/// have exactly 2 +1 and 0.
2143///
2144/// The plot range and the number of plot bins is determined by the parameters
2145/// of the plot variable of the frame (RooAbsReal::setPlotRange(), RooAbsReal::setPlotBins())
2146///
2147/// The optional cut string expression can be used to select the events to be plotted.
2148/// The cut specification may refer to any variable contained in the data set
2149///
2150/// The drawOptions are passed to the TH1::Draw() method
2151
2153{
2154 if(0 == frame) {
2155 coutE(Plotting) << ClassName() << "::" << GetName() << ":plotEffOn: frame is null" << endl;
2156 return 0;
2157 }
2159 if(0 == var) {
2160 coutE(Plotting) << ClassName() << "::" << GetName()
2161 << ":plotEffOn: frame does not specify a plot variable" << endl;
2162 return 0;
2163 }
2164
2165 // create and fill temporary histograms of this variable for each state
2166 TString hist1Name(GetName()),hist2Name(GetName());
2167 hist1Name.Append("_plot1");
2168 std::unique_ptr<TH1> hist1, hist2;
2169 hist2Name.Append("_plot2");
2170
2171 if (o.bins) {
2172 hist1.reset( var->createHistogram(hist1Name.Data(), "Events", *o.bins) );
2173 hist2.reset( var->createHistogram(hist2Name.Data(), "Events", *o.bins) );
2174 } else {
2175 hist1.reset( var->createHistogram(hist1Name.Data(), "Events",
2176 frame->GetXaxis()->GetXmin(), frame->GetXaxis()->GetXmax(),
2177 frame->GetNbinsX()) );
2178 hist2.reset( var->createHistogram(hist2Name.Data(), "Events",
2179 frame->GetXaxis()->GetXmin(), frame->GetXaxis()->GetXmax(),
2180 frame->GetNbinsX()) );
2181 }
2182
2183 assert(hist1 && hist2);
2184
2185 TString cuts1,cuts2 ;
2186 if (o.cuts && strlen(o.cuts)) {
2187 cuts1 = Form("(%s)&&(%s==1)",o.cuts,effCat.GetName());
2188 cuts2 = Form("(%s)&&(%s==0)",o.cuts,effCat.GetName());
2189 } else {
2190 cuts1 = Form("(%s==1)",effCat.GetName());
2191 cuts2 = Form("(%s==0)",effCat.GetName());
2192 }
2193
2194 if(! fillHistogram(hist1.get(), RooArgList(*var),cuts1.Data(),o.cutRange) ||
2195 ! fillHistogram(hist2.get(), RooArgList(*var),cuts2.Data(),o.cutRange)) {
2196 coutE(Plotting) << ClassName() << "::" << GetName()
2197 << ":plotEffOn: createHistogram() failed" << endl;
2198 return 0;
2199 }
2200
2201 // convert this histogram to a RooHist object on the heap
2202 RooHist *graph= new RooHist(*hist1,*hist2,0,1,o.etype,o.xErrorSize,kTRUE);
2203 graph->setYAxisLabel(Form("Efficiency of %s=%s", effCat.GetName(), effCat.lookupName(1).c_str()));
2204
2205 // initialize the frame's normalization setup, if necessary
2206 frame->updateNormVars(_vars);
2207
2208 // Rename graph if requested
2209 if (o.histName) {
2210 graph->SetName(o.histName) ;
2211 } else {
2212 TString hname(Form("h_%s_Eff[%s]",GetName(),effCat.GetName())) ;
2213 if (o.cutRange && strlen(o.cutRange)>0) {
2214 hname.Append(Form("_CutRange[%s]",o.cutRange)) ;
2215 }
2216 if (o.cuts && strlen(o.cuts)>0) {
2217 hname.Append(Form("_Cut[%s]",o.cuts)) ;
2218 }
2219 graph->SetName(hname.Data()) ;
2220 }
2221
2222 // add the RooHist to the specified plot
2224
2225 return frame;
2226}
2227
2228////////////////////////////////////////////////////////////////////////////////
2229/// Create and fill a 1-dimensional table for given category column
2230/// This functions is the equivalent of plotOn() for category dimensions.
2231///
2232/// The optional cut string expression can be used to select the events to be tabulated
2233/// The cut specification may refer to any variable contained in the data set
2234///
2235/// The option string is currently not used
2236
2237Roo1DTable* RooAbsData::table(const RooAbsCategory& cat, const char* cuts, const char* /*opts*/) const
2238{
2239 // First see if var is in data set
2240 RooAbsCategory* tableVar = (RooAbsCategory*) _vars.find(cat.GetName()) ;
2241 RooArgSet *tableSet = 0;
2242 Bool_t ownPlotVar(kFALSE) ;
2243 if (!tableVar) {
2244 if (!cat.dependsOn(_vars)) {
2245 coutE(Plotting) << "RooTreeData::Table(" << GetName() << "): Argument " << cat.GetName()
2246 << " is not in dataset and is also not dependent on data set" << endl ;
2247 return 0 ;
2248 }
2249
2250 // Clone derived variable
2251 tableSet = (RooArgSet*) RooArgSet(cat).snapshot(kTRUE) ;
2252 if (!tableSet) {
2253 coutE(Plotting) << "RooTreeData::table(" << GetName() << ") Couldn't deep-clone table category, abort." << endl ;
2254 return 0 ;
2255 }
2256 tableVar = (RooAbsCategory*) tableSet->find(cat.GetName()) ;
2257 ownPlotVar = kTRUE ;
2258
2259 //Redirect servers of derived clone to internal ArgSet representing the data in this set
2260 tableVar->recursiveRedirectServers(_vars) ;
2261 }
2262
2263 TString tableName(GetName()) ;
2264 if (cuts && strlen(cuts)) {
2265 tableName.Append("(") ;
2266 tableName.Append(cuts) ;
2267 tableName.Append(")") ;
2268 }
2269 Roo1DTable* table2 = tableVar->createTable(tableName) ;
2270
2271 // Make cut selector if cut is specified
2272 RooFormulaVar* cutVar = 0;
2273 if (cuts && strlen(cuts)) {
2274 cutVar = new RooFormulaVar("cutVar",cuts,_vars) ;
2275 }
2276
2277 // Dump contents
2278 Int_t nevent= numEntries() ;
2279 for(Int_t i=0; i < nevent; ++i) {
2280 get(i);
2281
2282 if (cutVar && cutVar->getVal()==0) continue ;
2283
2284 table2->fill(*tableVar,weight()) ;
2285 }
2286
2287 if (ownPlotVar) delete tableSet ;
2288 if (cutVar) delete cutVar ;
2289
2290 return table2 ;
2291}
2292
2293////////////////////////////////////////////////////////////////////////////////
2294/// Fill Doubles 'lowest' and 'highest' with the lowest and highest value of
2295/// observable 'var' in this dataset. If the return value is kTRUE and error
2296/// occurred
2297
2298Bool_t RooAbsData::getRange(const RooAbsRealLValue& var, Double_t& lowest, Double_t& highest, Double_t marginFrac, Bool_t symMode) const
2299{
2300 // Lookup variable in dataset
2301 const auto arg = _vars.find(var.GetName());
2302 if (!arg) {
2303 coutE(InputArguments) << "RooDataSet::getRange(" << GetName() << ") ERROR: unknown variable: " << var.GetName() << endl ;
2304 return kTRUE;
2305 }
2306
2307 auto varPtr = dynamic_cast<const RooRealVar*>(arg);
2308 // Check if found variable is of type RooRealVar
2309 if (!varPtr) {
2310 coutE(InputArguments) << "RooDataSet::getRange(" << GetName() << ") ERROR: variable " << var.GetName() << " is not of type RooRealVar" << endl ;
2311 return kTRUE;
2312 }
2313
2314 // Check if dataset is not empty
2315 if(sumEntries() == 0.) {
2316 coutE(InputArguments) << "RooDataSet::getRange(" << GetName() << ") WARNING: empty dataset" << endl ;
2317 return kTRUE;
2318 }
2319
2320 // Look for highest and lowest value
2321 lowest = RooNumber::infinity() ;
2322 highest = -RooNumber::infinity() ;
2323 for (Int_t i=0 ; i<numEntries() ; i++) {
2324 get(i) ;
2325 if (varPtr->getVal()<lowest) {
2326 lowest = varPtr->getVal() ;
2327 }
2328 if (varPtr->getVal()>highest) {
2329 highest = varPtr->getVal() ;
2330 }
2331 }
2332
2333 if (marginFrac>0) {
2334 if (symMode==kFALSE) {
2335
2336 Double_t margin = marginFrac*(highest-lowest) ;
2337 lowest -= margin ;
2338 highest += margin ;
2339 if (lowest<var.getMin()) lowest = var.getMin() ;
2340 if (highest>var.getMax()) highest = var.getMax() ;
2341
2342 } else {
2343
2344 Double_t mom1 = moment(*varPtr,1) ;
2345 Double_t delta = ((highest-mom1)>(mom1-lowest)?(highest-mom1):(mom1-lowest))*(1+marginFrac) ;
2346 lowest = mom1-delta ;
2347 highest = mom1+delta ;
2348 if (lowest<var.getMin()) lowest = var.getMin() ;
2349 if (highest>var.getMax()) highest = var.getMax() ;
2350
2351 }
2352 }
2353
2354 return kFALSE ;
2355}
2356
2357////////////////////////////////////////////////////////////////////////////////
2358/// Prepare dataset for use with cached constant terms listed in
2359/// 'cacheList' of expression 'arg'. Deactivate tree branches
2360/// for any dataset observable that is either not used at all,
2361/// or is used exclusively by cached branch nodes.
2362
2363void RooAbsData::optimizeReadingWithCaching(RooAbsArg& arg, const RooArgSet& cacheList, const RooArgSet& keepObsList)
2364{
2365 RooArgSet pruneSet ;
2366
2367 // Add unused observables in this dataset to pruneSet
2368 pruneSet.add(*get()) ;
2369 RooArgSet* usedObs = arg.getObservables(*this) ;
2370 pruneSet.remove(*usedObs,kTRUE,kTRUE) ;
2371
2372 // Add observables exclusively used to calculate cached observables to pruneSet
2373 TIterator* vIter = get()->createIterator() ;
2374 RooAbsArg *var ;
2375 while ((var=(RooAbsArg*) vIter->Next())) {
2376 if (allClientsCached(var,cacheList)) {
2377 pruneSet.add(*var) ;
2378 }
2379 }
2380 delete vIter ;
2381
2382
2383 if (pruneSet.getSize()!=0) {
2384
2385 // Go over all used observables and check if any of them have parameterized
2386 // ranges in terms of pruned observables. If so, remove those observable
2387 // from the pruning list
2388 for(auto const* rrv : dynamic_range_cast<RooRealVar*>(*usedObs)) {
2389 if (rrv && !rrv->getBinning().isShareable()) {
2390 RooArgSet depObs ;
2391 RooAbsReal* loFunc = rrv->getBinning().lowBoundFunc() ;
2392 RooAbsReal* hiFunc = rrv->getBinning().highBoundFunc() ;
2393 if (loFunc) {
2394 loFunc->leafNodeServerList(&depObs,0,kTRUE) ;
2395 }
2396 if (hiFunc) {
2397 hiFunc->leafNodeServerList(&depObs,0,kTRUE) ;
2398 }
2399 if (depObs.getSize()>0) {
2400 pruneSet.remove(depObs,kTRUE,kTRUE) ;
2401 }
2402 }
2403 }
2404 }
2405
2406
2407 // Remove all observables in keep list from prune list
2408 pruneSet.remove(keepObsList,kTRUE,kTRUE) ;
2409
2410 if (pruneSet.getSize()!=0) {
2411
2412 // Deactivate tree branches here
2413 cxcoutI(Optimization) << "RooTreeData::optimizeReadingForTestStatistic(" << GetName() << "): Observables " << pruneSet
2414 << " in dataset are either not used at all, orserving exclusively p.d.f nodes that are now cached, disabling reading of these observables for TTree" << endl ;
2415 setArgStatus(pruneSet,kFALSE) ;
2416 }
2417
2418 delete usedObs ;
2419
2420}
2421
2422////////////////////////////////////////////////////////////////////////////////
2423/// Utility function that determines if all clients of object 'var'
2424/// appear in given list of cached nodes.
2425
2427{
2428 Bool_t ret(kTRUE), anyClient(kFALSE) ;
2429
2430 for (const auto client : var->valueClients()) {
2431 anyClient = kTRUE ;
2432 if (!cacheList.find(client->GetName())) {
2433 // If client is not cached recurse
2434 ret &= allClientsCached(client,cacheList) ;
2435 }
2436 }
2437
2438 return anyClient?ret:kFALSE ;
2439}
2440
2441////////////////////////////////////////////////////////////////////////////////
2442
2444{
2445 _dstore->attachBuffers(extObs) ;
2446}
2447
2448////////////////////////////////////////////////////////////////////////////////
2449
2451{
2453}
2454
2455////////////////////////////////////////////////////////////////////////////////
2456
2458{
2459 if (_ownedComponents.size()>0) {
2460 return kTRUE ;
2461 }
2462 return kFALSE ;
2463}
2464
2465////////////////////////////////////////////////////////////////////////////////
2466
2468{
2469 map<string,RooAbsData*>::iterator i = _ownedComponents.find(name) ;
2470 if (i==_ownedComponents.end()) return 0 ;
2471 return i->second ;
2472}
2473
2474////////////////////////////////////////////////////////////////////////////////
2475
2476void RooAbsData::addOwnedComponent(const char* idxlabel, RooAbsData& data)
2477{
2478 _ownedComponents[idxlabel]= &data ;
2479}
2480
2481////////////////////////////////////////////////////////////////////////////////
2482/// Stream an object of class RooAbsData.
2483
2484void RooAbsData::Streamer(TBuffer &R__b)
2485{
2486 if (R__b.IsReading()) {
2487 R__b.ReadClassBuffer(RooAbsData::Class(),this);
2489
2490 // Convert on the fly to vector storage if that the current working default
2493 }
2494
2495 } else {
2496 R__b.WriteClassBuffer(RooAbsData::Class(),this);
2497 }
2498}
2499
2500////////////////////////////////////////////////////////////////////////////////
2501
2503{
2504 _dstore->checkInit() ;
2505}
2506
2507////////////////////////////////////////////////////////////////////////////////
2508/// Forward draw command to data store
2509
2511{
2512 if (_dstore) _dstore->Draw(option) ;
2513}
2514
2515////////////////////////////////////////////////////////////////////////////////
2516
2518{
2519 return _dstore->hasFilledCache() ;
2520}
2521
2522////////////////////////////////////////////////////////////////////////////////
2523/// Return a pointer to the TTree which stores the data. Returns a nullpointer
2524/// if vector-based storage is used. The RooAbsData remains owner of the tree.
2525/// GetClonedTree() can be used to get a tree even if the internal storage does not use one.
2526
2528{
2530 return _dstore->tree();
2531 } else {
2532 coutW(InputArguments) << "RooAbsData::tree(" << GetName() << ") WARNING: is not of StorageType::Tree. "
2533 << "Use GetClonedTree() instead or convert to tree storage." << endl;
2534 return (TTree *)nullptr;
2535 }
2536}
2537
2538////////////////////////////////////////////////////////////////////////////////
2539/// Return a clone of the TTree which stores the data or create such a tree
2540/// if vector storage is used. The user is responsible for deleting the tree
2541
2543{
2545 auto tmp = const_cast<TTree *>(_dstore->tree());
2546 return tmp->CloneTree();
2547 } else {
2548 RooTreeDataStore buffer(GetName(), GetTitle(), *get(), *_dstore);
2549 return buffer.tree().CloneTree();
2550 }
2551}
2552
2553////////////////////////////////////////////////////////////////////////////////
2554/// Convert vector-based storage to tree-based storage
2555
2557{
2560 delete _dstore;
2561 _dstore = newStore;
2563 }
2564}
2565
2566////////////////////////////////////////////////////////////////////////////////
2567/// If one of the TObject we have a referenced to is deleted, remove the
2568/// reference.
2569
2571{
2572 for(auto &iter : _ownedComponents) {
2573 if (iter.second == obj) {
2574 iter.second = nullptr;
2575 }
2576 }
2577}
2578
2579
2580////////////////////////////////////////////////////////////////////////////////
2581/// Sets the global observables stored in this data. A snapshot of the
2582/// observables will be saved.
2583/// \param[in] globalObservables The set of global observables to take a snapshot of.
2584
2585void RooAbsData::setGlobalObservables(RooArgSet const& globalObservables) {
2586 if(_globalObservables == nullptr) _globalObservables = std::make_unique<RooArgSet>();
2587 else _globalObservables->clear();
2588 globalObservables.snapshot(*_globalObservables);
2589 for(auto * arg : *_globalObservables) {
2590 arg->setAttribute("global",true);
2591 // Global observables are also always constant in fits
2592 if(auto lval = dynamic_cast<RooAbsRealLValue*>(arg)) lval->setConstant(true);
2593 if(auto lval = dynamic_cast<RooAbsCategoryLValue*>(arg)) lval->setConstant(true);
2594 }
2595}
2596
2597
2598////////////////////////////////////////////////////////////////////////////////
2599
2600void RooAbsData::SetName(const char* name)
2601{
2603 auto newPtr = RooNameReg::instance().constPtr(GetName()) ;
2604 if (newPtr != _namePtr) {
2605 //cout << "Rename '" << _namePtr->GetName() << "' to '" << name << "' (set flag in new name)" << endl;
2606 _namePtr = newPtr;
2609 }
2610}
2611
2612
2613
2614
2615////////////////////////////////////////////////////////////////////////////////
2616
2617void RooAbsData::SetNameTitle(const char *name, const char *title)
2618{
2619 TNamed::SetTitle(title) ;
2620 SetName(name);
2621}
2622
2623
2624
2625////////////////////////////////////////////////////////////////////////////////
2626/// Return sum of squared weights of this data.
2627
2629 const RooSpan<const double> eventWeights = getWeightBatch(0, numEntries(), /*sumW2=*/true);
2630 if (eventWeights.empty()) {
2631 return numEntries() * weightSquared();
2632 }
2633
2635 for (std::size_t i = 0; i < eventWeights.size(); ++i) {
2636 kahanWeight.AddIndexed(eventWeights[i], i);
2637 }
2638 return kahanWeight.Sum();
2639}
2640
2641
2642////////////////////////////////////////////////////////////////////////////////
2643/// Write information to retrieve data columns into `evalData.spans`.
2644/// All spans belonging to variables of this dataset are overwritten. Spans to other
2645/// variables remain intact.
2646/// \param[out] evalData Store references to all data batches in this struct's `spans`.
2647/// The key to retrieve an item is the pointer of the variable that owns the data.
2648/// \param first Index of first event that ends up in the batch.
2649/// \param len Number of events in each batch.
2650RooAbsData::RealSpans RooAbsData::getBatches(std::size_t begin, std::size_t len) const {
2651 return store()->getBatches(begin, len);
2652}
2653
2654
2656 return store()->getCategoryBatches(first, len);
2657}
static std::map< RooAbsData *, int > _dcc
TMatrixTSym< Double_t > TMatrixDSym
Definition RooAbsData.h:48
#define coutI(a)
#define cxcoutI(a)
#define coutW(a)
#define coutE(a)
int Int_t
Definition RtypesCore.h:45
float Size_t
Definition RtypesCore.h:96
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
static void indent(ostringstream &buf, int indent_level)
#define N
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
char * Form(const char *fmt,...)
The Kahan summation is a compensated summation algorithm, which significantly reduces numerical error...
Definition Util.h:122
T Sum() const
Definition Util.h:240
void AddIndexed(T input, std::size_t index)
Add input to the sum.
Definition Util.h:231
Roo1DTable implements a one-dimensional table.
Definition Roo1DTable.h:23
virtual void fill(RooAbsCategory &cat, Double_t weight=1.0)
Increment the counter of the table slot with the name corresponding to that of the current category s...
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Given a set of possible observables, return the observables that this PDF depends on.
Definition RooAbsArg.h:309
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.
void attachArgs(const RooAbsCollection &set)
Bind this node to objects in set.
void leafNodeServerList(RooAbsCollection *list, const RooAbsArg *arg=0, Bool_t recurseNonDerived=kFALSE) const
Fill supplied list with all leaf nodes of the arg tree, starting with ourself as top node.
const RefCountList_t & valueClients() const
List of all value clients of this object. Value clients receive value updates.
Definition RooAbsArg.h:186
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
virtual Bool_t isDerived() const
Does value or shape of this arg depend on any other arg?
Definition RooAbsArg.h:89
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Recursively replace all servers with the new servers in newSet.
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
void SetName(const char *name)
Set the name of the TNamed.
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
virtual Bool_t isUniform() const
virtual Double_t averageBinWidth() const =0
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
RooAbsCategory is the base class for objects that represent a discrete value with a finite number of ...
virtual const char * getCurrentLabel() const
Return label string of current state.
const std::string & lookupName(value_type index) const
Get the name corresponding to the given index.
Roo1DTable * createTable(const char *label) const
Create a table matching the shape of this category.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Int_t getSize() const
virtual RooAbsArg * addClone(const RooAbsArg &var, Bool_t silent=kFALSE)
Add a clone of the specified argument to list.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
Bool_t allInRange(const char *rangeSpec) const
Return true if all contained object report to have their value inside the specified range.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsDataStore is the abstract base class for data collection that use a TTree as internal storage m...
virtual void reset()=0
virtual void attachBuffers(const RooArgSet &extObs)=0
virtual Bool_t hasFilledCache() const
virtual RooAbsData::CategorySpans getCategoryBatches(std::size_t, std::size_t) const
virtual Bool_t changeObservableName(const char *from, const char *to)=0
virtual const RooArgSet * get(Int_t index) const =0
virtual void resetBuffers()=0
virtual const TTree * tree() const
virtual void cacheArgs(const RooAbsArg *cacheOwner, RooArgSet &varSet, const RooArgSet *nset=0, Bool_t skipZeroWeights=kFALSE)=0
virtual void checkInit() const
virtual RooAbsData::RealSpans getBatches(std::size_t first, std::size_t len) const =0
Retrieve batches for all observables in this data store.
virtual Int_t fill()=0
virtual RooAbsDataStore * clone(const char *newname=0) const =0
virtual void setDirtyProp(Bool_t flag)
virtual void resetCache()=0
void printMultiline(std::ostream &os, Int_t content, Bool_t verbose, TString indent) const override
Detailed printing interface.
virtual void attachCache(const RooAbsArg *newOwner, const RooArgSet &cachedVars)=0
virtual void setArgStatus(const RooArgSet &set, Bool_t active)=0
virtual Int_t numEntries() const =0
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:82
RooRealVar * meanVar(const RooRealVar &var, const char *cutSpec=0, const char *cutRange=0) const
Create a RooRealVar containing the mean of observable 'var' in this dataset.
virtual const RooArgSet * get() const
Definition RooAbsData.h:128
void SetNameTitle(const char *name, const char *title)
Set all the TNamed parameters (name and title).
const TNamed * _namePtr
Definition RooAbsData.h:381
void SetName(const char *name)
Set the name of the TNamed.
virtual void RecursiveRemove(TObject *obj)
If one of the TObject we have a referenced to is deleted, remove the reference.
RooAbsData()
Default constructor.
static void setDefaultStorageType(StorageType s)
virtual Bool_t changeObservableName(const char *from, const char *to)
CategorySpans getCategoryBatches(std::size_t first=0, std::size_t len=std::numeric_limits< std::size_t >::max()) const
RooRealVar * dataRealVar(const char *methodname, const RooRealVar &extVar) const
Internal method to check if given RooRealVar maps to a RooRealVar in this dataset.
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
void setGlobalObservables(RooArgSet const &globalObservables)
Sets the global observables stored in this data.
RooAbsDataStore * store()
Definition RooAbsData.h:104
std::map< RooFit::Detail::DataKey, RooSpan< const double > > RealSpans
Definition RooAbsData.h:159
virtual void reset()
static Bool_t releaseVars(RooAbsData *)
If return value is true variables can be deleted.
TMatrixDSym * corrcovMatrix(const RooArgList &vars, const char *cutSpec, const char *cutRange, Bool_t corr) const
Return covariance matrix from data for given list of observables.
virtual Bool_t isNonPoissonWeighted() const
Definition RooAbsData.h:184
virtual TList * split(const RooAbsCategory &splitCat, Bool_t createEmptyDataSets=kFALSE) const
Split dataset into subsets based on states of given splitCat in this dataset.
virtual Double_t sumEntries() const =0
Return effective number of entries in dataset, i.e., sum all weights.
virtual double weightError(ErrorType=Poisson) const
Return the symmetric error on the current weight.
Definition RooAbsData.h:140
void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Interface for detailed printing of object.
virtual TH1 * fillHistogram(TH1 *hist, const RooArgList &plotVars, const char *cuts="", const char *cutRange=0) const
Loop over columns of our tree data and fill the input histogram.
virtual RooPlot * statOn(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())
Add a box with statistics information to the specified frame.
void checkInit() const
virtual RooPlot * plotEffOn(RooPlot *frame, const RooAbsCategoryLValue &effCat, PlotOpt o) const
Create and fill a histogram with the efficiency N[1] / ( N[1] + N[0] ), where N(1/0) is the number of...
RealSpans getBatches(std::size_t first=0, std::size_t len=std::numeric_limits< std::size_t >::max()) const
Write information to retrieve data columns into evalData.spans.
virtual void optimizeReadingWithCaching(RooAbsArg &arg, const RooArgSet &cacheList, const RooArgSet &keepObsList)
Prepare dataset for use with cached constant terms listed in 'cacheList' of expression 'arg'.
static void claimVars(RooAbsData *)
virtual Double_t weight() const =0
static StorageType defaultStorageType
Definition RooAbsData.h:334
virtual Double_t weightSquared() const =0
Double_t standMoment(const RooRealVar &var, Double_t order, const char *cutSpec=0, const char *cutRange=0) const
Calculate standardized moment.
virtual Bool_t isWeighted() const
Definition RooAbsData.h:180
void addOwnedComponent(const char *idxlabel, RooAbsData &data)
virtual void fill()
RooArgSet _vars
Definition RooAbsData.h:372
virtual void printName(std::ostream &os) const
Print name of dataset.
virtual RooPlot * plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue &asymCat, PlotOpt o) const
Create and fill a histogram with the asymmetry N[+] - N[-] / ( N[+] + N[-] ), where N(+/-) is the num...
RooAbsData * getSimData(const char *idxstate)
virtual void Draw(Option_t *option="")
Forward draw command to data store.
virtual RooAbsData * reduceEng(const RooArgSet &varSubset, const RooFormulaVar *cutVar, const char *cutRange=0, std::size_t nStart=0, std::size_t=std::numeric_limits< std::size_t >::max(), Bool_t copyCache=kTRUE)=0
void copyGlobalObservables(const RooAbsData &other)
De-duplicated name pointer. This will be equal for all objects with the same name.
RooRealVar * rmsVar(const RooRealVar &var, const char *cutSpec=0, const char *cutRange=0) const
Create a RooRealVar containing the RMS of observable 'var' in this dataset.
double sumEntriesW2() const
Return sum of squared weights of this data.
virtual void attachCache(const RooAbsArg *newOwner, const RooArgSet &cachedVars)
Internal method – Attach dataset copied with cache contents to copied instances of functions.
std::map< RooFit::Detail::DataKey, RooSpan< const RooAbsCategory::value_type > > CategorySpans
Definition RooAbsData.h:160
Double_t corrcov(const RooRealVar &x, const RooRealVar &y, const char *cutSpec, const char *cutRange, Bool_t corr) const
Internal method to calculate single correlation and covariance elements.
void convertToVectorStore()
Convert tree-based storage to vector-based storage.
virtual RooSpan< const double > getWeightBatch(std::size_t first, std::size_t len, bool sumW2=false) const =0
Return event weights of all events in range [first, first+len).
RooArgSet _cachedVars
Definition RooAbsData.h:373
virtual void add(const RooArgSet &row, Double_t weight=1, Double_t weightError=0)=0
Bool_t getRange(const RooAbsRealLValue &var, Double_t &lowest, Double_t &highest, Double_t marginFrac=0, Bool_t symMode=kFALSE) const
Fill Doubles 'lowest' and 'highest' with the lowest and highest value of observable 'var' in this dat...
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
virtual void convertToTreeStore()
Convert vector-based storage to tree-based storage.
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 Int_t defaultPrintContents(Option_t *opt) const
Define default print options, for a given print style.
StorageType storageType
Definition RooAbsData.h:336
RooAbsData & operator=(const RooAbsData &other)
Double_t moment(const RooRealVar &var, Double_t order, const char *cutSpec=0, const char *cutRange=0) const
Calculate moment of requested order.
virtual void resetCache()
Internal method – Remove cached function values.
std::unique_ptr< RooArgSet > _globalObservables
Definition RooAbsData.h:379
Bool_t hasFilledCache() const
TTree * GetClonedTree() const
Return a clone of the TTree which stores the data or create such a tree if vector storage is used.
RooAbsData * reduce(const RooCmdArg &arg1, 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())
Create a reduced copy of this dataset.
virtual void cacheArgs(const RooAbsArg *owner, RooArgSet &varSet, const RooArgSet *nset=0, Bool_t skipZeroWeights=kFALSE)
Internal method – Cache given set of functions with data.
void attachBuffers(const RooArgSet &extObs)
virtual RooAbsData * emptyClone(const char *newName=0, const char *newTitle=0, const RooArgSet *vars=0, const char *wgtVarName=0) const =0
virtual void setArgStatus(const RooArgSet &set, Bool_t active)
virtual void printClassName(std::ostream &os) const
Print class name of dataset.
std::map< std::string, RooAbsData * > _ownedComponents
Definition RooAbsData.h:377
TH1 * createHistogram(const char *name, const RooAbsRealLValue &xvar, 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 createHistogram(const char *name, const RooAbsRealLValue& xvar, const RooLinkedList& argList) c...
static StorageType getDefaultStorageType()
void resetBuffers()
Bool_t canSplitFast() const
virtual void printTitle(std::ostream &os) const
Print title of dataset.
RooAbsDataStore * _dstore
External variables cached with this data set.
Definition RooAbsData.h:375
const TTree * tree() const
Return a pointer to the TTree which stores the data.
Bool_t allClientsCached(RooAbsArg *, const RooArgSet &)
Utility function that determines if all clients of object 'var' appear in given list of cached nodes.
void setDirtyProp(Bool_t flag)
Control propagation of dirty flags from observables in dataset.
virtual ~RooAbsData()
Destructor.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const =0
Retrive binning configuration with given name or default binning.
void setConstant(Bool_t value=kTRUE)
Bool_t hasMin(const char *name=0) const
Check if variable has a lower bound.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
TH1 * createHistogram(const char *name, 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
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
const char * getPlotLabel() const
Get the label associated with the variable.
void setPlotLabel(const char *label)
Set the label associated with this variable.
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:35
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:158
RooCategory is an object to represent discrete states.
Definition RooCategory.h:27
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition RooCmdArg.h:27
Double_t getDouble(Int_t idx) const
Definition RooCmdArg.h:84
Int_t getInt(Int_t idx) const
Definition RooCmdArg.h:80
virtual TObject * Clone(const char *newName=0) const
Make a clone of an object using the Streamer facility.
Definition RooCmdArg.h:51
Class RooCmdConfig is a configurable parser for RooCmdArg named arguments.
TObject * getObject(const char *name, TObject *obj=0)
Return TObject property registered with name 'name'.
Bool_t defineInt(const char *name, const char *argName, Int_t intNum, Int_t defValue=0)
Define integer property name 'name' mapped to integer in slot 'intNum' in RooCmdArg with name argName...
void defineMutex(const char *argName1, const char *argName2)
Define arguments named argName1 and argName2 mutually exclusive.
Bool_t defineObject(const char *name, const char *argName, Int_t setNum, const TObject *obj=0, Bool_t isArray=kFALSE)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
const char * getString(const char *name, const char *defaultValue="", Bool_t convEmptyToNull=kFALSE)
Return string property registered with name 'name'.
Int_t getInt(const char *name, Int_t defaultValue=0)
Return integer property registered with name 'name'.
Bool_t defineDouble(const char *name, const char *argName, Int_t doubleNum, Double_t defValue=0.)
Define Double_t property name 'name' mapped to Double_t in slot 'doubleNum' in RooCmdArg with name ar...
Double_t getDouble(const char *name, Double_t defaultValue=0)
Return Double_t property registered with name 'name'.
void allowUndefined(Bool_t flag=kTRUE)
void stripCmdList(RooLinkedList &cmdList, const char *cmdsToPurge)
Utility function that strips command names listed (comma separated) in cmdsToPurge from cmdList.
Bool_t defineString(const char *name, const char *argName, Int_t stringNum, const char *defValue="", Bool_t appendMode=kFALSE)
Define Double_t property name 'name' mapped to Double_t in slot 'stringNum' in RooCmdArg with name ar...
Bool_t ok(Bool_t verbose) const
Return true of parsing was successful.
Bool_t process(const RooCmdArg &arg)
Process given RooCmdArg.
Bool_t hasProcessed(const char *cmdName) const
Return true if RooCmdArg with name 'cmdName' has been processed.
RooCompositeDataStore combines several disjunct datasets into one.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
RooFormula internally uses ROOT's TFormula to compute user-defined expressions of RooAbsArgs.
Definition RooFormula.h:33
Double_t eval(const RooArgSet *nset=0) const
Evalute all parameters/observables, and then evaluate formula.
A RooHist is a graphical representation of binned data based on the TGraphAsymmErrors class.
Definition RooHist.h:27
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
TObject * FindObject(const char *name) const
Return pointer to obejct with given name.
void Delete(Option_t *o=0)
Remove all elements in collection and delete all elements NB: Collection does not own elements,...
TObject * find(const char *name) const
Return pointer to object with given name in collection.
virtual void Add(TObject *arg)
Bool_t Replace(const TObject *oldArg, const TObject *newArg)
Replace object 'oldArg' in collection with new object 'newArg'.
RooMultiCategory connects several RooAbsCategory objects into a single category.
RooNameReg is a registry for const char* names.
Definition RooNameReg.h:25
const TNamed * constPtr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
static RooNameReg & instance()
Return reference to singleton instance.
static void incrementRenameCounter()
The renaming counter has to be incremented every time a RooAbsArg is renamed.
static Double_t infinity()
Return internal infinity representation.
Definition RooNumber.cxx:49
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:44
TAttMarker * getAttMarker(const char *name=0) const
Return a pointer to the marker attributes of the named object in this plot, or zero if the named obje...
Definition RooPlot.cxx:856
void addObject(TObject *obj, Option_t *drawOptions="", Bool_t invisible=kFALSE)
Add a generic object to this plot.
Definition RooPlot.cxx:392
Double_t getFitRangeNEvt() const
Return the number of events in the fit range.
Definition RooPlot.h:138
Double_t getFitRangeBinW() const
Return the bin width that is being used to normalise the PDF.
Definition RooPlot.h:141
TAttFill * getAttFill(const char *name=0) const
Return a pointer to the fill attributes of the named object in this plot, or zero if the named object...
Definition RooPlot.cxx:846
RooAbsRealLValue * getPlotVar() const
Definition RooPlot.h:136
TAttLine * getAttLine(const char *name=0) const
Return a pointer to the line attributes of the named object in this plot, or zero if the named object...
Definition RooPlot.cxx:836
TAxis * GetXaxis() const
Definition RooPlot.cxx:1276
void updateNormVars(const RooArgSet &vars)
Install the given set of observables are reference normalization variables for this frame.
Definition RooPlot.cxx:380
void addPlotable(RooPlotable *plotable, Option_t *drawOptions="", Bool_t invisible=kFALSE, Bool_t refreshNorm=kFALSE)
Add the specified plotable object to our plot.
Definition RooPlot.cxx:540
Int_t GetNbinsX() const
Definition RooPlot.cxx:1280
TObject * findObject(const char *name, const TClass *clas=0) const
Find the named object in our list of items and return a pointer to it.
Definition RooPlot.cxx:952
RooPlotable is a 'mix-in' base class that define the standard RooFit plotting and printing methods.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
void setError(Double_t value)
Definition RooRealVar.h:64
TString * format(const RooCmdArg &formatArg) const
Format contents of RooRealVar for pretty printing on RooPlot parameter boxes.
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.
A simple container to hold a batch of data values.
Definition RooSpan.h:34
constexpr std::span< T >::index_type size() const noexcept
Definition RooSpan.h:121
constexpr bool empty() const noexcept
Definition RooSpan.h:125
static void destroy(const TObject *obj)
Register deletion of object 'obj'.
Definition RooTrace.cxx:81
static void create(const TObject *obj)
Register creation of object 'obj'.
Definition RooTrace.cxx:68
RooTreeDataStore is a TTree-backed data storage.
RooUniformBinning is an implementation of RooAbsBinning that provides a uniform binning in 'n' bins b...
RooVectorDataStore uses std::vectors to store data columns.
Int_t fN
Definition TArray.h:38
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:39
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:42
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:41
Double_t GetXmax() const
Definition TAxis.h:134
Double_t GetXmin() const
Definition TAxis.h:133
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:8893
virtual Int_t GetDimension() const
Definition TH1.h:282
virtual void SetBinError(Int_t bin, Double_t error)
Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empt...
Definition TH1.cxx:9036
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3351
virtual TArrayD * GetSumw2()
Definition TH1.h:312
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
Return Global bin number corresponding to x,y,z.
Definition TH1.cxx:3681
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition TH1.cxx:8850
Service class for 2-D histogram classes.
Definition TH2.h:30
The 3-D histogram classes derived from the 1-D histogram classes.
Definition TH3.h:31
Iterator abstract base class.
Definition TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
A doubly linked list.
Definition TList.h:38
virtual void Add(TObject *obj)
Definition TList.h:81
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition TList.cxx:578
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
virtual const char * GetTitle() const
Returns title of object.
Definition TNamed.h:48
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition TNamed.cxx:74
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Clear(Option_t *="")
Definition TObject.h:119
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:515
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:267
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
void ToUpper()
Change string to upper case.
Definition TString.cxx:1163
TString & Append(const char *cs)
Definition TString.h:564
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual TTree * CloneTree(Long64_t nentries=-1, Option_t *option="")
Create a clone of this tree and copy nentries.
Definition TTree.cxx:3133
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg ZVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg AxisLabel(const char *name)
RooCmdArg AutoBinning(Int_t nbins=100, Double_t marginFactor=0.1)
RooCmdArg Binning(const RooAbsBinning &binning)
const Double_t sigma
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:685
Definition first.py:1
Definition graph.py:1
static const char * what
Definition stlLoader.cc:6
const char * cuts
Definition RooAbsData.h:210
const char * cutRange
Definition RooAbsData.h:214
const char * histName
Definition RooAbsData.h:215
const char * addToHistName
Definition RooAbsData.h:217
RooAbsData::ErrorType etype
Definition RooAbsData.h:213
RooAbsBinning * bins
Definition RooAbsData.h:212
Option_t * drawOptions
Definition RooAbsData.h:211
auto * l
Definition textangle.C:4
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345