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