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