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