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 if(auto dHist = dynamic_cast<RooDataHist*>(*hiter)) {
331 dmap[token] = dHist;
332 }
333 if(auto hHist = dynamic_cast<TH1*>(*hiter)) {
334 hmap[token] = hHist;
335 }
336 ++hiter;
337 }
338 if(!dmap.empty() && !hmap.empty()) {
339 std::stringstream errorMsgStream;
340 errorMsgStream << "RooDataHist::ctor(" << GetName() << ") ERROR: you can't import mix of TH1 and RooDataHist";
341 const std::string errorMsg = errorMsgStream.str();
342 coutE(InputArguments) << errorMsg << std::endl;
343 throw std::invalid_argument(errorMsg);
344 }
345 if (!dmap.empty()) {
346 importDHistSet(vars,*indexCat,dmap,initWgt);
347 }
348 if (!hmap.empty()) {
349 importTH1Set(vars,*indexCat,hmap,initWgt,false);
350 }
351
352
353 } else {
354
355 // Initialize empty
356 initialize();
357 }
358
361
362}
363
364
365
366
367////////////////////////////////////////////////////////////////////////////////
368/// Import data from given TH1/2/3 into this RooDataHist
369
370void RooDataHist::importTH1(const RooArgList& vars, const TH1& histo, double wgt, bool doDensityCorrection)
371{
372 // Adjust binning of internal observables to match that of input THx
373 Int_t offset[3]{0, 0, 0};
374 adjustBinning(vars, histo, offset) ;
375
376 // Initialize internal data structure
377 initialize();
378
379 // Define x,y,z as 1st, 2nd and 3rd observable
380 RooRealVar* xvar = static_cast<RooRealVar*>(_vars.find(vars.at(0)->GetName())) ;
381 RooRealVar* yvar = static_cast<RooRealVar*>(vars.at(1) ? _vars.find(vars.at(1)->GetName()) : nullptr ) ;
382 RooRealVar* zvar = static_cast<RooRealVar*>(vars.at(2) ? _vars.find(vars.at(2)->GetName()) : nullptr ) ;
383
384 // Transfer contents
385 Int_t xmin(0);
386 Int_t ymin(0);
387 Int_t zmin(0);
389 xmin = offset[0] ;
390 if (yvar) {
391 vset.add(*yvar) ;
392 ymin = offset[1] ;
393 }
394 if (zvar) {
395 vset.add(*zvar) ;
396 zmin = offset[2] ;
397 }
398
399 Int_t iX(0);
400 Int_t iY(0);
401 Int_t iz(0);
402 for (iX=0 ; iX < xvar->getBins() ; iX++) {
403 xvar->setBin(iX) ;
404 if (yvar) {
405 for (iY=0 ; iY < yvar->getBins() ; iY++) {
406 yvar->setBin(iY) ;
407 if (zvar) {
408 for (iz=0 ; iz < zvar->getBins() ; iz++) {
409 zvar->setBin(iz) ;
410 double bv = doDensityCorrection ? binVolume(vset) : 1;
411 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)) ;
412 }
413 } else {
414 double bv = doDensityCorrection ? binVolume(vset) : 1;
415 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)) ;
416 }
417 }
418 } else {
419 double bv = doDensityCorrection ? binVolume(vset) : 1 ;
420 add(vset,bv*histo.GetBinContent(iX+1+xmin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin)*wgt,2)) ;
421 }
422 }
423
424}
425
426namespace {
427bool checkConsistentAxes(const TH1* first, const TH1* second) {
428 return first->GetDimension() == second->GetDimension()
429 && first->GetNbinsX() == second->GetNbinsX()
430 && first->GetNbinsY() == second->GetNbinsY()
431 && first->GetNbinsZ() == second->GetNbinsZ()
432 && first->GetXaxis()->GetXmin() == second->GetXaxis()->GetXmin()
433 && first->GetXaxis()->GetXmax() == second->GetXaxis()->GetXmax()
434 && (first->GetNbinsY() == 1 || (first->GetYaxis()->GetXmin() == second->GetYaxis()->GetXmin()
435 && first->GetYaxis()->GetXmax() == second->GetYaxis()->GetXmax() ) )
436 && (first->GetNbinsZ() == 1 || (first->GetZaxis()->GetXmin() == second->GetZaxis()->GetXmin()
437 && first->GetZaxis()->GetXmax() == second->GetZaxis()->GetXmax() ) );
438}
439
441
442 // Relative tolerance for bin boundary comparison
443 constexpr double tolerance = 1e-6;
444
445 auto const& vars1 = *h1.get();
446 auto const& vars2 = *h2.get();
447
448 // Check if number of variables and names is consistent
449 if(!vars1.hasSameLayout(vars2)) {
450 return false;
451 }
452
453 for(std::size_t iVar = 0; iVar < vars1.size(); ++iVar) {
454 auto * var1 = dynamic_cast<RooRealVar*>(vars1[iVar]);
455 auto * var2 = dynamic_cast<RooRealVar*>(vars2[iVar]);
456
457 // Check if variables are consistently real-valued
458 if((!var1 && var2) || (var1 && !var2)) return false;
459
460 // Not a real-valued variable
461 if(!var1) continue;
462
463 // Now check the binning
464 auto const& bng1 = var1->getBinning();
465 auto const& bng2 = var2->getBinning();
466
467 // Compare bin numbers
468 if(bng1.numBins() != bng2.numBins()) return false;
469
470 std::size_t nBins = bng1.numBins();
471
472 // Compare bin boundaries
473 for(std::size_t iBin = 0; iBin < nBins; ++iBin) {
474 double v1 = bng1.binLow(iBin);
475 double v2 = bng2.binLow(iBin);
476 if(std::abs((v1 - v2) / v1) > tolerance) return false;
477 }
478 double v1 = bng1.binHigh(nBins - 1);
479 double v2 = bng2.binHigh(nBins - 1);
480 if(std::abs((v1 - v2) / v1) > tolerance) return false;
481 }
482 return true;
483}
484}
485
486
487////////////////////////////////////////////////////////////////////////////////
488/// Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
489/// in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
490/// and the import source
491
492void RooDataHist::importTH1Set(const RooArgList& vars, RooCategory& indexCat, std::map<string,TH1*> hmap, double wgt, bool doDensityCorrection)
493{
494 RooCategory* icat = static_cast<RooCategory*>(_vars.find(indexCat.GetName())) ;
495
496 TH1* histo(nullptr) ;
497 bool init(false) ;
498 for (const auto& hiter : hmap) {
499 // Store pointer to first histogram from which binning specification will be taken
500 if (!histo) {
501 histo = hiter.second;
502 } else {
503 if (!checkConsistentAxes(histo, hiter.second)) {
504 coutE(InputArguments) << "Axes of histogram " << hiter.second->GetName() << " are not consistent with first processed "
505 << "histogram " << histo->GetName() << std::endl;
506 throw std::invalid_argument("Axes of inputs for RooDataHist are inconsistent");
507 }
508 }
509 // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
510 if (!indexCat.hasLabel(hiter.first)) {
511 indexCat.defineType(hiter.first) ;
512 coutI(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << ") defining state \"" << hiter.first << "\" in index category " << indexCat.GetName() << std::endl ;
513 }
514 if (!icat->hasLabel(hiter.first)) {
515 icat->defineType(hiter.first) ;
516 }
517 }
518
519 // Check consistency in number of dimensions
520 if (histo && int(vars.size()) != histo->GetDimension()) {
521 coutE(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << "): dimension of input histogram must match "
522 << "number of continuous variables" << std::endl ;
523 throw std::invalid_argument("Inputs histograms for RooDataHist are not compatible with dimensions of variables.");
524 }
525
526 // Copy bins and ranges from THx to dimension observables
527 Int_t offset[3] ;
528 adjustBinning(vars,*histo,offset) ;
529
530 // Initialize internal data structure
531 if (!init) {
532 initialize();
533 init = true;
534 }
535
536 // Define x,y,z as 1st, 2nd and 3rd observable
537 RooRealVar* xvar = static_cast<RooRealVar*>(_vars.find(vars.at(0)->GetName())) ;
538 RooRealVar* yvar = static_cast<RooRealVar*>(vars.at(1) ? _vars.find(vars.at(1)->GetName()) : nullptr ) ;
539 RooRealVar* zvar = static_cast<RooRealVar*>(vars.at(2) ? _vars.find(vars.at(2)->GetName()) : nullptr ) ;
540
541 // Transfer contents
542 Int_t xmin(0);
543 Int_t ymin(0);
544 Int_t zmin(0);
546 double volume = xvar->getMax()-xvar->getMin() ;
547 xmin = offset[0] ;
548 if (yvar) {
549 vset.add(*yvar) ;
550 ymin = offset[1] ;
551 volume *= (yvar->getMax()-yvar->getMin()) ;
552 }
553 if (zvar) {
554 vset.add(*zvar) ;
555 zmin = offset[2] ;
556 volume *= (zvar->getMax()-zvar->getMin()) ;
557 }
558 double avgBV = volume / numEntries() ;
559
560 Int_t ic(0);
561 Int_t iX(0);
562 Int_t iY(0);
563 Int_t iz(0);
564 for (ic=0 ; ic < icat->numBins(nullptr) ; ic++) {
565 icat->setBin(ic) ;
566 histo = hmap[icat->getCurrentLabel()] ;
567 for (iX=0 ; iX < xvar->getBins() ; iX++) {
568 xvar->setBin(iX) ;
569 if (yvar) {
570 for (iY=0 ; iY < yvar->getBins() ; iY++) {
571 yvar->setBin(iY) ;
572 if (zvar) {
573 for (iz=0 ; iz < zvar->getBins() ; iz++) {
574 zvar->setBin(iz) ;
575 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
576 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)) ;
577 }
578 } else {
579 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
580 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)) ;
581 }
582 }
583 } else {
584 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
585 add(vset,bv*histo->GetBinContent(iX+1+xmin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin)*wgt,2)) ;
586 }
587 }
588 }
589
590}
591
592
593
594////////////////////////////////////////////////////////////////////////////////
595/// Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
596/// in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
597/// and the import source
598
599void RooDataHist::importDHistSet(const RooArgList & /*vars*/, RooCategory &indexCat,
600 std::map<std::string, RooDataHist *> dmap, double initWgt)
601{
602 auto *icat = static_cast<RooCategory *>(_vars.find(indexCat.GetName()));
603
604 RooDataHist *dhistForBinning = nullptr;
605
606 for (const auto &diter : dmap) {
607
608 std::string const &label = diter.first;
609 RooDataHist *dhist = diter.second;
610
611 if (!dhistForBinning) {
613 } else {
615 coutE(InputArguments) << "Layout or binning of histogram " << dhist->GetName()
616 << " is not consistent with first processed "
617 << "histogram " << dhistForBinning->GetName() << std::endl;
618 throw std::invalid_argument("Layout or binning of inputs for RooDataHist is inconsistent");
619 }
620 }
621
622 // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
623 if (!indexCat.hasLabel(label)) {
624 indexCat.defineType(label);
625 coutI(InputArguments) << "RooDataHist::importDHistSet(" << GetName() << ") defining state \"" << label
626 << "\" in index category " << indexCat.GetName() << std::endl;
627 }
628 if (!icat->hasLabel(label)) {
629 icat->defineType(label);
630 }
631 }
632
633 // adjust the binning of the created histogram
635 auto *ourVar = dynamic_cast<RooRealVar *>(_vars.find(theirVar->GetName()));
636 if (!theirVar || !ourVar)
637 continue;
638 ourVar->setBinning(theirVar->getBinning());
639 }
640
641 initialize();
642
643 for (const auto &diter : dmap) {
644 std::string const &label = diter.first;
645 RooDataHist *dhist = diter.second;
646
647 icat->setLabel(label.c_str());
648
649 // Transfer contents
650 for (Int_t i = 0; i < dhist->numEntries(); i++) {
651 _vars.assign(*dhist->get(i));
652 add(_vars, dhist->weight() * initWgt, pow(dhist->weightError(SumW2), 2));
653 }
654 }
655}
656
657////////////////////////////////////////////////////////////////////////////////
658/// Helper doing the actual work of adjustBinning().
659
662{
663 const std::string ourVarName(ourVar->GetName() ? ourVar->GetName() : "");
664 const std::string ownName(GetName() ? GetName() : "");
665 // RooRealVar is derived from RooAbsRealLValue which is itself
666 // derived from RooAbsReal and a virtual class RooAbsLValue
667 // supplying setter functions, check if ourVar is indeed derived
668 // as real
669 if (!dynamic_cast<RooAbsReal *>(ourVar)) {
670 coutE(InputArguments) << "RooDataHist::adjustBinning(" << ownName << ") ERROR: dimension " << ourVarName
671 << " must be real\n";
672 throw std::logic_error("Incorrect type object (" + ourVarName +
673 ") passed as argument to RooDataHist::_adjustBinning. Please report this issue.");
674 }
675
676 const double xlo = theirVar.getMin();
677 const double xhi = theirVar.getMax();
678
679 const bool isUniform = !axis.GetXbins()->GetArray();
680 std::unique_ptr<RooAbsBinning> xbins;
681
682 if (!isUniform) {
683 xbins = std::make_unique<RooBinning>(axis.GetNbins(), axis.GetXbins()->GetArray());
684 } else {
685 xbins = std::make_unique<RooUniformBinning>(axis.GetXmin(), axis.GetXmax(), axis.GetNbins());
686 }
687
688 const double tolerance = 1e-6 * xbins->averageBinWidth();
689
690 // Adjust xlo/xhi to nearest boundary
691 const int iBinLo = xbins->binNumber(xlo + tolerance);
692 const int iBinHi = xbins->binNumber(xhi - tolerance);
693 const int nBinsAdj = iBinHi - iBinLo + 1;
694 const double xloAdj = xbins->binLow(iBinLo);
695 const double xhiAdj = xbins->binHigh(iBinHi);
696
697 if (isUniform) {
698 xbins = std::make_unique<RooUniformBinning>(xloAdj, xhiAdj, nBinsAdj);
699 theirVar.setRange(xloAdj, xhiAdj);
700 } else {
701 xbins->setRange(xloAdj, xhiAdj);
702 theirVar.setBinning(*xbins);
703 }
704
705 if (std::abs(xloAdj - xlo) > tolerance || std::abs(xhiAdj - xhi) > tolerance) {
706 coutI(DataHandling) << "RooDataHist::adjustBinning(" << ownName << "): fit range of variable " << ourVarName
707 << " expanded to nearest bin boundaries: [" << xlo << "," << xhi << "] --> [" << xloAdj << ","
708 << xhiAdj << "]"
709 << "\n";
710 }
711
712 ourVar->setBinning(*xbins);
713
714 // The offset is the bin number of the adjusted lower limit of the RooFit
715 // variable in the original TH1 histogram, starting from zero.
716 if (offset) {
717 *offset = axis.FindFixBin(xloAdj + tolerance) - 1;
718 }
719}
720
721////////////////////////////////////////////////////////////////////////////////
722/// Adjust binning specification on first and optionally second and third
723/// observable to binning in given reference TH1. Used by constructors
724/// that import data from an external TH1.
725/// Both the variables in vars and in this RooDataHist are adjusted.
726/// @param vars List with variables that are supposed to have their binning adjusted.
727/// @param href Reference histogram that dictates the binning
728/// @param offset If not nullptr, a possible bin count offset for the axes x,y,z is saved here as Int_t[3]
729
731{
732 auto xvar = static_cast<RooRealVar*>(_vars.find(*vars.at(0)) );
733 _adjustBinning(*static_cast<RooRealVar*>(vars.at(0)), *href.GetXaxis(), xvar, offset ? &offset[0] : nullptr);
734
735 if (vars.at(1)) {
736 auto yvar = static_cast<RooRealVar*>(_vars.find(*vars.at(1)));
737 if (yvar)
738 _adjustBinning(*static_cast<RooRealVar*>(vars.at(1)), *href.GetYaxis(), yvar, offset ? &offset[1] : nullptr);
739 }
740
741 if (vars.at(2)) {
742 auto zvar = static_cast<RooRealVar*>(_vars.find(*vars.at(2)));
743 if (zvar)
744 _adjustBinning(*static_cast<RooRealVar*>(vars.at(2)), *href.GetZaxis(), zvar, offset ? &offset[2] : nullptr);
745 }
746
747}
748
749
750namespace {
751/// Clone external weight arrays, unless the external array is nullptr.
752void cloneArray(double*& ours, const double* theirs, std::size_t n) {
753 if (ours) delete[] ours;
754 ours = nullptr;
755 if (!theirs) return;
756 ours = new double[n];
757 std::copy(theirs, theirs+n, ours);
758}
759
760/// Allocate and initialise an array with desired size and values.
761void initArray(double*& arr, std::size_t n, double val) {
762 if (arr) delete[] arr;
763 arr = nullptr;
764 if (n == 0) return;
765 arr = new double[n];
766 std::fill(arr, arr+n, val);
767}
768}
769
770
771////////////////////////////////////////////////////////////////////////////////
772/// Initialization procedure: allocate weights array, calculate
773/// multipliers needed for N-space to 1-dim array jump table,
774/// and fill the internal tree with all bin center coordinates
775
776void RooDataHist::initialize(const char* binningName, bool fillTree)
777{
778 _lvvars.clear();
779 _lvbins.clear();
780
781 // Fill array of LValue pointers to variables
782 for (unsigned int i = 0; i < _vars.size(); ++i) {
783 if (binningName) {
784 RooRealVar* rrv = dynamic_cast<RooRealVar*>(_vars[i]);
785 if (rrv) {
786 rrv->setBinning(rrv->getBinning(binningName));
787 }
788 }
789
790 auto lvarg = dynamic_cast<RooAbsLValue*>(_vars[i]);
791 assert(lvarg);
792 _lvvars.push_back(lvarg);
793
794 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
795 _lvbins.emplace_back(binning ? binning->clone() : nullptr);
796 }
797
798
799 // Allocate coefficients array
800 _idxMult.resize(_vars.size()) ;
801
802 _arrSize = 1 ;
803 unsigned int n = 0u;
804 for (const auto var : _vars) {
805 auto arg = dynamic_cast<const RooAbsLValue*>(var);
806 assert(arg);
807
808 // Calculate sub-index multipliers for master index
809 for (unsigned int i = 0u; i<n; i++) {
810 _idxMult[i] *= arg->numBins() ;
811 }
812 _idxMult[n++] = 1 ;
813
814 // Calculate dimension of weight array
815 _arrSize *= arg->numBins() ;
816 }
817
818 // Allocate and initialize weight array if necessary
819 if (!_wgt) {
820 initArray(_wgt, _arrSize, 0.);
821 delete[] _errLo; _errLo = nullptr;
822 delete[] _errHi; _errHi = nullptr;
823 delete[] _sumw2; _sumw2 = nullptr;
825
826 // Refill array pointers in data store when reading
827 // from Streamer
828 if (!fillTree) {
830 }
831 }
832
833 if (!fillTree) return ;
834
835 // Fill TTree with bin center coordinates
836 // Calculate plot bins of components from master index
837
838 for (Int_t ibin=0 ; ibin < _arrSize ; ibin++) {
839 Int_t j(0);
840 Int_t idx(0);
841 Int_t tmp(ibin);
842 double theBinVolume(1) ;
843 for (auto arg2 : _lvvars) {
844 idx = tmp / _idxMult[j] ;
845 tmp -= idx*_idxMult[j++] ;
846 arg2->setBin(idx) ;
847 theBinVolume *= arg2->getBinWidth(idx) ;
848 }
850
851 fill() ;
852 }
853
854
855}
856
857
858////////////////////////////////////////////////////////////////////////////////
859
861{
862 if (!_binbounds.empty()) return;
863 for (auto& it : _lvbins) {
864 _binbounds.push_back(std::vector<double>());
865 if (it) {
866 std::vector<double>& bounds = _binbounds.back();
867 bounds.reserve(2 * it->numBins());
868 for (Int_t i = 0; i < it->numBins(); ++i) {
869 bounds.push_back(it->binLow(i));
870 bounds.push_back(it->binHigh(i));
871 }
872 }
873 }
874}
875
876
877////////////////////////////////////////////////////////////////////////////////
878/// Copy constructor
879
881 RooAbsData(other,newname), RooDirItem(), _arrSize(other._arrSize), _idxMult(other._idxMult), _pbinvCache(other._pbinvCache)
882{
883 // Allocate and initialize weight array
884 assert(_arrSize == other._arrSize);
885 cloneArray(_wgt, other._wgt, other._arrSize);
886 cloneArray(_errLo, other._errLo, other._arrSize);
887 cloneArray(_errHi, other._errHi, other._arrSize);
888 cloneArray(_binv, other._binv, other._arrSize);
889 cloneArray(_sumw2, other._sumw2, other._arrSize);
890
891 // Fill array of LValue pointers to variables
892 for (const auto rvarg : _vars) {
893 auto lvarg = dynamic_cast<RooAbsLValue*>(rvarg);
894 assert(lvarg);
895 _lvvars.push_back(lvarg);
896 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
897 _lvbins.emplace_back(binning ? binning->clone() : nullptr) ;
898 }
899
901}
902
903
904////////////////////////////////////////////////////////////////////////////////
905/// Implementation of RooAbsData virtual method that drives the RooAbsData::reduce() methods
906
907std::unique_ptr<RooAbsData> RooDataHist::reduceEng(const RooArgSet& varSubset, const RooFormulaVar* cutVar, const char* cutRange,
908 std::size_t nStart, std::size_t nStop) const
909{
910 checkInit() ;
913 auto rdh = std::make_unique<RooDataHist>(GetName(), GetTitle(), myVarSubset);
914
915 RooFormulaVar* cloneVar = nullptr;
916 std::unique_ptr<RooArgSet> tmp;
917 if (cutVar) {
918 tmp = std::make_unique<RooArgSet>();
919 // Deep clone cutVar and attach clone to this dataset
920 if (RooArgSet(*cutVar).snapshot(*tmp)) {
921 coutE(DataHandling) << "RooDataHist::reduceEng(" << GetName() << ") Couldn't deep-clone cut variable, abort," << std::endl ;
922 return nullptr;
923 }
924 cloneVar = static_cast<RooFormulaVar*>(tmp->find(*cutVar));
925 cloneVar->attachDataSet(*this) ;
926 }
927
928 double lo;
929 double hi;
930 const std::size_t nevt = nStop < static_cast<std::size_t>(numEntries()) ? nStop : static_cast<std::size_t>(numEntries());
931 for (auto i=nStart; i<nevt ; i++) {
932 const RooArgSet* row = get(i) ;
933
934 bool doSelect(true) ;
935 if (cutRange) {
936 for (const auto arg : *row) {
937 if (!arg->inRange(cutRange)) {
938 doSelect = false ;
939 break ;
940 }
941 }
942 }
943 if (!doSelect) continue ;
944
945 if (!cloneVar || cloneVar->getVal()) {
946 weightError(lo,hi,SumW2) ;
947 rdh->add(*row,weight(),lo*lo) ;
948 }
949 }
950
951 return rdh ;
952}
953
954
955
956////////////////////////////////////////////////////////////////////////////////
957/// Destructor
958
960{
961 delete[] _wgt;
962 delete[] _errLo;
963 delete[] _errHi;
964 delete[] _sumw2;
965 delete[] _binv;
966
967 removeFromDir(this) ;
969}
970
971
972
973
974////////////////////////////////////////////////////////////////////////////////
975/// Calculate bin number of the given coordinates. If only a subset of the internal
976/// coordinates are passed, the missing coordinates are taken at their current value.
977/// \param[in] coord Variables that are representing the coordinates.
978/// \param[in] fast If the variables in `coord` and the ones of the data hist have the
979/// same size and layout, `fast` can be set to skip checking that all variables are
980/// present in `coord`.
982 checkInit() ;
983 return calcTreeIndex(coord, fast);
984}
985
987 bool correctForBinSize) const
988{
989 std::vector<double> vals(_arrSize);
990 for (std::size_t i = 0; i < vals.size(); ++i) {
991 vals[i] = correctForBinSize ? _wgt[i] / _binv[i] : _wgt[i];
992 }
993 return ctx.buildArg(vals);
994}
995
997 const RooAbsCollection &coords, bool reverse) const
998{
999 assert(coords.size() == _vars.size());
1000
1001 std::string code;
1002 int idxMult = 1;
1003
1004 for (std::size_t i = 0; i < _vars.size(); ++i) {
1005
1006 std::size_t iVar = reverse ? _vars.size() - 1 - i : i;
1007 const RooAbsArg *internalVar = _vars[iVar];
1008 const RooAbsArg *theVar = coords[iVar];
1009
1010 const RooAbsBinning *binning = _lvbins[iVar].get();
1011 if (!binning) {
1012 coutE(InputArguments) << "RooHistPdf::weight(" << GetName()
1013 << ") ERROR: Code Squashing currently does not support category values." << std::endl;
1014 return "";
1015 }
1016
1017 code += " + " + binning->translateBinNumber(ctx, *theVar, idxMult);
1018
1019 // Use RooAbsLValue here because it also generalized to categories, which
1020 // is useful in the future. dynamic_cast because it's a cross-cast.
1021 idxMult *= dynamic_cast<RooAbsLValue const *>(internalVar)->numBins();
1022 }
1023
1024 return "(" + code + ")";
1025}
1026
1027////////////////////////////////////////////////////////////////////////////////
1028/// Calculate the bin index corresponding to the coordinates passed as argument.
1029/// \param[in] coords Coordinates. If `fast == false`, these can be partial.
1030/// \param[in] fast Promise that the coordinates in `coords` have the same order
1031/// as the internal coordinates. In this case, values are looked up only by index.
1032std::size_t RooDataHist::calcTreeIndex(const RooAbsCollection& coords, bool fast) const
1033{
1034 // With fast, caller promises that layout of `coords` is identical to our internal `vars`.
1035 // Previously, this was verified with an assert in debug mode like this:
1036 //
1037 // assert(!fast || coords.hasSameLayout(_vars));
1038 //
1039 // However, there are usecases where the externally provided `coords` have
1040 // different names than the internal variables, even though they correspond
1041 // to each other. For example, if the observables in the computation graph
1042 // are renamed with `redirectServers`. Hence, we can't do a meaningful assert
1043 // here.
1044
1045 if (&_vars == &coords)
1046 fast = true;
1047
1048 std::size_t masterIdx = 0;
1049
1050 for (unsigned int i=0; i < _vars.size(); ++i) {
1051 const RooAbsArg* internalVar = _vars[i];
1052 const RooAbsBinning* binning = _lvbins[i].get();
1053
1054 // Find the variable that we need values from.
1055 // That's either the variable directly from the external coordinates
1056 // or we find the external one that has the same name as "internalVar".
1057 const RooAbsArg* theVar = fast ? coords[i] : coords.find(*internalVar);
1058 if (!theVar) {
1059 // Variable is not in external coordinates. Use current internal value.
1061 }
1062 // If fast is on, users promise that the sets have the same layout:
1063 //
1064 // assert(!fast || strcmp(internalVar->GetName(), theVar->GetName()) == 0);
1065 //
1066 // This assert is commented out for the same reasons that applied to the
1067 // other assert explained above.
1068
1069 if (binning) {
1070 assert(dynamic_cast<const RooAbsReal*>(theVar));
1071 const double val = static_cast<const RooAbsReal*>(theVar)->getVal();
1072 masterIdx += _idxMult[i] * binning->binNumber(val);
1073 } else {
1074 // We are a category. No binning.
1075 assert(dynamic_cast<const RooAbsCategoryLValue*>(theVar));
1076 auto cat = static_cast<const RooAbsCategoryLValue*>(theVar);
1077 masterIdx += _idxMult[i] * cat->getBin(static_cast<const char*>(nullptr));
1078 }
1079 }
1080
1081 return masterIdx ;
1082}
1083
1084
1085////////////////////////////////////////////////////////////////////////////////
1086/// Back end function to plotting functionality. Plot RooDataHist on given
1087/// frame in mode specified by plot options 'o'. The main purpose of
1088/// this function is to match the specified binning on 'o' to the
1089/// internal binning of the plot observable in this RooDataHist.
1090/// \note see RooAbsData::plotOnImpl() for plotting options.
1092{
1093 checkInit() ;
1094 if (o.bins) return RooAbsData::plotOnImpl(frame,o) ;
1095
1096 if(!frame) {
1097 coutE(InputArguments) << ClassName() << "::" << GetName() << ":plotOn: frame is null" << std::endl;
1098 return nullptr;
1099 }
1100 auto var= static_cast<RooAbsRealLValue*>(frame->getPlotVar());
1101 if(!var) {
1102 coutE(InputArguments) << ClassName() << "::" << GetName()
1103 << ":plotOn: frame does not specify a plot variable" << std::endl;
1104 return nullptr;
1105 }
1106
1107 auto dataVar = static_cast<RooRealVar*>(_vars.find(*var));
1108 if (!dataVar) {
1109 coutE(InputArguments) << ClassName() << "::" << GetName()
1110 << ":plotOn: dataset doesn't contain plot frame variable" << std::endl;
1111 return nullptr;
1112 }
1113
1114 o.bins = &dataVar->getBinning() ;
1115 return RooAbsData::plotOnImpl(frame,o) ;
1116}
1117
1118
1119////////////////////////////////////////////////////////////////////////////////
1120/// A vectorized version of interpolateDim for boundary safe quadratic
1121/// interpolation of one dimensional histograms.
1122///
1123/// \param[out] output An array of interpolated weights corresponding to the
1124/// values in xVals.
1125/// \param[in] xVals An array of event coordinates for which the weights should be
1126/// calculated.
1127/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1128/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1129/// Underflow bins are assumed to have weight zero and
1130/// overflow bins have weight one. Otherwise, the
1131/// histogram is mirrored at the boundaries for the
1132/// interpolation.
1133
1134void RooDataHist::interpolateQuadratic(double* output, std::span<const double> xVals,
1136{
1137 const std::size_t nBins = numEntries();
1138 const std::size_t nEvents = xVals.size();
1139
1140 RooAbsBinning const& binning = *_lvbins[0];
1141 // Reuse the output buffer for bin indices and zero-initialize it
1142 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1143 std::fill(binIndices, binIndices + nEvents, 0);
1144 binning.binNumbers(xVals.data(), binIndices, nEvents);
1145
1146 // Extend coordinates and weights with one extra point before the first bin
1147 // and one extra point after the last bin. This means the original histogram
1148 // bins span elements 1 to nBins in coordsExt and weightsExt
1149 std::vector<double> coordsExt(nBins+3);
1150 double* binCoords = coordsExt.data() + 2;
1151 binCoords[0] = binning.lowBound() + 0.5*_binv[0];
1152 for (std::size_t binIdx = 1; binIdx < nBins ; ++binIdx) {
1153 if (binning.isUniform()) {
1154 double binWidth = _binv[0];
1155 binCoords[binIdx] = binIdx*binWidth + binCoords[0];
1156 }
1157 else {
1158 double binCentDiff = 0.5*_binv[binIdx-1] + 0.5*_binv[binIdx];
1160 }
1161 }
1162
1163 std::vector<double> weightsExt(nBins+3);
1164 // Fill weights for bins that are inside histogram boundaries
1165 for (std::size_t binIdx = 0; binIdx < nBins; ++binIdx) {
1167 }
1168
1169 if (cdfBoundaries) {
1170 coordsExt[0] = - 1e-10 + binning.lowBound();
1171 weightsExt[0] = 0.;
1172
1173 coordsExt[1] = binning.lowBound();
1174 weightsExt[1] = 0.;
1175
1176 coordsExt[nBins+2] = binning.highBound();
1177 weightsExt[nBins+2] = 1.;
1178 }
1179 else {
1180 // Mirror first two bins and last bin
1181 coordsExt[0] = binCoords[1] - 2*_binv[0] - _binv[1];
1182 weightsExt[0] = weightsExt[3];
1183
1184 coordsExt[1] = binCoords[0] - _binv[0];
1185 weightsExt[1] = weightsExt[2];
1186
1187 coordsExt[nBins+2] = binCoords[nBins-1] + _binv[nBins-1];
1188 weightsExt[nBins+2] = weightsExt[nBins+1];
1189 }
1190
1191 // We use the current bin center and two bin centers on the left for
1192 // interpolation if xVal is to the left of the current bin center
1193 for (std::size_t i = 0; i < nEvents ; ++i) {
1194 double xVal = xVals[i];
1195 std::size_t binIdx = binIndices[i] + 2;
1196
1197 // If xVal is to the right of the current bin center, shift all bin
1198 // coordinates one step to the right and use that for the interpolation
1199 if (xVal > coordsExt[binIdx]) {
1200 binIdx += 1;
1201 }
1202
1203 double x1 = coordsExt[binIdx-2];
1204 double y1 = weightsExt[binIdx-2];
1205
1206 double x2 = coordsExt[binIdx-1];
1207 double y2 = weightsExt[binIdx-1];
1208
1209 double x3 = coordsExt[binIdx];
1210 double y3 = weightsExt[binIdx];
1211
1212 // Evaluate a few repeated factors
1213 double quotient = (x3-x1) / (x2-x1);
1214 double x1Sqrd = x1*x1;
1215 double x3Sqrd = x3*x3;
1216 // Solve coefficients in system of three quadratic equations!
1217 double secondCoeff = (y3 - y1 - (y2-y1) * quotient) / (x3Sqrd - x1Sqrd - (x2*x2 - x1Sqrd) * quotient);
1218 double firstCoeff = (y3 - y1 - secondCoeff*(x3Sqrd - x1Sqrd)) / (x3-x1);
1220 // Get the interpolated weight using the equation of a second degree polynomial
1222 }
1223}
1224
1225
1226////////////////////////////////////////////////////////////////////////////////
1227/// A vectorized version of interpolateDim for boundary safe linear
1228/// interpolation of one dimensional histograms.
1229///
1230/// \param[out] output An array of interpolated weights corresponding to the
1231/// values in xVals.
1232/// \param[in] xVals An array of event coordinates for which the weights should be
1233/// calculated.
1234/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1235/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1236/// Underflow bins are assumed to have weight zero and
1237/// overflow bins have weight one. Otherwise, the
1238/// histogram is mirrored at the boundaries for the
1239/// interpolation.
1240
1241void RooDataHist::interpolateLinear(double* output, std::span<const double> xVals,
1243{
1244 const std::size_t nBins = numEntries();
1245 const std::size_t nEvents = xVals.size();
1246
1247 RooAbsBinning const& binning = *_lvbins[0];
1248 // Reuse the output buffer for bin indices and zero-initialize it
1249 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1250 std::fill(binIndices, binIndices + nEvents, 0);
1251 binning.binNumbers(xVals.data(), binIndices, nEvents);
1252
1253 // Extend coordinates and weights with one extra point before the first bin
1254 // and one extra point after the last bin. This means the original histogram
1255 // bins span elements 1 to nBins in coordsExt and weightsExt
1256 std::vector<double> coordsExt(nBins+2);
1257 double* binCoords = coordsExt.data() + 1;
1258 binCoords[0] = binning.lowBound() + 0.5*_binv[0];
1259 for (std::size_t binIdx = 1; binIdx < nBins ; ++binIdx) {
1260 if (binning.isUniform()) {
1261 double binWidth = _binv[0];
1262 binCoords[binIdx] = binIdx*binWidth + binCoords[0];
1263 }
1264 else {
1265 double binCentDiff = 0.5*_binv[binIdx-1] + 0.5*_binv[binIdx];
1267 }
1268 }
1269
1270 std::vector<double> weightsExt(nBins+2);
1271 // Fill weights for bins that are inside histogram boundaries
1272 for (std::size_t binIdx = 0; binIdx < nBins; ++binIdx) {
1274 }
1275
1276 // Fill weights for bins that are outside histogram boundaries
1277 if (cdfBoundaries) {
1278 coordsExt[0] = binning.lowBound();
1279 weightsExt[0] = 0.;
1280 coordsExt[nBins+1] = binning.highBound();
1281 weightsExt[nBins+1] = 1.;
1282 }
1283 else {
1284 // Mirror first and last bins
1285 coordsExt[0] = binCoords[0] - _binv[0];
1286 weightsExt[0] = weightsExt[1];
1287 coordsExt[nBins+1] = binCoords[nBins-1] + _binv[nBins-1];
1288 weightsExt[nBins+1] = weightsExt[nBins];
1289 }
1290
1291 // Interpolate between current bin center and one bin center to the left
1292 // if xVal is to the left of the current bin center
1293 for (std::size_t i = 0; i < nEvents ; ++i) {
1294 double xVal = xVals[i];
1295 std::size_t binIdx = binIndices[i] + 1;
1296
1297 // If xVal is to the right of the current bin center, interpolate between
1298 // current bin center and one bin center to the right instead
1299 if (xVal > coordsExt[binIdx]) { binIdx += 1; }
1300
1301 double x1 = coordsExt[binIdx-1];
1302 double y1 = weightsExt[binIdx-1];
1303 double x2 = coordsExt[binIdx];
1304 double y2 = weightsExt[binIdx];
1305
1306 // Find coefficients by solving a system of two linear equations
1307 double firstCoeff = (y2-y1) / (x2-x1);
1308 double zerothCoeff = y1 - firstCoeff * x1;
1309 // Get the interpolated weight using the equation of a straight line
1311 }
1312}
1313
1314
1315////////////////////////////////////////////////////////////////////////////////
1316/// A vectorized version of RooDataHist::weight() for one dimensional histograms
1317/// with up to one dimensional interpolation.
1318/// \param[out] output An array of weights corresponding the values in xVals.
1319/// \param[in] xVals An array of coordinates for which the weights should be
1320/// calculated.
1321/// \param[in] intOrder Interpolation order; 0th and 1st order are supported.
1322/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1323/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1324/// Underflow bins are assumed to have weight zero and
1325/// overflow bins have weight one. Otherwise, the
1326/// histogram is mirrored at the boundaries for the
1327/// interpolation.
1328
1329void RooDataHist::weights(double* output, std::span<double const> xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
1330{
1331 auto const nEvents = xVals.size();
1332
1333 if (intOrder == 0) {
1334 RooAbsBinning const& binning = *_lvbins[0];
1335
1336 // Reuse the output buffer for bin indices and zero-initialize it
1337 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1338 std::fill(binIndices, binIndices + nEvents, 0);
1339 binning.binNumbers(xVals.data(), binIndices, nEvents);
1340
1341 for (std::size_t i=0; i < nEvents; ++i) {
1342 auto binIdx = binIndices[i];
1344 }
1345 }
1346 else if (intOrder == 1) {
1348 }
1349 else if (intOrder == 2) {
1351 }
1352 else {
1353 // Higher dimensional scenarios not yet implemented
1354 coutE(InputArguments) << "RooDataHist::weights(" << GetName() << ") interpolation in "
1355 << intOrder << " dimensions not yet implemented" << std::endl ;
1356 // Fall back to 1st order interpolation
1358 }
1359}
1360
1361
1362////////////////////////////////////////////////////////////////////////////////
1363/// A faster version of RooDataHist::weight that assumes the passed arguments
1364/// are aligned with the histogram variables.
1365/// \param[in] bin Coordinates for which the weight should be calculated.
1366/// Has to be aligned with the internal histogram variables.
1367/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1368/// used for the interpolation. If zero, the bare weight for
1369/// the bin enclosing the coordinatesis returned.
1370/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1371/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1372/// underflow bins are assumed to have weight zero and
1373/// overflow bins have weight one. Otherwise, the
1374/// histogram is mirrored at the boundaries for the
1375/// interpolation.
1376
1378{
1379 checkInit() ;
1380
1381 // Handle illegal intOrder values
1382 if (intOrder<0) {
1383 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << std::endl ;
1384 return 0 ;
1385 }
1386
1387 // Handle no-interpolation case
1388 if (intOrder==0) {
1389 const auto idx = calcTreeIndex(bin, true);
1390 return correctForBinSize ? _wgt[idx] / _binv[idx] : _wgt[idx];
1391 }
1392
1393 // Handle all interpolation cases
1395}
1396
1397
1398////////////////////////////////////////////////////////////////////////////////
1399/// Return the weight at given coordinates with optional interpolation.
1400/// \param[in] bin Coordinates for which the weight should be calculated.
1401/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1402/// used for the interpolation. If zero, the bare weight for
1403/// the bin enclosing the coordinatesis returned.
1404/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1405/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1406/// underflow bins are assumed to have weight zero and
1407/// overflow bins have weight one. Otherwise, the
1408/// histogram is mirrored at the boundaries for the
1409/// interpolation.
1410/// \param[in] oneSafe Ignored.
1411
1412double RooDataHist::weight(const RooArgSet& bin, Int_t intOrder, bool correctForBinSize, bool cdfBoundaries, bool /*oneSafe*/)
1413{
1414 checkInit() ;
1415
1416 // Handle illegal intOrder values
1417 if (intOrder<0) {
1418 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << std::endl ;
1419 return 0 ;
1420 }
1421
1422 // Handle no-interpolation case
1423 if (intOrder==0) {
1424 const auto idx = calcTreeIndex(bin, false);
1425 return correctForBinSize ? _wgt[idx] / _binv[idx] : _wgt[idx];
1426 }
1427
1428 // Handle all interpolation cases
1429 _vars.assignValueOnly(bin) ;
1430
1432}
1433
1434
1435////////////////////////////////////////////////////////////////////////////////
1436/// Return the weight at given coordinates with interpolation.
1437/// \param[in] bin Coordinates for which the weight should be calculated.
1438/// Has to be aligned with the internal histogram variables.
1439/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1440/// used for the interpolation.
1441/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1442/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1443/// underflow bins are assumed to have weight zero and
1444/// overflow bins have weight one. Otherwise, the
1445/// histogram is mirrored at the boundaries for the
1446/// interpolation.
1447
1449 VarInfo const& varInfo = getVarInfo();
1450
1451 const auto centralIdx = calcTreeIndex(bin, true);
1452
1453 double wInt{0} ;
1454 if (varInfo.nRealVars == 1) {
1455
1456 // buffer needs to be 2 x (interpolation order + 1), with the factor 2 for x and y.
1457 _interpolationBuffer.resize(2 * intOrder + 2);
1458
1459 // 1-dimensional interpolation
1460 auto const& realX = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx1]);
1462
1463 } else if (varInfo.nRealVars == 2) {
1464
1465 // buffer needs to be 2 x 2 x (interpolation order + 1), with one factor 2
1466 // for x and y, and the other for the number of dimensions.
1467 _interpolationBuffer.resize(4 * intOrder + 4);
1468
1469 // 2-dimensional interpolation
1470 auto const& realX = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx1]);
1471 auto const& realY = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx2]);
1472 double xval = realX.getVal() ;
1473 double yval = realY.getVal() ;
1474
1475 RooAbsBinning const& binningY = realY.getBinning();
1476
1477 int ybinC = binningY.binNumber(yval) ;
1478 int ybinLo = ybinC-intOrder/2 - ((yval<binningY.binCenter(ybinC))?1:0) ;
1479 int ybinM = binningY.numBins() ;
1480
1481 auto idxMultY = _idxMult[varInfo.realVarIdx2];
1483
1484 // Use a class-member buffer to avoid repeated heap allocations.
1485 double * yarr = _interpolationBuffer.data() + 2 * intOrder + 2; // add offset to skip part reserved for other dim
1486 double * xarr = yarr + intOrder + 1;
1487 for (int i=ybinLo ; i<=intOrder+ybinLo ; i++) {
1488 int ibin ;
1489 if (i>=0 && i<ybinM) {
1490 // In range
1491 ibin = i ;
1492 xarr[i-ybinLo] = binningY.binCenter(ibin) ;
1493 } else if (i>=ybinM) {
1494 // Overflow: mirror
1495 ibin = 2*ybinM-i-1 ;
1496 xarr[i-ybinLo] = 2*binningY.highBound()-binningY.binCenter(ibin) ;
1497 } else {
1498 // Underflow: mirror
1499 ibin = -i -1;
1500 xarr[i-ybinLo] = 2*binningY.lowBound()-binningY.binCenter(ibin) ;
1501 }
1504 }
1505
1506 if (gDebug>7) {
1507 std::cout << "RooDataHist interpolating data is" << std::endl ;
1508 std::cout << "xarr = " ;
1509 for (int q=0; q<=intOrder ; q++) std::cout << xarr[q] << " " ;
1510 std::cout << " yarr = " ;
1511 for (int q=0; q<=intOrder ; q++) std::cout << yarr[q] << " " ;
1512 std::cout << std::endl ;
1513 }
1515
1516 } else {
1517
1518 // Higher dimensional scenarios not yet implemented
1519 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") interpolation in "
1520 << varInfo.nRealVars << " dimensions not yet implemented" << std::endl ;
1522
1523 }
1524
1525 return wInt ;
1526}
1527
1528
1530 if (!_errLo || !_errHi) {
1531 initArray(_errLo, _arrSize, -1.);
1532 initArray(_errHi, _arrSize, -1.);
1534 }
1535}
1536
1537
1538////////////////////////////////////////////////////////////////////////////////
1539/// Return the asymmetric errors on the current weight.
1540/// \note see weightError(ErrorType) const for symmetric error.
1541/// \param[out] lo Low error.
1542/// \param[out] hi High error.
1543/// \param[in] etype Type of error to compute. May throw if not supported.
1544/// Supported errors are
1545/// - `Poisson` Default. Asymmetric Poisson errors (68% CL).
1546/// - `SumW2` The square root of the sum of weights. (Symmetric).
1547/// - `None` Return zero.
1548void RooDataHist::weightError(double& lo, double& hi, ErrorType etype) const
1549{
1550 checkInit() ;
1551
1552 switch (etype) {
1553
1554 case Auto:
1555 throw std::invalid_argument("RooDataHist::weightError(" + std::string(GetName()) + ") error type Auto not allowed here");
1556 break ;
1557
1558 case Expected:
1559 throw std::invalid_argument("RooDataHist::weightError(" + std::string(GetName()) + ") error type Expected not allowed here");
1560 break ;
1561
1562 case Poisson:
1563 if (_errLo && _errLo[_curIndex] >= 0.0) {
1564 // Weight is preset or precalculated
1565 lo = _errLo[_curIndex];
1566 hi = _errHi[_curIndex];
1567 return ;
1568 }
1569
1570 // We didn't track asymmetric errors so far, so now we need to allocate
1572
1573 // Calculate poisson errors
1574 double ym;
1575 double yp;
1576 RooHistError::instance().getPoissonInterval(Int_t(weight()+0.5),ym,yp,1) ;
1579 lo = _errLo[_curIndex];
1580 hi = _errHi[_curIndex];
1581 return ;
1582
1583 case SumW2:
1584 lo = std::sqrt(weightSquared(_curIndex));
1585 hi = lo;
1586 return ;
1587
1588 case None:
1589 lo = 0 ;
1590 hi = 0 ;
1591 return ;
1592 }
1593}
1594
1595
1596// wve adjust for variable bin sizes
1597
1598////////////////////////////////////////////////////////////////////////////////
1599/// Perform boundary safe 'intOrder'-th interpolation of weights in dimension 'dim'
1600/// at current value 'xval'
1601
1602/// \param[in] iDim Index of the histogram dimension along which to interpolate.
1603/// \param[in] xval Value of histogram variable at dimension `iDim` for which
1604/// we want to interpolate the histogram weight.
1605/// \param[in] centralIdx Index of the bin that the point at which we
1606/// interpolate the histogram weight falls into
1607/// (can be obtained with `RooDataHist::calcTreeIndex`).
1608/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1609/// used for the interpolation.
1610/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1611/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1612/// underflow bins are assumed to have weight zero and
1613/// overflow bins have weight one. Otherwise, the
1614/// histogram is mirrored at the boundaries for the
1615/// interpolation.
1617{
1618 auto const& binning = static_cast<RooRealVar&>(*_vars[iDim]).getBinning();
1619
1620 // Fill workspace arrays spanning interpolation area
1621 int fbinC = binning.binNumber(xval) ;
1622 int fbinLo = fbinC-intOrder/2 - ((xval<binning.binCenter(fbinC))?1:0) ;
1623 int fbinM = binning.numBins() ;
1624
1625 auto idxMult = _idxMult[iDim];
1626 auto offsetIdx = centralIdx - idxMult * fbinC;
1627
1628 // Use a class-member buffer to avoid repeated heap allocations.
1629 double * yarr = _interpolationBuffer.data();
1630 double * xarr = yarr + intOrder + 1;
1631
1632 for (int i=fbinLo ; i<=intOrder+fbinLo ; i++) {
1633 int ibin ;
1634 if (i>=0 && i<fbinM) {
1635 // In range
1636 ibin = i ;
1637 xarr[i-fbinLo] = binning.binCenter(ibin) ;
1638 auto idx = offsetIdx + idxMult * ibin;
1639 yarr[i - fbinLo] = _wgt[idx];
1640 if (correctForBinSize) yarr[i-fbinLo] /= _binv[idx] ;
1641 } else if (i>=fbinM) {
1642 // Overflow: mirror
1643 ibin = 2*fbinM-i-1 ;
1644 if (cdfBoundaries) {
1645 xarr[i-fbinLo] = binning.highBound()+1e-10*(i-fbinM+1) ;
1646 yarr[i-fbinLo] = 1.0 ;
1647 } else {
1648 auto idx = offsetIdx + idxMult * ibin;
1649 xarr[i-fbinLo] = 2*binning.highBound()-binning.binCenter(ibin) ;
1650 yarr[i - fbinLo] = _wgt[idx];
1652 yarr[i - fbinLo] /= _binv[idx];
1653 }
1654 } else {
1655 // Underflow: mirror
1656 ibin = -i - 1 ;
1657 if (cdfBoundaries) {
1658 xarr[i-fbinLo] = binning.lowBound()-ibin*(1e-10) ;
1659 yarr[i-fbinLo] = 0.0 ;
1660 } else {
1661 auto idx = offsetIdx + idxMult * ibin;
1662 xarr[i-fbinLo] = 2*binning.lowBound()-binning.binCenter(ibin) ;
1663 yarr[i - fbinLo] = _wgt[idx];
1665 yarr[i - fbinLo] /= _binv[idx];
1666 }
1667 }
1668 }
1670}
1671
1672
1673
1674
1675////////////////////////////////////////////////////////////////////////////////
1676/// Increment the bin content of the bin enclosing the given coordinates.
1677///
1678/// \param[in] row Coordinates of the bin.
1679/// \param[in] wgt Increment by this weight.
1680/// \param[in] sumw2 Optionally, track the sum of squared weights. If a value > 0 or
1681/// a weight != 1. is passed for the first time, a vector for the squared weights will be allocated.
1682void RooDataHist::add(const RooArgSet& row, double wgt, double sumw2)
1683{
1684 checkInit() ;
1685
1686 if ((sumw2 > 0. || wgt != 1.) && !_sumw2) {
1687 // Receiving a weighted entry. SumW2 != sumw from now on.
1688 _sumw2 = new double[_arrSize];
1689 std::copy(_wgt, _wgt+_arrSize, _sumw2);
1690
1692 }
1693
1694 const auto idx = calcTreeIndex(row, false);
1695
1696 _wgt[idx] += wgt ;
1697 if (_sumw2) _sumw2[idx] += (sumw2 > 0 ? sumw2 : wgt*wgt);
1698
1699 _cache_sum_valid = false;
1700}
1701
1702
1703
1704////////////////////////////////////////////////////////////////////////////////
1705/// Set a bin content.
1706/// \param[in] row Coordinates of the bin to be set.
1707/// \param[in] wgt New bin content.
1708/// \param[in] wgtErrLo Low error of the bin content.
1709/// \param[in] wgtErrHi High error of the bin content.
1710void RooDataHist::set(const RooArgSet& row, double wgt, double wgtErrLo, double wgtErrHi)
1711{
1712 checkInit() ;
1713
1715
1716 const auto idx = calcTreeIndex(row, false);
1717
1718 _wgt[idx] = wgt ;
1719 _errLo[idx] = wgtErrLo ;
1720 _errHi[idx] = wgtErrHi ;
1721
1722 _cache_sum_valid = false;
1723}
1724
1725
1726
1727////////////////////////////////////////////////////////////////////////////////
1728/// Set bin content of bin that was last loaded with get(std::size_t).
1729/// \param[in] binNumber Optional bin number to set. If empty, currently active bin is set.
1730/// \param[in] wgt New bin content.
1731/// \param[in] wgtErr Error of the new bin content. If the weight need not have an error, use 0. or a negative number.
1732void RooDataHist::set(std::size_t binNumber, double wgt, double wgtErr) {
1733 checkInit() ;
1734
1735 if (wgtErr > 0. && !_sumw2) {
1736 // Receiving a weighted entry. Need to track sumw2 from now on:
1738
1740 }
1741
1742 _wgt[binNumber] = wgt ;
1743 if (_errLo) _errLo[binNumber] = wgtErr;
1744 if (_errHi) _errHi[binNumber] = wgtErr;
1745 if (_sumw2) _sumw2[binNumber] = wgtErr*wgtErr;
1746
1748}
1749
1750
1751////////////////////////////////////////////////////////////////////////////////
1752/// Set bin content of bin that was last loaded with get(std::size_t).
1753/// \param[in] wgt New bin content.
1754/// \param[in] wgtErr Optional error of the bin content.
1755void RooDataHist::set(double wgt, double wgtErr) {
1756 if (_curIndex == std::numeric_limits<std::size_t>::max()) {
1757 _curIndex = calcTreeIndex(_vars, true) ;
1758 }
1759
1761}
1762
1763
1764////////////////////////////////////////////////////////////////////////////////
1765/// Set a bin content.
1766/// \param[in] row Coordinates to compute the bin from.
1767/// \param[in] wgt New bin content.
1768/// \param[in] wgtErr Optional error of the bin content.
1769void RooDataHist::set(const RooArgSet& row, double wgt, double wgtErr) {
1770 set(calcTreeIndex(row, false), wgt, wgtErr);
1771}
1772
1773
1774
1775////////////////////////////////////////////////////////////////////////////////
1776/// Add all data points contained in 'dset' to this data set with given weight.
1777/// Optional cut string expression selects the data points to be added and can
1778/// reference any variable contained in this data set
1779
1780void RooDataHist::add(const RooAbsData& dset, const char* cut, double wgt)
1781{
1782 RooFormulaVar cutVar("select",cut,*dset.get()) ;
1783 add(dset,&cutVar,wgt) ;
1784}
1785
1786
1787
1788////////////////////////////////////////////////////////////////////////////////
1789/// Add all data points contained in 'dset' to this data set with given weight.
1790/// Optional RooFormulaVar pointer selects the data points to be added.
1791
1793{
1794 checkInit() ;
1795
1796 RooFormulaVar* cloneVar = nullptr;
1797 std::unique_ptr<RooArgSet> tmp;
1798 if (cutVar) {
1799 // Deep clone cutVar and attach clone to this dataset
1800 tmp = std::make_unique<RooArgSet>();
1801 if(RooArgSet(*cutVar).snapshot(*tmp)) {
1802 coutE(DataHandling) << "RooDataHist::add(" << GetName() << ") Couldn't deep-clone cut variable, abort," << std::endl ;
1803 return ;
1804 }
1805
1806 cloneVar = static_cast<RooFormulaVar*>(tmp->find(*cutVar)) ;
1807 cloneVar->attachDataSet(dset) ;
1808 }
1809
1810
1811 Int_t i ;
1812 for (i=0 ; i<dset.numEntries() ; i++) {
1813 const RooArgSet* row = dset.get(i) ;
1814 if (!cloneVar || cloneVar->getVal()) {
1815 add(*row,wgt*dset.weight(), wgt*wgt*dset.weightSquared()) ;
1816 }
1817 }
1818
1820}
1821
1822
1823
1824////////////////////////////////////////////////////////////////////////////////
1825/// Return the sum of the weights of all bins in the histogram.
1826///
1827/// \param[in] correctForBinSize Multiply the sum of weights in each bin
1828/// with the N-dimensional bin volume, making the return value
1829/// the integral over the function represented by this histogram.
1830/// \param[in] inverseBinCor Divide by the N-dimensional bin volume.
1832{
1833 checkInit() ;
1834
1835 // Check if result was cached
1837 if (_cache_sum_valid == static_cast<Int_t>(cache_code)) {
1838 return _cache_sum ;
1839 }
1840
1842 for (Int_t i=0; i < _arrSize; i++) {
1843 const double theBinVolume = correctForBinSize ? (inverseBinCor ? 1/_binv[i] : _binv[i]) : 1.0 ;
1844 kahanSum += _wgt[i] * theBinVolume;
1845 }
1846
1847 // Store result in cache
1849 _cache_sum = kahanSum.Sum();
1850
1851 return kahanSum.Sum();
1852}
1853
1854
1855
1856////////////////////////////////////////////////////////////////////////////////
1857/// Return the sum of the weights of a multi-dimensional slice of the histogram
1858/// by summing only over the dimensions specified in sumSet.
1859///
1860/// The coordinates of all other dimensions are fixed to those given in sliceSet
1861///
1862/// If correctForBinSize is specified, the sum of weights
1863/// is multiplied by the M-dimensional bin volume, (M = N(sumSet)),
1864/// making the return value the integral over the function
1865/// represented by this histogram
1866
1868{
1869 checkInit() ;
1870
1872 varSave.addClone(_vars) ;
1873
1875 sliceOnlySet.remove(sumSet,true,true) ;
1876
1878 std::vector<double> const * pbinv = nullptr;
1879
1882 } else if(correctForBinSize && !inverseBinCor) {
1884 }
1885
1886 // Calculate mask and reference plot bins for non-iterating variables
1887 std::vector<bool> mask(_vars.size());
1888 std::vector<int> refBin(_vars.size());
1889
1890 for (unsigned int i = 0; i < _vars.size(); ++i) {
1891 const RooAbsArg* arg = _vars[i];
1892 const RooAbsLValue* argLv = _lvvars[i]; // Same as above, but cross-cast
1893
1894 if (sumSet.find(*arg)) {
1895 mask[i] = false ;
1896 } else {
1897 mask[i] = true ;
1898 refBin[i] = argLv->getBin();
1899 }
1900 }
1901
1902 // Loop over entire data set, skipping masked entries
1904 for (Int_t ibin=0; ibin < _arrSize; ++ibin) {
1905
1906 std::size_t tmpibin = ibin;
1907 bool skip(false) ;
1908
1909 // Check if this bin belongs in selected slice
1910 for (unsigned int ivar = 0; !skip && ivar < _vars.size(); ++ivar) {
1911 const Int_t idx = tmpibin / _idxMult[ivar] ;
1912 tmpibin -= idx*_idxMult[ivar] ;
1913 if (mask[ivar] && idx!=refBin[ivar])
1914 skip = true ;
1915 }
1916
1917 if (!skip) {
1918 const double theBinVolume = correctForBinSize ? (inverseBinCor ? 1/(*pbinv)[ibin] : (*pbinv)[ibin] ) : 1.0 ;
1920 }
1921 }
1922
1924
1925 return total.Sum();
1926}
1927
1928////////////////////////////////////////////////////////////////////////////////
1929/// Return the sum of the weights of a multi-dimensional slice of the histogram
1930/// by summing only over the dimensions specified in sumSet.
1931///
1932/// The coordinates of all other dimensions are fixed to those given in sliceSet
1933///
1934/// If correctForBinSize is specified, the sum of weights
1935/// is multiplied by the M-dimensional bin volume, (M = N(sumSet)),
1936/// or the fraction of it that falls inside the range rangeName,
1937/// making the return value the integral over the function
1938/// represented by this histogram.
1939///
1940/// If correctForBinSize is not specified, the weights are multiplied by the
1941/// fraction of the bin volume that falls inside the range, i.e. a factor of
1942/// binVolumeInRange/totalBinVolume.
1943
1946 const std::map<const RooAbsArg*, std::pair<double, double> >& ranges,
1947 std::function<double(int)> getBinScale)
1948{
1949 checkInit();
1952 varSave.addClone(_vars);
1953 {
1955 sliceOnlySet.remove(sumSet, true, true);
1957 }
1958
1959 // Calculate mask and reference plot bins for non-iterating variables,
1960 // and get ranges for iterating variables
1961 std::vector<bool> mask(_vars.size());
1962 std::vector<int> refBin(_vars.size());
1963 std::vector<double> rangeLo(_vars.size(), -std::numeric_limits<double>::infinity());
1964 std::vector<double> rangeHi(_vars.size(), +std::numeric_limits<double>::infinity());
1965
1966 for (std::size_t i = 0; i < _vars.size(); ++i) {
1967 const RooAbsArg* arg = _vars[i];
1968 const RooAbsLValue* argLV = _lvvars[i]; // Same object as above, but cross cast
1969
1970 RooAbsArg* sumsetv = sumSet.find(*arg);
1971 RooAbsArg* slicesetv = sliceSet.find(*arg);
1972 mask[i] = !sumsetv;
1973 if (mask[i]) {
1974 assert(argLV);
1975 refBin[i] = argLV->getBin();
1976 }
1977
1978 auto it = ranges.find(sumsetv ? sumsetv : slicesetv);
1979 if (ranges.end() != it) {
1980 rangeLo[i] = it->second.first;
1981 rangeHi[i] = it->second.second;
1982 }
1983 }
1984
1985 // Loop over entire data set, skipping masked entries
1987 for (Int_t ibin = 0; ibin < _arrSize; ++ibin) {
1988 // Check if this bin belongs in selected slice
1989 bool skip{false};
1990 for (int ivar = 0, tmp = ibin; !skip && ivar < int(_vars.size()); ++ivar) {
1991 const Int_t idx = tmp / _idxMult[ivar];
1992 tmp -= idx*_idxMult[ivar];
1993 if (mask[ivar] && idx!=refBin[ivar]) skip = true;
1994 }
1995
1996 if (skip) continue;
1997
1998 // Work out bin volume
1999 // It's not necessary to figure out the bin volume for the slice-only set explicitly here.
2000 // We need to loop over the sumSet anyway to get the partial bin containment correction,
2001 // so we can get the slice-only set volume later by dividing _binv[ibin] / binVolumeSumSetFull.
2002 double binVolumeSumSetFull = 1.;
2003 double binVolumeSumSetInRange = 1.;
2004 for (Int_t ivar = 0, tmp = ibin; ivar < (int)_vars.size(); ++ivar) {
2005 const Int_t idx = tmp / _idxMult[ivar];
2006 tmp -= idx*_idxMult[ivar];
2007
2008 // If the current variable is not in the sumSet, it should not be considered for the bin volume
2009 const auto arg = _vars[ivar];
2010 if (!sumSet.find(*arg)) {
2011 continue;
2012 }
2013
2014 if (_binbounds[ivar].empty()) continue;
2015 const double binLo = _binbounds[ivar][2 * idx];
2016 const double binHi = _binbounds[ivar][2 * idx + 1];
2017 if (binHi < rangeLo[ivar] || binLo > rangeHi[ivar]) {
2018 // bin is outside of allowed range - effective bin volume is zero
2020 break;
2021 }
2022
2024 binVolumeSumSetInRange *= std::min(rangeHi[ivar], binHi) - std::max(rangeLo[ivar], binLo);
2025 }
2027 if (0. == corrPartial) continue;
2029 total += getBinScale(ibin)*(_wgt[ibin] * corr * corrPartial);
2030 }
2031
2033
2034 return total.Sum();
2035}
2036
2037
2038
2039////////////////////////////////////////////////////////////////////////////////
2040/// Fill the transient cache with partial bin volumes with up-to-date
2041/// values for the partial volume specified by observables 'dimSet'
2042
2043const std::vector<double>& RooDataHist::calculatePartialBinVolume(const RooArgSet& dimSet) const
2044{
2045 // The code bitset has all bits set to one whose position corresponds to arguments in dimSet.
2046 // It is used as the key for the bin volume caching hash map.
2047 int code{0};
2048 {
2049 int i{0} ;
2050 for (auto const& v : _vars) {
2051 code += ((dimSet.find(*v) ? 1 : 0) << i) ;
2052 ++i;
2053 }
2054 }
2055
2056 auto& pbinv = _pbinvCache[code];
2057 if(!pbinv.empty()) {
2058 return pbinv;
2059 }
2060 pbinv.resize(_arrSize);
2061
2062 // Calculate plot bins of components from master index
2063 std::vector<bool> selDim(_vars.size());
2064 for (std::size_t i = 0; i < selDim.size(); ++i) {
2065 selDim[i] = (code >> i) & 1 ;
2066 }
2067
2068 // Recalculate partial bin volume cache
2069 for (Int_t ibin=0; ibin < _arrSize ;ibin++) {
2070 Int_t idx(0);
2071 Int_t tmp(ibin);
2072 double theBinVolume(1) ;
2073 for (unsigned int j=0; j < _lvvars.size(); ++j) {
2074 const RooAbsLValue* arg = _lvvars[j];
2075 assert(arg);
2076
2077 idx = tmp / _idxMult[j];
2078 tmp -= idx*_idxMult[j];
2079 if (selDim[j]) {
2080 theBinVolume *= arg->getBinWidth(idx) ;
2081 }
2082 }
2084 }
2085
2086 return pbinv;
2087}
2088
2089
2090////////////////////////////////////////////////////////////////////////////////
2091/// Sum the weights of all bins.
2095
2096
2097
2098////////////////////////////////////////////////////////////////////////////////
2099/// Return the sum of weights in all entries matching cutSpec (if specified)
2100/// and in named range cutRange (if specified)
2101/// Return the
2102
2103double RooDataHist::sumEntries(const char* cutSpec, const char* cutRange) const
2104{
2105 checkInit() ;
2106
2107 if (cutSpec==nullptr && cutRange==nullptr) {
2108 return sumEntries();
2109 } else {
2110
2111 // Setup RooFormulaVar for cutSpec if it is present
2112 std::unique_ptr<RooFormula> select;
2113 if (cutSpec) {
2114 select = std::make_unique<RooFormula>("select",cutSpec,*get());
2115 }
2116
2117 // Otherwise sum the weights in the event
2118 ROOT::Math::KahanSum<> kahanSum;
2119 for (Int_t i=0; i < _arrSize; i++) {
2120 get(i) ;
2121 if ((select && select->eval() == 0.) || (cutRange && !_vars.allInRange(cutRange)))
2122 continue;
2123
2124 kahanSum += weight(i);
2125 }
2126
2127 return kahanSum.Sum();
2128 }
2129}
2130
2131
2132
2133////////////////////////////////////////////////////////////////////////////////
2134/// Reset all bin weights to zero
2135
2137{
2138 // WVE DO NOT CALL RooTreeData::reset() for binned
2139 // datasets as this will delete the bin definitions
2140
2141 std::fill(_wgt, _wgt + _arrSize, 0.);
2142 delete[] _errLo; _errLo = nullptr;
2143 delete[] _errHi; _errHi = nullptr;
2144 delete[] _sumw2; _sumw2 = nullptr;
2145
2147
2148 _cache_sum_valid = false;
2149}
2150
2151
2152
2153////////////////////////////////////////////////////////////////////////////////
2154/// Load bin `binNumber`, and return an argset with the coordinates of the bin centre.
2155/// \note The argset is owned by this data hist, and this function has a side effect, because
2156/// it alters the currently active bin.
2157const RooArgSet* RooDataHist::get(Int_t binNumber) const
2158{
2159 checkInit() ;
2160 _curIndex = binNumber;
2161
2162 return RooAbsData::get(_curIndex);
2163}
2164
2165
2166
2167////////////////////////////////////////////////////////////////////////////////
2168/// Return a RooArgSet with whose coordinates denote the bin centre of the bin
2169/// enclosing the point in `coord`.
2170/// \note The argset is owned by this data hist, and this function has a side effect, because
2171/// it alters the currently active bin.
2173 return get(calcTreeIndex(coord, false));
2174}
2175
2176
2177
2178////////////////////////////////////////////////////////////////////////////////
2179/// Return the volume of the bin enclosing coordinates 'coord'.
2181 checkInit() ;
2182 return _binv[calcTreeIndex(coord, false)] ;
2183}
2184
2185
2186////////////////////////////////////////////////////////////////////////////////
2187/// Create an iterator over all bins in a slice defined by the subset of observables
2188/// listed in sliceArg. The position of the slice is given by otherArgs
2189
2191{
2192 // Update to current position
2194 _curIndex = calcTreeIndex(_vars, true);
2195
2197 if (!intArg) {
2198 coutE(InputArguments) << "RooDataHist::sliceIterator() variable " << sliceArg.GetName() << " is not part of this RooDataHist" << std::endl ;
2199 return nullptr ;
2200 }
2201 return new RooDataHistSliceIter(*this,*intArg) ;
2202}
2203
2204
2205////////////////////////////////////////////////////////////////////////////////
2206/// Change the name of the RooDataHist
2207
2209{
2210 if (_dir) _dir->GetList()->Remove(this);
2211 // We need to use the function from RooAbsData, because it already overrides TNamed::SetName
2213 if (_dir) _dir->GetList()->Add(this);
2214}
2215
2216
2217////////////////////////////////////////////////////////////////////////////////
2218/// Change the title of this RooDataHist
2219
2220void RooDataHist::SetNameTitle(const char *name, const char* title)
2221{
2222 SetName(name);
2223 SetTitle(title);
2224}
2225
2226
2227////////////////////////////////////////////////////////////////////////////////
2228/// Print value of the dataset, i.e. the sum of weights contained in the dataset
2229
2230void RooDataHist::printValue(ostream& os) const
2231{
2232 os << numEntries() << " bins (" << sumEntries() << " weights)" ;
2233}
2234
2235
2236
2237
2238////////////////////////////////////////////////////////////////////////////////
2239/// Print argument of dataset, i.e. the observable names
2240
2241void RooDataHist::printArgs(ostream& os) const
2242{
2243 os << "[" ;
2244 bool first(true) ;
2245 for (const auto arg : _vars) {
2246 if (first) {
2247 first=false ;
2248 } else {
2249 os << "," ;
2250 }
2251 os << arg->GetName() ;
2252 }
2253 os << "]" ;
2254}
2255
2256
2257
2258////////////////////////////////////////////////////////////////////////////////
2259/// Returns true if dataset contains entries with a non-integer weight.
2260
2262{
2263 for (Int_t i=0; i < _arrSize; ++i) {
2264 const double wgt = _wgt[i];
2265 double intpart;
2266 if (std::abs(std::modf(wgt, &intpart)) > 1.E-10)
2267 return true;
2268 }
2269
2270 return false;
2271}
2272
2273
2274////////////////////////////////////////////////////////////////////////////////
2275/// Print the details on the dataset contents
2276
2277void RooDataHist::printMultiline(ostream& os, Int_t content, bool verbose, TString indent) const
2278{
2280
2281 os << indent << "Binned Dataset " << GetName() << " (" << GetTitle() << ")" << std::endl ;
2282 os << indent << " Contains " << numEntries() << " bins with a total weight of " << sumEntries() << std::endl;
2283
2284 if (!verbose) {
2285 os << indent << " Observables " << _vars << std::endl ;
2286 } else {
2287 os << indent << " Observables: " ;
2289 }
2290
2291 if(verbose) {
2292 if (!_cachedVars.empty()) {
2293 os << indent << " Caches " << _cachedVars << std::endl ;
2294 }
2295 }
2296}
2297
2298/**
2299 * \brief Prints the contents of the RooDataHist to the specified output stream.
2300 *
2301 * This function iterates through all bins of the histogram and prints the
2302 * coordinates of each bin, along with its weight and statistical error.
2303 * It is designed to be robust, handling empty or invalid datasets,
2304 * and works for histograms of any dimension.
2305 *
2306 * \param os The output stream (e.g., std::cout) to write the contents to.
2307 */
2308void RooDataHist::printContents(std::ostream& os) const
2309{
2310 os << "Contents of RooDataHist \"" << GetName() << "\"" << std::endl;
2311
2312 if (numEntries() == 0) {
2313 os << "(dataset is empty)" << std::endl;
2314 return;
2315 }
2316
2317 for (int i = 0; i < numEntries(); ++i) {
2318 const RooArgSet* obs = get(i); // load i-th bin
2319 os << " Bin " << i << ": ";
2320
2321 bool first = true;
2322 for (const auto* var : *obs) {
2323 if (!first) os << ", ";
2324 first = false;
2325
2326 os << var->GetName() << "=";
2327 if (auto realVar = dynamic_cast<const RooRealVar*>(var)) {
2328 os << realVar->getVal();
2329 } else if (auto catVar = dynamic_cast<const RooCategory*>(var)) {
2330 os << catVar->getCurrentLabel();
2331 } else {
2332 os << "(unsupported type)"; //added as a precaution
2333 }
2334 }
2335
2336 double lo, hi;
2338 os << ", weight=" << weight() << " +/- [" << lo << "," << hi << "]"
2339 << std::endl;
2340 }
2341}
2342
2343
2344////////////////////////////////////////////////////////////////////////////////
2345/// Stream an object of class RooDataHist.
2347 if (R__b.IsReading()) {
2348
2349 UInt_t R__s;
2350 UInt_t R__c;
2351 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2352
2353 if (R__v > 2) {
2354 R__b.ReadClassBuffer(RooDataHist::Class(),this,R__v,R__s,R__c);
2355 R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());
2356 initialize(nullptr, false);
2357 } else {
2358
2359 // Legacy dataset conversion happens here. Legacy RooDataHist inherits from RooTreeData
2360 // which in turn inherits from RooAbsData. Manually stream RooTreeData contents on
2361 // file here and convert it into a RooTreeDataStore which is installed in the
2362 // new-style RooAbsData base class
2363
2364 // --- This is the contents of the streamer code of RooTreeData version 2 ---
2365 UInt_t R__s1;
2366 UInt_t R__c1;
2367 Version_t R__v1 = R__b.ReadVersion(&R__s1, &R__c1); if (R__v1) { }
2368
2370 TTree* X_tree(nullptr) ; R__b >> X_tree;
2371 RooArgSet X_truth ; X_truth.Streamer(R__b);
2373 R__b.CheckByteCount(R__s1, R__c1, TClass::GetClass("RooTreeData"));
2374 // --- End of RooTreeData-v1 streamer
2375
2376 // Construct RooTreeDataStore from X_tree and complete initialization of new-style RooAbsData
2377 _dstore = std::make_unique<RooTreeDataStore>(X_tree,_vars);
2378 _dstore->SetName(GetName()) ;
2379 _dstore->SetTitle(GetTitle()) ;
2380 _dstore->checkInit() ;
2381
2383 R__b >> _arrSize;
2384 delete [] _wgt;
2385 _wgt = new double[_arrSize];
2386 R__b.ReadFastArray(_wgt,_arrSize);
2387 delete [] _errLo;
2388 _errLo = new double[_arrSize];
2389 R__b.ReadFastArray(_errLo,_arrSize);
2390 delete [] _errHi;
2391 _errHi = new double[_arrSize];
2392 R__b.ReadFastArray(_errHi,_arrSize);
2393 delete [] _sumw2;
2394 _sumw2 = new double[_arrSize];
2395 R__b.ReadFastArray(_sumw2,_arrSize);
2396 delete [] _binv;
2397 _binv = new double[_arrSize];
2399 tmpSet.Streamer(R__b);
2400 double tmp;
2401 R__b >> tmp; //_curWeight;
2402 R__b >> tmp; //_curWgtErrLo;
2403 R__b >> tmp; //_curWgtErrHi;
2404 R__b >> tmp; //_curSumW2;
2405 R__b >> tmp; //_curVolume;
2406 R__b >> _curIndex;
2407 R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());
2408 }
2409
2410 } else {
2411
2412 R__b.WriteClassBuffer(RooDataHist::Class(),this);
2413 }
2414}
2415
2416
2417////////////////////////////////////////////////////////////////////////////////
2418/// Return event weights of all events in range [first, first+len).
2419/// If cacheValidEntries() has been called, out-of-range events will have a weight of 0.
2420std::span<const double> RooDataHist::getWeightBatch(std::size_t first, std::size_t len, bool sumW2 /*=false*/) const {
2421 return {(sumW2 && _sumw2 ? _sumw2 : _wgt) + first, len};
2422}
2423
2424
2425////////////////////////////////////////////////////////////////////////////////
2426/// Hand over pointers to our weight arrays to the data store implementation.
2428 _dstore->setExternalWeightArray(_wgt, _errLo, _errHi, _sumw2);
2429}
2430
2431
2432////////////////////////////////////////////////////////////////////////////////
2433/// Return reference to VarInfo struct with cached histogram variable
2434/// information that is frequently used for histogram weights retrieval.
2435///
2436/// If the `_varInfo` struct was not initialized yet, it will be initialized in
2437/// this function.
2439
2440 if(_varInfo.initialized) return _varInfo;
2441
2442 auto& info = _varInfo;
2443
2444 {
2445 // count the number of real vars and get their indices
2446 info.nRealVars = 0;
2447 size_t iVar = 0;
2448 for (const auto real : _vars) {
2449 if (dynamic_cast<RooRealVar*>(real)) {
2450 if(info.nRealVars == 0) info.realVarIdx1 = iVar;
2451 if(info.nRealVars == 1) info.realVarIdx2 = iVar;
2452 ++info.nRealVars;
2453 }
2454 ++iVar;
2455 }
2456 }
2457
2458 {
2459 // assert that the variables are either real values or categories
2460 for (unsigned int i=0; i < _vars.size(); ++i) {
2461 if (_lvbins[i].get()) {
2462 assert(dynamic_cast<const RooAbsReal*>(_vars[i]));
2463 } else {
2464 assert(dynamic_cast<const RooAbsCategoryLValue*>(_vars[i]));
2465 }
2466 }
2467 }
2468
2469 info.initialized = true;
2470
2471 return info;
2472}
#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:110
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:2973
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:9142
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:5103
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:226
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()