Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooDataHist.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 RooDataHist.cxx
19\class RooDataHist
20\ingroup Roofitcore
21
22Container class to hold N-dimensional binned data. Each bin's central
23coordinates in N-dimensional space are represented by a RooArgSet containing RooRealVar, RooCategory
24or RooStringVar objects, thus data can be binned in real and/or discrete dimensions.
25
26There is an unbinned equivalent, RooDataSet.
27
28### Inspecting a datahist
29Inspect a datahist using Print() to get the coordinates and `weight()` to get the bin contents:
30```
31datahist->Print("V");
32datahist->get(0)->Print("V"); std::cout << "w=" << datahist->weight(0) << std::endl;
33datahist->get(1)->Print("V"); std::cout << "w=" << datahist->weight(1) << std::endl;
34...
35```
36
37### Plotting data.
38See RooAbsData::plotOn().
39
40### Creating a datahist using RDataFrame
41See RooAbsDataHelper, rf408_RDataFrameToRooFit.C
42
43**/
44
45#include "RooDataHist.h"
46
47#include "Riostream.h"
48#include "RooMsgService.h"
50#include "RooAbsLValue.h"
51#include "RooArgList.h"
52#include "RooRealVar.h"
53#include "RooMath.h"
54#include "RooBinning.h"
55#include "RooPlot.h"
56#include "RooHistError.h"
57#include "RooCategory.h"
58#include "RooCmdConfig.h"
59#include "RooLinkedListIter.h"
60#include "RooTreeDataStore.h"
61#include "RooVectorDataStore.h"
62#include "RooTrace.h"
63#include "RooFormulaVar.h"
64#include "RooFormula.h"
65#include "RooUniformBinning.h"
66
67#include "RooFitImplHelpers.h"
68
69#include <ROOT/RSpan.hxx>
70#include <ROOT/StringUtils.hxx>
71
72#include "TAxis.h"
73#include "TH1.h"
74#include "TTree.h"
75#include "TBuffer.h"
76#include "TMath.h"
77#include "Math/Util.h"
78
79using std::string, std::ostream;
80
81
82
83////////////////////////////////////////////////////////////////////////////////
84/// Default constructor
85
90
91
92std::unique_ptr<RooAbsDataStore>
94{
96 ? static_cast<std::unique_ptr<RooAbsDataStore>>(std::make_unique<RooTreeDataStore>(name, title, vars))
97 : static_cast<std::unique_ptr<RooAbsDataStore>>(std::make_unique<RooVectorDataStore>(name, title, vars));
98}
99
100
101////////////////////////////////////////////////////////////////////////////////
102/// Constructor of an empty data hist from a RooArgSet defining the dimensions
103/// of the data space. The range and number of bins in each dimensions are taken
104/// from getMin()getMax(),getBins() of each RooAbsArg representing that
105/// dimension.
106///
107/// For real dimensions, the fit range and number of bins can be set independently
108/// of the plot range and number of bins, but it is advisable to keep the
109/// ratio of the plot bin width and the fit bin width an integer value.
110/// For category dimensions, the fit ranges always comprises all defined states
111/// and each state is always has its individual bin
112///
113/// To effectively bin real dimensions with variable bin sizes,
114/// construct a RooThresholdCategory of the real dimension to be binned variably.
115/// Set the thresholds at the desired bin boundaries, and construct the
116/// data hist as a function of the threshold category instead of the real variable.
117RooDataHist::RooDataHist(RooStringView name, RooStringView title, const RooArgSet& vars, const char* binningName) :
118 RooAbsData(name,title,vars)
119{
120 // Initialize datastore
122
123 initialize(binningName) ;
124
126
128}
129
130
131
132////////////////////////////////////////////////////////////////////////////////
133/// Constructor of a data hist from an existing data collection (binned or unbinned)
134/// The RooArgSet 'vars' defines the dimensions of the histogram.
135/// The range and number of bins in each dimensions are taken
136/// from getMin(), getMax(), getBins() of each argument passed.
137///
138/// For real dimensions, the fit range and number of bins can be set independently
139/// of the plot range and number of bins, but it is advisable to keep the
140/// ratio of the plot bin width and the fit bin width an integer value.
141/// For category dimensions, the fit ranges always comprises all defined states
142/// and each state is always has its individual bin
143///
144/// To effectively bin real dimensions with variable bin sizes,
145/// construct a RooThresholdCategory of the real dimension to be binned variably.
146/// Set the thresholds at the desired bin boundaries, and construct the
147/// data hist as a function of the threshold category instead of the real variable.
148///
149/// If the constructed data hist has less dimensions that in source data collection,
150/// all missing dimensions will be projected.
151
153 RooDataHist(name,title,vars)
154{
155 add(data,static_cast<const RooFormulaVar*>(nullptr),wgt);
156}
157
158
159
160////////////////////////////////////////////////////////////////////////////////
161/// Constructor of a data hist from a map of TH1,TH2 or TH3 that are collated into a x+1 dimensional
162/// RooDataHist where the added dimension is a category that labels the input source as defined
163/// in the histMap argument. The state names used in histMap must correspond to predefined states
164/// 'indexCat'
165///
166/// The RooArgList 'vars' defines the dimensions of the histogram.
167/// The ranges and number of bins are taken from the input histogram and must be the same in all histograms
168
170 std::map<string,TH1*> histMap, double wgt) :
171 RooAbsData(name,title,RooArgSet(vars,&indexCat))
172{
173 // Initialize datastore
175
176 importTH1Set(vars, indexCat, histMap, wgt, false) ;
177
180}
181
182
183
184////////////////////////////////////////////////////////////////////////////////
185/// Constructor of a data hist from a map of RooDataHists that are collated into a x+1 dimensional
186/// RooDataHist where the added dimension is a category that labels the input source as defined
187/// in the histMap argument. The state names used in histMap must correspond to predefined states
188/// 'indexCat'
189///
190/// The RooArgList 'vars' defines the dimensions of the histogram.
191/// The ranges and number of bins are taken from the input histogram and must be the same in all histograms
192
194 std::map<string,RooDataHist*> dhistMap, double wgt) :
195 RooAbsData(name,title,RooArgSet(vars,&indexCat))
196{
197 // Initialize datastore
199
200 importDHistSet(vars, indexCat, dhistMap, wgt) ;
201
204}
205
206
207
208////////////////////////////////////////////////////////////////////////////////
209/// Constructor of a data hist from an TH1,TH2 or TH3
210/// The RooArgSet 'vars' defines the dimensions of the histogram. The ranges
211/// and number of bins are taken from the input histogram, and the corresponding
212/// values are set accordingly on the arguments in 'vars'
213
214RooDataHist::RooDataHist(RooStringView name, RooStringView title, const RooArgList& vars, const TH1* hist, double wgt) :
215 RooAbsData(name,title,vars)
216{
217 // Initialize datastore
219
220 // Check consistency in number of dimensions
221 if (int(vars.size()) != hist->GetDimension()) {
222 std::stringstream errorMsgStream;
223 errorMsgStream << "RooDataHist::ctor(" << GetName() << ") ERROR: dimension of input histogram must match "
224 << "number of dimension variables";
225 const std::string errorMsg = errorMsgStream.str();
226 coutE(InputArguments) << errorMsg << std::endl;
227 throw std::invalid_argument(errorMsg);
228 }
229
230 importTH1(vars,*hist,wgt, false) ;
231
234}
235
236
237
238////////////////////////////////////////////////////////////////////////////////
239/// Constructor of a binned dataset from a RooArgSet defining the dimensions
240/// of the data space. The range and number of bins in each dimensions are taken
241/// from getMin() getMax(),getBins() of each RooAbsArg representing that
242/// dimension.
243///
244/// <table>
245/// <tr><th> Optional Argument <th> Effect
246/// <tr><td> Import(TH1&, bool impDens) <td> Import contents of the given TH1/2/3 into this binned dataset. The
247/// ranges and binning of the binned dataset are automatically adjusted to
248/// match those of the imported histogram.
249///
250/// Please note: for TH1& with unequal binning _only_,
251/// you should decide if you want to import the absolute bin content,
252/// or the bin content expressed as density. The latter is default and will
253/// result in the same histogram as the original TH1. For certain types of
254/// bin contents (containing efficiencies, asymmetries, or ratio is general)
255/// you should import the absolute value and set impDens to false
256///
257///
258/// <tr><td> Weight(double) <td> Apply given weight factor when importing histograms
259///
260/// <tr><td> Index(RooCategory&) <td> Prepare import of multiple TH1/1/2/3 into a N+1 dimensional RooDataHist
261/// where the extra discrete dimension labels the source of the imported histogram
262/// If the index category defines states for which no histogram is be imported
263/// the corresponding bins will be left empty.
264///
265/// <tr><td> Import(const char*, TH1&) <td> Import a THx to be associated with the given state name of the index category
266/// specified in Index(). If the given state name is not yet defined in the index
267/// category it will be added on the fly. The import command can be specified
268/// multiple times.
269/// <tr><td> Import(map<string,TH1*>&) <td> As above, but allows specification of many imports in a single operation
270/// <tr><td> `GlobalObservables(const RooArgSet&)` <td> Define the set of global observables to be stored in this RooDataHist.
271/// A snapshot of the passed RooArgSet is stored, meaning the values wont't change unexpectedly.
272/// </table>
273///
274
276 const RooCmdArg& arg4,const RooCmdArg& arg5,const RooCmdArg& arg6,const RooCmdArg& arg7,const RooCmdArg& arg8) :
277 RooAbsData(name,title,RooArgSet(vars,static_cast<RooAbsArg*>(RooCmdConfig::decodeObjOnTheFly("RooDataHist::RooDataHist", "IndexCat",0,nullptr,arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8))))
278{
279 // Initialize datastore
281
282 // Define configuration for this method
283 RooCmdConfig pc("RooDataHist::ctor(" + std::string(GetName()) + ")");
284 pc.defineObject("impHist","ImportHisto",0) ;
285 pc.defineInt("impDens","ImportHisto",0) ;
286 pc.defineObject("indexCat","IndexCat",0) ;
287 pc.defineObject("impSliceData","ImportDataSlice",0,nullptr,true) ; // array
288 pc.defineString("impSliceState","ImportDataSlice",0,"",true) ; // array
289 pc.defineDouble("weight","Weight",0,1) ;
290 pc.defineObject("dummy1","ImportDataSliceMany",0) ;
291 pc.defineSet("glObs","GlobalObservables",0,nullptr) ;
292 pc.defineMutex("ImportHisto","ImportDataSlice");
293 pc.defineDependency("ImportDataSlice","IndexCat") ;
294
296 l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
297 l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
298 l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
299 l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
300
301 // Process & check varargs
302 pc.process(l) ;
303 if (!pc.ok(true)) {
304 throw std::invalid_argument("Invalid command arguments passed to RooDataHist constructor!");
305 }
306
307 if(pc.getSet("glObs")) setGlobalObservables(*pc.getSet("glObs"));
308
309 TH1* impHist = static_cast<TH1*>(pc.getObject("impHist")) ;
310 bool impDens = pc.getInt("impDens") ;
311 double initWgt = pc.getDouble("weight") ;
312 RooCategory* indexCat = static_cast<RooCategory*>(pc.getObject("indexCat")) ;
313 const char* impSliceNames = pc.getString("impSliceState","",true) ;
314 const RooLinkedList& impSliceHistos = pc.getObjectList("impSliceData") ;
315
316
317 if (impHist) {
318
319 // Initialize importing contents from TH1
321
322 } else if (indexCat) {
323
324
325 // Initialize importing mapped set of RooDataHists and TH1s
326 std::map<std::string,RooDataHist*> dmap ;
327 std::map<std::string,TH1*> hmap ;
328 auto hiter = impSliceHistos.begin() ;
329 for (const auto& token : ROOT::Split(impSliceNames, ",", /*skipEmpty=*/true)) {
330
331 if (!indexCat->hasLabel(token)) {
332 std::stringstream errorMsgStream;
333 errorMsgStream << "RooDataHist::RooDataHist(\"" << GetName() << "\") "
334 << "you are providing import data for the category state \"" << token
335 << "\", but the index category \"" << indexCat->GetName() << "\" has no such state!";
336 const std::string errorMsg = errorMsgStream.str();
337 coutE(InputArguments) << errorMsg << std::endl;
338 throw std::invalid_argument(errorMsg);
339 }
340
341 if(auto dHist = dynamic_cast<RooDataHist*>(*hiter)) {
342 dmap[token] = dHist;
343 }
344 if(auto hHist = dynamic_cast<TH1*>(*hiter)) {
345 hmap[token] = hHist;
346 }
347 ++hiter;
348 }
349 if(!dmap.empty() && !hmap.empty()) {
350 std::stringstream errorMsgStream;
351 errorMsgStream << "RooDataHist::ctor(" << GetName() << ") ERROR: you can't import mix of TH1 and RooDataHist";
352 const std::string errorMsg = errorMsgStream.str();
353 coutE(InputArguments) << errorMsg << std::endl;
354 throw std::invalid_argument(errorMsg);
355 }
356 if (!dmap.empty()) {
357 importDHistSet(vars,*indexCat,dmap,initWgt);
358 }
359 if (!hmap.empty()) {
360 importTH1Set(vars,*indexCat,hmap,initWgt,false);
361 }
362
363
364 } else {
365
366 // Initialize empty
367 initialize();
368 }
369
372
373}
374
375
376
377
378////////////////////////////////////////////////////////////////////////////////
379/// Import data from given TH1/2/3 into this RooDataHist
380
381void RooDataHist::importTH1(const RooArgList& vars, const TH1& histo, double wgt, bool doDensityCorrection)
382{
383 // Adjust binning of internal observables to match that of input THx
384 Int_t offset[3]{0, 0, 0};
385 adjustBinning(vars, histo, offset) ;
386
387 // Initialize internal data structure
388 initialize();
389
390 // Define x,y,z as 1st, 2nd and 3rd observable
391 RooRealVar* xvar = static_cast<RooRealVar*>(_vars.find(vars.at(0)->GetName())) ;
392 RooRealVar* yvar = static_cast<RooRealVar*>(vars.at(1) ? _vars.find(vars.at(1)->GetName()) : nullptr ) ;
393 RooRealVar* zvar = static_cast<RooRealVar*>(vars.at(2) ? _vars.find(vars.at(2)->GetName()) : nullptr ) ;
394
395 // Transfer contents
396 Int_t xmin(0);
397 Int_t ymin(0);
398 Int_t zmin(0);
400 xmin = offset[0] ;
401 if (yvar) {
402 vset.add(*yvar) ;
403 ymin = offset[1] ;
404 }
405 if (zvar) {
406 vset.add(*zvar) ;
407 zmin = offset[2] ;
408 }
409
410 Int_t iX(0);
411 Int_t iY(0);
412 Int_t iz(0);
413 for (iX=0 ; iX < xvar->getBins() ; iX++) {
414 xvar->setBin(iX) ;
415 if (yvar) {
416 for (iY=0 ; iY < yvar->getBins() ; iY++) {
417 yvar->setBin(iY) ;
418 if (zvar) {
419 for (iz=0 ; iz < zvar->getBins() ; iz++) {
420 zvar->setBin(iz) ;
421 double bv = doDensityCorrection ? binVolume(vset) : 1;
422 add(vset,bv*histo.GetBinContent(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,2)) ;
423 }
424 } else {
425 double bv = doDensityCorrection ? binVolume(vset) : 1;
426 add(vset,bv*histo.GetBinContent(iX+1+xmin,iY+1+ymin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin,iY+1+ymin)*wgt,2)) ;
427 }
428 }
429 } else {
430 double bv = doDensityCorrection ? binVolume(vset) : 1 ;
431 add(vset,bv*histo.GetBinContent(iX+1+xmin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin)*wgt,2)) ;
432 }
433 }
434
435}
436
437namespace {
438bool checkConsistentAxes(const TH1* first, const TH1* second) {
439 return first->GetDimension() == second->GetDimension()
440 && first->GetNbinsX() == second->GetNbinsX()
441 && first->GetNbinsY() == second->GetNbinsY()
442 && first->GetNbinsZ() == second->GetNbinsZ()
443 && first->GetXaxis()->GetXmin() == second->GetXaxis()->GetXmin()
444 && first->GetXaxis()->GetXmax() == second->GetXaxis()->GetXmax()
445 && (first->GetNbinsY() == 1 || (first->GetYaxis()->GetXmin() == second->GetYaxis()->GetXmin()
446 && first->GetYaxis()->GetXmax() == second->GetYaxis()->GetXmax() ) )
447 && (first->GetNbinsZ() == 1 || (first->GetZaxis()->GetXmin() == second->GetZaxis()->GetXmin()
448 && first->GetZaxis()->GetXmax() == second->GetZaxis()->GetXmax() ) );
449}
450
452
453 // Relative tolerance for bin boundary comparison
454 constexpr double tolerance = 1e-6;
455
456 auto const& vars1 = *h1.get();
457 auto const& vars2 = *h2.get();
458
459 // Check if number of variables and names is consistent
460 if(!vars1.hasSameLayout(vars2)) {
461 return false;
462 }
463
464 for(std::size_t iVar = 0; iVar < vars1.size(); ++iVar) {
465 auto * var1 = dynamic_cast<RooRealVar*>(vars1[iVar]);
466 auto * var2 = dynamic_cast<RooRealVar*>(vars2[iVar]);
467
468 // Check if variables are consistently real-valued
469 if((!var1 && var2) || (var1 && !var2)) return false;
470
471 // Not a real-valued variable
472 if(!var1) continue;
473
474 // Now check the binning
475 auto const& bng1 = var1->getBinning();
476 auto const& bng2 = var2->getBinning();
477
478 // Compare bin numbers
479 if(bng1.numBins() != bng2.numBins()) return false;
480
481 std::size_t nBins = bng1.numBins();
482
483 // Compare bin boundaries
484 for(std::size_t iBin = 0; iBin < nBins; ++iBin) {
485 double v1 = bng1.binLow(iBin);
486 double v2 = bng2.binLow(iBin);
487 if(std::abs((v1 - v2) / v1) > tolerance) return false;
488 }
489 double v1 = bng1.binHigh(nBins - 1);
490 double v2 = bng2.binHigh(nBins - 1);
491 if(std::abs((v1 - v2) / v1) > tolerance) return false;
492 }
493 return true;
494}
495}
496
497
498////////////////////////////////////////////////////////////////////////////////
499/// Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
500/// in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
501/// and the import source
502
503void RooDataHist::importTH1Set(const RooArgList& vars, RooCategory& indexCat, std::map<string,TH1*> hmap, double wgt, bool doDensityCorrection)
504{
505 RooCategory* icat = static_cast<RooCategory*>(_vars.find(indexCat.GetName())) ;
506
507 TH1* histo(nullptr) ;
508 bool init(false) ;
509 for (const auto& hiter : hmap) {
510 // Store pointer to first histogram from which binning specification will be taken
511 if (!histo) {
512 histo = hiter.second;
513 } else {
514 if (!checkConsistentAxes(histo, hiter.second)) {
515 coutE(InputArguments) << "Axes of histogram " << hiter.second->GetName() << " are not consistent with first processed "
516 << "histogram " << histo->GetName() << std::endl;
517 throw std::invalid_argument("Axes of inputs for RooDataHist are inconsistent");
518 }
519 }
520 // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
521 if (!indexCat.hasLabel(hiter.first)) {
522 indexCat.defineType(hiter.first) ;
523 coutI(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << ") defining state \"" << hiter.first << "\" in index category " << indexCat.GetName() << std::endl ;
524 }
525 if (!icat->hasLabel(hiter.first)) {
526 icat->defineType(hiter.first) ;
527 }
528 }
529
530 // Check consistency in number of dimensions
531 if (histo && int(vars.size()) != histo->GetDimension()) {
532 coutE(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << "): dimension of input histogram must match "
533 << "number of continuous variables" << std::endl ;
534 throw std::invalid_argument("Inputs histograms for RooDataHist are not compatible with dimensions of variables.");
535 }
536
537 // Copy bins and ranges from THx to dimension observables
538 Int_t offset[3] ;
539 adjustBinning(vars,*histo,offset) ;
540
541 // Initialize internal data structure
542 if (!init) {
543 initialize();
544 init = true;
545 }
546
547 // Define x,y,z as 1st, 2nd and 3rd observable
548 RooRealVar* xvar = static_cast<RooRealVar*>(_vars.find(vars.at(0)->GetName())) ;
549 RooRealVar* yvar = static_cast<RooRealVar*>(vars.at(1) ? _vars.find(vars.at(1)->GetName()) : nullptr ) ;
550 RooRealVar* zvar = static_cast<RooRealVar*>(vars.at(2) ? _vars.find(vars.at(2)->GetName()) : nullptr ) ;
551
552 // Transfer contents
553 Int_t xmin(0);
554 Int_t ymin(0);
555 Int_t zmin(0);
557 double volume = xvar->getMax()-xvar->getMin() ;
558 xmin = offset[0] ;
559 if (yvar) {
560 vset.add(*yvar) ;
561 ymin = offset[1] ;
562 volume *= (yvar->getMax()-yvar->getMin()) ;
563 }
564 if (zvar) {
565 vset.add(*zvar) ;
566 zmin = offset[2] ;
567 volume *= (zvar->getMax()-zvar->getMin()) ;
568 }
569 double avgBV = volume / numEntries() ;
570
571 Int_t ic(0);
572 Int_t iX(0);
573 Int_t iY(0);
574 Int_t iz(0);
575 for (ic=0 ; ic < icat->numBins(nullptr) ; ic++) {
576 icat->setBin(ic) ;
577 histo = hmap[icat->getCurrentLabel()] ;
578 for (iX=0 ; iX < xvar->getBins() ; iX++) {
579 xvar->setBin(iX) ;
580 if (yvar) {
581 for (iY=0 ; iY < yvar->getBins() ; iY++) {
582 yvar->setBin(iY) ;
583 if (zvar) {
584 for (iz=0 ; iz < zvar->getBins() ; iz++) {
585 zvar->setBin(iz) ;
586 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
587 add(vset,bv*histo->GetBinContent(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,2)) ;
588 }
589 } else {
590 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
591 add(vset,bv*histo->GetBinContent(iX+1+xmin,iY+1+ymin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin,iY+1+ymin)*wgt,2)) ;
592 }
593 }
594 } else {
595 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
596 add(vset,bv*histo->GetBinContent(iX+1+xmin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin)*wgt,2)) ;
597 }
598 }
599 }
600
601}
602
603
604
605////////////////////////////////////////////////////////////////////////////////
606/// Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
607/// in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
608/// and the import source
609
610void RooDataHist::importDHistSet(const RooArgList & /*vars*/, RooCategory &indexCat,
611 std::map<std::string, RooDataHist *> dmap, double initWgt)
612{
613 auto *icat = static_cast<RooCategory *>(_vars.find(indexCat.GetName()));
614
615 RooDataHist *dhistForBinning = nullptr;
616
617 for (const auto &diter : dmap) {
618
619 std::string const &label = diter.first;
620 RooDataHist *dhist = diter.second;
621
622 if (!dhistForBinning) {
624 } else {
626 coutE(InputArguments) << "Layout or binning of histogram " << dhist->GetName()
627 << " is not consistent with first processed "
628 << "histogram " << dhistForBinning->GetName() << std::endl;
629 throw std::invalid_argument("Layout or binning of inputs for RooDataHist is inconsistent");
630 }
631 }
632
633 // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
634 if (!indexCat.hasLabel(label)) {
635 indexCat.defineType(label);
636 coutI(InputArguments) << "RooDataHist::importDHistSet(" << GetName() << ") defining state \"" << label
637 << "\" in index category " << indexCat.GetName() << std::endl;
638 }
639 if (!icat->hasLabel(label)) {
640 icat->defineType(label);
641 }
642 }
643
644 // adjust the binning of the created histogram
646 auto *ourVar = dynamic_cast<RooRealVar *>(_vars.find(theirVar->GetName()));
647 if (!theirVar || !ourVar)
648 continue;
649 ourVar->setBinning(theirVar->getBinning());
650 }
651
652 initialize();
653
654 for (const auto &diter : dmap) {
655 std::string const &label = diter.first;
656 RooDataHist *dhist = diter.second;
657
658 icat->setLabel(label.c_str());
659
660 // Transfer contents
661 for (Int_t i = 0; i < dhist->numEntries(); i++) {
662 _vars.assign(*dhist->get(i));
663 add(_vars, dhist->weight() * initWgt, pow(dhist->weightError(SumW2), 2));
664 }
665 }
666}
667
668////////////////////////////////////////////////////////////////////////////////
669/// Helper doing the actual work of adjustBinning().
670
673{
674 const std::string ourVarName(ourVar->GetName() ? ourVar->GetName() : "");
675 const std::string ownName(GetName() ? GetName() : "");
676 // RooRealVar is derived from RooAbsRealLValue which is itself
677 // derived from RooAbsReal and a virtual class RooAbsLValue
678 // supplying setter functions, check if ourVar is indeed derived
679 // as real
680 if (!dynamic_cast<RooAbsReal *>(ourVar)) {
681 coutE(InputArguments) << "RooDataHist::adjustBinning(" << ownName << ") ERROR: dimension " << ourVarName
682 << " must be real\n";
683 throw std::logic_error("Incorrect type object (" + ourVarName +
684 ") passed as argument to RooDataHist::_adjustBinning. Please report this issue.");
685 }
686
687 const double xlo = theirVar.getMin();
688 const double xhi = theirVar.getMax();
689
690 const bool isUniform = !axis.GetXbins()->GetArray();
691 std::unique_ptr<RooAbsBinning> xbins;
692
693 if (!isUniform) {
694 xbins = std::make_unique<RooBinning>(axis.GetNbins(), axis.GetXbins()->GetArray());
695 } else {
696 xbins = std::make_unique<RooUniformBinning>(axis.GetXmin(), axis.GetXmax(), axis.GetNbins());
697 }
698
699 const double tolerance = 1e-6 * xbins->averageBinWidth();
700
701 // Adjust xlo/xhi to nearest boundary
702 const int iBinLo = xbins->binNumber(xlo + tolerance);
703 const int iBinHi = xbins->binNumber(xhi - tolerance);
704 const int nBinsAdj = iBinHi - iBinLo + 1;
705 const double xloAdj = xbins->binLow(iBinLo);
706 const double xhiAdj = xbins->binHigh(iBinHi);
707
708 if (isUniform) {
709 xbins = std::make_unique<RooUniformBinning>(xloAdj, xhiAdj, nBinsAdj);
710 theirVar.setRange(xloAdj, xhiAdj);
711 } else {
712 xbins->setRange(xloAdj, xhiAdj);
713 theirVar.setBinning(*xbins);
714 }
715
716 if (std::abs(xloAdj - xlo) > tolerance || std::abs(xhiAdj - xhi) > tolerance) {
717 coutI(DataHandling) << "RooDataHist::adjustBinning(" << ownName << "): fit range of variable " << ourVarName
718 << " expanded to nearest bin boundaries: [" << xlo << "," << xhi << "] --> [" << xloAdj << ","
719 << xhiAdj << "]"
720 << "\n";
721 }
722
723 ourVar->setBinning(*xbins);
724
725 // The offset is the bin number of the adjusted lower limit of the RooFit
726 // variable in the original TH1 histogram, starting from zero.
727 if (offset) {
728 *offset = axis.FindFixBin(xloAdj + tolerance) - 1;
729 }
730}
731
732////////////////////////////////////////////////////////////////////////////////
733/// Adjust binning specification on first and optionally second and third
734/// observable to binning in given reference TH1. Used by constructors
735/// that import data from an external TH1.
736/// Both the variables in vars and in this RooDataHist are adjusted.
737/// @param vars List with variables that are supposed to have their binning adjusted.
738/// @param href Reference histogram that dictates the binning
739/// @param offset If not nullptr, a possible bin count offset for the axes x,y,z is saved here as Int_t[3]
740
742{
743 auto xvar = static_cast<RooRealVar*>(_vars.find(*vars.at(0)) );
744 _adjustBinning(*static_cast<RooRealVar*>(vars.at(0)), *href.GetXaxis(), xvar, offset ? &offset[0] : nullptr);
745
746 if (vars.at(1)) {
747 auto yvar = static_cast<RooRealVar*>(_vars.find(*vars.at(1)));
748 if (yvar)
749 _adjustBinning(*static_cast<RooRealVar*>(vars.at(1)), *href.GetYaxis(), yvar, offset ? &offset[1] : nullptr);
750 }
751
752 if (vars.at(2)) {
753 auto zvar = static_cast<RooRealVar*>(_vars.find(*vars.at(2)));
754 if (zvar)
755 _adjustBinning(*static_cast<RooRealVar*>(vars.at(2)), *href.GetZaxis(), zvar, offset ? &offset[2] : nullptr);
756 }
757
758}
759
760
761namespace {
762/// Clone external weight arrays, unless the external array is nullptr.
763void cloneArray(double*& ours, const double* theirs, std::size_t n) {
764 if (ours) delete[] ours;
765 ours = nullptr;
766 if (!theirs) return;
767 ours = new double[n];
768 std::copy(theirs, theirs+n, ours);
769}
770
771/// Allocate and initialise an array with desired size and values.
772void initArray(double*& arr, std::size_t n, double val) {
773 if (arr) delete[] arr;
774 arr = nullptr;
775 if (n == 0) return;
776 arr = new double[n];
777 std::fill(arr, arr+n, val);
778}
779}
780
781
782////////////////////////////////////////////////////////////////////////////////
783/// Initialization procedure: allocate weights array, calculate
784/// multipliers needed for N-space to 1-dim array jump table,
785/// and fill the internal tree with all bin center coordinates
786
787void RooDataHist::initialize(const char* binningName, bool fillTree)
788{
789 _lvvars.clear();
790 _lvbins.clear();
791
792 // Fill array of LValue pointers to variables
793 for (unsigned int i = 0; i < _vars.size(); ++i) {
794 if (binningName) {
795 RooRealVar* rrv = dynamic_cast<RooRealVar*>(_vars[i]);
796 if (rrv) {
797 rrv->setBinning(rrv->getBinning(binningName));
798 }
799 }
800
801 auto lvarg = dynamic_cast<RooAbsLValue*>(_vars[i]);
802 assert(lvarg);
803 _lvvars.push_back(lvarg);
804
805 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
806 _lvbins.emplace_back(binning ? binning->clone() : nullptr);
807 }
808
809
810 // Allocate coefficients array
811 _idxMult.resize(_vars.size()) ;
812
813 _arrSize = 1 ;
814 unsigned int n = 0u;
815 for (const auto var : _vars) {
816 auto arg = dynamic_cast<const RooAbsLValue*>(var);
817 assert(arg);
818
819 // Calculate sub-index multipliers for master index
820 for (unsigned int i = 0u; i<n; i++) {
821 _idxMult[i] *= arg->numBins() ;
822 }
823 _idxMult[n++] = 1 ;
824
825 // Calculate dimension of weight array
826 _arrSize *= arg->numBins() ;
827 }
828
829 // Allocate and initialize weight array if necessary
830 if (!_wgt) {
831 initArray(_wgt, _arrSize, 0.);
832 delete[] _errLo; _errLo = nullptr;
833 delete[] _errHi; _errHi = nullptr;
834 delete[] _sumw2; _sumw2 = nullptr;
836
837 // Refill array pointers in data store when reading
838 // from Streamer
839 if (!fillTree) {
841 }
842 }
843
844 if (!fillTree) return ;
845
846 // Fill TTree with bin center coordinates
847 // Calculate plot bins of components from master index
848
849 for (Int_t ibin=0 ; ibin < _arrSize ; ibin++) {
850 Int_t j(0);
851 Int_t idx(0);
852 Int_t tmp(ibin);
853 double theBinVolume(1) ;
854 for (auto arg2 : _lvvars) {
855 idx = tmp / _idxMult[j] ;
856 tmp -= idx*_idxMult[j++] ;
857 arg2->setBin(idx) ;
858 theBinVolume *= arg2->getBinWidth(idx) ;
859 }
861
862 fill() ;
863 }
864
865
866}
867
868
869////////////////////////////////////////////////////////////////////////////////
870
872{
873 if (!_binbounds.empty()) return;
874 for (auto& it : _lvbins) {
875 _binbounds.push_back(std::vector<double>());
876 if (it) {
877 std::vector<double>& bounds = _binbounds.back();
878 bounds.reserve(2 * it->numBins());
879 for (Int_t i = 0; i < it->numBins(); ++i) {
880 bounds.push_back(it->binLow(i));
881 bounds.push_back(it->binHigh(i));
882 }
883 }
884 }
885}
886
887
888////////////////////////////////////////////////////////////////////////////////
889/// Copy constructor
890
892 RooAbsData(other,newname), RooDirItem(), _arrSize(other._arrSize), _idxMult(other._idxMult), _pbinvCache(other._pbinvCache)
893{
894 // Allocate and initialize weight array
895 assert(_arrSize == other._arrSize);
896 cloneArray(_wgt, other._wgt, other._arrSize);
897 cloneArray(_errLo, other._errLo, other._arrSize);
898 cloneArray(_errHi, other._errHi, other._arrSize);
899 cloneArray(_binv, other._binv, other._arrSize);
900 cloneArray(_sumw2, other._sumw2, other._arrSize);
901
902 // Fill array of LValue pointers to variables
903 for (const auto rvarg : _vars) {
904 auto lvarg = dynamic_cast<RooAbsLValue*>(rvarg);
905 assert(lvarg);
906 _lvvars.push_back(lvarg);
907 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
908 _lvbins.emplace_back(binning ? binning->clone() : nullptr) ;
909 }
910
912}
913
914
915////////////////////////////////////////////////////////////////////////////////
916/// Implementation of RooAbsData virtual method that drives the RooAbsData::reduce() methods
917
918std::unique_ptr<RooAbsData> RooDataHist::reduceEng(const RooArgSet& varSubset, const RooFormulaVar* cutVar, const char* cutRange,
919 std::size_t nStart, std::size_t nStop) const
920{
921 checkInit() ;
924 auto rdh = std::make_unique<RooDataHist>(GetName(), GetTitle(), myVarSubset);
925
926 RooFormulaVar* cloneVar = nullptr;
927 std::unique_ptr<RooArgSet> tmp;
928 if (cutVar) {
929 tmp = std::make_unique<RooArgSet>();
930 // Deep clone cutVar and attach clone to this dataset
931 if (RooArgSet(*cutVar).snapshot(*tmp)) {
932 coutE(DataHandling) << "RooDataHist::reduceEng(" << GetName() << ") Couldn't deep-clone cut variable, abort," << std::endl ;
933 return nullptr;
934 }
935 cloneVar = static_cast<RooFormulaVar*>(tmp->find(*cutVar));
936 cloneVar->attachDataSet(*this) ;
937 }
938
939 double lo;
940 double hi;
941 const std::size_t nevt = nStop < static_cast<std::size_t>(numEntries()) ? nStop : static_cast<std::size_t>(numEntries());
942 for (auto i=nStart; i<nevt ; i++) {
943 const RooArgSet* row = get(i) ;
944
945 bool doSelect(true) ;
946 if (cutRange) {
947 for (const auto arg : *row) {
948 if (!arg->inRange(cutRange)) {
949 doSelect = false ;
950 break ;
951 }
952 }
953 }
954 if (!doSelect) continue ;
955
956 if (!cloneVar || cloneVar->getVal()) {
957 weightError(lo,hi,SumW2) ;
958 rdh->add(*row,weight(),lo*lo) ;
959 }
960 }
961
962 return rdh ;
963}
964
965
966
967////////////////////////////////////////////////////////////////////////////////
968/// Destructor
969
971{
972 delete[] _wgt;
973 delete[] _errLo;
974 delete[] _errHi;
975 delete[] _sumw2;
976 delete[] _binv;
977
978 removeFromDir(this) ;
980}
981
982
983
984
985////////////////////////////////////////////////////////////////////////////////
986/// Calculate bin number of the given coordinates. If only a subset of the internal
987/// coordinates are passed, the missing coordinates are taken at their current value.
988/// \param[in] coord Variables that are representing the coordinates.
989/// \param[in] fast If the variables in `coord` and the ones of the data hist have the
990/// same size and layout, `fast` can be set to skip checking that all variables are
991/// present in `coord`.
993 checkInit() ;
994 return calcTreeIndex(coord, fast);
995}
996
998 bool correctForBinSize) const
999{
1000 std::vector<double> vals(_arrSize);
1001 for (std::size_t i = 0; i < vals.size(); ++i) {
1002 vals[i] = correctForBinSize ? _wgt[i] / _binv[i] : _wgt[i];
1003 }
1004 return ctx.buildArg(vals);
1005}
1006
1008 const RooAbsCollection &coords, bool reverse) const
1009{
1010 assert(coords.size() == _vars.size());
1011
1012 std::string code;
1013 int idxMult = 1;
1014
1015 for (std::size_t i = 0; i < _vars.size(); ++i) {
1016
1017 std::size_t iVar = reverse ? _vars.size() - 1 - i : i;
1018 const RooAbsArg *internalVar = _vars[iVar];
1019 const RooAbsArg *theVar = coords[iVar];
1020
1021 const RooAbsBinning *binning = _lvbins[iVar].get();
1022 if (!binning) {
1023 coutE(InputArguments) << "RooHistPdf::weight(" << GetName()
1024 << ") ERROR: Code Squashing currently does not support category values." << std::endl;
1025 return "";
1026 }
1027
1028 if (i > 0)
1029 code += " + ";
1030 code += binning->translateBinNumber(ctx, *theVar, idxMult);
1031
1032 // Use RooAbsLValue here because it also generalized to categories, which
1033 // is useful in the future. dynamic_cast because it's a cross-cast.
1034 idxMult *= dynamic_cast<RooAbsLValue const *>(internalVar)->numBins();
1035 }
1036
1037 return _vars.size() == 1 ? code : "(" + code + ")";
1038}
1039
1040////////////////////////////////////////////////////////////////////////////////
1041/// Calculate the bin index corresponding to the coordinates passed as argument.
1042/// \param[in] coords Coordinates. If `fast == false`, these can be partial.
1043/// \param[in] fast Promise that the coordinates in `coords` have the same order
1044/// as the internal coordinates. In this case, values are looked up only by index.
1045std::size_t RooDataHist::calcTreeIndex(const RooAbsCollection& coords, bool fast) const
1046{
1047 // With fast, caller promises that layout of `coords` is identical to our internal `vars`.
1048 // Previously, this was verified with an assert in debug mode like this:
1049 //
1050 // assert(!fast || coords.hasSameLayout(_vars));
1051 //
1052 // However, there are usecases where the externally provided `coords` have
1053 // different names than the internal variables, even though they correspond
1054 // to each other. For example, if the observables in the computation graph
1055 // are renamed with `redirectServers`. Hence, we can't do a meaningful assert
1056 // here.
1057
1058 if (&_vars == &coords)
1059 fast = true;
1060
1061 std::size_t masterIdx = 0;
1062
1063 for (unsigned int i=0; i < _vars.size(); ++i) {
1064 const RooAbsArg* internalVar = _vars[i];
1065 const RooAbsBinning* binning = _lvbins[i].get();
1066
1067 // Find the variable that we need values from.
1068 // That's either the variable directly from the external coordinates
1069 // or we find the external one that has the same name as "internalVar".
1070 const RooAbsArg* theVar = fast ? coords[i] : coords.find(*internalVar);
1071 if (!theVar) {
1072 // Variable is not in external coordinates. Use current internal value.
1074 }
1075 // If fast is on, users promise that the sets have the same layout:
1076 //
1077 // assert(!fast || strcmp(internalVar->GetName(), theVar->GetName()) == 0);
1078 //
1079 // This assert is commented out for the same reasons that applied to the
1080 // other assert explained above.
1081
1082 if (binning) {
1083 assert(dynamic_cast<const RooAbsReal*>(theVar));
1084 const double val = static_cast<const RooAbsReal*>(theVar)->getVal();
1085 masterIdx += _idxMult[i] * binning->binNumber(val);
1086 } else {
1087 // We are a category. No binning.
1088 assert(dynamic_cast<const RooAbsCategoryLValue*>(theVar));
1089 auto cat = static_cast<const RooAbsCategoryLValue*>(theVar);
1090 masterIdx += _idxMult[i] * cat->getBin(static_cast<const char*>(nullptr));
1091 }
1092 }
1093
1094 return masterIdx ;
1095}
1096
1097
1098////////////////////////////////////////////////////////////////////////////////
1099/// Back end function to plotting functionality. Plot RooDataHist on given
1100/// frame in mode specified by plot options 'o'. The main purpose of
1101/// this function is to match the specified binning on 'o' to the
1102/// internal binning of the plot observable in this RooDataHist.
1103/// \note see RooAbsData::plotOnImpl() for plotting options.
1105{
1106 checkInit() ;
1107 if (o.bins) return RooAbsData::plotOnImpl(frame,o) ;
1108
1109 if(!frame) {
1110 coutE(InputArguments) << ClassName() << "::" << GetName() << ":plotOn: frame is null" << std::endl;
1111 return nullptr;
1112 }
1113 auto var= static_cast<RooAbsRealLValue*>(frame->getPlotVar());
1114 if(!var) {
1115 coutE(InputArguments) << ClassName() << "::" << GetName()
1116 << ":plotOn: frame does not specify a plot variable" << std::endl;
1117 return nullptr;
1118 }
1119
1120 auto dataVar = static_cast<RooRealVar*>(_vars.find(*var));
1121 if (!dataVar) {
1122 coutE(InputArguments) << ClassName() << "::" << GetName()
1123 << ":plotOn: dataset doesn't contain plot frame variable" << std::endl;
1124 return nullptr;
1125 }
1126
1127 o.bins = &dataVar->getBinning() ;
1128 return RooAbsData::plotOnImpl(frame,o) ;
1129}
1130
1131
1132////////////////////////////////////////////////////////////////////////////////
1133/// A vectorized version of interpolateDim for boundary safe quadratic
1134/// interpolation of one dimensional histograms.
1135///
1136/// \param[out] output An array of interpolated weights corresponding to the
1137/// values in xVals.
1138/// \param[in] xVals An array of event coordinates for which the weights should be
1139/// calculated.
1140/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1141/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1142/// Underflow bins are assumed to have weight zero and
1143/// overflow bins have weight one. Otherwise, the
1144/// histogram is mirrored at the boundaries for the
1145/// interpolation.
1146
1147void RooDataHist::interpolateQuadratic(double* output, std::span<const double> xVals,
1149{
1150 const std::size_t nBins = numEntries();
1151 const std::size_t nEvents = xVals.size();
1152
1153 RooAbsBinning const& binning = *_lvbins[0];
1154 // Reuse the output buffer for bin indices and zero-initialize it
1155 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1156 std::fill(binIndices, binIndices + nEvents, 0);
1157 binning.binNumbers(xVals.data(), binIndices, nEvents);
1158
1159 // Extend coordinates and weights with one extra point before the first bin
1160 // and one extra point after the last bin. This means the original histogram
1161 // bins span elements 1 to nBins in coordsExt and weightsExt
1162 std::vector<double> coordsExt(nBins+3);
1163 double* binCoords = coordsExt.data() + 2;
1164 binCoords[0] = binning.lowBound() + 0.5*_binv[0];
1165 for (std::size_t binIdx = 1; binIdx < nBins ; ++binIdx) {
1166 if (binning.isUniform()) {
1167 double binWidth = _binv[0];
1168 binCoords[binIdx] = binIdx*binWidth + binCoords[0];
1169 }
1170 else {
1171 double binCentDiff = 0.5*_binv[binIdx-1] + 0.5*_binv[binIdx];
1173 }
1174 }
1175
1176 std::vector<double> weightsExt(nBins+3);
1177 // Fill weights for bins that are inside histogram boundaries
1178 for (std::size_t binIdx = 0; binIdx < nBins; ++binIdx) {
1180 }
1181
1182 if (cdfBoundaries) {
1183 coordsExt[0] = - 1e-10 + binning.lowBound();
1184 weightsExt[0] = 0.;
1185
1186 coordsExt[1] = binning.lowBound();
1187 weightsExt[1] = 0.;
1188
1189 coordsExt[nBins+2] = binning.highBound();
1190 weightsExt[nBins+2] = 1.;
1191 }
1192 else {
1193 // Mirror first two bins and last bin
1194 coordsExt[0] = binCoords[1] - 2*_binv[0] - _binv[1];
1195 weightsExt[0] = weightsExt[3];
1196
1197 coordsExt[1] = binCoords[0] - _binv[0];
1198 weightsExt[1] = weightsExt[2];
1199
1200 coordsExt[nBins+2] = binCoords[nBins-1] + _binv[nBins-1];
1201 weightsExt[nBins+2] = weightsExt[nBins+1];
1202 }
1203
1204 // We use the current bin center and two bin centers on the left for
1205 // interpolation if xVal is to the left of the current bin center
1206 for (std::size_t i = 0; i < nEvents ; ++i) {
1207 double xVal = xVals[i];
1208 std::size_t binIdx = binIndices[i] + 2;
1209
1210 // If xVal is to the right of the current bin center, shift all bin
1211 // coordinates one step to the right and use that for the interpolation
1212 if (xVal > coordsExt[binIdx]) {
1213 binIdx += 1;
1214 }
1215
1216 double x1 = coordsExt[binIdx-2];
1217 double y1 = weightsExt[binIdx-2];
1218
1219 double x2 = coordsExt[binIdx-1];
1220 double y2 = weightsExt[binIdx-1];
1221
1222 double x3 = coordsExt[binIdx];
1223 double y3 = weightsExt[binIdx];
1224
1225 // Evaluate a few repeated factors
1226 double quotient = (x3-x1) / (x2-x1);
1227 double x1Sqrd = x1*x1;
1228 double x3Sqrd = x3*x3;
1229 // Solve coefficients in system of three quadratic equations!
1230 double secondCoeff = (y3 - y1 - (y2-y1) * quotient) / (x3Sqrd - x1Sqrd - (x2*x2 - x1Sqrd) * quotient);
1231 double firstCoeff = (y3 - y1 - secondCoeff*(x3Sqrd - x1Sqrd)) / (x3-x1);
1233 // Get the interpolated weight using the equation of a second degree polynomial
1235 }
1236}
1237
1238
1239////////////////////////////////////////////////////////////////////////////////
1240/// A vectorized version of interpolateDim for boundary safe linear
1241/// interpolation of one dimensional histograms.
1242///
1243/// \param[out] output An array of interpolated weights corresponding to the
1244/// values in xVals.
1245/// \param[in] xVals An array of event coordinates for which the weights should be
1246/// calculated.
1247/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1248/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1249/// Underflow bins are assumed to have weight zero and
1250/// overflow bins have weight one. Otherwise, the
1251/// histogram is mirrored at the boundaries for the
1252/// interpolation.
1253
1254void RooDataHist::interpolateLinear(double* output, std::span<const double> xVals,
1256{
1257 const std::size_t nBins = numEntries();
1258 const std::size_t nEvents = xVals.size();
1259
1260 RooAbsBinning const& binning = *_lvbins[0];
1261 // Reuse the output buffer for bin indices and zero-initialize it
1262 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1263 std::fill(binIndices, binIndices + nEvents, 0);
1264 binning.binNumbers(xVals.data(), binIndices, nEvents);
1265
1266 // Extend coordinates and weights with one extra point before the first bin
1267 // and one extra point after the last bin. This means the original histogram
1268 // bins span elements 1 to nBins in coordsExt and weightsExt
1269 std::vector<double> coordsExt(nBins+2);
1270 double* binCoords = coordsExt.data() + 1;
1271 binCoords[0] = binning.lowBound() + 0.5*_binv[0];
1272 for (std::size_t binIdx = 1; binIdx < nBins ; ++binIdx) {
1273 if (binning.isUniform()) {
1274 double binWidth = _binv[0];
1275 binCoords[binIdx] = binIdx*binWidth + binCoords[0];
1276 }
1277 else {
1278 double binCentDiff = 0.5*_binv[binIdx-1] + 0.5*_binv[binIdx];
1280 }
1281 }
1282
1283 std::vector<double> weightsExt(nBins+2);
1284 // Fill weights for bins that are inside histogram boundaries
1285 for (std::size_t binIdx = 0; binIdx < nBins; ++binIdx) {
1287 }
1288
1289 // Fill weights for bins that are outside histogram boundaries
1290 if (cdfBoundaries) {
1291 coordsExt[0] = binning.lowBound();
1292 weightsExt[0] = 0.;
1293 coordsExt[nBins+1] = binning.highBound();
1294 weightsExt[nBins+1] = 1.;
1295 }
1296 else {
1297 // Mirror first and last bins
1298 coordsExt[0] = binCoords[0] - _binv[0];
1299 weightsExt[0] = weightsExt[1];
1300 coordsExt[nBins+1] = binCoords[nBins-1] + _binv[nBins-1];
1301 weightsExt[nBins+1] = weightsExt[nBins];
1302 }
1303
1304 // Interpolate between current bin center and one bin center to the left
1305 // if xVal is to the left of the current bin center
1306 for (std::size_t i = 0; i < nEvents ; ++i) {
1307 double xVal = xVals[i];
1308 std::size_t binIdx = binIndices[i] + 1;
1309
1310 // If xVal is to the right of the current bin center, interpolate between
1311 // current bin center and one bin center to the right instead
1312 if (xVal > coordsExt[binIdx]) { binIdx += 1; }
1313
1314 double x1 = coordsExt[binIdx-1];
1315 double y1 = weightsExt[binIdx-1];
1316 double x2 = coordsExt[binIdx];
1317 double y2 = weightsExt[binIdx];
1318
1319 // Find coefficients by solving a system of two linear equations
1320 double firstCoeff = (y2-y1) / (x2-x1);
1321 double zerothCoeff = y1 - firstCoeff * x1;
1322 // Get the interpolated weight using the equation of a straight line
1324 }
1325}
1326
1327
1328////////////////////////////////////////////////////////////////////////////////
1329/// A vectorized version of RooDataHist::weight() for one dimensional histograms
1330/// with up to one dimensional interpolation.
1331/// \param[out] output An array of weights corresponding the values in xVals.
1332/// \param[in] xVals An array of coordinates for which the weights should be
1333/// calculated.
1334/// \param[in] intOrder Interpolation order; 0th and 1st order are supported.
1335/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1336/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1337/// Underflow bins are assumed to have weight zero and
1338/// overflow bins have weight one. Otherwise, the
1339/// histogram is mirrored at the boundaries for the
1340/// interpolation.
1341
1342void RooDataHist::weights(double* output, std::span<double const> xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
1343{
1344 auto const nEvents = xVals.size();
1345
1346 if (intOrder == 0) {
1347 RooAbsBinning const& binning = *_lvbins[0];
1348
1349 // Reuse the output buffer for bin indices and zero-initialize it
1350 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1351 std::fill(binIndices, binIndices + nEvents, 0);
1352 binning.binNumbers(xVals.data(), binIndices, nEvents);
1353
1354 for (std::size_t i=0; i < nEvents; ++i) {
1355 auto binIdx = binIndices[i];
1357 }
1358 }
1359 else if (intOrder == 1) {
1361 }
1362 else if (intOrder == 2) {
1364 }
1365 else {
1366 // Higher dimensional scenarios not yet implemented
1367 coutE(InputArguments) << "RooDataHist::weights(" << GetName() << ") interpolation in "
1368 << intOrder << " dimensions not yet implemented" << std::endl ;
1369 // Fall back to 1st order interpolation
1371 }
1372}
1373
1374
1375////////////////////////////////////////////////////////////////////////////////
1376/// A faster version of RooDataHist::weight that assumes the passed arguments
1377/// are aligned with the histogram variables.
1378/// \param[in] bin Coordinates for which the weight should be calculated.
1379/// Has to be aligned with the internal histogram variables.
1380/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1381/// used for the interpolation. If zero, the bare weight for
1382/// the bin enclosing the coordinatesis returned.
1383/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1384/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1385/// underflow bins are assumed to have weight zero and
1386/// overflow bins have weight one. Otherwise, the
1387/// histogram is mirrored at the boundaries for the
1388/// interpolation.
1389
1391{
1392 checkInit() ;
1393
1394 // Handle illegal intOrder values
1395 if (intOrder<0) {
1396 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << std::endl ;
1397 return 0 ;
1398 }
1399
1400 // Handle no-interpolation case
1401 if (intOrder==0) {
1402 const auto idx = calcTreeIndex(bin, true);
1403 return correctForBinSize ? _wgt[idx] / _binv[idx] : _wgt[idx];
1404 }
1405
1406 // Handle all interpolation cases
1408}
1409
1410
1411////////////////////////////////////////////////////////////////////////////////
1412/// Return the weight at given coordinates with optional interpolation.
1413/// \param[in] bin Coordinates for which the weight should be calculated.
1414/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1415/// used for the interpolation. If zero, the bare weight for
1416/// the bin enclosing the coordinatesis returned.
1417/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1418/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1419/// underflow bins are assumed to have weight zero and
1420/// overflow bins have weight one. Otherwise, the
1421/// histogram is mirrored at the boundaries for the
1422/// interpolation.
1423/// \param[in] oneSafe Ignored.
1424
1425double RooDataHist::weight(const RooArgSet& bin, Int_t intOrder, bool correctForBinSize, bool cdfBoundaries, bool /*oneSafe*/)
1426{
1427 checkInit() ;
1428
1429 // Handle illegal intOrder values
1430 if (intOrder<0) {
1431 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << std::endl ;
1432 return 0 ;
1433 }
1434
1435 // Handle no-interpolation case
1436 if (intOrder==0) {
1437 const auto idx = calcTreeIndex(bin, false);
1438 return correctForBinSize ? _wgt[idx] / _binv[idx] : _wgt[idx];
1439 }
1440
1441 // Handle all interpolation cases
1442 _vars.assignValueOnly(bin) ;
1443
1445}
1446
1447
1448////////////////////////////////////////////////////////////////////////////////
1449/// Return the weight at given coordinates with interpolation.
1450/// \param[in] bin Coordinates for which the weight should be calculated.
1451/// Has to be aligned with the internal histogram variables.
1452/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1453/// used for the interpolation.
1454/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1455/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1456/// underflow bins are assumed to have weight zero and
1457/// overflow bins have weight one. Otherwise, the
1458/// histogram is mirrored at the boundaries for the
1459/// interpolation.
1460
1462 VarInfo const& varInfo = getVarInfo();
1463
1464 const auto centralIdx = calcTreeIndex(bin, true);
1465
1466 double wInt{0} ;
1467 if (varInfo.nRealVars == 1) {
1468
1469 // buffer needs to be 2 x (interpolation order + 1), with the factor 2 for x and y.
1470 _interpolationBuffer.resize(2 * intOrder + 2);
1471
1472 // 1-dimensional interpolation
1473 auto const& realX = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx1]);
1475
1476 } else if (varInfo.nRealVars == 2) {
1477
1478 // buffer needs to be 2 x 2 x (interpolation order + 1), with one factor 2
1479 // for x and y, and the other for the number of dimensions.
1480 _interpolationBuffer.resize(4 * intOrder + 4);
1481
1482 // 2-dimensional interpolation
1483 auto const& realX = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx1]);
1484 auto const& realY = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx2]);
1485 double xval = realX.getVal() ;
1486 double yval = realY.getVal() ;
1487
1488 RooAbsBinning const& binningY = realY.getBinning();
1489
1490 int ybinC = binningY.binNumber(yval) ;
1491 int ybinLo = ybinC-intOrder/2 - ((yval<binningY.binCenter(ybinC))?1:0) ;
1492 int ybinM = binningY.numBins() ;
1493
1494 auto idxMultY = _idxMult[varInfo.realVarIdx2];
1496
1497 // Use a class-member buffer to avoid repeated heap allocations.
1498 double * yarr = _interpolationBuffer.data() + 2 * intOrder + 2; // add offset to skip part reserved for other dim
1499 double * xarr = yarr + intOrder + 1;
1500 for (int i=ybinLo ; i<=intOrder+ybinLo ; i++) {
1501 int ibin ;
1502 if (i>=0 && i<ybinM) {
1503 // In range
1504 ibin = i ;
1505 xarr[i-ybinLo] = binningY.binCenter(ibin) ;
1506 } else if (i>=ybinM) {
1507 // Overflow: mirror
1508 ibin = 2*ybinM-i-1 ;
1509 xarr[i-ybinLo] = 2*binningY.highBound()-binningY.binCenter(ibin) ;
1510 } else {
1511 // Underflow: mirror
1512 ibin = -i -1;
1513 xarr[i-ybinLo] = 2*binningY.lowBound()-binningY.binCenter(ibin) ;
1514 }
1517 }
1518
1519 if (gDebug>7) {
1520 std::cout << "RooDataHist interpolating data is" << std::endl ;
1521 std::cout << "xarr = " ;
1522 for (int q=0; q<=intOrder ; q++) std::cout << xarr[q] << " " ;
1523 std::cout << " yarr = " ;
1524 for (int q=0; q<=intOrder ; q++) std::cout << yarr[q] << " " ;
1525 std::cout << std::endl ;
1526 }
1528
1529 } else {
1530
1531 // Higher dimensional scenarios not yet implemented
1532 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") interpolation in "
1533 << varInfo.nRealVars << " dimensions not yet implemented" << std::endl ;
1535
1536 }
1537
1538 return wInt ;
1539}
1540
1541
1543 if (!_errLo || !_errHi) {
1544 initArray(_errLo, _arrSize, -1.);
1545 initArray(_errHi, _arrSize, -1.);
1547 }
1548}
1549
1550
1551////////////////////////////////////////////////////////////////////////////////
1552/// Return the asymmetric errors on the current weight.
1553/// \note see weightError(ErrorType) const for symmetric error.
1554/// \param[out] lo Low error.
1555/// \param[out] hi High error.
1556/// \param[in] etype Type of error to compute. May throw if not supported.
1557/// Supported errors are
1558/// - `Poisson` Default. Asymmetric Poisson errors (68% CL).
1559/// - `SumW2` The square root of the sum of weights. (Symmetric).
1560/// - `None` Return zero.
1561void RooDataHist::weightError(double& lo, double& hi, ErrorType etype) const
1562{
1563 checkInit() ;
1564
1565 switch (etype) {
1566
1567 case Auto:
1568 throw std::invalid_argument("RooDataHist::weightError(" + std::string(GetName()) + ") error type Auto not allowed here");
1569 break ;
1570
1571 case Expected:
1572 throw std::invalid_argument("RooDataHist::weightError(" + std::string(GetName()) + ") error type Expected not allowed here");
1573 break ;
1574
1575 case Poisson:
1576 if (_errLo && _errLo[_curIndex] >= 0.0) {
1577 // Weight is preset or precalculated
1578 lo = _errLo[_curIndex];
1579 hi = _errHi[_curIndex];
1580 return ;
1581 }
1582
1583 // We didn't track asymmetric errors so far, so now we need to allocate
1585
1586 // Calculate poisson errors
1587 double ym;
1588 double yp;
1589 RooHistError::instance().getPoissonInterval(Int_t(weight()+0.5),ym,yp,1) ;
1592 lo = _errLo[_curIndex];
1593 hi = _errHi[_curIndex];
1594 return ;
1595
1596 case SumW2:
1597 lo = std::sqrt(weightSquared(_curIndex));
1598 hi = lo;
1599 return ;
1600
1601 case None:
1602 lo = 0 ;
1603 hi = 0 ;
1604 return ;
1605 }
1606}
1607
1608
1609// wve adjust for variable bin sizes
1610
1611////////////////////////////////////////////////////////////////////////////////
1612/// Perform boundary safe 'intOrder'-th interpolation of weights in dimension 'dim'
1613/// at current value 'xval'
1614
1615/// \param[in] iDim Index of the histogram dimension along which to interpolate.
1616/// \param[in] xval Value of histogram variable at dimension `iDim` for which
1617/// we want to interpolate the histogram weight.
1618/// \param[in] centralIdx Index of the bin that the point at which we
1619/// interpolate the histogram weight falls into
1620/// (can be obtained with `RooDataHist::calcTreeIndex`).
1621/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1622/// used for the interpolation.
1623/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1624/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1625/// underflow bins are assumed to have weight zero and
1626/// overflow bins have weight one. Otherwise, the
1627/// histogram is mirrored at the boundaries for the
1628/// interpolation.
1630{
1631 auto const& binning = static_cast<RooRealVar&>(*_vars[iDim]).getBinning();
1632
1633 // Fill workspace arrays spanning interpolation area
1634 int fbinC = binning.binNumber(xval) ;
1635 int fbinLo = fbinC-intOrder/2 - ((xval<binning.binCenter(fbinC))?1:0) ;
1636 int fbinM = binning.numBins() ;
1637
1638 auto idxMult = _idxMult[iDim];
1639 auto offsetIdx = centralIdx - idxMult * fbinC;
1640
1641 // Use a class-member buffer to avoid repeated heap allocations.
1642 double * yarr = _interpolationBuffer.data();
1643 double * xarr = yarr + intOrder + 1;
1644
1645 for (int i=fbinLo ; i<=intOrder+fbinLo ; i++) {
1646 int ibin ;
1647 if (i>=0 && i<fbinM) {
1648 // In range
1649 ibin = i ;
1650 xarr[i-fbinLo] = binning.binCenter(ibin) ;
1651 auto idx = offsetIdx + idxMult * ibin;
1652 yarr[i - fbinLo] = _wgt[idx];
1653 if (correctForBinSize) yarr[i-fbinLo] /= _binv[idx] ;
1654 } else if (i>=fbinM) {
1655 // Overflow: mirror
1656 ibin = 2*fbinM-i-1 ;
1657 if (cdfBoundaries) {
1658 xarr[i-fbinLo] = binning.highBound()+1e-10*(i-fbinM+1) ;
1659 yarr[i-fbinLo] = 1.0 ;
1660 } else {
1661 auto idx = offsetIdx + idxMult * ibin;
1662 xarr[i-fbinLo] = 2*binning.highBound()-binning.binCenter(ibin) ;
1663 yarr[i - fbinLo] = _wgt[idx];
1665 yarr[i - fbinLo] /= _binv[idx];
1666 }
1667 } else {
1668 // Underflow: mirror
1669 ibin = -i - 1 ;
1670 if (cdfBoundaries) {
1671 xarr[i-fbinLo] = binning.lowBound()-ibin*(1e-10) ;
1672 yarr[i-fbinLo] = 0.0 ;
1673 } else {
1674 auto idx = offsetIdx + idxMult * ibin;
1675 xarr[i-fbinLo] = 2*binning.lowBound()-binning.binCenter(ibin) ;
1676 yarr[i - fbinLo] = _wgt[idx];
1678 yarr[i - fbinLo] /= _binv[idx];
1679 }
1680 }
1681 }
1683}
1684
1685
1686
1687
1688////////////////////////////////////////////////////////////////////////////////
1689/// Increment the bin content of the bin enclosing the given coordinates.
1690///
1691/// \param[in] row Coordinates of the bin.
1692/// \param[in] wgt Increment by this weight.
1693/// \param[in] sumw2 Optionally, track the sum of squared weights. If a value > 0 or
1694/// a weight != 1. is passed for the first time, a vector for the squared weights will be allocated.
1695void RooDataHist::add(const RooArgSet& row, double wgt, double sumw2)
1696{
1697 checkInit() ;
1698
1699 if ((sumw2 > 0. || wgt != 1.) && !_sumw2) {
1700 // Receiving a weighted entry. SumW2 != sumw from now on.
1701 _sumw2 = new double[_arrSize];
1702 std::copy(_wgt, _wgt+_arrSize, _sumw2);
1703
1705 }
1706
1707 const auto idx = calcTreeIndex(row, false);
1708
1709 _wgt[idx] += wgt ;
1710 if (_sumw2) _sumw2[idx] += (sumw2 > 0 ? sumw2 : wgt*wgt);
1711
1712 _cache_sum_valid = false;
1713}
1714
1715
1716
1717////////////////////////////////////////////////////////////////////////////////
1718/// Set a bin content.
1719/// \param[in] row Coordinates of the bin to be set.
1720/// \param[in] wgt New bin content.
1721/// \param[in] wgtErrLo Low error of the bin content.
1722/// \param[in] wgtErrHi High error of the bin content.
1723void RooDataHist::set(const RooArgSet& row, double wgt, double wgtErrLo, double wgtErrHi)
1724{
1725 checkInit() ;
1726
1728
1729 const auto idx = calcTreeIndex(row, false);
1730
1731 _wgt[idx] = wgt ;
1732 _errLo[idx] = wgtErrLo ;
1733 _errHi[idx] = wgtErrHi ;
1734
1735 _cache_sum_valid = false;
1736}
1737
1738
1739
1740////////////////////////////////////////////////////////////////////////////////
1741/// Set bin content of bin that was last loaded with get(std::size_t).
1742/// \param[in] binNumber Optional bin number to set. If empty, currently active bin is set.
1743/// \param[in] wgt New bin content.
1744/// \param[in] wgtErr Error of the new bin content. If the weight need not have an error, use 0. or a negative number.
1745void RooDataHist::set(std::size_t binNumber, double wgt, double wgtErr) {
1746 checkInit() ;
1747
1748 if (wgtErr > 0. && !_sumw2) {
1749 // Receiving a weighted entry. Need to track sumw2 from now on:
1751
1753 }
1754
1755 _wgt[binNumber] = wgt ;
1756 if (_errLo) _errLo[binNumber] = wgtErr;
1757 if (_errHi) _errHi[binNumber] = wgtErr;
1758 if (_sumw2) _sumw2[binNumber] = wgtErr*wgtErr;
1759
1761}
1762
1763
1764////////////////////////////////////////////////////////////////////////////////
1765/// Set bin content of bin that was last loaded with get(std::size_t).
1766/// \param[in] wgt New bin content.
1767/// \param[in] wgtErr Optional error of the bin content.
1768void RooDataHist::set(double wgt, double wgtErr) {
1769 if (_curIndex == std::numeric_limits<std::size_t>::max()) {
1770 _curIndex = calcTreeIndex(_vars, true) ;
1771 }
1772
1774}
1775
1776
1777////////////////////////////////////////////////////////////////////////////////
1778/// Set a bin content.
1779/// \param[in] row Coordinates to compute the bin from.
1780/// \param[in] wgt New bin content.
1781/// \param[in] wgtErr Optional error of the bin content.
1782void RooDataHist::set(const RooArgSet& row, double wgt, double wgtErr) {
1783 set(calcTreeIndex(row, false), wgt, wgtErr);
1784}
1785
1786
1787
1788////////////////////////////////////////////////////////////////////////////////
1789/// Add all data points contained in 'dset' to this data set with given weight.
1790/// Optional cut string expression selects the data points to be added and can
1791/// reference any variable contained in this data set
1792
1793void RooDataHist::add(const RooAbsData& dset, const char* cut, double wgt)
1794{
1795 RooFormulaVar cutVar("select",cut,*dset.get()) ;
1796 add(dset,&cutVar,wgt) ;
1797}
1798
1799
1800
1801////////////////////////////////////////////////////////////////////////////////
1802/// Add all data points contained in 'dset' to this data set with given weight.
1803/// Optional RooFormulaVar pointer selects the data points to be added.
1804
1806{
1807 checkInit() ;
1808
1809 RooFormulaVar* cloneVar = nullptr;
1810 std::unique_ptr<RooArgSet> tmp;
1811 if (cutVar) {
1812 // Deep clone cutVar and attach clone to this dataset
1813 tmp = std::make_unique<RooArgSet>();
1814 if(RooArgSet(*cutVar).snapshot(*tmp)) {
1815 coutE(DataHandling) << "RooDataHist::add(" << GetName() << ") Couldn't deep-clone cut variable, abort," << std::endl ;
1816 return ;
1817 }
1818
1819 cloneVar = static_cast<RooFormulaVar*>(tmp->find(*cutVar)) ;
1820 cloneVar->attachDataSet(dset) ;
1821 }
1822
1823
1824 Int_t i ;
1825 for (i=0 ; i<dset.numEntries() ; i++) {
1826 const RooArgSet* row = dset.get(i) ;
1827 if (!cloneVar || cloneVar->getVal()) {
1828 add(*row,wgt*dset.weight(), wgt*wgt*dset.weightSquared()) ;
1829 }
1830 }
1831
1833}
1834
1835
1836
1837////////////////////////////////////////////////////////////////////////////////
1838/// Return the sum of the weights of all bins in the histogram.
1839///
1840/// \param[in] correctForBinSize Multiply the sum of weights in each bin
1841/// with the N-dimensional bin volume, making the return value
1842/// the integral over the function represented by this histogram.
1843/// \param[in] inverseBinCor Divide by the N-dimensional bin volume.
1845{
1846 checkInit() ;
1847
1848 // Check if result was cached
1850 if (_cache_sum_valid == static_cast<Int_t>(cache_code)) {
1851 return _cache_sum ;
1852 }
1853
1855 for (Int_t i=0; i < _arrSize; i++) {
1856 const double theBinVolume = correctForBinSize ? (inverseBinCor ? 1/_binv[i] : _binv[i]) : 1.0 ;
1857 kahanSum += _wgt[i] * theBinVolume;
1858 }
1859
1860 // Store result in cache
1862 _cache_sum = kahanSum.Sum();
1863
1864 return kahanSum.Sum();
1865}
1866
1867
1868
1869////////////////////////////////////////////////////////////////////////////////
1870/// Return the sum of the weights of a multi-dimensional slice of the histogram
1871/// by summing only over the dimensions specified in sumSet.
1872///
1873/// The coordinates of all other dimensions are fixed to those given in sliceSet
1874///
1875/// If correctForBinSize is specified, the sum of weights
1876/// is multiplied by the M-dimensional bin volume, (M = N(sumSet)),
1877/// making the return value the integral over the function
1878/// represented by this histogram
1879
1881{
1882 checkInit() ;
1883
1885 varSave.addClone(_vars) ;
1886
1888 sliceOnlySet.remove(sumSet,true,true) ;
1889
1891 std::vector<double> const * pbinv = nullptr;
1892
1895 } else if(correctForBinSize && !inverseBinCor) {
1897 }
1898
1899 // Calculate mask and reference plot bins for non-iterating variables
1900 std::vector<bool> mask(_vars.size());
1901 std::vector<int> refBin(_vars.size());
1902
1903 for (unsigned int i = 0; i < _vars.size(); ++i) {
1904 const RooAbsArg* arg = _vars[i];
1905 const RooAbsLValue* argLv = _lvvars[i]; // Same as above, but cross-cast
1906
1907 if (sumSet.find(*arg)) {
1908 mask[i] = false ;
1909 } else {
1910 mask[i] = true ;
1911 refBin[i] = argLv->getBin();
1912 }
1913 }
1914
1915 // Loop over entire data set, skipping masked entries
1917 for (Int_t ibin=0; ibin < _arrSize; ++ibin) {
1918
1919 std::size_t tmpibin = ibin;
1920 bool skip(false) ;
1921
1922 // Check if this bin belongs in selected slice
1923 for (unsigned int ivar = 0; !skip && ivar < _vars.size(); ++ivar) {
1924 const Int_t idx = tmpibin / _idxMult[ivar] ;
1925 tmpibin -= idx*_idxMult[ivar] ;
1926 if (mask[ivar] && idx!=refBin[ivar])
1927 skip = true ;
1928 }
1929
1930 if (!skip) {
1931 const double theBinVolume = correctForBinSize ? (inverseBinCor ? 1/(*pbinv)[ibin] : (*pbinv)[ibin] ) : 1.0 ;
1933 }
1934 }
1935
1937
1938 return total.Sum();
1939}
1940
1941////////////////////////////////////////////////////////////////////////////////
1942/// Return the sum of the weights of a multi-dimensional slice of the histogram
1943/// by summing only over the dimensions specified in sumSet.
1944///
1945/// The coordinates of all other dimensions are fixed to those given in sliceSet
1946///
1947/// If correctForBinSize is specified, the sum of weights
1948/// is multiplied by the M-dimensional bin volume, (M = N(sumSet)),
1949/// or the fraction of it that falls inside the range rangeName,
1950/// making the return value the integral over the function
1951/// represented by this histogram.
1952///
1953/// If correctForBinSize is not specified, the weights are multiplied by the
1954/// fraction of the bin volume that falls inside the range, i.e. a factor of
1955/// binVolumeInRange/totalBinVolume.
1956
1959 const std::map<const RooAbsArg*, std::pair<double, double> >& ranges,
1960 std::function<double(int)> getBinScale)
1961{
1962 checkInit();
1965 varSave.addClone(_vars);
1966 {
1968 sliceOnlySet.remove(sumSet, true, true);
1970 }
1971
1972 // Calculate mask and reference plot bins for non-iterating variables,
1973 // and get ranges for iterating variables
1974 std::vector<bool> mask(_vars.size());
1975 std::vector<int> refBin(_vars.size());
1976 std::vector<double> rangeLo(_vars.size(), -std::numeric_limits<double>::infinity());
1977 std::vector<double> rangeHi(_vars.size(), +std::numeric_limits<double>::infinity());
1978
1979 for (std::size_t i = 0; i < _vars.size(); ++i) {
1980 const RooAbsArg* arg = _vars[i];
1981 const RooAbsLValue* argLV = _lvvars[i]; // Same object as above, but cross cast
1982
1983 RooAbsArg* sumsetv = sumSet.find(*arg);
1984 RooAbsArg* slicesetv = sliceSet.find(*arg);
1985 mask[i] = !sumsetv;
1986 if (mask[i]) {
1987 assert(argLV);
1988 refBin[i] = argLV->getBin();
1989 }
1990
1991 auto it = ranges.find(sumsetv ? sumsetv : slicesetv);
1992 if (ranges.end() != it) {
1993 rangeLo[i] = it->second.first;
1994 rangeHi[i] = it->second.second;
1995 }
1996 }
1997
1998 // Loop over entire data set, skipping masked entries
2000 for (Int_t ibin = 0; ibin < _arrSize; ++ibin) {
2001 // Check if this bin belongs in selected slice
2002 bool skip{false};
2003 for (int ivar = 0, tmp = ibin; !skip && ivar < int(_vars.size()); ++ivar) {
2004 const Int_t idx = tmp / _idxMult[ivar];
2005 tmp -= idx*_idxMult[ivar];
2006 if (mask[ivar] && idx!=refBin[ivar]) skip = true;
2007 }
2008
2009 if (skip) continue;
2010
2011 // Work out bin volume
2012 // It's not necessary to figure out the bin volume for the slice-only set explicitly here.
2013 // We need to loop over the sumSet anyway to get the partial bin containment correction,
2014 // so we can get the slice-only set volume later by dividing _binv[ibin] / binVolumeSumSetFull.
2015 double binVolumeSumSetFull = 1.;
2016 double binVolumeSumSetInRange = 1.;
2017 for (Int_t ivar = 0, tmp = ibin; ivar < (int)_vars.size(); ++ivar) {
2018 const Int_t idx = tmp / _idxMult[ivar];
2019 tmp -= idx*_idxMult[ivar];
2020
2021 // If the current variable is not in the sumSet, it should not be considered for the bin volume
2022 const auto arg = _vars[ivar];
2023 if (!sumSet.find(*arg)) {
2024 continue;
2025 }
2026
2027 if (_binbounds[ivar].empty()) continue;
2028 const double binLo = _binbounds[ivar][2 * idx];
2029 const double binHi = _binbounds[ivar][2 * idx + 1];
2030 if (binHi < rangeLo[ivar] || binLo > rangeHi[ivar]) {
2031 // bin is outside of allowed range - effective bin volume is zero
2033 break;
2034 }
2035
2037 binVolumeSumSetInRange *= std::min(rangeHi[ivar], binHi) - std::max(rangeLo[ivar], binLo);
2038 }
2040 if (0. == corrPartial) continue;
2042 total += getBinScale(ibin)*(_wgt[ibin] * corr * corrPartial);
2043 }
2044
2046
2047 return total.Sum();
2048}
2049
2050
2051
2052////////////////////////////////////////////////////////////////////////////////
2053/// Fill the transient cache with partial bin volumes with up-to-date
2054/// values for the partial volume specified by observables 'dimSet'
2055
2056const std::vector<double>& RooDataHist::calculatePartialBinVolume(const RooArgSet& dimSet) const
2057{
2058 // The code bitset has all bits set to one whose position corresponds to arguments in dimSet.
2059 // It is used as the key for the bin volume caching hash map.
2060 int code{0};
2061 {
2062 int i{0} ;
2063 for (auto const& v : _vars) {
2064 code += ((dimSet.find(*v) ? 1 : 0) << i) ;
2065 ++i;
2066 }
2067 }
2068
2069 auto& pbinv = _pbinvCache[code];
2070 if(!pbinv.empty()) {
2071 return pbinv;
2072 }
2073 pbinv.resize(_arrSize);
2074
2075 // Calculate plot bins of components from master index
2076 std::vector<bool> selDim(_vars.size());
2077 for (std::size_t i = 0; i < selDim.size(); ++i) {
2078 selDim[i] = (code >> i) & 1 ;
2079 }
2080
2081 // Recalculate partial bin volume cache
2082 for (Int_t ibin=0; ibin < _arrSize ;ibin++) {
2083 Int_t idx(0);
2084 Int_t tmp(ibin);
2085 double theBinVolume(1) ;
2086 for (unsigned int j=0; j < _lvvars.size(); ++j) {
2087 const RooAbsLValue* arg = _lvvars[j];
2088 assert(arg);
2089
2090 idx = tmp / _idxMult[j];
2091 tmp -= idx*_idxMult[j];
2092 if (selDim[j]) {
2093 theBinVolume *= arg->getBinWidth(idx) ;
2094 }
2095 }
2097 }
2098
2099 return pbinv;
2100}
2101
2102
2103////////////////////////////////////////////////////////////////////////////////
2104/// Sum the weights of all bins.
2108
2109
2110
2111////////////////////////////////////////////////////////////////////////////////
2112/// Return the sum of weights in all entries matching cutSpec (if specified)
2113/// and in named range cutRange (if specified)
2114/// Return the
2115
2116double RooDataHist::sumEntries(const char* cutSpec, const char* cutRange) const
2117{
2118 checkInit() ;
2119
2120 if (cutSpec==nullptr && cutRange==nullptr) {
2121 return sumEntries();
2122 } else {
2123
2124 // Setup RooFormulaVar for cutSpec if it is present
2125 std::unique_ptr<RooFormula> select;
2126 if (cutSpec) {
2127 select = std::make_unique<RooFormula>("select",cutSpec,*get());
2128 }
2129
2130 // Otherwise sum the weights in the event
2131 ROOT::Math::KahanSum<> kahanSum;
2132 for (Int_t i=0; i < _arrSize; i++) {
2133 get(i) ;
2134 if ((select && select->eval() == 0.) || (cutRange && !_vars.allInRange(cutRange)))
2135 continue;
2136
2137 kahanSum += weight(i);
2138 }
2139
2140 return kahanSum.Sum();
2141 }
2142}
2143
2144
2145
2146////////////////////////////////////////////////////////////////////////////////
2147/// Reset all bin weights to zero
2148
2150{
2151 // WVE DO NOT CALL RooTreeData::reset() for binned
2152 // datasets as this will delete the bin definitions
2153
2154 std::fill(_wgt, _wgt + _arrSize, 0.);
2155 delete[] _errLo; _errLo = nullptr;
2156 delete[] _errHi; _errHi = nullptr;
2157 delete[] _sumw2; _sumw2 = nullptr;
2158
2160
2161 _cache_sum_valid = false;
2162}
2163
2164
2165
2166////////////////////////////////////////////////////////////////////////////////
2167/// Load bin `binNumber`, and return an argset with the coordinates of the bin centre.
2168/// \note The argset is owned by this data hist, and this function has a side effect, because
2169/// it alters the currently active bin.
2170const RooArgSet* RooDataHist::get(Int_t binNumber) const
2171{
2172 checkInit() ;
2173 _curIndex = binNumber;
2174
2175 return RooAbsData::get(_curIndex);
2176}
2177
2178
2179
2180////////////////////////////////////////////////////////////////////////////////
2181/// Return a RooArgSet with whose coordinates denote the bin centre of the bin
2182/// enclosing the point in `coord`.
2183/// \note The argset is owned by this data hist, and this function has a side effect, because
2184/// it alters the currently active bin.
2186 return get(calcTreeIndex(coord, false));
2187}
2188
2189
2190
2191////////////////////////////////////////////////////////////////////////////////
2192/// Return the volume of the bin enclosing coordinates 'coord'.
2194 checkInit() ;
2195 return _binv[calcTreeIndex(coord, false)] ;
2196}
2197
2198
2199////////////////////////////////////////////////////////////////////////////////
2200/// Create an iterator over all bins in a slice defined by the subset of observables
2201/// listed in sliceArg. The position of the slice is given by otherArgs
2202
2204{
2205 // Update to current position
2207 _curIndex = calcTreeIndex(_vars, true);
2208
2210 if (!intArg) {
2211 coutE(InputArguments) << "RooDataHist::sliceIterator() variable " << sliceArg.GetName() << " is not part of this RooDataHist" << std::endl ;
2212 return nullptr ;
2213 }
2214 return new RooDataHistSliceIter(*this,*intArg) ;
2215}
2216
2217
2218////////////////////////////////////////////////////////////////////////////////
2219/// Change the name of the RooDataHist
2220
2222{
2223 if (_dir) _dir->GetList()->Remove(this);
2224 // We need to use the function from RooAbsData, because it already overrides TNamed::SetName
2226 if (_dir) _dir->GetList()->Add(this);
2227}
2228
2229
2230////////////////////////////////////////////////////////////////////////////////
2231/// Change the title of this RooDataHist
2232
2233void RooDataHist::SetNameTitle(const char *name, const char* title)
2234{
2235 SetName(name);
2236 SetTitle(title);
2237}
2238
2239
2240////////////////////////////////////////////////////////////////////////////////
2241/// Print value of the dataset, i.e. the sum of weights contained in the dataset
2242
2243void RooDataHist::printValue(ostream& os) const
2244{
2245 os << numEntries() << " bins (" << sumEntries() << " weights)" ;
2246}
2247
2248
2249
2250
2251////////////////////////////////////////////////////////////////////////////////
2252/// Print argument of dataset, i.e. the observable names
2253
2254void RooDataHist::printArgs(ostream& os) const
2255{
2256 os << "[" ;
2257 bool first(true) ;
2258 for (const auto arg : _vars) {
2259 if (first) {
2260 first=false ;
2261 } else {
2262 os << "," ;
2263 }
2264 os << arg->GetName() ;
2265 }
2266 os << "]" ;
2267}
2268
2269
2270
2271////////////////////////////////////////////////////////////////////////////////
2272/// Returns true if dataset contains entries with a non-integer weight.
2273
2275{
2276 for (Int_t i=0; i < _arrSize; ++i) {
2277 const double wgt = _wgt[i];
2278 double intpart;
2279 if (std::abs(std::modf(wgt, &intpart)) > 1.E-10)
2280 return true;
2281 }
2282
2283 return false;
2284}
2285
2286
2287////////////////////////////////////////////////////////////////////////////////
2288/// Print the details on the dataset contents
2289
2290void RooDataHist::printMultiline(ostream& os, Int_t content, bool verbose, TString indent) const
2291{
2293
2294 os << indent << "Binned Dataset " << GetName() << " (" << GetTitle() << ")" << std::endl ;
2295 os << indent << " Contains " << numEntries() << " bins with a total weight of " << sumEntries() << std::endl;
2296
2297 if (!verbose) {
2298 os << indent << " Observables " << _vars << std::endl ;
2299 } else {
2300 os << indent << " Observables: " ;
2302 }
2303
2304 if(verbose) {
2305 if (!_cachedVars.empty()) {
2306 os << indent << " Caches " << _cachedVars << std::endl ;
2307 }
2308 }
2309}
2310
2311/**
2312 * \brief Prints the contents of the RooDataHist to the specified output stream.
2313 *
2314 * This function iterates through all bins of the histogram and prints the
2315 * coordinates of each bin, along with its weight and statistical error.
2316 * It is designed to be robust, handling empty or invalid datasets,
2317 * and works for histograms of any dimension.
2318 *
2319 * \param os The output stream (e.g., std::cout) to write the contents to.
2320 */
2321void RooDataHist::printContents(std::ostream& os) const
2322{
2323 os << "Contents of RooDataHist \"" << GetName() << "\"" << std::endl;
2324
2325 if (numEntries() == 0) {
2326 os << "(dataset is empty)" << std::endl;
2327 return;
2328 }
2329
2330 for (int i = 0; i < numEntries(); ++i) {
2331 const RooArgSet* obs = get(i); // load i-th bin
2332 os << " Bin " << i << ": ";
2333
2334 bool first = true;
2335 for (const auto* var : *obs) {
2336 if (!first) os << ", ";
2337 first = false;
2338
2339 os << var->GetName() << "=";
2340 if (auto realVar = dynamic_cast<const RooRealVar*>(var)) {
2341 os << realVar->getVal();
2342 } else if (auto catVar = dynamic_cast<const RooCategory*>(var)) {
2343 os << catVar->getCurrentLabel();
2344 } else {
2345 os << "(unsupported type)"; //added as a precaution
2346 }
2347 }
2348
2349 double lo, hi;
2351 os << ", weight=" << weight() << " +/- [" << lo << "," << hi << "]"
2352 << std::endl;
2353 }
2354}
2355
2356
2357////////////////////////////////////////////////////////////////////////////////
2358/// Stream an object of class RooDataHist.
2360 if (R__b.IsReading()) {
2361
2362 UInt_t R__s;
2363 UInt_t R__c;
2364 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2365
2366 if (R__v > 2) {
2367 R__b.ReadClassBuffer(RooDataHist::Class(),this,R__v,R__s,R__c);
2368 R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());
2369 initialize(nullptr, false);
2370 } else {
2371
2372 // Legacy dataset conversion happens here. Legacy RooDataHist inherits from RooTreeData
2373 // which in turn inherits from RooAbsData. Manually stream RooTreeData contents on
2374 // file here and convert it into a RooTreeDataStore which is installed in the
2375 // new-style RooAbsData base class
2376
2377 // --- This is the contents of the streamer code of RooTreeData version 2 ---
2378 UInt_t R__s1;
2379 UInt_t R__c1;
2380 Version_t R__v1 = R__b.ReadVersion(&R__s1, &R__c1); if (R__v1) { }
2381
2383 TTree* X_tree(nullptr) ; R__b >> X_tree;
2384 RooArgSet X_truth ; X_truth.Streamer(R__b);
2386 R__b.CheckByteCount(R__s1, R__c1, TClass::GetClass("RooTreeData"));
2387 // --- End of RooTreeData-v1 streamer
2388
2389 // Construct RooTreeDataStore from X_tree and complete initialization of new-style RooAbsData
2390 _dstore = std::make_unique<RooTreeDataStore>(X_tree,_vars);
2391 _dstore->SetName(GetName()) ;
2392 _dstore->SetTitle(GetTitle()) ;
2393 _dstore->checkInit() ;
2394
2396 R__b >> _arrSize;
2397 delete [] _wgt;
2398 _wgt = new double[_arrSize];
2399 R__b.ReadFastArray(_wgt,_arrSize);
2400 delete [] _errLo;
2401 _errLo = new double[_arrSize];
2402 R__b.ReadFastArray(_errLo,_arrSize);
2403 delete [] _errHi;
2404 _errHi = new double[_arrSize];
2405 R__b.ReadFastArray(_errHi,_arrSize);
2406 delete [] _sumw2;
2407 _sumw2 = new double[_arrSize];
2408 R__b.ReadFastArray(_sumw2,_arrSize);
2409 delete [] _binv;
2410 _binv = new double[_arrSize];
2412 tmpSet.Streamer(R__b);
2413 double tmp;
2414 R__b >> tmp; //_curWeight;
2415 R__b >> tmp; //_curWgtErrLo;
2416 R__b >> tmp; //_curWgtErrHi;
2417 R__b >> tmp; //_curSumW2;
2418 R__b >> tmp; //_curVolume;
2419 R__b >> _curIndex;
2420 R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());
2421 }
2422
2423 } else {
2424
2425 R__b.WriteClassBuffer(RooDataHist::Class(),this);
2426 }
2427}
2428
2429
2430////////////////////////////////////////////////////////////////////////////////
2431/// Return event weights of all events in range [first, first+len).
2432/// If cacheValidEntries() has been called, out-of-range events will have a weight of 0.
2433std::span<const double> RooDataHist::getWeightBatch(std::size_t first, std::size_t len, bool sumW2 /*=false*/) const {
2434 return {(sumW2 && _sumw2 ? _sumw2 : _wgt) + first, len};
2435}
2436
2437
2438////////////////////////////////////////////////////////////////////////////////
2439/// Hand over pointers to our weight arrays to the data store implementation.
2441 _dstore->setExternalWeightArray(_wgt, _errLo, _errHi, _sumw2);
2442}
2443
2444
2445////////////////////////////////////////////////////////////////////////////////
2446/// Return reference to VarInfo struct with cached histogram variable
2447/// information that is frequently used for histogram weights retrieval.
2448///
2449/// If the `_varInfo` struct was not initialized yet, it will be initialized in
2450/// this function.
2452
2453 if(_varInfo.initialized) return _varInfo;
2454
2455 auto& info = _varInfo;
2456
2457 {
2458 // count the number of real vars and get their indices
2459 info.nRealVars = 0;
2460 size_t iVar = 0;
2461 for (const auto real : _vars) {
2462 if (dynamic_cast<RooRealVar*>(real)) {
2463 if(info.nRealVars == 0) info.realVarIdx1 = iVar;
2464 if(info.nRealVars == 1) info.realVarIdx2 = iVar;
2465 ++info.nRealVars;
2466 }
2467 ++iVar;
2468 }
2469 }
2470
2471 {
2472 // assert that the variables are either real values or categories
2473 for (unsigned int i=0; i < _vars.size(); ++i) {
2474 if (_lvbins[i].get()) {
2475 assert(dynamic_cast<const RooAbsReal*>(_vars[i]));
2476 } else {
2477 assert(dynamic_cast<const RooAbsCategoryLValue*>(_vars[i]));
2478 }
2479 }
2480 }
2481
2482 info.initialized = true;
2483
2484 return info;
2485}
#define e(i)
Definition RSha256.hxx:103
#define coutI(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
short Version_t
Class version identifier (short)
Definition RtypesCore.h:79
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
static unsigned int total
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 mask
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 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
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:157
float xmin
#define hi
float * q
float ymin
Int_t gDebug
Global variable setting the debug level. Set to 0 to disable, increase it in steps of 1 to increase t...
Definition TROOT.cxx:627
The Kahan summation is a compensated summation algorithm, which significantly reduces numerical error...
Definition Util.h:141
static KahanSum< T, N > Accumulate(Iterator begin, Iterator end, T initialValue=T{})
Iterate over a range and return an instance of a KahanSum.
Definition Util.h:230
const_iterator begin() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
void attachDataSet(const RooAbsData &set)
Replace server nodes with names matching the dataset variable names with those data set variables,...
Abstract base class for RooRealVar binning definitions.
int binNumber(double x) const
Returns the bin number corresponding to the value x.
virtual void binNumbers(double const *x, int *bins, std::size_t n, int coef=1) const =0
Compute the bin indices for multiple values of x.
virtual bool isUniform() const
virtual double highBound() const =0
virtual double lowBound() const =0
virtual std::string translateBinNumber(RooFit::Experimental::CodegenContext &ctx, RooAbsArg const &var, int coef) const
virtual RooAbsBinning * clone(const char *name=nullptr) const =0
Abstract base class for objects that represent a discrete value that can be set from the outside,...
bool hasLabel(const std::string &label) const
Check if a state with name label exists.
Abstract container object that can hold multiple RooAbsArg objects.
RooAbsCollection & assignValueOnly(const RooAbsCollection &other, bool forceIfSizeOne=false)
Sets the value of any argument in our set that also appears in the other set.
bool allInRange(const char *rangeSpec) const
Return true if all contained object report to have their value inside the specified range.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:56
virtual const RooArgSet * get() const
Definition RooAbsData.h:100
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Interface for detailed printing of object.
void SetName(const char *name) override
Set the name of the TNamed.
void setGlobalObservables(RooArgSet const &globalObservables)
Sets the global observables stored in this data.
void checkInit() const
static StorageType defaultStorageType
Definition RooAbsData.h:298
std::unique_ptr< RooAbsDataStore > _dstore
Data storage implementation.
Definition RooAbsData.h:358
virtual void fill()
RooArgSet _vars
Dimensions of this data set.
Definition RooAbsData.h:355
RooArgSet _cachedVars
! External variables cached with this data set
Definition RooAbsData.h:356
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
void Streamer(TBuffer &) override
Stream an object of class RooAbsData.
virtual RooPlot * plotOnImpl(RooPlot *frame, PlotOpt o) const
Create and fill a histogram of the frame's variable and append it to the frame.
Abstract base class for objects that are lvalues, i.e.
virtual double getBinWidth(Int_t i, const char *rangeName=nullptr) const =0
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
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:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
Object to represent discrete states.
Definition RooCategory.h:28
bool defineType(const std::string &label)
Define a state with given name.
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Configurable parser for RooCmdArg named arguments.
void defineMutex(const char *head, Args_t &&... tail)
Define arguments where any pair is mutually exclusive.
bool process(const RooCmdArg &arg)
Process given RooCmdArg.
double getDouble(const char *name, double defaultValue=0.0) const
Return double property registered with name 'name'.
void defineDependency(const char *refArgName, const char *neededArgName)
Define that processing argument name refArgName requires processing of argument named neededArgName t...
bool defineDouble(const char *name, const char *argName, int doubleNum, double defValue=0.0)
Define double property name 'name' mapped to double in slot 'doubleNum' in RooCmdArg with name argNam...
RooArgSet * getSet(const char *name, RooArgSet *set=nullptr) const
Return RooArgSet property registered with name 'name'.
bool defineSet(const char *name, const char *argName, int setNum, const RooArgSet *set=nullptr)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
bool ok(bool verbose) const
Return true of parsing was successful.
bool defineObject(const char *name, const char *argName, int setNum, const TObject *obj=nullptr, bool isArray=false)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
const char * getString(const char *name, const char *defaultValue="", bool convEmptyToNull=false) const
Return string property registered with name 'name'.
bool defineString(const char *name, const char *argName, int stringNum, const char *defValue="", bool appendMode=false)
Define double property name 'name' mapped to double in slot 'stringNum' in RooCmdArg with name argNam...
const RooLinkedList & getObjectList(const char *name) const
Return list of objects registered with name 'name'.
bool defineInt(const char *name, const char *argName, int intNum, int defValue=0)
Define integer property name 'name' mapped to integer in slot 'intNum' in RooCmdArg with name argName...
int getInt(const char *name, int defaultValue=0) const
Return integer property registered with name 'name'.
TObject * getObject(const char *name, TObject *obj=nullptr) const
Return TObject property registered with name 'name'.
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
std::span< const double > getWeightBatch(std::size_t first, std::size_t len, bool sumW2=false) const override
Return event weights of all events in range [first, first+len).
void interpolateQuadratic(double *output, std::span< const double > xVals, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of interpolateDim for boundary safe quadratic interpolation of one dimensional h...
double sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
void weights(double *output, std::span< double const > xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of RooDataHist::weight() for one dimensional histograms with up to one dimension...
Int_t _cache_sum_valid
! Is cache sum valid? Needs to be Int_t instead of CacheSumState_t for subclasses.
void printContents(std::ostream &os=std::cout) const override
Print the contents of the dataset to the specified output stream.
double interpolateDim(int iDim, double xval, size_t centralIdx, int intOrder, bool correctForBinSize, bool cdfBoundaries)
Perform boundary safe 'intOrder'-th interpolation of weights in dimension 'dim' at current value 'xva...
double weightSquared() const override
Return squared weight of last bin that was requested with get().
friend class RooDataHistSliceIter
void importTH1(const RooArgList &vars, const TH1 &histo, double initWgt, bool doDensityCorrection)
Import data from given TH1/2/3 into this RooDataHist.
static TClass * Class()
TClass * IsA() const override
void SetNameTitle(const char *name, const char *title) override
Change the title of this RooDataHist.
double _cache_sum
! Cache for sum of entries ;
void initialize(const char *binningName=nullptr, bool fillTree=true)
Initialization procedure: allocate weights array, calculate multipliers needed for N-space to 1-dim a...
VarInfo _varInfo
!
std::string declWeightArrayForCodeSquash(RooFit::Experimental::CodegenContext &ctx, bool correctForBinSize) const
Int_t getIndex(const RooAbsCollection &coord, bool fast=false) const
Calculate bin number of the given coordinates.
void add(const RooArgSet &row, double wgt=1.0) override
Add wgt to the bin content enclosed by the coordinates passed in row.
Definition RooDataHist.h:72
const std::vector< double > & calculatePartialBinVolume(const RooArgSet &dimSet) const
Fill the transient cache with partial bin volumes with up-to-date values for the partial volume speci...
static std::unique_ptr< RooAbsDataStore > makeDefaultDataStore(RooStringView name, RooStringView title, RooArgSet const &vars)
double weightInterpolated(const RooArgSet &bin, int intOrder, bool correctForBinSize, bool cdfBoundaries)
Return the weight at given coordinates with interpolation.
std::unordered_map< int, std::vector< double > > _pbinvCache
! Cache for arrays of partial bin volumes
void checkBinBounds() const
void initializeAsymErrArrays() const
void set(std::size_t binNumber, double weight, double wgtErr)
Set bin content of bin that was last loaded with get(std::size_t).
void weightError(double &lo, double &hi, ErrorType etype=Poisson) const override
Return the asymmetric errors on the current weight.
double * _errHi
[_arrSize] High-side error on weight array
void importTH1Set(const RooArgList &vars, RooCategory &indexCat, std::map< std::string, TH1 * > hmap, double initWgt, bool doDensityCorrection)
Import data from given set of TH1/2/3 into this RooDataHist.
void adjustBinning(const RooArgList &vars, const TH1 &href, Int_t *offset=nullptr)
Adjust binning specification on first and optionally second and third observable to binning in given ...
double * _binv
[_arrSize] Bin volume array
RooDataHist()
Default constructor.
ULong64_t _curIndex
Current index.
std::string calculateTreeIndexForCodeSquash(RooFit::Experimental::CodegenContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
double weightFast(const RooArgSet &bin, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A faster version of RooDataHist::weight that assumes the passed arguments are aligned with the histog...
double weight() const override
Return weight of last bin that was requested with get().
std::vector< std::vector< double > > _binbounds
! list of bin bounds per dimension
void printArgs(std::ostream &os) const override
Print argument of dataset, i.e. the observable names.
void importDHistSet(const RooArgList &vars, RooCategory &indexCat, std::map< std::string, RooDataHist * > dmap, double initWgt)
Import data from given set of TH1/2/3 into this RooDataHist.
void _adjustBinning(RooRealVar &theirVar, const TAxis &axis, RooRealVar *ourVar, Int_t *offset)
Helper doing the actual work of adjustBinning().
void printMultiline(std::ostream &os, Int_t content, bool verbose=false, TString indent="") const override
Print the details on the dataset contents.
double * _sumw2
[_arrSize] Sum of weights^2
TIterator * sliceIterator(RooAbsArg &sliceArg, const RooArgSet &otherArgs)
Create an iterator over all bins in a slice defined by the subset of observables listed in sliceArg.
Int_t calcTreeIndex() const
Legacy overload to calculate the tree index from the current value of _vars.
~RooDataHist() override
Destructor.
bool isNonPoissonWeighted() const override
Returns true if dataset contains entries with a non-integer weight.
std::vector< RooAbsLValue * > _lvvars
! List of observables casted as RooAbsLValue
void SetName(const char *name) override
Change the name of the RooDataHist.
std::vector< std::unique_ptr< const RooAbsBinning > > _lvbins
! List of used binnings associated with lvalues
void Streamer(TBuffer &) override
Stream an object of class RooDataHist.
std::vector< double > _interpolationBuffer
! Buffer to contain values used for weight interpolation
std::vector< Int_t > _idxMult
void registerWeightArraysToDataStore() const
Hand over pointers to our weight arrays to the data store implementation.
void reset() override
Reset all bin weights to zero.
double * _errLo
[_arrSize] Low-side error on weight array
double * _wgt
[_arrSize] Weight array
RooPlot * plotOnImpl(RooPlot *frame, PlotOpt o) const override
Back end function to plotting functionality.
void printValue(std::ostream &os) const override
Print value of the dataset, i.e. the sum of weights contained in the dataset.
VarInfo const & getVarInfo()
Return reference to VarInfo struct with cached histogram variable information that is frequently used...
std::unique_ptr< RooAbsData > reduceEng(const RooArgSet &varSubset, const RooFormulaVar *cutVar, const char *cutRange=nullptr, std::size_t nStart=0, std::size_t nStop=std::numeric_limits< std::size_t >::max()) const override
Implementation of RooAbsData virtual method that drives the RooAbsData::reduce() methods.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
void interpolateLinear(double *output, std::span< const double > xVals, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of interpolateDim for boundary safe linear interpolation of one dimensional hist...
double binVolume() const
Return volume of current bin.
double sumEntries() const override
Sum the weights of all bins.
Utility base class for RooFit objects that are to be attached to ROOT directories.
Definition RooDirItem.h:22
virtual void Streamer(TBuffer &)
void removeFromDir(TObject *obj)
Remove object from directory it was added to.
TDirectory * _dir
! Associated directory
Definition RooDirItem.h:33
A class to maintain the context for squashing of RooFit models into code.
std::string buildArg(RooAbsCollection const &x, std::string const &arrayType="double")
Function to save a RooListProxy as an array in the squashed code.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
static const RooHistError & instance()
Return a reference to a singleton object that is created the first time this method is called.
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
static double interpolate(double yArr[], Int_t nOrder, double x)
Definition RooMath.cxx:78
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
RooAbsRealLValue * getPlotVar() const
Definition RooPlot.h:137
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
The RooStringView is a wrapper around a C-style string that can also be constructed from a std::strin...
Class to manage histogram axis.
Definition TAxis.h:32
const TArrayD * GetXbins() const
Definition TAxis.h:138
Double_t GetXmax() const
Definition TAxis.h:142
virtual Int_t FindFixBin(Double_t x) const
Find bin number corresponding to abscissa x
Definition TAxis.cxx:422
Double_t GetXmin() const
Definition TAxis.h:141
Int_t GetNbins() const
Definition TAxis.h:127
Buffer base class used for serializing objects.
Definition TBuffer.h:43
static TClass * GetClass(const char *name, Bool_t load=kTRUE, Bool_t silent=kFALSE)
Static method returning pointer to TClass of the specified class name.
Definition TClass.cxx:2979
virtual TList * GetList() const
Definition TDirectory.h:223
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
TAxis * GetZaxis()
Definition TH1.h:573
virtual Int_t GetNbinsY() const
Definition TH1.h:542
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:9176
virtual Int_t GetNbinsZ() const
Definition TH1.h:543
virtual Int_t GetDimension() const
Definition TH1.h:527
TAxis * GetXaxis()
Definition TH1.h:571
virtual Int_t GetNbinsX() const
Definition TH1.h:541
TAxis * GetYaxis()
Definition TH1.h:572
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5122
Iterator abstract base class.
Definition TIterator.h:30
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
Mother of all ROOT objects.
Definition TObject.h:42
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:224
Basic string class.
Definition TString.h:138
A TTree represents a columnar dataset.
Definition TTree.h:89
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5
std::vector< std::string > Split(std::string_view str, std::string_view delims, bool skipEmpty=false)
Splits a string at each character in delims.
RooAbsBinning * bins
Definition RooAbsData.h:313
Structure to cache information on the histogram variable that is frequently used for histogram weight...
TLine l
Definition textangle.C:4
static void output()