Logo ROOT  
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#include "RooFit.h"
30
31#include <iostream>
32
33#include "TBuffer.h"
34#include "TClass.h"
35#include "TMath.h"
36#include "TTree.h"
37
38#include "RooFormulaVar.h"
39#include "RooCmdConfig.h"
40#include "RooAbsRealLValue.h"
41#include "RooMsgService.h"
42#include "RooMultiCategory.h"
43#include "Roo1DTable.h"
44#include "RooAbsDataStore.h"
45#include "RooVectorDataStore.h"
46#include "RooTreeDataStore.h"
47#include "RooDataHist.h"
49#include "RooCategory.h"
50#include "RooTrace.h"
51#include "RooUniformBinning.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) ;
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// Default constructor
116
118{
119 claimVars(this) ;
120 _dstore = 0 ;
122
123 RooTrace::create(this) ;
124}
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)) {
134 storageType = RooAbsData::Tree;
135 } else if (dynamic_cast<RooVectorDataStore *>(dstore)) {
136 storageType = RooAbsData::Vector;
137 } else {
138 storageType = RooAbsData::Composite;
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;
156
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
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,std::numeric_limits<int>::max()) ;
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",std::numeric_limits<int>::max()) ;
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 for (const auto arg : varSubset) {
417 if (!_vars.find(arg->GetName())) {
418 coutW(InputArguments) << "RooAbsData::reduce(" << GetName() << ") WARNING: variable "
419 << arg->GetName() << " not in dataset, ignored" << endl ;
420 varSubset.remove(*arg) ;
421 }
422 }
423 } else {
424 varSubset.add(*get()) ;
425 }
426
427 RooAbsData* ret = 0 ;
428 if (cutSpec) {
429
430 RooFormulaVar cutVarTmp(cutSpec,cutSpec,*get()) ;
431 ret = reduceEng(varSubset,&cutVarTmp,cutRange,nStart,nStop,kFALSE) ;
432
433 } else if (cutVar) {
434
435 ret = reduceEng(varSubset,cutVar,cutRange,nStart,nStop,kFALSE) ;
436
437 } else {
438
439 ret = reduceEng(varSubset,0,cutRange,nStart,nStop,kFALSE) ;
440
441 }
442
443 if (!ret) return 0 ;
444
445 if (name) {
446 ret->SetName(name) ;
447 }
448 if (title) {
449 ret->SetTitle(title) ;
450 }
451
452 return ret ;
453}
454
455////////////////////////////////////////////////////////////////////////////////
456/// Create a subset of the data set by applying the given cut on the data points.
457/// The cut expression can refer to any variable in the data set. For cuts involving
458/// other variables, such as intermediate formula objects, use the equivalent
459/// reduce method specifying the as a RooFormulVar reference.
460
462{
463 RooFormulaVar cutVar(cut,cut,*get()) ;
464 return reduceEng(*get(),&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
465}
466
467////////////////////////////////////////////////////////////////////////////////
468/// Create a subset of the data set by applying the given cut on the data points.
469/// The 'cutVar' formula variable is used to select the subset of data points to be
470/// retained in the reduced data collection.
471
473{
474 return reduceEng(*get(),&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
475}
476
477////////////////////////////////////////////////////////////////////////////////
478/// Create a subset of the data set by applying the given cut on the data points
479/// and reducing the dimensions to the specified set.
480///
481/// The cut expression can refer to any variable in the data set. For cuts involving
482/// other variables, such as intermediate formula objects, use the equivalent
483/// reduce method specifying the as a RooFormulVar reference.
484
485RooAbsData* RooAbsData::reduce(const RooArgSet& varSubset, const char* cut)
486{
487 // Make sure varSubset doesn't contain any variable not in this dataset
488 RooArgSet varSubset2(varSubset) ;
489 for (const auto arg : varSubset) {
490 if (!_vars.find(arg->GetName())) {
491 coutW(InputArguments) << "RooAbsData::reduce(" << GetName() << ") WARNING: variable "
492 << arg->GetName() << " not in dataset, ignored" << endl ;
493 varSubset2.remove(*arg) ;
494 }
495 }
496
497 if (cut && strlen(cut)>0) {
498 RooFormulaVar cutVar(cut, cut, *get(), false);
499 return reduceEng(varSubset2,&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
500 }
501 return reduceEng(varSubset2,0,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
502}
503
504////////////////////////////////////////////////////////////////////////////////
505/// Create a subset of the data set by applying the given cut on the data points
506/// and reducing the dimensions to the specified set.
507///
508/// The 'cutVar' formula variable is used to select the subset of data points to be
509/// retained in the reduced data collection.
510
511RooAbsData* RooAbsData::reduce(const RooArgSet& varSubset, const RooFormulaVar& cutVar)
512{
513 // Make sure varSubset doesn't contain any variable not in this dataset
514 RooArgSet varSubset2(varSubset) ;
515 TIterator* iter = varSubset.createIterator() ;
516 RooAbsArg* arg ;
517 while((arg=(RooAbsArg*)iter->Next())) {
518 if (!_vars.find(arg->GetName())) {
519 coutW(InputArguments) << "RooAbsData::reduce(" << GetName() << ") WARNING: variable "
520 << arg->GetName() << " not in dataset, ignored" << endl ;
521 varSubset2.remove(*arg) ;
522 }
523 }
524 delete iter ;
525
526 return reduceEng(varSubset2,&cutVar,0,0,std::numeric_limits<std::size_t>::max(),kFALSE) ;
527}
528
529////////////////////////////////////////////////////////////////////////////////
530/// Return error on current weight (dummy implementation returning zero)
531
533{
534 return 0 ;
535}
536
537////////////////////////////////////////////////////////////////////////////////
538/// Return asymmetric error on weight. (Dummy implementation returning zero)
539
541{
542 lo=0 ; hi=0 ;
543}
544
545
546RooPlot* RooAbsData::plotOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2,
547 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
548 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
549{
551 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
552 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
553 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
554 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
555 return plotOn(frame,l) ;
556}
557
558////////////////////////////////////////////////////////////////////////////////
559/// Create and fill a ROOT histogram TH1,TH2 or TH3 with the values of this dataset for the variables with given names
560/// The range of each observable that is histogrammed is always automatically calculated from the distribution in
561/// the dataset. The number of bins can be controlled using the [xyz]bins parameters. For a greater degree of control
562/// use the createHistogram() method below with named arguments
563///
564/// The caller takes ownership of the returned histogram
565
566TH1 *RooAbsData::createHistogram(const char* varNameList, Int_t xbins, Int_t ybins, Int_t zbins) const
567{
568 // Parse list of variable names
569 char buf[1024] ;
570 strlcpy(buf,varNameList,1024) ;
571 char* varName = strtok(buf,",:") ;
572
573 RooRealVar* xvar = (RooRealVar*) get()->find(varName) ;
574 if (!xvar) {
575 coutE(InputArguments) << "RooAbsData::createHistogram(" << GetName() << ") ERROR: dataset does not contain an observable named " << varName << endl ;
576 return 0 ;
577 }
578 varName = strtok(0,",") ;
579 RooRealVar* yvar = varName ? (RooRealVar*) get()->find(varName) : 0 ;
580 if (varName && !yvar) {
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* zvar = varName ? (RooRealVar*) get()->find(varName) : 0 ;
586 if (varName && !zvar) {
587 coutE(InputArguments) << "RooAbsData::createHistogram(" << GetName() << ") ERROR: dataset does not contain an observable named " << varName << endl ;
588 return 0 ;
589 }
590
591 // Construct list of named arguments to pass to the implementation version of createHistogram()
592
593 RooLinkedList argList ;
594 if (xbins<=0 || !xvar->hasMax() || !xvar->hasMin() ) {
595 argList.Add(RooFit::AutoBinning(xbins==0?xvar->numBins():abs(xbins)).Clone()) ;
596 } else {
597 argList.Add(RooFit::Binning(xbins).Clone()) ;
598 }
599
600 if (yvar) {
601 if (ybins<=0 || !yvar->hasMax() || !yvar->hasMin() ) {
602 argList.Add(RooFit::YVar(*yvar,RooFit::AutoBinning(ybins==0?yvar->numBins():abs(ybins))).Clone()) ;
603 } else {
604 argList.Add(RooFit::YVar(*yvar,RooFit::Binning(ybins)).Clone()) ;
605 }
606 }
607
608 if (zvar) {
609 if (zbins<=0 || !zvar->hasMax() || !zvar->hasMin() ) {
610 argList.Add(RooFit::ZVar(*zvar,RooFit::AutoBinning(zbins==0?zvar->numBins():abs(zbins))).Clone()) ;
611 } else {
612 argList.Add(RooFit::ZVar(*zvar,RooFit::Binning(zbins)).Clone()) ;
613 }
614 }
615
616
617
618 // Call implementation function
619 TH1* result = createHistogram(GetName(),*xvar,argList) ;
620
621 // Delete temporary list of RooCmdArgs
622 argList.Delete() ;
623
624 return result ;
625}
626
627
629 const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
630 const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
631{
633 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
634 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
635 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
636 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
637
638 return createHistogram(name,xvar,l) ;
639}
640
641////////////////////////////////////////////////////////////////////////////////
642///
643/// This function accepts the following arguments
644///
645/// \param[in] name Name of the ROOT histogram
646/// \param[in] xvar Observable to be mapped on x axis of ROOT histogram
647/// \return Histogram now owned by user.
648///
649/// <table>
650/// <tr><td> `AutoBinning(Int_t nbins, Double_y margin)` <td> Automatically calculate range with given added fractional margin, set binning to nbins
651/// <tr><td> `AutoSymBinning(Int_t nbins, Double_y margin)` <td> Automatically calculate range with given added fractional margin,
652/// with additional constraint that mean of data is in center of range, set binning to nbins
653/// <tr><td> `Binning(const char* name)` <td> Apply binning with given name to x axis of histogram
654/// <tr><td> `Binning(RooAbsBinning& binning)` <td> Apply specified binning to x axis of histogram
655/// <tr><td> `Binning(int nbins, double lo, double hi)` <td> Apply specified binning to x axis of histogram
656///
657/// <tr><td> `YVar(const RooAbsRealLValue& var,...)` <td> Observable to be mapped on y axis of ROOT histogram
658/// <tr><td> `ZVar(const RooAbsRealLValue& var,...)` <td> Observable to be mapped on z axis of ROOT histogram
659/// </table>
660///
661/// 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.
662/// ```
663/// createHistogram("histo",x,Binning(-1,1,20), YVar(y,Binning(-1,1,30)), ZVar(z,Binning("zbinning")))
664/// ```
665///
666/// The caller takes ownership of the returned histogram
667
668TH1 *RooAbsData::createHistogram(const char *name, const RooAbsRealLValue& xvar, const RooLinkedList& argListIn) const
669{
670 RooLinkedList argList(argListIn) ;
671
672 // Define configuration for this method
673 RooCmdConfig pc(Form("RooAbsData::createHistogram(%s)",GetName())) ;
674 pc.defineString("cutRange","CutRange",0,"",kTRUE) ;
675 pc.defineString("cutString","CutSpec",0,"") ;
676 pc.defineObject("yvar","YVar",0,0) ;
677 pc.defineObject("zvar","ZVar",0,0) ;
678 pc.allowUndefined() ;
679
680 // Process & check varargs
681 pc.process(argList) ;
682 if (!pc.ok(kTRUE)) {
683 return 0 ;
684 }
685
686 const char* cutSpec = pc.getString("cutString",0,kTRUE) ;
687 const char* cutRange = pc.getString("cutRange",0,kTRUE) ;
688
689 RooArgList vars(xvar) ;
690 RooAbsArg* yvar = static_cast<RooAbsArg*>(pc.getObject("yvar")) ;
691 if (yvar) {
692 vars.add(*yvar) ;
693 }
694 RooAbsArg* zvar = static_cast<RooAbsArg*>(pc.getObject("zvar")) ;
695 if (zvar) {
696 vars.add(*zvar) ;
697 }
698
699 pc.stripCmdList(argList,"CutRange,CutSpec") ;
700
701 // Swap Auto(Sym)RangeData with a Binning command
702 RooLinkedList ownedCmds ;
703 RooCmdArg* autoRD = (RooCmdArg*) argList.find("AutoRangeData") ;
704 if (autoRD) {
706 getRange((RooRealVar&)xvar,xmin,xmax,autoRD->getDouble(0),autoRD->getInt(0)) ;
707 RooCmdArg* bincmd = (RooCmdArg*) RooFit::Binning(autoRD->getInt(1),xmin,xmax).Clone() ;
708 ownedCmds.Add(bincmd) ;
709 argList.Replace(autoRD,bincmd) ;
710 }
711
712 if (yvar) {
713 RooCmdArg* autoRDY = (RooCmdArg*) ((RooCmdArg*)argList.find("YVar"))->subArgs().find("AutoRangeData") ;
714 if (autoRDY) {
716 getRange((RooRealVar&)(*yvar),ymin,ymax,autoRDY->getDouble(0),autoRDY->getInt(0)) ;
717 RooCmdArg* bincmd = (RooCmdArg*) RooFit::Binning(autoRDY->getInt(1),ymin,ymax).Clone() ;
718 //ownedCmds.Add(bincmd) ;
719 ((RooCmdArg*)argList.find("YVar"))->subArgs().Replace(autoRDY,bincmd) ;
720 delete autoRDY ;
721 }
722 }
723
724 if (zvar) {
725 RooCmdArg* autoRDZ = (RooCmdArg*) ((RooCmdArg*)argList.find("ZVar"))->subArgs().find("AutoRangeData") ;
726 if (autoRDZ) {
727 Double_t zmin,zmax ;
728 getRange((RooRealVar&)(*zvar),zmin,zmax,autoRDZ->getDouble(0),autoRDZ->getInt(0)) ;
729 RooCmdArg* bincmd = (RooCmdArg*) RooFit::Binning(autoRDZ->getInt(1),zmin,zmax).Clone() ;
730 //ownedCmds.Add(bincmd) ;
731 ((RooCmdArg*)argList.find("ZVar"))->subArgs().Replace(autoRDZ,bincmd) ;
732 delete autoRDZ ;
733 }
734 }
735
736
737 TH1* histo = xvar.createHistogram(name,argList) ;
738 fillHistogram(histo,vars,cutSpec,cutRange) ;
739
740 ownedCmds.Delete() ;
741
742 return histo ;
743}
744
745////////////////////////////////////////////////////////////////////////////////
746/// Construct table for product of categories in catSet
747
748Roo1DTable* RooAbsData::table(const RooArgSet& catSet, const char* cuts, const char* opts) const
749{
750 RooArgSet catSet2 ;
751
752 string prodName("(") ;
753 TIterator* iter = catSet.createIterator() ;
754 RooAbsArg* arg ;
755 while((arg=(RooAbsArg*)iter->Next())) {
756 if (dynamic_cast<RooAbsCategory*>(arg)) {
757 RooAbsCategory* varsArg = dynamic_cast<RooAbsCategory*>(_vars.find(arg->GetName())) ;
758 if (varsArg != 0) catSet2.add(*varsArg) ;
759 else catSet2.add(*arg) ;
760 if (prodName.length()>1) {
761 prodName += " x " ;
762 }
763 prodName += arg->GetName() ;
764 } else {
765 coutW(InputArguments) << "RooAbsData::table(" << GetName() << ") non-RooAbsCategory input argument " << arg->GetName() << " ignored" << endl ;
766 }
767 }
768 prodName += ")" ;
769 delete iter ;
770
771 RooMultiCategory tmp(prodName.c_str(),prodName.c_str(),catSet2) ;
772 return table(tmp,cuts,opts) ;
773}
774
775////////////////////////////////////////////////////////////////////////////////
776/// Print name of dataset
777
778void RooAbsData::printName(ostream& os) const
779{
780 os << GetName() ;
781}
782
783////////////////////////////////////////////////////////////////////////////////
784/// Print title of dataset
785
786void RooAbsData::printTitle(ostream& os) const
787{
788 os << GetTitle() ;
789}
790
791////////////////////////////////////////////////////////////////////////////////
792/// Print class name of dataset
793
794void RooAbsData::printClassName(ostream& os) const
795{
796 os << IsA()->GetName() ;
797}
798
799////////////////////////////////////////////////////////////////////////////////
800
801void RooAbsData::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
802{
803 _dstore->printMultiline(os,contents,verbose,indent) ;
804}
805
806////////////////////////////////////////////////////////////////////////////////
807/// Define default print options, for a given print style
808
810{
812}
813
814////////////////////////////////////////////////////////////////////////////////
815/// Calculate standardized moment.
816///
817/// \param[in] var Variable to be used for calculating the moment.
818/// \param[in] order Order of the moment.
819/// \param[in] cutSpec If specified, the moment is calculated on the subset of the data which pass the C++ cut specification expression 'cutSpec'
820/// \param[in] cutRange If specified, calculate inside the range named 'cutRange' (also applies cut spec)
821/// \return \f$ \frac{\left< \left( X - \left< X \right> \right)^n \right>}{\sigma^n} \f$, where n = order.
822
823Double_t RooAbsData::standMoment(const RooRealVar &var, Double_t order, const char* cutSpec, const char* cutRange) const
824{
825 // Hardwire invariant answer for first and second moment
826 if (order==1) return 0 ;
827 if (order==2) return 1 ;
828
829 return moment(var,order,cutSpec,cutRange) / TMath::Power(sigma(var,cutSpec,cutRange),order) ;
830}
831
832////////////////////////////////////////////////////////////////////////////////
833/// Calculate moment of requested order.
834///
835/// \param[in] var Variable to be used for calculating the moment.
836/// \param[in] order Order of the moment.
837/// \param[in] cutSpec If specified, the moment is calculated on the subset of the data which pass the C++ cut specification expression 'cutSpec'
838/// \param[in] cutRange If specified, calculate inside the range named 'cutRange' (also applies cut spec)
839/// \return \f$ \left< \left( X - \left< X \right> \right)^n \right> \f$ of order \f$n\f$.
840///
841
842Double_t RooAbsData::moment(const RooRealVar& var, Double_t order, const char* cutSpec, const char* cutRange) const
843{
844 Double_t offset = order>1 ? moment(var,1,cutSpec,cutRange) : 0 ;
845 return moment(var,order,offset,cutSpec,cutRange) ;
846
847}
848
849////////////////////////////////////////////////////////////////////////////////
850/// Return the 'order'-ed moment of observable 'var' in this dataset. If offset is non-zero it is subtracted
851/// from the values of 'var' prior to the moment calculation. If cutSpec and/or cutRange are specified
852/// the moment is calculated on the subset of the data which pass the C++ cut specification expression 'cutSpec'
853/// and/or are inside the range named 'cutRange'
854
855Double_t RooAbsData::moment(const RooRealVar& var, Double_t order, Double_t offset, const char* cutSpec, const char* cutRange) const
856{
857 // Lookup variable in dataset
858 auto arg = _vars.find(var.GetName());
859 if (!arg) {
860 coutE(InputArguments) << "RooDataSet::moment(" << GetName() << ") ERROR: unknown variable: " << var.GetName() << endl ;
861 return 0;
862 }
863
864 auto varPtr = dynamic_cast<const RooRealVar*>(arg);
865 // Check if found variable is of type RooRealVar
866 if (!varPtr) {
867 coutE(InputArguments) << "RooDataSet::moment(" << GetName() << ") ERROR: variable " << var.GetName() << " is not of type RooRealVar" << endl ;
868 return 0;
869 }
870
871 // Check if dataset is not empty
872 if(sumEntries(cutSpec, cutRange) == 0.) {
873 coutE(InputArguments) << "RooDataSet::moment(" << GetName() << ") WARNING: empty dataset" << endl ;
874 return 0;
875 }
876
877 // Setup RooFormulaVar for cutSpec if it is present
878 std::unique_ptr<RooFormula> select;
879 if (cutSpec) {
880 select.reset(new RooFormula("select",cutSpec,*get()));
881 }
882
883
884 // Calculate requested moment
885 Double_t sum(0);
886 const RooArgSet* vars ;
887 for(Int_t index= 0; index < numEntries(); index++) {
888 vars = get(index) ;
889 if (select && select->eval()==0) continue ;
890 if (cutRange && vars->allInRange(cutRange)) continue ;
891
892 sum+= weight() * TMath::Power(varPtr->getVal() - offset,order);
893 }
894 return sum/sumEntries(cutSpec, cutRange);
895}
896
897////////////////////////////////////////////////////////////////////////////////
898/// Internal method to check if given RooRealVar maps to a RooRealVar in this dataset
899
900RooRealVar* RooAbsData::dataRealVar(const char* methodname, const RooRealVar& extVar) const
901{
902 // Lookup variable in dataset
903 RooRealVar *xdata = (RooRealVar*) _vars.find(extVar.GetName());
904 if(!xdata) {
905 coutE(InputArguments) << "RooDataSet::" << methodname << "(" << GetName() << ") ERROR: variable : " << extVar.GetName() << " is not in data" << endl ;
906 return 0;
907 }
908 // Check if found variable is of type RooRealVar
909 if (!dynamic_cast<RooRealVar*>(xdata)) {
910 coutE(InputArguments) << "RooDataSet::" << methodname << "(" << GetName() << ") ERROR: variable : " << extVar.GetName() << " is not of type RooRealVar in data" << endl ;
911 return 0;
912 }
913 return xdata;
914}
915
916////////////////////////////////////////////////////////////////////////////////
917/// Internal method to calculate single correlation and covariance elements
918
919Double_t RooAbsData::corrcov(const RooRealVar &x, const RooRealVar &y, const char* cutSpec, const char* cutRange, Bool_t corr) const
920{
921 // Lookup variable in dataset
922 RooRealVar *xdata = dataRealVar(corr?"correlation":"covariance",x) ;
923 RooRealVar *ydata = dataRealVar(corr?"correlation":"covariance",y) ;
924 if (!xdata||!ydata) return 0 ;
925
926 // Check if dataset is not empty
927 if(sumEntries(cutSpec, cutRange) == 0.) {
928 coutW(InputArguments) << "RooDataSet::" << (corr?"correlation":"covariance") << "(" << GetName() << ") WARNING: empty dataset, returning zero" << endl ;
929 return 0;
930 }
931
932 // Setup RooFormulaVar for cutSpec if it is present
933 RooFormula* select = cutSpec ? new RooFormula("select",cutSpec,*get()) : 0 ;
934
935 // Calculate requested moment
936 Double_t xysum(0),xsum(0),ysum(0),x2sum(0),y2sum(0);
937 const RooArgSet* vars ;
938 for(Int_t index= 0; index < numEntries(); index++) {
939 vars = get(index) ;
940 if (select && select->eval()==0) continue ;
941 if (cutRange && vars->allInRange(cutRange)) continue ;
942
943 xysum += weight()*xdata->getVal()*ydata->getVal() ;
944 xsum += weight()*xdata->getVal() ;
945 ysum += weight()*ydata->getVal() ;
946 if (corr) {
947 x2sum += weight()*xdata->getVal()*xdata->getVal() ;
948 y2sum += weight()*ydata->getVal()*ydata->getVal() ;
949 }
950 }
951
952 // Normalize entries
953 xysum/=sumEntries(cutSpec, cutRange) ;
954 xsum/=sumEntries(cutSpec, cutRange) ;
955 ysum/=sumEntries(cutSpec, cutRange) ;
956 if (corr) {
957 x2sum/=sumEntries(cutSpec, cutRange) ;
958 y2sum/=sumEntries(cutSpec, cutRange) ;
959 }
960
961 // Cleanup
962 if (select) delete select ;
963
964 // Return covariance or correlation as requested
965 if (corr) {
966 return (xysum-xsum*ysum)/(sqrt(x2sum-(xsum*xsum))*sqrt(y2sum-(ysum*ysum))) ;
967 } else {
968 return (xysum-xsum*ysum);
969 }
970}
971
972////////////////////////////////////////////////////////////////////////////////
973/// Return covariance matrix from data for given list of observables
974
975TMatrixDSym* RooAbsData::corrcovMatrix(const RooArgList& vars, const char* cutSpec, const char* cutRange, Bool_t corr) const
976{
977 RooArgList varList ;
978 TIterator* iter = vars.createIterator() ;
979 RooRealVar* var ;
980 while((var=(RooRealVar*)iter->Next())) {
981 RooRealVar* datavar = dataRealVar("covarianceMatrix",*var) ;
982 if (!datavar) {
983 delete iter ;
984 return 0 ;
985 }
986 varList.add(*datavar) ;
987 }
988 delete iter ;
989
990
991 // Check if dataset is not empty
992 if(sumEntries(cutSpec, cutRange) == 0.) {
993 coutW(InputArguments) << "RooDataSet::covariance(" << GetName() << ") WARNING: empty dataset, returning zero" << endl ;
994 return 0;
995 }
996
997 // Setup RooFormulaVar for cutSpec if it is present
998 RooFormula* select = cutSpec ? new RooFormula("select",cutSpec,*get()) : 0 ;
999
1000 iter = varList.createIterator() ;
1001 TIterator* iter2 = varList.createIterator() ;
1002
1003 TMatrixDSym xysum(varList.getSize()) ;
1004 vector<double> xsum(varList.getSize()) ;
1005 vector<double> x2sum(varList.getSize()) ;
1006
1007 // Calculate <x_i> and <x_i y_j>
1008 for(Int_t index= 0; index < numEntries(); index++) {
1009 const RooArgSet* dvars = get(index) ;
1010 if (select && select->eval()==0) continue ;
1011 if (cutRange && dvars->allInRange(cutRange)) continue ;
1012
1013 RooRealVar* varx, *vary ;
1014 iter->Reset() ;
1015 Int_t ix=0,iy=0 ;
1016 while((varx=(RooRealVar*)iter->Next())) {
1017 xsum[ix] += weight()*varx->getVal() ;
1018 if (corr) {
1019 x2sum[ix] += weight()*varx->getVal()*varx->getVal() ;
1020 }
1021
1022 *iter2=*iter ; iy=ix ;
1023 vary=varx ;
1024 while(vary) {
1025 xysum(ix,iy) += weight()*varx->getVal()*vary->getVal() ;
1026 xysum(iy,ix) = xysum(ix,iy) ;
1027 iy++ ;
1028 vary=(RooRealVar*)iter2->Next() ;
1029 }
1030 ix++ ;
1031 }
1032
1033 }
1034
1035 // Normalize sums
1036 for (Int_t ix=0 ; ix<varList.getSize() ; ix++) {
1037 xsum[ix] /= sumEntries(cutSpec, cutRange) ;
1038 if (corr) {
1039 x2sum[ix] /= sumEntries(cutSpec, cutRange) ;
1040 }
1041 for (Int_t iy=0 ; iy<varList.getSize() ; iy++) {
1042 xysum(ix,iy) /= sumEntries(cutSpec, cutRange) ;
1043 }
1044 }
1045
1046 // Calculate covariance matrix
1047 TMatrixDSym* C = new TMatrixDSym(varList.getSize()) ;
1048 for (Int_t ix=0 ; ix<varList.getSize() ; ix++) {
1049 for (Int_t iy=0 ; iy<varList.getSize() ; iy++) {
1050 (*C)(ix,iy) = xysum(ix,iy)-xsum[ix]*xsum[iy] ;
1051 if (corr) {
1052 (*C)(ix,iy) /= sqrt((x2sum[ix]-(xsum[ix]*xsum[ix]))*(x2sum[iy]-(xsum[iy]*xsum[iy]))) ;
1053 }
1054 }
1055 }
1056
1057 if (select) delete select ;
1058 delete iter ;
1059 delete iter2 ;
1060
1061 return C ;
1062}
1063
1064////////////////////////////////////////////////////////////////////////////////
1065/// Create a RooRealVar containing the mean of observable 'var' in
1066/// this dataset. If cutSpec and/or cutRange are specified the
1067/// moment is calculated on the subset of the data which pass the C++
1068/// cut specification expression 'cutSpec' and/or are inside the
1069/// range named 'cutRange'
1070
1071RooRealVar* RooAbsData::meanVar(const RooRealVar &var, const char* cutSpec, const char* cutRange) const
1072{
1073 // Create a new variable with appropriate strings. The error is calculated as
1074 // RMS/Sqrt(N) which is generally valid.
1075
1076 // Create holder variable for mean
1077 TString name(var.GetName()),title("Mean of ") ;
1078 name.Append("Mean");
1079 title.Append(var.GetTitle());
1080 RooRealVar *meanv= new RooRealVar(name,title,0) ;
1081 meanv->setConstant(kFALSE) ;
1082
1083 // Adjust plot label
1084 TString label("<") ;
1085 label.Append(var.getPlotLabel());
1086 label.Append(">");
1087 meanv->setPlotLabel(label.Data());
1088
1089 // fill in this variable's value and error
1090 Double_t meanVal=moment(var,1,0,cutSpec,cutRange) ;
1091 Double_t N(sumEntries(cutSpec,cutRange)) ;
1092
1093 Double_t rmsVal= sqrt(moment(var,2,meanVal,cutSpec,cutRange)*N/(N-1));
1094 meanv->setVal(meanVal) ;
1095 meanv->setError(N > 0 ? rmsVal/sqrt(N) : 0);
1096
1097 return meanv;
1098}
1099
1100////////////////////////////////////////////////////////////////////////////////
1101/// Create a RooRealVar containing the RMS of observable 'var' in
1102/// this dataset. If cutSpec and/or cutRange are specified the
1103/// moment is calculated on the subset of the data which pass the C++
1104/// cut specification expression 'cutSpec' and/or are inside the
1105/// range named 'cutRange'
1106
1107RooRealVar* RooAbsData::rmsVar(const RooRealVar &var, const char* cutSpec, const char* cutRange) const
1108{
1109 // Create a new variable with appropriate strings. The error is calculated as
1110 // RMS/(2*Sqrt(N)) which is only valid if the variable has a Gaussian distribution.
1111
1112 // Create RMS value holder
1113 TString name(var.GetName()),title("RMS of ") ;
1114 name.Append("RMS");
1115 title.Append(var.GetTitle());
1116 RooRealVar *rms= new RooRealVar(name,title,0) ;
1117 rms->setConstant(kFALSE) ;
1118
1119 // Adjust plot label
1120 TString label(var.getPlotLabel());
1121 label.Append("_{RMS}");
1122 rms->setPlotLabel(label);
1123
1124 // Fill in this variable's value and error
1125 Double_t meanVal(moment(var,1,0,cutSpec,cutRange)) ;
1126 Double_t N(sumEntries(cutSpec, cutRange));
1127 Double_t rmsVal= sqrt(moment(var,2,meanVal,cutSpec,cutRange)*N/(N-1));
1128 rms->setVal(rmsVal) ;
1129 rms->setError(rmsVal/sqrt(2*N));
1130
1131 return rms;
1132}
1133
1134////////////////////////////////////////////////////////////////////////////////
1135/// Add a box with statistics information to the specified frame. By default a box with the
1136/// event count, mean and rms of the plotted variable is added.
1137///
1138/// The following optional named arguments are accepted
1139/// <table>
1140/// <tr><td> `What(const char* whatstr)` <td> Controls what is printed: "N" = count, "M" is mean, "R" is RMS.
1141/// <tr><td> `Format(const char* optStr)` <td> \deprecated Classing parameter formatting options, provided for backward compatibility
1142///
1143/// <tr><td> `Format(const char* what,...)` <td> Parameter formatting options.
1144/// <table>
1145/// <tr><td> const char* what <td> Controls what is shown:
1146/// - "N" adds name
1147/// - "E" adds error
1148/// - "A" shows asymmetric error
1149/// - "U" shows unit
1150/// - "H" hides the value
1151/// <tr><td> `FixedPrecision(int n)` <td> Controls precision, set fixed number of digits
1152/// <tr><td> `AutoPrecision(int n)` <td> Controls precision. Number of shown digits is calculated from error + n specified additional digits (1 is sensible default)
1153/// <tr><td> `VerbatimName(Bool_t flag)` <td> Put variable name in a \\verb+ + clause.
1154/// </table>
1155/// <tr><td> `Label(const chat* label)` <td> Add header label to parameter box
1156/// <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
1157/// bottom of box is calculated automatically from number lines in box
1158/// <tr><td> `Cut(const char* expression)` <td> Apply given cut expression to data when calculating statistics
1159/// <tr><td> `CutRange(const char* rangeName)` <td> Only consider events within given range when calculating statistics. Multiple
1160/// CutRange() argument may be specified to combine ranges.
1161///
1162/// </table>
1163
1164RooPlot* RooAbsData::statOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2,
1165 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
1166 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
1167{
1168 // Stuff all arguments in a list
1169 RooLinkedList cmdList;
1170 cmdList.Add(const_cast<RooCmdArg*>(&arg1)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg2)) ;
1171 cmdList.Add(const_cast<RooCmdArg*>(&arg3)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg4)) ;
1172 cmdList.Add(const_cast<RooCmdArg*>(&arg5)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg6)) ;
1173 cmdList.Add(const_cast<RooCmdArg*>(&arg7)) ; cmdList.Add(const_cast<RooCmdArg*>(&arg8)) ;
1174
1175 // Select the pdf-specific commands
1176 RooCmdConfig pc(Form("RooTreeData::statOn(%s)",GetName())) ;
1177 pc.defineString("what","What",0,"MNR") ;
1178 pc.defineString("label","Label",0,"") ;
1179 pc.defineDouble("xmin","Layout",0,0.65) ;
1180 pc.defineDouble("xmax","Layout",1,0.99) ;
1181 pc.defineInt("ymaxi","Layout",0,Int_t(0.95*10000)) ;
1182 pc.defineString("formatStr","Format",0,"NELU") ;
1183 pc.defineInt("sigDigit","Format",0,2) ;
1184 pc.defineInt("dummy","FormatArgs",0,0) ;
1185 pc.defineString("cutRange","CutRange",0,"",kTRUE) ;
1186 pc.defineString("cutString","CutSpec",0,"") ;
1187 pc.defineMutex("Format","FormatArgs") ;
1188
1189 // Process and check varargs
1190 pc.process(cmdList) ;
1191 if (!pc.ok(kTRUE)) {
1192 return frame ;
1193 }
1194
1195 const char* label = pc.getString("label") ;
1196 Double_t xmin = pc.getDouble("xmin") ;
1197 Double_t xmax = pc.getDouble("xmax") ;
1198 Double_t ymax = pc.getInt("ymaxi") / 10000. ;
1199 const char* formatStr = pc.getString("formatStr") ;
1200 Int_t sigDigit = pc.getInt("sigDigit") ;
1201 const char* what = pc.getString("what") ;
1202
1203 const char* cutSpec = pc.getString("cutString",0,kTRUE) ;
1204 const char* cutRange = pc.getString("cutRange",0,kTRUE) ;
1205
1206 if (pc.hasProcessed("FormatArgs")) {
1207 RooCmdArg* formatCmd = static_cast<RooCmdArg*>(cmdList.FindObject("FormatArgs")) ;
1208 return statOn(frame,what,label,0,0,xmin,xmax,ymax,cutSpec,cutRange,formatCmd) ;
1209 } else {
1210 return statOn(frame,what,label,sigDigit,formatStr,xmin,xmax,ymax,cutSpec,cutRange) ;
1211 }
1212}
1213
1214////////////////////////////////////////////////////////////////////////////////
1215/// Implementation back-end of statOn() method with named arguments
1216
1217RooPlot* RooAbsData::statOn(RooPlot* frame, const char* what, const char *label, Int_t sigDigits,
1219 const char* cutSpec, const char* cutRange, const RooCmdArg* formatCmd)
1220{
1221 Bool_t showLabel= (label != 0 && strlen(label) > 0);
1222
1223 TString whatStr(what) ;
1224 whatStr.ToUpper() ;
1225 Bool_t showN = whatStr.Contains("N") ;
1226 Bool_t showR = whatStr.Contains("R") ;
1227 Bool_t showM = whatStr.Contains("M") ;
1228 Int_t nPar= 0;
1229 if (showN) nPar++ ;
1230 if (showR) nPar++ ;
1231 if (showM) nPar++ ;
1232
1233 // calculate the box's size
1234 Double_t dy(0.06), ymin(ymax-nPar*dy);
1235 if(showLabel) ymin-= dy;
1236
1237 // create the box and set its options
1238 TPaveText *box= new TPaveText(xmin,ymax,xmax,ymin,"BRNDC");
1239 if(!box) return 0;
1240 box->SetName(Form("%s_statBox",GetName())) ;
1241 box->SetFillColor(0);
1242 box->SetBorderSize(1);
1243 box->SetTextAlign(12);
1244 box->SetTextSize(0.04F);
1245 box->SetFillStyle(1001);
1246
1247 // add formatted text for each statistic
1248 RooRealVar N("N","Number of Events",sumEntries(cutSpec,cutRange));
1249 N.setPlotLabel("Entries") ;
1250 RooRealVar *meanv= meanVar(*(RooRealVar*)frame->getPlotVar(),cutSpec,cutRange);
1251 meanv->setPlotLabel("Mean") ;
1252 RooRealVar *rms= rmsVar(*(RooRealVar*)frame->getPlotVar(),cutSpec,cutRange);
1253 rms->setPlotLabel("RMS") ;
1254 TString *rmsText, *meanText, *NText ;
1255 if (options) {
1256 rmsText= rms->format(sigDigits,options);
1257 meanText= meanv->format(sigDigits,options);
1258 NText= N.format(sigDigits,options);
1259 } else {
1260 rmsText= rms->format(*formatCmd);
1261 meanText= meanv->format(*formatCmd);
1262 NText= N.format(*formatCmd);
1263 }
1264 if (showR) box->AddText(rmsText->Data());
1265 if (showM) box->AddText(meanText->Data());
1266 if (showN) box->AddText(NText->Data());
1267
1268 // cleanup heap memory
1269 delete NText;
1270 delete meanText;
1271 delete rmsText;
1272 delete meanv;
1273 delete rms;
1274
1275 // add the optional label if specified
1276 if(showLabel) box->AddText(label);
1277
1278 frame->addObject(box) ;
1279 return frame ;
1280}
1281
1282////////////////////////////////////////////////////////////////////////////////
1283/// Loop over columns of our tree data and fill the input histogram. Returns a pointer to the
1284/// input histogram, or zero in case of an error. The input histogram can be any TH1 subclass, and
1285/// therefore of arbitrary dimension. Variables are matched with the (x,y,...) dimensions of the input
1286/// histogram according to the order in which they appear in the input plotVars list.
1287
1288TH1 *RooAbsData::fillHistogram(TH1 *hist, const RooArgList &plotVars, const char *cuts, const char* cutRange) const
1289{
1290 // Do we have a valid histogram to use?
1291 if(0 == hist) {
1292 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: no valid histogram to fill" << endl;
1293 return 0;
1294 }
1295
1296 // Check that the number of plotVars matches the input histogram's dimension
1297 Int_t hdim= hist->GetDimension();
1298 if(hdim != plotVars.getSize()) {
1299 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: plotVars has the wrong dimension" << endl;
1300 return 0;
1301 }
1302
1303 // Check that the plot variables are all actually RooAbsReal's and print a warning if we do not
1304 // explicitly depend on one of them. Clone any variables that we do not contain directly and
1305 // redirect them to use our event data.
1306 RooArgSet plotClones,localVars;
1307 for(Int_t index= 0; index < plotVars.getSize(); index++) {
1308 const RooAbsArg *var= plotVars.at(index);
1309 const RooAbsReal *realVar= dynamic_cast<const RooAbsReal*>(var);
1310 if(0 == realVar) {
1311 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: cannot plot variable \"" << var->GetName()
1312 << "\" of type " << var->ClassName() << endl;
1313 return 0;
1314 }
1315 RooAbsArg *found= _vars.find(realVar->GetName());
1316 if(!found) {
1317 RooAbsArg *clone= plotClones.addClone(*realVar,kTRUE); // do not complain about duplicates
1318 assert(0 != clone);
1319 if(!clone->dependsOn(_vars)) {
1320 coutE(InputArguments) << ClassName() << "::" << GetName()
1321 << ":fillHistogram: Data does not contain the variable '" << realVar->GetName() << "'." << endl;
1322 return nullptr;
1323 }
1324 else {
1326 }
1327 localVars.add(*clone);
1328 }
1329 else {
1330 localVars.add(*found);
1331 }
1332 }
1333
1334 // Create selection formula if selection cuts are specified
1335 std::unique_ptr<RooFormula> select;
1336 if (cuts != nullptr && strlen(cuts) > 0) {
1337 select.reset(new RooFormula(cuts, cuts, _vars, false));
1338 if (!select || !select->ok()) {
1339 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: invalid cuts \"" << cuts << "\"" << endl;
1340 return 0 ;
1341 }
1342 }
1343
1344 // Lookup each of the variables we are binning in our tree variables
1345 const RooAbsReal *xvar = 0;
1346 const RooAbsReal *yvar = 0;
1347 const RooAbsReal *zvar = 0;
1348 switch(hdim) {
1349 case 3:
1350 zvar= dynamic_cast<RooAbsReal*>(localVars.find(plotVars.at(2)->GetName()));
1351 assert(0 != zvar);
1352 // fall through to next case...
1353 case 2:
1354 yvar= dynamic_cast<RooAbsReal*>(localVars.find(plotVars.at(1)->GetName()));
1355 assert(0 != yvar);
1356 // fall through to next case...
1357 case 1:
1358 xvar= dynamic_cast<RooAbsReal*>(localVars.find(plotVars.at(0)->GetName()));
1359 assert(0 != xvar);
1360 break;
1361 default:
1362 coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: cannot fill histogram with "
1363 << hdim << " dimensions" << endl;
1364 break;
1365 }
1366
1367 // Parse cutRange specification
1368 vector<string> cutVec ;
1369 if (cutRange && strlen(cutRange)>0) {
1370 if (strchr(cutRange,',')==0) {
1371 cutVec.push_back(cutRange) ;
1372 } else {
1373 const size_t bufSize = strlen(cutRange)+1;
1374 char* buf = new char[bufSize] ;
1375 strlcpy(buf,cutRange,bufSize) ;
1376 const char* oneRange = strtok(buf,",") ;
1377 while(oneRange) {
1378 cutVec.push_back(oneRange) ;
1379 oneRange = strtok(0,",") ;
1380 }
1381 delete[] buf ;
1382 }
1383 }
1384
1385 // Loop over events and fill the histogram
1386 if (hist->GetSumw2()->fN==0) {
1387 hist->Sumw2() ;
1388 }
1389 Int_t nevent= numEntries() ; //(Int_t)_tree->GetEntries();
1390 for(Int_t i=0; i < nevent; ++i) {
1391
1392 //Int_t entryNumber= _tree->GetEntryNumber(i);
1393 //if (entryNumber<0) break;
1394 get(i);
1395
1396 // Apply expression based selection criteria
1397 if (select && select->eval()==0) {
1398 continue ;
1399 }
1400
1401
1402 // Apply range based selection criteria
1403 Bool_t selectByRange = kTRUE ;
1404 if (cutRange) {
1405 for (const auto arg : _vars) {
1406 Bool_t selectThisArg = kFALSE ;
1407 UInt_t icut ;
1408 for (icut=0 ; icut<cutVec.size() ; icut++) {
1409 if (arg->inRange(cutVec[icut].c_str())) {
1410 selectThisArg = kTRUE ;
1411 break ;
1412 }
1413 }
1414 if (!selectThisArg) {
1415 selectByRange = kFALSE ;
1416 break ;
1417 }
1418 }
1419 }
1420
1421 if (!selectByRange) {
1422 // Go to next event in loop over events
1423 continue ;
1424 }
1425
1426 Int_t bin(0);
1427 switch(hdim) {
1428 case 1:
1429 bin= hist->FindBin(xvar->getVal());
1430 hist->Fill(xvar->getVal(),weight()) ;
1431 break;
1432 case 2:
1433 bin= hist->FindBin(xvar->getVal(),yvar->getVal());
1434 static_cast<TH2*>(hist)->Fill(xvar->getVal(),yvar->getVal(),weight()) ;
1435 break;
1436 case 3:
1437 bin= hist->FindBin(xvar->getVal(),yvar->getVal(),zvar->getVal());
1438 static_cast<TH3*>(hist)->Fill(xvar->getVal(),yvar->getVal(),zvar->getVal(),weight()) ;
1439 break;
1440 default:
1441 assert(hdim < 3);
1442 break;
1443 }
1444
1445
1446 Double_t error2 = TMath::Power(hist->GetBinError(bin),2)-TMath::Power(weight(),2) ;
1448 if (we==0) we = weight() ;
1449 error2 += TMath::Power(we,2) ;
1450
1451
1452// Double_t we = weightError(RooAbsData::SumW2) ;
1453// Double_t error2(0) ;
1454// if (we==0) {
1455// we = weight() ; //sqrt(weight()) ;
1456// error2 = TMath::Power(hist->GetBinError(bin),2)-TMath::Power(weight(),2) + TMath::Power(we,2) ;
1457// } else {
1458// error2 = TMath::Power(hist->GetBinError(bin),2)-TMath::Power(weight(),2) + TMath::Power(we,2) ;
1459// }
1460 //hist->AddBinContent(bin,weight());
1461 hist->SetBinError(bin,sqrt(error2)) ;
1462
1463 //cout << "RooTreeData::fillHistogram() bin = " << bin << " weight() = " << weight() << " we = " << we << endl ;
1464
1465 }
1466
1467 return hist;
1468}
1469
1470////////////////////////////////////////////////////////////////////////////////
1471/// Split dataset into subsets based on states of given splitCat in this dataset.
1472/// A TList of RooDataSets is returned in which each RooDataSet is named
1473/// after the state name of splitCat of which it contains the dataset subset.
1474/// The observables splitCat itself is no longer present in the sub datasets.
1475/// If createEmptyDataSets is kFALSE (default) this method only creates datasets for states
1476/// which have at least one entry The caller takes ownership of the returned list and its contents
1477
1478TList* RooAbsData::split(const RooAbsCategory& splitCat, Bool_t createEmptyDataSets) const
1479{
1480 // Sanity check
1481 if (!splitCat.dependsOn(*get())) {
1482 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") ERROR category " << splitCat.GetName()
1483 << " doesn't depend on any variable in this dataset" << endl ;
1484 return 0 ;
1485 }
1486
1487 // Clone splitting category and attach to self
1488 RooAbsCategory* cloneCat =0;
1489 RooArgSet* cloneSet = 0;
1490 if (splitCat.isDerived()) {
1491 cloneSet = (RooArgSet*) RooArgSet(splitCat).snapshot(kTRUE) ;
1492 if (!cloneSet) {
1493 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") Couldn't deep-clone splitting category, abort." << endl ;
1494 return 0 ;
1495 }
1496 cloneCat = (RooAbsCategory*) cloneSet->find(splitCat.GetName()) ;
1497 cloneCat->attachDataSet(*this) ;
1498 } else {
1499 cloneCat = dynamic_cast<RooAbsCategory*>(get()->find(splitCat.GetName())) ;
1500 if (!cloneCat) {
1501 coutE(InputArguments) << "RooTreeData::split(" << GetName() << ") ERROR category " << splitCat.GetName()
1502 << " is fundamental and does not appear in this dataset" << endl ;
1503 return 0 ;
1504 }
1505 }
1506
1507 // Split a dataset in a series of subsets, each corresponding
1508 // to a state of splitCat
1509 TList* dsetList = new TList ;
1510
1511 // Construct set of variables to be included in split sets = full set - split category
1512 RooArgSet subsetVars(*get()) ;
1513 if (splitCat.isDerived()) {
1514 RooArgSet* vars = splitCat.getVariables() ;
1515 subsetVars.remove(*vars,kTRUE,kTRUE) ;
1516 delete vars ;
1517 } else {
1518 subsetVars.remove(splitCat,kTRUE,kTRUE) ;
1519 }
1520
1521 // Add weight variable explicitly if dataset has weights, but no top-level weight
1522 // variable exists (can happen with composite datastores)
1523 Bool_t addWV(kFALSE) ;
1524 RooRealVar newweight("weight","weight",-1e9,1e9) ;
1525 if (isWeighted() && !IsA()->InheritsFrom(RooDataHist::Class())) {
1526 subsetVars.add(newweight) ;
1527 addWV = kTRUE ;
1528 }
1529
1530 // If createEmptyDataSets is true, prepopulate with empty sets corresponding to all states
1531 if (createEmptyDataSets) {
1532 for (const auto& nameIdx : *cloneCat) {
1533 RooAbsData* subset = emptyClone(nameIdx.first.c_str(), nameIdx.first.c_str(), &subsetVars,(addWV?"weight":0)) ;
1534 dsetList->Add((RooAbsArg*)subset) ;
1535 }
1536 }
1537
1538
1539 // Loop over dataset and copy event to matching subset
1540 const bool propWeightSquared = isWeighted();
1541 for (Int_t i = 0; i < numEntries(); ++i) {
1542 const RooArgSet* row = get(i);
1543 RooAbsData* subset = (RooAbsData*) dsetList->FindObject(cloneCat->getCurrentLabel());
1544 if (!subset) {
1545 subset = emptyClone(cloneCat->getCurrentLabel(),cloneCat->getCurrentLabel(),&subsetVars,(addWV?"weight":0));
1546 dsetList->Add((RooAbsArg*)subset);
1547 }
1548 if (!propWeightSquared) {
1549 subset->add(*row, weight());
1550 } else {
1551 subset->add(*row, weight(), weightSquared());
1552 }
1553 }
1554
1555 delete cloneSet;
1556 return dsetList;
1557}
1558
1559////////////////////////////////////////////////////////////////////////////////
1560/// Plot dataset on specified frame. By default an unbinned dataset will use the default binning of
1561/// the target frame. A binned dataset will by default retain its intrinsic binning.
1562///
1563/// The following optional named arguments can be used to modify the default behavior
1564///
1565/// <table>
1566/// <tr><th> <th> Data representation options
1567/// <tr><td> `Asymmetry(const RooCategory& c)` <td> Show the asymmetry of the data in given two-state category [F(+)-F(-)] / [F(+)+F(-)].
1568/// Category must have two states with indices -1 and +1 or three states with indices -1,0 and +1.
1569/// <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
1570/// <tr><td> `DataError(RooAbsData::EType)` <td> Select the type of error drawn:
1571/// - `Auto(default)` results in Poisson for unweighted data and SumW2 for weighted data
1572/// - `Poisson` draws asymmetric Poisson confidence intervals.
1573/// - `SumW2` draws symmetric sum-of-weights error ( sum(w)^2/sum(w^2) )
1574/// - `None` draws no error bars
1575/// <tr><td> `Binning(int nbins, double xlo, double xhi)` <td> Use specified binning to draw dataset
1576/// <tr><td> `Binning(const RooAbsBinning&)` <td> Use specified binning to draw dataset
1577/// <tr><td> `Binning(const char* name)` <td> Use binning with specified name to draw dataset
1578/// <tr><td> `RefreshNorm(Bool_t flag)` <td> Force refreshing for PDF normalization information in frame.
1579/// If set, any subsequent PDF will normalize to this dataset, even if it is
1580/// not the first one added to the frame. By default only the 1st dataset
1581/// added to a frame will update the normalization information
1582/// <tr><td> `Rescale(Double_t f)` <td> Rescale drawn histogram by given factor
1583/// <tr><td> `CutRange(const char*)` <td> Apply cuts to dataset.
1584/// \note This often requires passing the normalisation when plotting the PDF because RooFit does not save
1585/// how many events were being plotted (it will only work for cutting slices out of uniformly distributed variables).
1586/// ```
1587/// data->plotOn(frame01, CutRange("SB1"));
1588/// const double nData = data->sumEntries("", "SB1");
1589/// // Make clear that the target normalisation is nData. The enumerator NumEvent
1590/// // is needed to switch between relative and absolute scaling.
1591/// model.plotOn(frame01, Normalization(nData, RooAbsReal::NumEvent),
1592/// ProjectionRange("SB1"));
1593/// ```
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.lookupName(1).c_str()));
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(const RooAbsRealLValue& var, Double_t& lowest, Double_t& highest, Double_t marginFrac, Bool_t symMode) const
2115{
2116 // Lookup variable in dataset
2117 const auto arg = _vars.find(var.GetName());
2118 if (!arg) {
2119 coutE(InputArguments) << "RooDataSet::getRange(" << GetName() << ") ERROR: unknown variable: " << var.GetName() << endl ;
2120 return kTRUE;
2121 }
2122
2123 auto varPtr = dynamic_cast<const RooRealVar*>(arg);
2124 // Check if found variable is of type RooRealVar
2125 if (!varPtr) {
2126 coutE(InputArguments) << "RooDataSet::getRange(" << GetName() << ") ERROR: variable " << var.GetName() << " is not of type RooRealVar" << endl ;
2127 return kTRUE;
2128 }
2129
2130 // Check if dataset is not empty
2131 if(sumEntries() == 0.) {
2132 coutE(InputArguments) << "RooDataSet::getRange(" << GetName() << ") WARNING: empty dataset" << endl ;
2133 return kTRUE;
2134 }
2135
2136 // Look for highest and lowest value
2137 lowest = RooNumber::infinity() ;
2138 highest = -RooNumber::infinity() ;
2139 for (Int_t i=0 ; i<numEntries() ; i++) {
2140 get(i) ;
2141 if (varPtr->getVal()<lowest) {
2142 lowest = varPtr->getVal() ;
2143 }
2144 if (varPtr->getVal()>highest) {
2145 highest = varPtr->getVal() ;
2146 }
2147 }
2148
2149 if (marginFrac>0) {
2150 if (symMode==kFALSE) {
2151
2152 Double_t margin = marginFrac*(highest-lowest) ;
2153 lowest -= margin ;
2154 highest += margin ;
2155 if (lowest<var.getMin()) lowest = var.getMin() ;
2156 if (highest>var.getMax()) highest = var.getMax() ;
2157
2158 } else {
2159
2160 Double_t mom1 = moment(*varPtr,1) ;
2161 Double_t delta = ((highest-mom1)>(mom1-lowest)?(highest-mom1):(mom1-lowest))*(1+marginFrac) ;
2162 lowest = mom1-delta ;
2163 highest = mom1+delta ;
2164 if (lowest<var.getMin()) lowest = var.getMin() ;
2165 if (highest>var.getMax()) highest = var.getMax() ;
2166
2167 }
2168 }
2169
2170 return kFALSE ;
2171}
2172
2173////////////////////////////////////////////////////////////////////////////////
2174/// Prepare dataset for use with cached constant terms listed in
2175/// 'cacheList' of expression 'arg'. Deactivate tree branches
2176/// for any dataset observable that is either not used at all,
2177/// or is used exclusively by cached branch nodes.
2178
2179void RooAbsData::optimizeReadingWithCaching(RooAbsArg& arg, const RooArgSet& cacheList, const RooArgSet& keepObsList)
2180{
2181 RooArgSet pruneSet ;
2182
2183 // Add unused observables in this dataset to pruneSet
2184 pruneSet.add(*get()) ;
2185 RooArgSet* usedObs = arg.getObservables(*this) ;
2186 pruneSet.remove(*usedObs,kTRUE,kTRUE) ;
2187
2188 // Add observables exclusively used to calculate cached observables to pruneSet
2189 TIterator* vIter = get()->createIterator() ;
2190 RooAbsArg *var ;
2191 while ((var=(RooAbsArg*) vIter->Next())) {
2192 if (allClientsCached(var,cacheList)) {
2193 pruneSet.add(*var) ;
2194 }
2195 }
2196 delete vIter ;
2197
2198
2199 if (pruneSet.getSize()!=0) {
2200
2201 // Go over all used observables and check if any of them have parameterized
2202 // ranges in terms of pruned observables. If so, remove those observable
2203 // from the pruning list
2204 TIterator* uIter = usedObs->createIterator() ;
2205 RooAbsArg* obs ;
2206 while((obs=(RooAbsArg*)uIter->Next())) {
2207 RooRealVar* rrv = dynamic_cast<RooRealVar*>(obs) ;
2208 if (rrv && !rrv->getBinning().isShareable()) {
2209 RooArgSet depObs ;
2210 RooAbsReal* loFunc = rrv->getBinning().lowBoundFunc() ;
2211 RooAbsReal* hiFunc = rrv->getBinning().highBoundFunc() ;
2212 if (loFunc) {
2213 loFunc->leafNodeServerList(&depObs,0,kTRUE) ;
2214 }
2215 if (hiFunc) {
2216 hiFunc->leafNodeServerList(&depObs,0,kTRUE) ;
2217 }
2218 if (depObs.getSize()>0) {
2219 pruneSet.remove(depObs,kTRUE,kTRUE) ;
2220 }
2221 }
2222 }
2223 delete uIter ;
2224 }
2225
2226
2227 // Remove all observables in keep list from prune list
2228 pruneSet.remove(keepObsList,kTRUE,kTRUE) ;
2229
2230 if (pruneSet.getSize()!=0) {
2231
2232 // Deactivate tree branches here
2233 cxcoutI(Optimization) << "RooTreeData::optimizeReadingForTestStatistic(" << GetName() << "): Observables " << pruneSet
2234 << " 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 ;
2235 setArgStatus(pruneSet,kFALSE) ;
2236 }
2237
2238 delete usedObs ;
2239
2240}
2241
2242////////////////////////////////////////////////////////////////////////////////
2243/// Utility function that determines if all clients of object 'var'
2244/// appear in given list of cached nodes.
2245
2247{
2248 Bool_t ret(kTRUE), anyClient(kFALSE) ;
2249
2250 for (const auto client : var->valueClients()) {
2251 anyClient = kTRUE ;
2252 if (!cacheList.find(client->GetName())) {
2253 // If client is not cached recurse
2254 ret &= allClientsCached(client,cacheList) ;
2255 }
2256 }
2257
2258 return anyClient?ret:kFALSE ;
2259}
2260
2261////////////////////////////////////////////////////////////////////////////////
2262
2264{
2265 _dstore->attachBuffers(extObs) ;
2266}
2267
2268////////////////////////////////////////////////////////////////////////////////
2269
2271{
2273}
2274
2275////////////////////////////////////////////////////////////////////////////////
2276
2278{
2279 if (_ownedComponents.size()>0) {
2280 return kTRUE ;
2281 }
2282 return kFALSE ;
2283}
2284
2285////////////////////////////////////////////////////////////////////////////////
2286
2288{
2289 map<string,RooAbsData*>::iterator i = _ownedComponents.find(name) ;
2290 if (i==_ownedComponents.end()) return 0 ;
2291 return i->second ;
2292}
2293
2294////////////////////////////////////////////////////////////////////////////////
2295
2296void RooAbsData::addOwnedComponent(const char* idxlabel, RooAbsData& data)
2297{
2298 _ownedComponents[idxlabel]= &data ;
2299}
2300
2301////////////////////////////////////////////////////////////////////////////////
2302/// Stream an object of class RooAbsData.
2303
2304void RooAbsData::Streamer(TBuffer &R__b)
2305{
2306 if (R__b.IsReading()) {
2308
2309 // Convert on the fly to vector storage if that the current working default
2312 }
2313
2314 } else {
2316 }
2317}
2318
2319////////////////////////////////////////////////////////////////////////////////
2320
2322{
2323 _dstore->checkInit() ;
2324}
2325
2326////////////////////////////////////////////////////////////////////////////////
2327/// Forward draw command to data store
2328
2330{
2331 if (_dstore) _dstore->Draw(option) ;
2332}
2333
2334////////////////////////////////////////////////////////////////////////////////
2335
2337{
2338 return _dstore->hasFilledCache() ;
2339}
2340
2341////////////////////////////////////////////////////////////////////////////////
2342/// Return a pointer to the TTree which stores the data. Returns a nullpointer
2343/// if vector-based storage is used. The RooAbsData remains owner of the tree.
2344/// GetClonedTree() can be used to get a tree even if the internal storage does not use one.
2345
2347{
2349 return _dstore->tree();
2350 } else {
2351 coutW(InputArguments) << "RooAbsData::tree(" << GetName() << ") WARNING: is not of StorageType::Tree. "
2352 << "Use GetClonedTree() instead or convert to tree storage." << endl;
2353 return (TTree *)nullptr;
2354 }
2355}
2356
2357////////////////////////////////////////////////////////////////////////////////
2358/// Return a clone of the TTree which stores the data or create such a tree
2359/// if vector storage is used. The user is responsible for deleting the tree
2360
2362{
2364 auto tmp = const_cast<TTree *>(_dstore->tree());
2365 return tmp->CloneTree();
2366 } else {
2367 RooTreeDataStore buffer(GetName(), GetTitle(), *get(), *_dstore);
2368 return buffer.tree().CloneTree();
2369 }
2370}
2371
2372////////////////////////////////////////////////////////////////////////////////
2373/// Convert vector-based storage to tree-based storage
2374
2376{
2379 delete _dstore;
2380 _dstore = newStore;
2382 }
2383}
2384
2385////////////////////////////////////////////////////////////////////////////////
2386/// If one of the TObject we have a referenced to is deleted, remove the
2387/// reference.
2388
2390{
2391 for(auto &iter : _ownedComponents) {
2392 if (iter.second == obj) {
2393 iter.second = nullptr;
2394 }
2395 }
2396}
void Class()
Definition: Class.C:29
static std::map< RooAbsData *, int > _dcc
Definition: RooAbsData.cxx:69
#define coutI(a)
Definition: RooMsgService.h:30
#define cxcoutI(a)
Definition: RooMsgService.h:85
#define coutW(a)
Definition: RooMsgService.h:32
#define coutE(a)
Definition: RooMsgService.h:33
int Int_t
Definition: RtypesCore.h:43
float Size_t
Definition: RtypesCore.h:85
const Bool_t kFALSE
Definition: RtypesCore.h:90
const Bool_t kTRUE
Definition: RtypesCore.h:89
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
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:23
virtual void fill(RooAbsCategory &cat, Double_t weight=1.0)
Increment the counter of the table slot with the name corresponding to that of the current category s...
Definition: Roo1DTable.cxx:104
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:73
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Return the observables of this pdf given a set of observables.
Definition: RooAbsArg.h:276
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:730
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:475
const RefCountList_t & valueClients() const
List of all value clients of this object. Value clients receive value updates.
Definition: RooAbsArg.h:189
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:243
RooArgSet * getVariables(Bool_t stripDisconnected=kTRUE) const
Return RooArgSet with all variables (tree leaf nodes of expresssion tree)
Definition: RooAbsArg.cxx:1909
virtual Bool_t isDerived() const
Does value or shape of this arg depend on any other arg?
Definition: RooAbsArg.h:93
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Definition: RooAbsArg.cxx:1064
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
Definition: RooAbsArg.cxx:1467
void SetName(const char *name)
Set the name of the TNamed.
Definition: RooAbsArg.cxx:2216
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 base class for objects that represent a discrete value with a finite number of ...
virtual const char * getCurrentLabel() const
Return label string of current state.
const std::string & lookupName(value_type index) const
Get the name corresponding to the given index.
Roo1DTable * createTable(const char *label) const
Create a table matching the shape of this category.
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
Bool_t allInRange(const char *rangeSpec) const
Return true if all contained object report to have their value inside the specified range.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsDataStore is the abstract base class for data collection that use a TTree as internal storage m...
virtual void reset()=0
virtual void attachBuffers(const RooArgSet &extObs)=0
virtual Bool_t hasFilledCache() const
virtual 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:44
RooRealVar * meanVar(const RooRealVar &var, const char *cutSpec=0, const char *cutRange=0) const
Create a RooRealVar containing the mean of observable 'var' in this dataset.
virtual const RooArgSet * get() const
Definition: RooAbsData.h:87
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
RooRealVar * dataRealVar(const char *methodname, const RooRealVar &extVar) const
Internal method to check if given RooRealVar maps to a RooRealVar in this dataset.
Definition: RooAbsData.cxx:900
virtual Roo1DTable * table(const RooArgSet &catSet, const char *cuts="", const char *opts="") const
Construct table for product of categories in catSet.
Definition: RooAbsData.cxx:748
RooAbsDataStore * store()
Definition: RooAbsData.h:63
virtual void reset()
Definition: RooAbsData.cxx:313
Double_t sigma(const RooRealVar &var, const char *cutSpec=0, const char *cutRange=0) const
Definition: RooAbsData.h:194
virtual Double_t weightError(ErrorType etype=Poisson) const
Return error on current weight (dummy implementation returning zero)
Definition: RooAbsData.cxx:532
static Bool_t releaseVars(RooAbsData *)
If return value is true variables can be deleted.
Definition: RooAbsData.cxx:104
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:975
virtual Bool_t isNonPoissonWeighted() const
Definition: RooAbsData.h:112
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.
void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Interface for detailed printing of object.
Definition: RooAbsData.cxx:801
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.
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:628
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:238
virtual Double_t weightSquared() const =0
Double_t standMoment(const RooRealVar &var, Double_t order, const char *cutSpec=0, const char *cutRange=0) const
Calculate standardized moment.
Definition: RooAbsData.cxx:823
virtual Bool_t isWeighted() const
Definition: RooAbsData.h:108
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)
virtual void fill()
Definition: RooAbsData.cxx:299
RooArgSet _vars
Definition: RooAbsData.h:273
virtual void printName(std::ostream &os) const
Print name of dataset.
Definition: RooAbsData.cxx:778
virtual RooPlot * plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue &asymCat, PlotOpt o) const
Create and fill a histogram with the asymmetry N[+] - N[-] / ( N[+] + N[-] ), where N(+/-) is the num...
RooAbsData * getSimData(const char *idxstate)
virtual void Draw(Option_t *option="")
Forward draw command to data store.
virtual RooAbsData * reduceEng(const RooArgSet &varSubset, const RooFormulaVar *cutVar, const char *cutRange=0, std::size_t nStart=0, std::size_t=std::numeric_limits< std::size_t >::max(), Bool_t copyCache=kTRUE)=0
RooRealVar * rmsVar(const RooRealVar &var, const char *cutSpec=0, const char *cutRange=0) const
Create a RooRealVar containing the RMS of observable 'var' in this dataset.
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
Double_t corrcov(const RooRealVar &x, const RooRealVar &y, const char *cutSpec, const char *cutRange, Bool_t corr) const
Internal method to calculate single correlation and covariance elements.
Definition: RooAbsData.cxx:919
void convertToVectorStore()
Convert tree-based storage to vector-based storage.
Definition: RooAbsData.cxx:274
RooArgSet _cachedVars
Definition: RooAbsData.h:274
virtual void add(const RooArgSet &row, Double_t weight=1, Double_t weightError=0)=0
Bool_t getRange(const RooAbsRealLValue &var, Double_t &lowest, Double_t &highest, Double_t marginFrac=0, Bool_t symMode=kFALSE) const
Fill Doubles 'lowest' and 'highest' with the lowest and highest value of observable 'var' in this dat...
virtual Int_t numEntries() const
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:809
StorageType storageType
Definition: RooAbsData.h:240
RooAbsData & operator=(const RooAbsData &other)
Definition: RooAbsData.cxx:210
Double_t moment(const RooRealVar &var, Double_t order, const char *cutSpec=0, const char *cutRange=0) const
Calculate moment of requested order.
Definition: RooAbsData.cxx:842
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:546
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:794
std::map< std::string, RooAbsData * > _ownedComponents
Definition: RooAbsData.h:278
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:786
RooAbsDataStore * _dstore
External variables cached with this data set.
Definition: RooAbsData.h:276
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:60
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
const char * getPlotLabel() const
Get the label associated with the variable.
Definition: RooAbsReal.cxx:426
void setPlotLabel(const char *label)
Set the label associated with this variable.
Definition: RooAbsReal.cxx:436
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition: RooArgSet.h:126
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
RooCategory is an object to represent discrete states.
Definition: RooCategory.h:23
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition: RooCmdArg.h:28
Double_t getDouble(Int_t idx) const
Definition: RooCmdArg.h:85
Int_t getInt(Int_t idx) const
Definition: RooCmdArg.h:81
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:29
RooFormula internally uses ROOT's TFormula to compute user-defined expressions of RooAbsArgs.
Definition: RooFormula.h:28
Double_t eval(const RooArgSet *nset=0) const
Evalute all parameters/observables, and then evaluate formula.
Definition: RooFormula.cxx:340
A RooHist is a graphical representation of binned data based on the TGraphAsymmErrors class.
Definition: RooHist.h:27
RooLinkedList is an collection class for internal use, storing a collection of RooAbsArg pointers in ...
Definition: RooLinkedList.h:35
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:62
Bool_t Replace(const TObject *oldArg, const TObject *newArg)
Replace object 'oldArg' in collection with new object 'newArg'.
RooMultiCategory connects several RooAbsCategory objects into a single category.
static Double_t infinity()
Return internal infinity representation.
Definition: RooNumber.cxx:49
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:44
TAttMarker * getAttMarker(const char *name=0) const
Return a pointer to the marker attributes of the named object in this plot, or zero if the named obje...
Definition: RooPlot.cxx:907
void addObject(TObject *obj, Option_t *drawOptions="", Bool_t invisible=kFALSE)
Add a generic object to this plot.
Definition: RooPlot.cxx:443
Double_t getFitRangeNEvt() const
Return the number of events in the fit range.
Definition: RooPlot.h:141
Double_t getFitRangeBinW() const
Return the bin width that is being used to normalise the PDF.
Definition: RooPlot.h:144
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:897
RooAbsRealLValue * getPlotVar() const
Definition: RooPlot.h:139
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:887
TAxis * GetXaxis() const
Definition: RooPlot.cxx:1276
void updateNormVars(const RooArgSet &vars)
Install the given set of observables are reference normalization variables for this frame.
Definition: RooPlot.cxx:401
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:591
Int_t GetNbinsX() const
Definition: RooPlot.cxx:1278
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:1003
RooPlotable is a 'mix-in' base class that define the standard RooFit plotting and printing methods.
Definition: RooPrintable.h:25
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
void setError(Double_t value)
Definition: RooRealVar.h:64
TString * format(const RooCmdArg &formatArg) const
Format contents of RooRealVar for pretty printing on RooPlot parameter boxes.
Definition: RooRealVar.cxx:862
const RooAbsBinning & getBinning(const char *name=0, Bool_t verbose=kTRUE, Bool_t createOnTheFly=kFALSE) const
Return binning definition with name.
Definition: RooRealVar.cxx:323
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:261
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:8519
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:8662
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3275
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:3596
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition: TH1.cxx:8476
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:577
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:115
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:78
virtual TTree * CloneTree(Long64_t nentries=-1, Option_t *option="")
Create a clone of this tree and copy nentries.
Definition: TTree.cxx:3094
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg ZVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg AxisLabel(const char *name)
RooCmdArg AutoBinning(Int_t nbins=100, Double_t marginFactor=0.1)
RooCmdArg Binning(const RooAbsBinning &binning)
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
static double C[]
@ Optimization
Definition: RooGlobalFunc.h:68
@ InputArguments
Definition: RooGlobalFunc.h:68
static constexpr double s
static constexpr double pc
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:725
Definition: graph.py:1
static const char * what
Definition: stlLoader.cc:6
Double_t scaleFactor
Definition: RooAbsData.h:151
const char * cuts
Definition: RooAbsData.h:137
const char * cutRange
Definition: RooAbsData.h:142
Double_t addToWgtOther
Definition: RooAbsData.h:147
const char * histName
Definition: RooAbsData.h:143
const char * addToHistName
Definition: RooAbsData.h:145
Double_t addToWgtSelf
Definition: RooAbsData.h:146
RooAbsData::ErrorType etype
Definition: RooAbsData.h:141
RooAbsBinning * bins
Definition: RooAbsData.h:140
Bool_t correctForBinWidth
Definition: RooAbsData.h:150
Option_t * drawOptions
Definition: RooAbsData.h:139
auto * l
Definition: textangle.C:4
static long int sum(long int i)
Definition: Factory.cxx:2275