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