Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RooDataHist.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooDataHist.cxx
19\class RooDataHist
20\ingroup Roofitcore
21
22Container class to hold N-dimensional binned data. Each bin's central
23coordinates in N-dimensional space are represented by a RooArgSet containing RooRealVar, RooCategory
24or RooStringVar objects, thus data can be binned in real and/or discrete dimensions.
25
26There is an unbinned equivalent, RooDataSet.
27
28### Inspecting a datahist
29Inspect a datahist using Print() to get the coordinates and `weight()` to get the bin contents:
30```
31datahist->Print("V");
32datahist->get(0)->Print("V"); std::cout << "w=" << datahist->weight(0) << std::endl;
33datahist->get(1)->Print("V"); std::cout << "w=" << datahist->weight(1) << std::endl;
34...
35```
36
37### Plotting data.
38See RooAbsData::plotOn().
39
40### Creating a datahist using RDataFrame
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 std::cout, std::endl, std::string, std::ostream;
80
82
83
84////////////////////////////////////////////////////////////////////////////////
85/// Default constructor
86
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 (int(vars.size()) != 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
278 const RooCmdArg& arg4,const RooCmdArg& arg5,const RooCmdArg& arg6,const RooCmdArg& arg7,const RooCmdArg& arg8) :
279 RooAbsData(name,title,RooArgSet(vars,static_cast<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
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 = static_cast<RooRealVar*>(_vars.find(vars.at(0)->GetName())) ;
386 RooRealVar* yvar = static_cast<RooRealVar*>(vars.at(1) ? _vars.find(vars.at(1)->GetName()) : nullptr ) ;
387 RooRealVar* zvar = static_cast<RooRealVar*>(vars.at(2) ? _vars.find(vars.at(2)->GetName()) : nullptr ) ;
388
389 // Transfer contents
390 Int_t xmin(0);
391 Int_t ymin(0);
392 Int_t zmin(0);
394 xmin = offset[0] ;
395 if (yvar) {
396 vset.add(*yvar) ;
397 ymin = offset[1] ;
398 }
399 if (zvar) {
400 vset.add(*zvar) ;
401 zmin = offset[2] ;
402 }
403
404 Int_t iX(0);
405 Int_t iY(0);
406 Int_t iz(0);
407 for (iX=0 ; iX < xvar->getBins() ; iX++) {
408 xvar->setBin(iX) ;
409 if (yvar) {
410 for (iY=0 ; iY < yvar->getBins() ; iY++) {
411 yvar->setBin(iY) ;
412 if (zvar) {
413 for (iz=0 ; iz < zvar->getBins() ; iz++) {
414 zvar->setBin(iz) ;
415 double bv = doDensityCorrection ? binVolume(vset) : 1;
416 add(vset,bv*histo.GetBinContent(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,2)) ;
417 }
418 } else {
419 double bv = doDensityCorrection ? binVolume(vset) : 1;
420 add(vset,bv*histo.GetBinContent(iX+1+xmin,iY+1+ymin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin,iY+1+ymin)*wgt,2)) ;
421 }
422 }
423 } else {
424 double bv = doDensityCorrection ? binVolume(vset) : 1 ;
425 add(vset,bv*histo.GetBinContent(iX+1+xmin)*wgt,bv*std::pow(histo.GetBinError(iX+1+xmin)*wgt,2)) ;
426 }
427 }
428
429}
430
431namespace {
432bool checkConsistentAxes(const TH1* first, const TH1* second) {
433 return first->GetDimension() == second->GetDimension()
434 && first->GetNbinsX() == second->GetNbinsX()
435 && first->GetNbinsY() == second->GetNbinsY()
436 && first->GetNbinsZ() == second->GetNbinsZ()
437 && first->GetXaxis()->GetXmin() == second->GetXaxis()->GetXmin()
438 && first->GetXaxis()->GetXmax() == second->GetXaxis()->GetXmax()
439 && (first->GetNbinsY() == 1 || (first->GetYaxis()->GetXmin() == second->GetYaxis()->GetXmin()
440 && first->GetYaxis()->GetXmax() == second->GetYaxis()->GetXmax() ) )
441 && (first->GetNbinsZ() == 1 || (first->GetZaxis()->GetXmin() == second->GetZaxis()->GetXmin()
442 && first->GetZaxis()->GetXmax() == second->GetZaxis()->GetXmax() ) );
443}
444
446
447 // Relative tolerance for bin boundary comparison
448 constexpr double tolerance = 1e-6;
449
450 auto const& vars1 = *h1.get();
451 auto const& vars2 = *h2.get();
452
453 // Check if number of variables and names is consistent
454 if(!vars1.hasSameLayout(vars2)) {
455 return false;
456 }
457
458 for(std::size_t iVar = 0; iVar < vars1.size(); ++iVar) {
459 auto * var1 = dynamic_cast<RooRealVar*>(vars1[iVar]);
460 auto * var2 = dynamic_cast<RooRealVar*>(vars2[iVar]);
461
462 // Check if variables are consistently real-valued
463 if((!var1 && var2) || (var1 && !var2)) return false;
464
465 // Not a real-valued variable
466 if(!var1) continue;
467
468 // Now check the binning
469 auto const& bng1 = var1->getBinning();
470 auto const& bng2 = var2->getBinning();
471
472 // Compare bin numbers
473 if(bng1.numBins() != bng2.numBins()) return false;
474
475 std::size_t nBins = bng1.numBins();
476
477 // Compare bin boundaries
478 for(std::size_t iBin = 0; iBin < nBins; ++iBin) {
479 double v1 = bng1.binLow(iBin);
480 double v2 = bng2.binLow(iBin);
481 if(std::abs((v1 - v2) / v1) > tolerance) return false;
482 }
483 double v1 = bng1.binHigh(nBins - 1);
484 double v2 = bng2.binHigh(nBins - 1);
485 if(std::abs((v1 - v2) / v1) > tolerance) return false;
486 }
487 return true;
488}
489}
490
491
492////////////////////////////////////////////////////////////////////////////////
493/// Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
494/// in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
495/// and the import source
496
497void RooDataHist::importTH1Set(const RooArgList& vars, RooCategory& indexCat, std::map<string,TH1*> hmap, double wgt, bool doDensityCorrection)
498{
499 RooCategory* icat = static_cast<RooCategory*>(_vars.find(indexCat.GetName())) ;
500
501 TH1* histo(nullptr) ;
502 bool init(false) ;
503 for (const auto& hiter : hmap) {
504 // Store pointer to first histogram from which binning specification will be taken
505 if (!histo) {
506 histo = hiter.second;
507 } else {
508 if (!checkConsistentAxes(histo, hiter.second)) {
509 coutE(InputArguments) << "Axes of histogram " << hiter.second->GetName() << " are not consistent with first processed "
510 << "histogram " << histo->GetName() << std::endl;
511 throw std::invalid_argument("Axes of inputs for RooDataHist are inconsistent");
512 }
513 }
514 // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
515 if (!indexCat.hasLabel(hiter.first)) {
516 indexCat.defineType(hiter.first) ;
517 coutI(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << ") defining state \"" << hiter.first << "\" in index category " << indexCat.GetName() << endl ;
518 }
519 if (!icat->hasLabel(hiter.first)) {
520 icat->defineType(hiter.first) ;
521 }
522 }
523
524 // Check consistency in number of dimensions
525 if (histo && int(vars.size()) != histo->GetDimension()) {
526 coutE(InputArguments) << "RooDataHist::importTH1Set(" << GetName() << "): dimension of input histogram must match "
527 << "number of continuous variables" << endl ;
528 throw std::invalid_argument("Inputs histograms for RooDataHist are not compatible with dimensions of variables.");
529 }
530
531 // Copy bins and ranges from THx to dimension observables
532 Int_t offset[3] ;
533 adjustBinning(vars,*histo,offset) ;
534
535 // Initialize internal data structure
536 if (!init) {
537 initialize() ;
538 appendToDir(this,true) ;
539 init = true ;
540 }
541
542 // Define x,y,z as 1st, 2nd and 3rd observable
543 RooRealVar* xvar = static_cast<RooRealVar*>(_vars.find(vars.at(0)->GetName())) ;
544 RooRealVar* yvar = static_cast<RooRealVar*>(vars.at(1) ? _vars.find(vars.at(1)->GetName()) : nullptr ) ;
545 RooRealVar* zvar = static_cast<RooRealVar*>(vars.at(2) ? _vars.find(vars.at(2)->GetName()) : nullptr ) ;
546
547 // Transfer contents
548 Int_t xmin(0);
549 Int_t ymin(0);
550 Int_t zmin(0);
552 double volume = xvar->getMax()-xvar->getMin() ;
553 xmin = offset[0] ;
554 if (yvar) {
555 vset.add(*yvar) ;
556 ymin = offset[1] ;
557 volume *= (yvar->getMax()-yvar->getMin()) ;
558 }
559 if (zvar) {
560 vset.add(*zvar) ;
561 zmin = offset[2] ;
562 volume *= (zvar->getMax()-zvar->getMin()) ;
563 }
564 double avgBV = volume / numEntries() ;
565
566 Int_t ic(0);
567 Int_t iX(0);
568 Int_t iY(0);
569 Int_t iz(0);
570 for (ic=0 ; ic < icat->numBins(nullptr) ; ic++) {
571 icat->setBin(ic) ;
572 histo = hmap[icat->getCurrentLabel()] ;
573 for (iX=0 ; iX < xvar->getBins() ; iX++) {
574 xvar->setBin(iX) ;
575 if (yvar) {
576 for (iY=0 ; iY < yvar->getBins() ; iY++) {
577 yvar->setBin(iY) ;
578 if (zvar) {
579 for (iz=0 ; iz < zvar->getBins() ; iz++) {
580 zvar->setBin(iz) ;
581 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
582 add(vset,bv*histo->GetBinContent(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin,iY+1+ymin,iz+1+zmin)*wgt,2)) ;
583 }
584 } else {
585 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
586 add(vset,bv*histo->GetBinContent(iX+1+xmin,iY+1+ymin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin,iY+1+ymin)*wgt,2)) ;
587 }
588 }
589 } else {
590 double bv = doDensityCorrection ? binVolume(vset)/avgBV : 1;
591 add(vset,bv*histo->GetBinContent(iX+1+xmin)*wgt,bv*std::pow(histo->GetBinError(iX+1+xmin)*wgt,2)) ;
592 }
593 }
594 }
595
596}
597
598
599
600////////////////////////////////////////////////////////////////////////////////
601/// Import data from given set of TH1/2/3 into this RooDataHist. The category indexCat labels the sources
602/// in the constructed RooDataHist. The stl map provides the mapping between the indexCat state labels
603/// and the import source
604
605void RooDataHist::importDHistSet(const RooArgList & /*vars*/, RooCategory &indexCat,
606 std::map<std::string, RooDataHist *> dmap, double initWgt)
607{
608 auto *icat = static_cast<RooCategory *>(_vars.find(indexCat.GetName()));
609
610 RooDataHist *dhistForBinning = nullptr;
611
612 for (const auto &diter : dmap) {
613
614 std::string const &label = diter.first;
615 RooDataHist *dhist = diter.second;
616
617 if (!dhistForBinning) {
619 } else {
621 coutE(InputArguments) << "Layout or binning of histogram " << dhist->GetName()
622 << " is not consistent with first processed "
623 << "histogram " << dhistForBinning->GetName() << std::endl;
624 throw std::invalid_argument("Layout or binning of inputs for RooDataHist is inconsistent");
625 }
626 }
627
628 // Define state labels in index category (both in provided indexCat and in internal copy in dataset)
629 if (!indexCat.hasLabel(label)) {
630 indexCat.defineType(label);
631 coutI(InputArguments) << "RooDataHist::importDHistSet(" << GetName() << ") defining state \"" << label
632 << "\" in index category " << indexCat.GetName() << endl;
633 }
634 if (!icat->hasLabel(label)) {
635 icat->defineType(label);
636 }
637 }
638
639 // adjust the binning of the created histogram
641 auto *ourVar = dynamic_cast<RooRealVar *>(_vars.find(theirVar->GetName()));
642 if (!theirVar || !ourVar)
643 continue;
644 ourVar->setBinning(theirVar->getBinning());
645 }
646
647 initialize();
648 appendToDir(this, true);
649
650 for (const auto &diter : dmap) {
651 std::string const &label = diter.first;
652 RooDataHist *dhist = diter.second;
653
654 icat->setLabel(label.c_str());
655
656 // Transfer contents
657 for (Int_t i = 0; i < dhist->numEntries(); i++) {
658 _vars.assign(*dhist->get(i));
659 add(_vars, dhist->weight() * initWgt, pow(dhist->weightError(SumW2), 2));
660 }
661 }
662}
663
664////////////////////////////////////////////////////////////////////////////////
665/// Helper doing the actual work of adjustBinning().
666
669{
670 const std::string ourVarName(ourVar->GetName() ? ourVar->GetName() : "");
671 const std::string ownName(GetName() ? GetName() : "");
672 // RooRealVar is derived from RooAbsRealLValue which is itself
673 // derived from RooAbsReal and a virtual class RooAbsLValue
674 // supplying setter functions, check if ourVar is indeed derived
675 // as real
676 if (!dynamic_cast<RooAbsReal *>(ourVar)) {
677 coutE(InputArguments) << "RooDataHist::adjustBinning(" << ownName << ") ERROR: dimension " << ourVarName
678 << " must be real\n";
679 throw std::logic_error("Incorrect type object (" + ourVarName +
680 ") passed as argument to RooDataHist::_adjustBinning. Please report this issue.");
681 }
682
683 const double xlo = theirVar.getMin();
684 const double xhi = theirVar.getMax();
685
686 const bool isUniform = !axis.GetXbins()->GetArray();
687 std::unique_ptr<RooAbsBinning> xbins;
688
689 if (!isUniform) {
690 xbins = std::make_unique<RooBinning>(axis.GetNbins(), axis.GetXbins()->GetArray());
691 } else {
692 xbins = std::make_unique<RooUniformBinning>(axis.GetXmin(), axis.GetXmax(), axis.GetNbins());
693 }
694
695 const double tolerance = 1e-6 * xbins->averageBinWidth();
696
697 // Adjust xlo/xhi to nearest boundary
698 const int iBinLo = xbins->binNumber(xlo + tolerance);
699 const int iBinHi = xbins->binNumber(xhi - tolerance);
700 const int nBinsAdj = iBinHi - iBinLo + 1;
701 const double xloAdj = xbins->binLow(iBinLo);
702 const double xhiAdj = xbins->binHigh(iBinHi);
703
704 if (isUniform) {
705 xbins = std::make_unique<RooUniformBinning>(xloAdj, xhiAdj, nBinsAdj);
706 theirVar.setRange(xloAdj, xhiAdj);
707 } else {
708 xbins->setRange(xloAdj, xhiAdj);
709 theirVar.setBinning(*xbins);
710 }
711
712 if (std::abs(xloAdj - xlo) > tolerance || std::abs(xhiAdj - xhi) > tolerance) {
713 coutI(DataHandling) << "RooDataHist::adjustBinning(" << ownName << "): fit range of variable " << ourVarName
714 << " expanded to nearest bin boundaries: [" << xlo << "," << xhi << "] --> [" << xloAdj << ","
715 << xhiAdj << "]"
716 << "\n";
717 }
718
719 ourVar->setBinning(*xbins);
720
721 // The offset is the bin number of the adjusted lower limit of the RooFit
722 // variable in the original TH1 histogram, starting from zero.
723 if (offset) {
724 *offset = axis.FindFixBin(xloAdj + tolerance) - 1;
725 }
726}
727
728////////////////////////////////////////////////////////////////////////////////
729/// Adjust binning specification on first and optionally second and third
730/// observable to binning in given reference TH1. Used by constructors
731/// that import data from an external TH1.
732/// Both the variables in vars and in this RooDataHist are adjusted.
733/// @param vars List with variables that are supposed to have their binning adjusted.
734/// @param href Reference histogram that dictates the binning
735/// @param offset If not nullptr, a possible bin count offset for the axes x,y,z is saved here as Int_t[3]
736
738{
739 auto xvar = static_cast<RooRealVar*>(_vars.find(*vars.at(0)) );
740 _adjustBinning(*static_cast<RooRealVar*>(vars.at(0)), *href.GetXaxis(), xvar, offset ? &offset[0] : nullptr);
741
742 if (vars.at(1)) {
743 auto yvar = static_cast<RooRealVar*>(_vars.find(*vars.at(1)));
744 if (yvar)
745 _adjustBinning(*static_cast<RooRealVar*>(vars.at(1)), *href.GetYaxis(), yvar, offset ? &offset[1] : nullptr);
746 }
747
748 if (vars.at(2)) {
749 auto zvar = static_cast<RooRealVar*>(_vars.find(*vars.at(2)));
750 if (zvar)
751 _adjustBinning(*static_cast<RooRealVar*>(vars.at(2)), *href.GetZaxis(), zvar, offset ? &offset[2] : nullptr);
752 }
753
754}
755
756
757namespace {
758/// Clone external weight arrays, unless the external array is nullptr.
759void cloneArray(double*& ours, const double* theirs, std::size_t n) {
760 if (ours) delete[] ours;
761 ours = nullptr;
762 if (!theirs) return;
763 ours = new double[n];
764 std::copy(theirs, theirs+n, ours);
765}
766
767/// Allocate and initialise an array with desired size and values.
768void initArray(double*& arr, std::size_t n, double val) {
769 if (arr) delete[] arr;
770 arr = nullptr;
771 if (n == 0) return;
772 arr = new double[n];
773 std::fill(arr, arr+n, val);
774}
775}
776
777
778////////////////////////////////////////////////////////////////////////////////
779/// Initialization procedure: allocate weights array, calculate
780/// multipliers needed for N-space to 1-dim array jump table,
781/// and fill the internal tree with all bin center coordinates
782
783void RooDataHist::initialize(const char* binningName, bool fillTree)
784{
785 _lvvars.clear();
786 _lvbins.clear();
787
788 // Fill array of LValue pointers to variables
789 for (unsigned int i = 0; i < _vars.size(); ++i) {
790 if (binningName) {
791 RooRealVar* rrv = dynamic_cast<RooRealVar*>(_vars[i]);
792 if (rrv) {
793 rrv->setBinning(rrv->getBinning(binningName));
794 }
795 }
796
797 auto lvarg = dynamic_cast<RooAbsLValue*>(_vars[i]);
798 assert(lvarg);
799 _lvvars.push_back(lvarg);
800
801 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
802 _lvbins.emplace_back(binning ? binning->clone() : nullptr);
803 }
804
805
806 // Allocate coefficients array
807 _idxMult.resize(_vars.size()) ;
808
809 _arrSize = 1 ;
810 unsigned int n = 0u;
811 for (const auto var : _vars) {
812 auto arg = dynamic_cast<const RooAbsLValue*>(var);
813 assert(arg);
814
815 // Calculate sub-index multipliers for master index
816 for (unsigned int i = 0u; i<n; i++) {
817 _idxMult[i] *= arg->numBins() ;
818 }
819 _idxMult[n++] = 1 ;
820
821 // Calculate dimension of weight array
822 _arrSize *= arg->numBins() ;
823 }
824
825 // Allocate and initialize weight array if necessary
826 if (!_wgt) {
827 initArray(_wgt, _arrSize, 0.);
828 delete[] _errLo; _errLo = nullptr;
829 delete[] _errHi; _errHi = nullptr;
830 delete[] _sumw2; _sumw2 = nullptr;
832
833 // Refill array pointers in data store when reading
834 // from Streamer
835 if (!fillTree) {
837 }
838 }
839
840 if (!fillTree) return ;
841
842 // Fill TTree with bin center coordinates
843 // Calculate plot bins of components from master index
844
845 for (Int_t ibin=0 ; ibin < _arrSize ; ibin++) {
846 Int_t j(0);
847 Int_t idx(0);
848 Int_t tmp(ibin);
849 double theBinVolume(1) ;
850 for (auto arg2 : _lvvars) {
851 idx = tmp / _idxMult[j] ;
852 tmp -= idx*_idxMult[j++] ;
853 arg2->setBin(idx) ;
854 theBinVolume *= arg2->getBinWidth(idx) ;
855 }
857
858 fill() ;
859 }
860
861
862}
863
864
865////////////////////////////////////////////////////////////////////////////////
866
868{
869 if (!_binbounds.empty()) return;
870 for (auto& it : _lvbins) {
871 _binbounds.push_back(std::vector<double>());
872 if (it) {
873 std::vector<double>& bounds = _binbounds.back();
874 bounds.reserve(2 * it->numBins());
875 for (Int_t i = 0; i < it->numBins(); ++i) {
876 bounds.push_back(it->binLow(i));
877 bounds.push_back(it->binHigh(i));
878 }
879 }
880 }
881}
882
883
884////////////////////////////////////////////////////////////////////////////////
885/// Copy constructor
886
888 RooAbsData(other,newname), RooDirItem(), _arrSize(other._arrSize), _idxMult(other._idxMult), _pbinvCache(other._pbinvCache)
889{
890 // Allocate and initialize weight array
891 assert(_arrSize == other._arrSize);
892 cloneArray(_wgt, other._wgt, other._arrSize);
893 cloneArray(_errLo, other._errLo, other._arrSize);
894 cloneArray(_errHi, other._errHi, other._arrSize);
895 cloneArray(_binv, other._binv, other._arrSize);
896 cloneArray(_sumw2, other._sumw2, other._arrSize);
897
898 // Fill array of LValue pointers to variables
899 for (const auto rvarg : _vars) {
900 auto lvarg = dynamic_cast<RooAbsLValue*>(rvarg);
901 assert(lvarg);
902 _lvvars.push_back(lvarg);
903 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
904 _lvbins.emplace_back(binning ? binning->clone() : nullptr) ;
905 }
906
908
909 appendToDir(this,true) ;
910}
911
912
913////////////////////////////////////////////////////////////////////////////////
914/// Implementation of RooAbsData virtual method that drives the RooAbsData::reduce() methods
915
916std::unique_ptr<RooAbsData> RooDataHist::reduceEng(const RooArgSet& varSubset, const RooFormulaVar* cutVar, const char* cutRange,
917 std::size_t nStart, std::size_t nStop) const
918{
919 checkInit() ;
922 auto rdh = std::make_unique<RooDataHist>(GetName(), GetTitle(), myVarSubset);
923
924 RooFormulaVar* cloneVar = nullptr;
925 std::unique_ptr<RooArgSet> tmp;
926 if (cutVar) {
927 tmp = std::make_unique<RooArgSet>();
928 // Deep clone cutVar and attach clone to this dataset
929 if (RooArgSet(*cutVar).snapshot(*tmp)) {
930 coutE(DataHandling) << "RooDataHist::reduceEng(" << GetName() << ") Couldn't deep-clone cut variable, abort," << endl ;
931 return nullptr;
932 }
933 cloneVar = static_cast<RooFormulaVar*>(tmp->find(*cutVar));
934 cloneVar->attachDataSet(*this) ;
935 }
936
937 double lo;
938 double hi;
939 const std::size_t nevt = nStop < static_cast<std::size_t>(numEntries()) ? nStop : static_cast<std::size_t>(numEntries());
940 for (auto i=nStart; i<nevt ; i++) {
941 const RooArgSet* row = get(i) ;
942
943 bool doSelect(true) ;
944 if (cutRange) {
945 for (const auto arg : *row) {
946 if (!arg->inRange(cutRange)) {
947 doSelect = false ;
948 break ;
949 }
950 }
951 }
952 if (!doSelect) continue ;
953
954 if (!cloneVar || cloneVar->getVal()) {
955 weightError(lo,hi,SumW2) ;
956 rdh->add(*row,weight(),lo*lo) ;
957 }
958 }
959
960 return rdh ;
961}
962
963
964
965////////////////////////////////////////////////////////////////////////////////
966/// Destructor
967
969{
970 delete[] _wgt;
971 delete[] _errLo;
972 delete[] _errHi;
973 delete[] _sumw2;
974 delete[] _binv;
975
976 removeFromDir(this) ;
978}
979
980
981
982
983////////////////////////////////////////////////////////////////////////////////
984/// Calculate bin number of the given coordinates. If only a subset of the internal
985/// coordinates are passed, the missing coordinates are taken at their current value.
986/// \param[in] coord Variables that are representing the coordinates.
987/// \param[in] fast If the variables in `coord` and the ones of the data hist have the
988/// same size and layout, `fast` can be set to skip checking that all variables are
989/// present in `coord`.
991 checkInit() ;
992 return calcTreeIndex(coord, fast);
993}
994
996 bool correctForBinSize) const
997{
998 std::vector<double> vals(_arrSize);
999 for (std::size_t i = 0; i < vals.size(); ++i) {
1000 vals[i] = correctForBinSize ? _wgt[i] / _binv[i] : _wgt[i];
1001 }
1002 return ctx.buildArg(vals);
1003}
1004
1006 const RooAbsCollection &coords, bool reverse) const
1007{
1008 assert(coords.size() == _vars.size());
1009
1010 std::string code;
1011 int idxMult = 1;
1012
1013 for (std::size_t i = 0; i < _vars.size(); ++i) {
1014
1015 std::size_t iVar = reverse ? _vars.size() - 1 - i : i;
1016 const RooAbsArg *internalVar = _vars[iVar];
1017 const RooAbsArg *theVar = coords[iVar];
1018
1019 const RooAbsBinning *binning = _lvbins[iVar].get();
1020 if (!binning) {
1021 coutE(InputArguments) << "RooHistPdf::weight(" << GetName()
1022 << ") ERROR: Code Squashing currently does not support category values." << std::endl;
1023 return "";
1024 }
1025
1026 code += " + " + binning->translateBinNumber(ctx, *theVar, idxMult);
1027
1028 // Use RooAbsLValue here because it also generalized to categories, which
1029 // is useful in the future. dynamic_cast because it's a cross-cast.
1030 idxMult *= dynamic_cast<RooAbsLValue const *>(internalVar)->numBins();
1031 }
1032
1033 return "(" + code + ")";
1034}
1035
1036////////////////////////////////////////////////////////////////////////////////
1037/// Calculate the bin index corresponding to the coordinates passed as argument.
1038/// \param[in] coords Coordinates. If `fast == false`, these can be partial.
1039/// \param[in] fast Promise that the coordinates in `coords` have the same order
1040/// as the internal coordinates. In this case, values are looked up only by index.
1041std::size_t RooDataHist::calcTreeIndex(const RooAbsCollection& coords, bool fast) const
1042{
1043 // With fast, caller promises that layout of `coords` is identical to our internal `vars`.
1044 // Previously, this was verified with an assert in debug mode like this:
1045 //
1046 // assert(!fast || coords.hasSameLayout(_vars));
1047 //
1048 // However, there are usecases where the externally provided `coords` have
1049 // different names than the internal variables, even though they correspond
1050 // to each other. For example, if the observables in the computation graph
1051 // are renamed with `redirectServers`. Hence, we can't do a meaningful assert
1052 // here.
1053
1054 if (&_vars == &coords)
1055 fast = true;
1056
1057 std::size_t masterIdx = 0;
1058
1059 for (unsigned int i=0; i < _vars.size(); ++i) {
1060 const RooAbsArg* internalVar = _vars[i];
1061 const RooAbsBinning* binning = _lvbins[i].get();
1062
1063 // Find the variable that we need values from.
1064 // That's either the variable directly from the external coordinates
1065 // or we find the external one that has the same name as "internalVar".
1066 const RooAbsArg* theVar = fast ? coords[i] : coords.find(*internalVar);
1067 if (!theVar) {
1068 // Variable is not in external coordinates. Use current internal value.
1070 }
1071 // If fast is on, users promise that the sets have the same layout:
1072 //
1073 // assert(!fast || strcmp(internalVar->GetName(), theVar->GetName()) == 0);
1074 //
1075 // This assert is commented out for the same reasons that applied to the
1076 // other assert explained above.
1077
1078 if (binning) {
1079 assert(dynamic_cast<const RooAbsReal*>(theVar));
1080 const double val = static_cast<const RooAbsReal*>(theVar)->getVal();
1081 masterIdx += _idxMult[i] * binning->binNumber(val);
1082 } else {
1083 // We are a category. No binning.
1084 assert(dynamic_cast<const RooAbsCategoryLValue*>(theVar));
1085 auto cat = static_cast<const RooAbsCategoryLValue*>(theVar);
1086 masterIdx += _idxMult[i] * cat->getBin(static_cast<const char*>(nullptr));
1087 }
1088 }
1089
1090 return masterIdx ;
1091}
1092
1093
1094////////////////////////////////////////////////////////////////////////////////
1095/// Back end function to plotting functionality. Plot RooDataHist on given
1096/// frame in mode specified by plot options 'o'. The main purpose of
1097/// this function is to match the specified binning on 'o' to the
1098/// internal binning of the plot observable in this RooDataHist.
1099/// \see RooAbsData::plotOn() for plotting options.
1101{
1102 checkInit() ;
1103 if (o.bins) return RooAbsData::plotOn(frame,o) ;
1104
1105 if(!frame) {
1106 coutE(InputArguments) << ClassName() << "::" << GetName() << ":plotOn: frame is null" << endl;
1107 return nullptr;
1108 }
1109 auto var= static_cast<RooAbsRealLValue*>(frame->getPlotVar());
1110 if(!var) {
1111 coutE(InputArguments) << ClassName() << "::" << GetName()
1112 << ":plotOn: frame does not specify a plot variable" << endl;
1113 return nullptr;
1114 }
1115
1116 auto dataVar = static_cast<RooRealVar*>(_vars.find(*var));
1117 if (!dataVar) {
1118 coutE(InputArguments) << ClassName() << "::" << GetName()
1119 << ":plotOn: dataset doesn't contain plot frame variable" << endl;
1120 return nullptr;
1121 }
1122
1123 o.bins = &dataVar->getBinning() ;
1124 return RooAbsData::plotOn(frame,o) ;
1125}
1126
1127
1128////////////////////////////////////////////////////////////////////////////////
1129/// A vectorized version of interpolateDim for boundary safe quadratic
1130/// interpolation of one dimensional histograms.
1131///
1132/// \param[out] output An array of interpolated weights corresponding to the
1133/// values in xVals.
1134/// \param[in] xVals An array of event coordinates for which the weights should be
1135/// calculated.
1136/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1137/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1138/// Underflow bins are assumed to have weight zero and
1139/// overflow bins have weight one. Otherwise, the
1140/// histogram is mirrored at the boundaries for the
1141/// interpolation.
1142
1143void RooDataHist::interpolateQuadratic(double* output, std::span<const double> xVals,
1145{
1146 const std::size_t nBins = numEntries();
1147 const std::size_t nEvents = xVals.size();
1148
1149 RooAbsBinning const& binning = *_lvbins[0];
1150 // Reuse the output buffer for bin indices and zero-initialize it
1151 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1152 std::fill(binIndices, binIndices + nEvents, 0);
1153 binning.binNumbers(xVals.data(), binIndices, nEvents);
1154
1155 // Extend coordinates and weights with one extra point before the first bin
1156 // and one extra point after the last bin. This means the original histogram
1157 // bins span elements 1 to nBins in coordsExt and weightsExt
1158 std::vector<double> coordsExt(nBins+3);
1159 double* binCoords = coordsExt.data() + 2;
1160 binCoords[0] = binning.lowBound() + 0.5*_binv[0];
1161 for (std::size_t binIdx = 1; binIdx < nBins ; ++binIdx) {
1162 if (binning.isUniform()) {
1163 double binWidth = _binv[0];
1164 binCoords[binIdx] = binIdx*binWidth + binCoords[0];
1165 }
1166 else {
1167 double binCentDiff = 0.5*_binv[binIdx-1] + 0.5*_binv[binIdx];
1169 }
1170 }
1171
1172 std::vector<double> weightsExt(nBins+3);
1173 // Fill weights for bins that are inside histogram boundaries
1174 for (std::size_t binIdx = 0; binIdx < nBins; ++binIdx) {
1176 }
1177
1178 if (cdfBoundaries) {
1179 coordsExt[0] = - 1e-10 + binning.lowBound();
1180 weightsExt[0] = 0.;
1181
1182 coordsExt[1] = binning.lowBound();
1183 weightsExt[1] = 0.;
1184
1185 coordsExt[nBins+2] = binning.highBound();
1186 weightsExt[nBins+2] = 1.;
1187 }
1188 else {
1189 // Mirror first two bins and last bin
1190 coordsExt[0] = binCoords[1] - 2*_binv[0] - _binv[1];
1191 weightsExt[0] = weightsExt[3];
1192
1193 coordsExt[1] = binCoords[0] - _binv[0];
1194 weightsExt[1] = weightsExt[2];
1195
1196 coordsExt[nBins+2] = binCoords[nBins-1] + _binv[nBins-1];
1197 weightsExt[nBins+2] = weightsExt[nBins+1];
1198 }
1199
1200 // We use the current bin center and two bin centers on the left for
1201 // interpolation if xVal is to the left of the current bin center
1202 for (std::size_t i = 0; i < nEvents ; ++i) {
1203 double xVal = xVals[i];
1204 std::size_t binIdx = binIndices[i] + 2;
1205
1206 // If xVal is to the right of the current bin center, shift all bin
1207 // coordinates one step to the right and use that for the interpolation
1208 if (xVal > coordsExt[binIdx]) {
1209 binIdx += 1;
1210 }
1211
1212 double x1 = coordsExt[binIdx-2];
1213 double y1 = weightsExt[binIdx-2];
1214
1215 double x2 = coordsExt[binIdx-1];
1216 double y2 = weightsExt[binIdx-1];
1217
1218 double x3 = coordsExt[binIdx];
1219 double y3 = weightsExt[binIdx];
1220
1221 // Evaluate a few repeated factors
1222 double quotient = (x3-x1) / (x2-x1);
1223 double x1Sqrd = x1*x1;
1224 double x3Sqrd = x3*x3;
1225 // Solve coefficients in system of three quadratic equations!
1226 double secondCoeff = (y3 - y1 - (y2-y1) * quotient) / (x3Sqrd - x1Sqrd - (x2*x2 - x1Sqrd) * quotient);
1227 double firstCoeff = (y3 - y1 - secondCoeff*(x3Sqrd - x1Sqrd)) / (x3-x1);
1229 // Get the interpolated weight using the equation of a second degree polynomial
1231 }
1232}
1233
1234
1235////////////////////////////////////////////////////////////////////////////////
1236/// A vectorized version of interpolateDim for boundary safe linear
1237/// interpolation of one dimensional histograms.
1238///
1239/// \param[out] output An array of interpolated weights corresponding to the
1240/// values in xVals.
1241/// \param[in] xVals An array of event coordinates for which the weights should be
1242/// calculated.
1243/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1244/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1245/// Underflow bins are assumed to have weight zero and
1246/// overflow bins have weight one. Otherwise, the
1247/// histogram is mirrored at the boundaries for the
1248/// interpolation.
1249
1250void RooDataHist::interpolateLinear(double* output, std::span<const double> xVals,
1252{
1253 const std::size_t nBins = numEntries();
1254 const std::size_t nEvents = xVals.size();
1255
1256 RooAbsBinning const& binning = *_lvbins[0];
1257 // Reuse the output buffer for bin indices and zero-initialize it
1258 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1259 std::fill(binIndices, binIndices + nEvents, 0);
1260 binning.binNumbers(xVals.data(), binIndices, nEvents);
1261
1262 // Extend coordinates and weights with one extra point before the first bin
1263 // and one extra point after the last bin. This means the original histogram
1264 // bins span elements 1 to nBins in coordsExt and weightsExt
1265 std::vector<double> coordsExt(nBins+2);
1266 double* binCoords = coordsExt.data() + 1;
1267 binCoords[0] = binning.lowBound() + 0.5*_binv[0];
1268 for (std::size_t binIdx = 1; binIdx < nBins ; ++binIdx) {
1269 if (binning.isUniform()) {
1270 double binWidth = _binv[0];
1271 binCoords[binIdx] = binIdx*binWidth + binCoords[0];
1272 }
1273 else {
1274 double binCentDiff = 0.5*_binv[binIdx-1] + 0.5*_binv[binIdx];
1276 }
1277 }
1278
1279 std::vector<double> weightsExt(nBins+2);
1280 // Fill weights for bins that are inside histogram boundaries
1281 for (std::size_t binIdx = 0; binIdx < nBins; ++binIdx) {
1283 }
1284
1285 // Fill weights for bins that are outside histogram boundaries
1286 if (cdfBoundaries) {
1287 coordsExt[0] = binning.lowBound();
1288 weightsExt[0] = 0.;
1289 coordsExt[nBins+1] = binning.highBound();
1290 weightsExt[nBins+1] = 1.;
1291 }
1292 else {
1293 // Mirror first and last bins
1294 coordsExt[0] = binCoords[0] - _binv[0];
1295 weightsExt[0] = weightsExt[1];
1296 coordsExt[nBins+1] = binCoords[nBins-1] + _binv[nBins-1];
1297 weightsExt[nBins+1] = weightsExt[nBins];
1298 }
1299
1300 // Interpolate between current bin center and one bin center to the left
1301 // if xVal is to the left of the current bin center
1302 for (std::size_t i = 0; i < nEvents ; ++i) {
1303 double xVal = xVals[i];
1304 std::size_t binIdx = binIndices[i] + 1;
1305
1306 // If xVal is to the right of the current bin center, interpolate between
1307 // current bin center and one bin center to the right instead
1308 if (xVal > coordsExt[binIdx]) { binIdx += 1; }
1309
1310 double x1 = coordsExt[binIdx-1];
1311 double y1 = weightsExt[binIdx-1];
1312 double x2 = coordsExt[binIdx];
1313 double y2 = weightsExt[binIdx];
1314
1315 // Find coefficients by solving a system of two linear equations
1316 double firstCoeff = (y2-y1) / (x2-x1);
1317 double zerothCoeff = y1 - firstCoeff * x1;
1318 // Get the interpolated weight using the equation of a straight line
1320 }
1321}
1322
1323
1324////////////////////////////////////////////////////////////////////////////////
1325/// A vectorized version of RooDataHist::weight() for one dimensional histograms
1326/// with up to one dimensional interpolation.
1327/// \param[out] output An array of weights corresponding the values in xVals.
1328/// \param[in] xVals An array of coordinates for which the weights should be
1329/// calculated.
1330/// \param[in] intOrder Interpolation order; 0th and 1st order are supported.
1331/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1332/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1333/// Underflow bins are assumed to have weight zero and
1334/// overflow bins have weight one. Otherwise, the
1335/// histogram is mirrored at the boundaries for the
1336/// interpolation.
1337
1338void RooDataHist::weights(double* output, std::span<double const> xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
1339{
1340 auto const nEvents = xVals.size();
1341
1342 if (intOrder == 0) {
1343 RooAbsBinning const& binning = *_lvbins[0];
1344
1345 // Reuse the output buffer for bin indices and zero-initialize it
1346 auto binIndices = reinterpret_cast<int*>(output + nEvents) - nEvents;
1347 std::fill(binIndices, binIndices + nEvents, 0);
1348 binning.binNumbers(xVals.data(), binIndices, nEvents);
1349
1350 for (std::size_t i=0; i < nEvents; ++i) {
1351 auto binIdx = binIndices[i];
1353 }
1354 }
1355 else if (intOrder == 1) {
1357 }
1358 else if (intOrder == 2) {
1360 }
1361 else {
1362 // Higher dimensional scenarios not yet implemented
1363 coutE(InputArguments) << "RooDataHist::weights(" << GetName() << ") interpolation in "
1364 << intOrder << " dimensions not yet implemented" << std::endl ;
1365 // Fall back to 1st order interpolation
1367 }
1368}
1369
1370
1371////////////////////////////////////////////////////////////////////////////////
1372/// A faster version of RooDataHist::weight that assumes the passed arguments
1373/// are aligned with the histogram variables.
1374/// \param[in] bin Coordinates for which the weight should be calculated.
1375/// Has to be aligned with the internal histogram variables.
1376/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1377/// used for the interpolation. If zero, the bare weight for
1378/// the bin enclosing the coordinatesis returned.
1379/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1380/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1381/// underflow bins are assumed to have weight zero and
1382/// overflow bins have weight one. Otherwise, the
1383/// histogram is mirrored at the boundaries for the
1384/// interpolation.
1385
1387{
1388 checkInit() ;
1389
1390 // Handle illegal intOrder values
1391 if (intOrder<0) {
1392 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << endl ;
1393 return 0 ;
1394 }
1395
1396 // Handle no-interpolation case
1397 if (intOrder==0) {
1398 const auto idx = calcTreeIndex(bin, true);
1399 return correctForBinSize ? _wgt[idx] / _binv[idx] : _wgt[idx];
1400 }
1401
1402 // Handle all interpolation cases
1404}
1405
1406
1407////////////////////////////////////////////////////////////////////////////////
1408/// Return the weight at given coordinates with optional interpolation.
1409/// \param[in] bin Coordinates for which the weight should be calculated.
1410/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1411/// used for the interpolation. If zero, the bare weight for
1412/// the bin enclosing the coordinatesis returned.
1413/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1414/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1415/// underflow bins are assumed to have weight zero and
1416/// overflow bins have weight one. Otherwise, the
1417/// histogram is mirrored at the boundaries for the
1418/// interpolation.
1419/// \param[in] oneSafe Ignored.
1420
1421double RooDataHist::weight(const RooArgSet& bin, Int_t intOrder, bool correctForBinSize, bool cdfBoundaries, bool /*oneSafe*/)
1422{
1423 checkInit() ;
1424
1425 // Handle illegal intOrder values
1426 if (intOrder<0) {
1427 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") ERROR: interpolation order must be positive" << endl ;
1428 return 0 ;
1429 }
1430
1431 // Handle no-interpolation case
1432 if (intOrder==0) {
1433 const auto idx = calcTreeIndex(bin, false);
1434 return correctForBinSize ? _wgt[idx] / _binv[idx] : _wgt[idx];
1435 }
1436
1437 // Handle all interpolation cases
1438 _vars.assignValueOnly(bin) ;
1439
1441}
1442
1443
1444////////////////////////////////////////////////////////////////////////////////
1445/// Return the weight at given coordinates with interpolation.
1446/// \param[in] bin Coordinates for which the weight should be calculated.
1447/// Has to be aligned with the internal histogram variables.
1448/// \param[in] intOrder Interpolation order, i.e. how many neighbouring bins are
1449/// used for the interpolation.
1450/// \param[in] correctForBinSize Enable the inverse bin volume correction factor.
1451/// \param[in] cdfBoundaries Enable the special boundary condition for a cdf:
1452/// underflow bins are assumed to have weight zero and
1453/// overflow bins have weight one. Otherwise, the
1454/// histogram is mirrored at the boundaries for the
1455/// interpolation.
1456
1458 VarInfo const& varInfo = getVarInfo();
1459
1460 const auto centralIdx = calcTreeIndex(bin, true);
1461
1462 double wInt{0} ;
1463 if (varInfo.nRealVars == 1) {
1464
1465 // buffer needs to be 2 x (interpolation order + 1), with the factor 2 for x and y.
1466 _interpolationBuffer.resize(2 * intOrder + 2);
1467
1468 // 1-dimensional interpolation
1469 auto const& realX = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx1]);
1471
1472 } else if (varInfo.nRealVars == 2) {
1473
1474 // buffer needs to be 2 x 2 x (interpolation order + 1), with one factor 2
1475 // for x and y, and the other for the number of dimensions.
1476 _interpolationBuffer.resize(4 * intOrder + 4);
1477
1478 // 2-dimensional interpolation
1479 auto const& realX = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx1]);
1480 auto const& realY = static_cast<RooRealVar const&>(*bin[varInfo.realVarIdx2]);
1481 double xval = realX.getVal() ;
1482 double yval = realY.getVal() ;
1483
1484 RooAbsBinning const& binningY = realY.getBinning();
1485
1486 int ybinC = binningY.binNumber(yval) ;
1487 int ybinLo = ybinC-intOrder/2 - ((yval<binningY.binCenter(ybinC))?1:0) ;
1488 int ybinM = binningY.numBins() ;
1489
1490 auto idxMultY = _idxMult[varInfo.realVarIdx2];
1492
1493 // Use a class-member buffer to avoid repeated heap allocations.
1494 double * yarr = _interpolationBuffer.data() + 2 * intOrder + 2; // add offset to skip part reserved for other dim
1495 double * xarr = yarr + intOrder + 1;
1496 for (int i=ybinLo ; i<=intOrder+ybinLo ; i++) {
1497 int ibin ;
1498 if (i>=0 && i<ybinM) {
1499 // In range
1500 ibin = i ;
1501 xarr[i-ybinLo] = binningY.binCenter(ibin) ;
1502 } else if (i>=ybinM) {
1503 // Overflow: mirror
1504 ibin = 2*ybinM-i-1 ;
1505 xarr[i-ybinLo] = 2*binningY.highBound()-binningY.binCenter(ibin) ;
1506 } else {
1507 // Underflow: mirror
1508 ibin = -i -1;
1509 xarr[i-ybinLo] = 2*binningY.lowBound()-binningY.binCenter(ibin) ;
1510 }
1513 }
1514
1515 if (gDebug>7) {
1516 cout << "RooDataHist interpolating data is" << endl ;
1517 cout << "xarr = " ;
1518 for (int q=0; q<=intOrder ; q++) cout << xarr[q] << " " ;
1519 cout << " yarr = " ;
1520 for (int q=0; q<=intOrder ; q++) cout << yarr[q] << " " ;
1521 cout << endl ;
1522 }
1524
1525 } else {
1526
1527 // Higher dimensional scenarios not yet implemented
1528 coutE(InputArguments) << "RooDataHist::weight(" << GetName() << ") interpolation in "
1529 << varInfo.nRealVars << " dimensions not yet implemented" << endl ;
1531
1532 }
1533
1534 return wInt ;
1535}
1536
1537
1539 if (!_errLo || !_errHi) {
1540 initArray(_errLo, _arrSize, -1.);
1541 initArray(_errHi, _arrSize, -1.);
1543 }
1544}
1545
1546
1547////////////////////////////////////////////////////////////////////////////////
1548/// Return the asymmetric errors on the current weight.
1549/// \see weightError(ErrorType) const for symmetric error.
1550/// \param[out] lo Low error.
1551/// \param[out] hi High error.
1552/// \param[in] etype Type of error to compute. May throw if not supported.
1553/// Supported errors are
1554/// - `Poisson` Default. Asymmetric Poisson errors (68% CL).
1555/// - `SumW2` The square root of the sum of weights. (Symmetric).
1556/// - `None` Return zero.
1557void RooDataHist::weightError(double& lo, double& hi, ErrorType etype) const
1558{
1559 checkInit() ;
1560
1561 switch (etype) {
1562
1563 case Auto:
1564 throw std::invalid_argument("RooDataHist::weightError(" + std::string(GetName()) + ") error type Auto not allowed here");
1565 break ;
1566
1567 case Expected:
1568 throw std::invalid_argument("RooDataHist::weightError(" + std::string(GetName()) + ") error type Expected not allowed here");
1569 break ;
1570
1571 case Poisson:
1572 if (_errLo && _errLo[_curIndex] >= 0.0) {
1573 // Weight is preset or precalculated
1574 lo = _errLo[_curIndex];
1575 hi = _errHi[_curIndex];
1576 return ;
1577 }
1578
1579 // We didn't track asymmetric errors so far, so now we need to allocate
1581
1582 // Calculate poisson errors
1583 double ym;
1584 double yp;
1585 RooHistError::instance().getPoissonInterval(Int_t(weight()+0.5),ym,yp,1) ;
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.
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];
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];
1674 yarr[i - fbinLo] /= _binv[idx];
1675 }
1676 }
1677 }
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:
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
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
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
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 = static_cast<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
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.
1841{
1842 checkInit() ;
1843
1844 // Check if result was cached
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
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
1877{
1878 checkInit() ;
1879
1881 varSave.addClone(_vars) ;
1882
1884 sliceOnlySet.remove(sumSet,true,true) ;
1885
1887 std::vector<double> const * pbinv = nullptr;
1888
1891 } else if(correctForBinSize && !inverseBinCor) {
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 ;
1929 }
1930 }
1931
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
1955 const std::map<const RooAbsArg*, std::pair<double, double> >& ranges,
1956 std::function<double(int)> getBinScale)
1957{
1958 checkInit();
1961 varSave.addClone(_vars);
1962 {
1964 sliceOnlySet.remove(sumSet, true, true);
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
2029 break;
2030 }
2031
2033 binVolumeSumSetInRange *= std::min(rangeHi[ivar], binHi) - std::max(rangeLo[ivar], binLo);
2034 }
2036 if (0. == corrPartial) continue;
2039 }
2040
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);
2080 Int_t tmp(ibin);
2081 double theBinVolume(1) ;
2082 for (unsigned int j=0; j < _lvvars.size(); ++j) {
2083 const RooAbsLValue* arg = _lvvars[j];
2084 assert(arg);
2085
2086 idx = tmp / _idxMult[j];
2087 tmp -= idx*_idxMult[j];
2088 if (selDim[j]) {
2089 theBinVolume *= arg->getBinWidth(idx) ;
2090 }
2091 }
2093 }
2094
2095 return pbinv;
2096}
2097
2098
2099////////////////////////////////////////////////////////////////////////////////
2100/// Sum the weights of all bins.
2104
2105
2106
2107////////////////////////////////////////////////////////////////////////////////
2108/// Return the sum of weights in all entries matching cutSpec (if specified)
2109/// and in named range cutRange (if specified)
2110/// Return the
2111
2112double RooDataHist::sumEntries(const char* cutSpec, const char* cutRange) const
2113{
2114 checkInit() ;
2115
2116 if (cutSpec==nullptr && cutRange==nullptr) {
2117 return sumEntries();
2118 } else {
2119
2120 // Setup RooFormulaVar for cutSpec if it is present
2121 std::unique_ptr<RooFormula> select;
2122 if (cutSpec) {
2123 select = std::make_unique<RooFormula>("select",cutSpec,*get());
2124 }
2125
2126 // Otherwise sum the weights in the event
2127 ROOT::Math::KahanSum<> kahanSum;
2128 for (Int_t i=0; i < _arrSize; i++) {
2129 get(i) ;
2130 if ((select && select->eval() == 0.) || (cutRange && !_vars.allInRange(cutRange)))
2131 continue;
2132
2133 kahanSum += weight(i);
2134 }
2135
2136 return kahanSum.Sum();
2137 }
2138}
2139
2140
2141
2142////////////////////////////////////////////////////////////////////////////////
2143/// Reset all bin weights to zero
2144
2146{
2147 // WVE DO NOT CALL RooTreeData::reset() for binned
2148 // datasets as this will delete the bin definitions
2149
2150 std::fill(_wgt, _wgt + _arrSize, 0.);
2151 delete[] _errLo; _errLo = nullptr;
2152 delete[] _errHi; _errHi = nullptr;
2153 delete[] _sumw2; _sumw2 = nullptr;
2154
2156
2157 _cache_sum_valid = false;
2158}
2159
2160
2161
2162////////////////////////////////////////////////////////////////////////////////
2163/// Load bin `binNumber`, and return an argset with the coordinates of the bin centre.
2164/// \note The argset is owned by this data hist, and this function has a side effect, because
2165/// it alters the currently active bin.
2166const RooArgSet* RooDataHist::get(Int_t binNumber) const
2167{
2168 checkInit() ;
2169 _curIndex = binNumber;
2170
2171 return RooAbsData::get(_curIndex);
2172}
2173
2174
2175
2176////////////////////////////////////////////////////////////////////////////////
2177/// Return a RooArgSet with whose coordinates denote the bin centre of the bin
2178/// enclosing the point in `coord`.
2179/// \note The argset is owned by this data hist, and this function has a side effect, because
2180/// it alters the currently active bin.
2182 return get(calcTreeIndex(coord, false));
2183}
2184
2185
2186
2187////////////////////////////////////////////////////////////////////////////////
2188/// Return the volume of the bin enclosing coordinates 'coord'.
2190 checkInit() ;
2191 return _binv[calcTreeIndex(coord, false)] ;
2192}
2193
2194
2195////////////////////////////////////////////////////////////////////////////////
2196/// Create an iterator over all bins in a slice defined by the subset of observables
2197/// listed in sliceArg. The position of the slice is given by otherArgs
2198
2200{
2201 // Update to current position
2203 _curIndex = calcTreeIndex(_vars, true);
2204
2206 if (!intArg) {
2207 coutE(InputArguments) << "RooDataHist::sliceIterator() variable " << sliceArg.GetName() << " is not part of this RooDataHist" << endl ;
2208 return nullptr ;
2209 }
2210 return new RooDataHistSliceIter(*this,*intArg) ;
2211}
2212
2213
2214////////////////////////////////////////////////////////////////////////////////
2215/// Change the name of the RooDataHist
2216
2218{
2219 if (_dir) _dir->GetList()->Remove(this);
2220 // We need to use the function from RooAbsData, because it already overrides TNamed::SetName
2222 if (_dir) _dir->GetList()->Add(this);
2223}
2224
2225
2226////////////////////////////////////////////////////////////////////////////////
2227/// Change the title of this RooDataHist
2228
2229void RooDataHist::SetNameTitle(const char *name, const char* title)
2230{
2231 SetName(name);
2232 SetTitle(title);
2233}
2234
2235
2236////////////////////////////////////////////////////////////////////////////////
2237/// Print value of the dataset, i.e. the sum of weights contained in the dataset
2238
2239void RooDataHist::printValue(ostream& os) const
2240{
2241 os << numEntries() << " bins (" << sumEntries() << " weights)" ;
2242}
2243
2244
2245
2246
2247////////////////////////////////////////////////////////////////////////////////
2248/// Print argument of dataset, i.e. the observable names
2249
2250void RooDataHist::printArgs(ostream& os) const
2251{
2252 os << "[" ;
2253 bool first(true) ;
2254 for (const auto arg : _vars) {
2255 if (first) {
2256 first=false ;
2257 } else {
2258 os << "," ;
2259 }
2260 os << arg->GetName() ;
2261 }
2262 os << "]" ;
2263}
2264
2265
2266
2267////////////////////////////////////////////////////////////////////////////////
2268/// Returns true if dataset contains entries with a non-integer weight.
2269
2271{
2272 for (Int_t i=0; i < _arrSize; ++i) {
2273 const double wgt = _wgt[i];
2274 double intpart;
2275 if (std::abs(std::modf(wgt, &intpart)) > 1.E-10)
2276 return true;
2277 }
2278
2279 return false;
2280}
2281
2282
2283////////////////////////////////////////////////////////////////////////////////
2284/// Print the details on the dataset contents
2285
2286void RooDataHist::printMultiline(ostream& os, Int_t content, bool verbose, TString indent) const
2287{
2289
2290 os << indent << "Binned Dataset " << GetName() << " (" << GetTitle() << ")" << endl ;
2291 os << indent << " Contains " << numEntries() << " bins with a total weight of " << sumEntries() << endl;
2292
2293 if (!verbose) {
2294 os << indent << " Observables " << _vars << endl ;
2295 } else {
2296 os << indent << " Observables: " ;
2298 }
2299
2300 if(verbose) {
2301 if (!_cachedVars.empty()) {
2302 os << indent << " Caches " << _cachedVars << endl ;
2303 }
2304 }
2305}
2306
2307void RooDataHist::printDataHistogram(ostream& os, RooRealVar* obs) const
2308{
2309 for(Int_t i=0; i<obs->getBins(); ++i){
2310 this->get(i);
2311 obs->setBin(i);
2312 os << this->weight() << " +/- " << this->weightSquared() << endl;
2313 }
2314}
2315
2316
2317////////////////////////////////////////////////////////////////////////////////
2318/// Stream an object of class RooDataHist.
2320 if (R__b.IsReading()) {
2321
2322 UInt_t R__s;
2323 UInt_t R__c;
2324 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2325
2326 if (R__v > 2) {
2327 R__b.ReadClassBuffer(RooDataHist::Class(),this,R__v,R__s,R__c);
2328 R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());
2329 initialize(nullptr, false);
2330 } else {
2331
2332 // Legacy dataset conversion happens here. Legacy RooDataHist inherits from RooTreeData
2333 // which in turn inherits from RooAbsData. Manually stream RooTreeData contents on
2334 // file here and convert it into a RooTreeDataStore which is installed in the
2335 // new-style RooAbsData base class
2336
2337 // --- This is the contents of the streamer code of RooTreeData version 2 ---
2338 UInt_t R__s1;
2339 UInt_t R__c1;
2340 Version_t R__v1 = R__b.ReadVersion(&R__s1, &R__c1); if (R__v1) { }
2341
2343 TTree* X_tree(nullptr) ; R__b >> X_tree;
2344 RooArgSet X_truth ; X_truth.Streamer(R__b);
2346 R__b.CheckByteCount(R__s1, R__c1, TClass::GetClass("RooTreeData"));
2347 // --- End of RooTreeData-v1 streamer
2348
2349 // Construct RooTreeDataStore from X_tree and complete initialization of new-style RooAbsData
2350 _dstore = std::make_unique<RooTreeDataStore>(X_tree,_vars);
2351 _dstore->SetName(GetName()) ;
2352 _dstore->SetTitle(GetTitle()) ;
2353 _dstore->checkInit() ;
2354
2356 R__b >> _arrSize;
2357 delete [] _wgt;
2358 _wgt = new double[_arrSize];
2359 R__b.ReadFastArray(_wgt,_arrSize);
2360 delete [] _errLo;
2361 _errLo = new double[_arrSize];
2362 R__b.ReadFastArray(_errLo,_arrSize);
2363 delete [] _errHi;
2364 _errHi = new double[_arrSize];
2365 R__b.ReadFastArray(_errHi,_arrSize);
2366 delete [] _sumw2;
2367 _sumw2 = new double[_arrSize];
2368 R__b.ReadFastArray(_sumw2,_arrSize);
2369 delete [] _binv;
2370 _binv = new double[_arrSize];
2372 tmpSet.Streamer(R__b);
2373 double tmp;
2374 R__b >> tmp; //_curWeight;
2375 R__b >> tmp; //_curWgtErrLo;
2376 R__b >> tmp; //_curWgtErrHi;
2377 R__b >> tmp; //_curSumW2;
2378 R__b >> tmp; //_curVolume;
2379 R__b >> _curIndex;
2380 R__b.CheckByteCount(R__s, R__c, RooDataHist::IsA());
2381 }
2382
2383 } else {
2384
2385 R__b.WriteClassBuffer(RooDataHist::Class(),this);
2386 }
2387}
2388
2389
2390////////////////////////////////////////////////////////////////////////////////
2391/// Return event weights of all events in range [first, first+len).
2392/// If cacheValidEntries() has been called, out-of-range events will have a weight of 0.
2393std::span<const double> RooDataHist::getWeightBatch(std::size_t first, std::size_t len, bool sumW2 /*=false*/) const {
2394 return {(sumW2 && _sumw2 ? _sumw2 : _wgt) + first, len};
2395}
2396
2397
2398////////////////////////////////////////////////////////////////////////////////
2399/// Hand over pointers to our weight arrays to the data store implementation.
2401 _dstore->setExternalWeightArray(_wgt, _errLo, _errHi, _sumw2);
2402}
2403
2404
2405////////////////////////////////////////////////////////////////////////////////
2406/// Return reference to VarInfo struct with cached histogram variable
2407/// information that is frequently used for histogram weights retrieval.
2408///
2409/// If the `_varInfo` struct was not initialized yet, it will be initialized in
2410/// this function.
2412
2413 if(_varInfo.initialized) return _varInfo;
2414
2415 auto& info = _varInfo;
2416
2417 {
2418 // count the number of real vars and get their indices
2419 info.nRealVars = 0;
2420 size_t iVar = 0;
2421 for (const auto real : _vars) {
2422 if (dynamic_cast<RooRealVar*>(real)) {
2423 if(info.nRealVars == 0) info.realVarIdx1 = iVar;
2424 if(info.nRealVars == 1) info.realVarIdx2 = iVar;
2425 ++info.nRealVars;
2426 }
2427 ++iVar;
2428 }
2429 }
2430
2431 {
2432 // assert that the variables are either real values or categories
2433 for (unsigned int i=0; i < _vars.size(); ++i) {
2434 if (_lvbins[i].get()) {
2435 assert(dynamic_cast<const RooAbsReal*>(_vars[i]));
2436 } else {
2437 assert(dynamic_cast<const RooAbsCategoryLValue*>(_vars[i]));
2438 }
2439 }
2440 }
2441
2442 info.initialized = true;
2443
2444 return info;
2445}
#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:382
static void indent(ostringstream &buf, int indent_level)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
static unsigned int total
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t mask
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
float xmin
#define hi
float * q
float ymin
Int_t gDebug
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
const_iterator begin() const
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,...
Abstract base class for RooRealVar binning definitions.
int binNumber(double x) const
Returns the bin number corresponding to the value x.
virtual void binNumbers(double const *x, int *bins, std::size_t n, int coef=1) const =0
Compute the bin indices for multiple values of x.
virtual bool isUniform() const
virtual double highBound() const =0
virtual double lowBound() const =0
virtual std::string translateBinNumber(RooFit::Detail::CodeSquashContext &ctx, RooAbsArg const &var, int coef) const
virtual RooAbsBinning * clone(const char *name=nullptr) const =0
Abstract base class for objects that represent a discrete value that can be set from the outside,...
bool hasLabel(const std::string &label) const
Check if a state with name label exists.
Abstract container object that can hold multiple RooAbsArg objects.
RooAbsCollection & assignValueOnly(const RooAbsCollection &other, bool forceIfSizeOne=false)
Sets the value of any argument in our set that also appears in the other set.
bool allInRange(const char *rangeSpec) const
Return true if all contained object report to have their value inside the specified range.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
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.
void Streamer(TBuffer &) override
Stream an object of class RooAbsData.
Abstract base class for objects that are lvalues, i.e.
virtual double getBinWidth(Int_t i, const char *rangeName=nullptr) const =0
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
virtual Int_t getBins(const char *name=nullptr) const
Get number of bins 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...
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:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:159
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
Object to represent discrete states.
Definition RooCategory.h:28
bool defineType(const std::string &label)
Define a state with given name.
Named container for two doubles, two integers two object points and three string pointers that can be...
Definition RooCmdArg.h:26
Configurable parser for RooCmdArg named arguments.
void defineMutex(const char *head, Args_t &&... tail)
Define arguments where any pair is mutually exclusive.
bool process(const RooCmdArg &arg)
Process given RooCmdArg.
double getDouble(const char *name, double defaultValue=0.0) const
Return double property registered with name 'name'.
void defineDependency(const char *refArgName, const char *neededArgName)
Define that processing argument name refArgName requires processing of argument named neededArgName t...
bool defineDouble(const char *name, const char *argName, int doubleNum, double defValue=0.0)
Define double property name 'name' mapped to double in slot 'doubleNum' in RooCmdArg with name argNam...
RooArgSet * getSet(const char *name, RooArgSet *set=nullptr) const
Return RooArgSet property registered with name 'name'.
bool defineSet(const char *name, const char *argName, int setNum, const RooArgSet *set=nullptr)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
bool ok(bool verbose) const
Return true of parsing was successful.
bool defineObject(const char *name, const char *argName, int setNum, const TObject *obj=nullptr, bool isArray=false)
Define TObject property name 'name' mapped to object in slot 'setNum' in RooCmdArg with name argName ...
const char * getString(const char *name, const char *defaultValue="", bool convEmptyToNull=false) const
Return string property registered with name 'name'.
bool defineString(const char *name, const char *argName, int stringNum, const char *defValue="", bool appendMode=false)
Define double property name 'name' mapped to double in slot 'stringNum' in RooCmdArg with name argNam...
const RooLinkedList & getObjectList(const char *name) const
Return list of objects registered with name 'name'.
bool defineInt(const char *name, const char *argName, int intNum, int defValue=0)
Define integer property name 'name' mapped to integer in slot 'intNum' in RooCmdArg with name argName...
int getInt(const char *name, int defaultValue=0) const
Return integer property registered with name 'name'.
TObject * getObject(const char *name, TObject *obj=nullptr) const
Return TObject property registered with name 'name'.
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
std::span< const double > getWeightBatch(std::size_t first, std::size_t len, bool sumW2=false) const override
Return event weights of all events in range [first, first+len).
void interpolateQuadratic(double *output, std::span< const double > xVals, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of interpolateDim for boundary safe quadratic interpolation of one dimensional h...
double sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
void weights(double *output, std::span< double const > xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of RooDataHist::weight() for one dimensional histograms with up to one dimension...
Int_t _cache_sum_valid
! Is cache sum valid? Needs to be Int_t instead of CacheSumState_t for subclasses.
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...
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:72
const std::vector< double > & calculatePartialBinVolume(const RooArgSet &dimSet) const
Fill the transient cache with partial bin volumes with up-to-date values for the partial volume speci...
static std::unique_ptr< RooAbsDataStore > makeDefaultDataStore(RooStringView name, RooStringView title, RooArgSet const &vars)
double weightInterpolated(const RooArgSet &bin, int intOrder, bool correctForBinSize, bool cdfBoundaries)
Return the weight at given coordinates with interpolation.
std::unordered_map< int, std::vector< double > > _pbinvCache
! Cache for arrays of partial bin volumes
void checkBinBounds() const
void initializeAsymErrArrays() const
void set(std::size_t binNumber, double weight, double wgtErr)
Set bin content of bin that was last loaded with get(std::size_t).
void weightError(double &lo, double &hi, ErrorType etype=Poisson) const override
Return the asymmetric errors on the current weight.
double * _errHi
[_arrSize] High-side error on weight array
void importTH1Set(const RooArgList &vars, RooCategory &indexCat, std::map< std::string, TH1 * > hmap, double initWgt, bool doDensityCorrection)
Import data from given set of TH1/2/3 into this RooDataHist.
void adjustBinning(const RooArgList &vars, const TH1 &href, Int_t *offset=nullptr)
Adjust binning specification on first and optionally second and third observable to binning in given ...
double * _binv
[_arrSize] Bin volume array
RooDataHist()
Default constructor.
ULong64_t _curIndex
Current index.
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.
std::string declWeightArrayForCodeSquash(RooFit::Detail::CodeSquashContext &ctx, bool correctForBinSize) const
bool isNonPoissonWeighted() const override
Returns true if dataset contains entries with a non-integer weight.
std::vector< RooAbsLValue * > _lvvars
! List of observables casted as RooAbsLValue
void SetName(const char *name) override
Change the name of the RooDataHist.
std::vector< std::unique_ptr< const RooAbsBinning > > _lvbins
! List of used binnings associated with lvalues
void Streamer(TBuffer &) override
Stream an object of class RooDataHist.
std::vector< double > _interpolationBuffer
! Buffer to contain values used for weight interpolation
std::vector< Int_t > _idxMult
void registerWeightArraysToDataStore() const
Hand over pointers to our weight arrays to the data store implementation.
void reset() override
Reset all bin weights to zero.
double * _errLo
[_arrSize] Low-side error on weight array
std::string calculateTreeIndexForCodeSquash(RooFit::Detail::CodeSquashContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
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()) const override
Implementation of RooAbsData virtual method that drives the RooAbsData::reduce() methods.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
void interpolateLinear(double *output, std::span< const double > xVals, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of interpolateDim for boundary safe linear interpolation of one dimensional hist...
double binVolume() const
Return volume of current bin.
double sumEntries() const override
Sum the weights of all bins.
Utility base class for RooFit objects that are to be attached to ROOT directories.
Definition RooDirItem.h:22
virtual void Streamer(TBuffer &)
void removeFromDir(TObject *obj)
Remove object from directory it was added to.
TDirectory * _dir
! Associated directory
Definition RooDirItem.h:33
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 buildArg(RooAbsCollection const &x)
Function to save a RooListProxy as an array in the squashed code.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
static const RooHistError & instance()
Return a reference to a singleton object that is created the first time this method is called.
Collection class for internal use, storing a collection of RooAbsArg pointers in a doubly linked list...
static double interpolate(double yArr[], Int_t nOrder, double x)
Definition RooMath.cxx:78
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:45
RooAbsRealLValue * getPlotVar() const
Definition RooPlot.h:143
virtual void printStream(std::ostream &os, Int_t contents, StyleOption style, TString indent="") const
Print description of object on ostream, printing contents set by contents integer,...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
The RooStringView is a wrapper around a C-style string that can also be constructed from a std::strin...
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
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:3037
virtual TList * GetList() const
Definition TDirectory.h:222
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
TAxis * GetZaxis()
Definition TH1.h:326
virtual Int_t GetNbinsY() const
Definition TH1.h:298
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:9084
virtual Int_t GetNbinsZ() const
Definition TH1.h:299
virtual Int_t GetDimension() const
Definition TH1.h:283
TAxis * GetXaxis()
Definition TH1.h:324
virtual Int_t GetNbinsX() const
Definition TH1.h:297
TAxis * GetYaxis()
Definition TH1.h:325
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5082
Iterator abstract base class.
Definition TIterator.h:30
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:213
Basic string class.
Definition TString.h:139
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.
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()