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