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