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