Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ParamHistFunc.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id: cranmer $
2// Author: Kyle Cranmer, George Lewis
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11////////////////////////////////////////////////////////////////////////////////
12
13/** \class ParamHistFunc
14 * \ingroup HistFactory
15 * A class which maps the current values of a RooRealVar
16 * (or a set of RooRealVars) to one of a number of RooAbsReal
17 * (nominally RooRealVar):
18 *
19 * `ParamHistFunc: {val1, val2, ...} -> {gamma (RooRealVar)}`
20 *
21 * The intended interpretation is that each parameter in the
22 * range represent the height of a bin over the domain
23 * space.
24 *
25 * The 'createParamSet' is an easy way to create these
26 * parameters from a set of observables. They are
27 * stored using the "TH1" ordering convention (as compared
28 * to the RooDataHist convention, which is used internally
29 * and one must map between the two).
30 *
31 * All indices include '0':<br>
32 * \f$ \gamma_{i,j} \f$ = `paramSet[ size(i)*j + i ]`
33 *
34 * ie assuming the dimensions are 5*5:<br>
35 * \f$ \gamma_{2,1} \f$ = `paramSet[ 5*1 + 2 ] = paramSet[7]`
36 */
37
38
40
41#include "RooConstVar.h"
42#include "RooBinning.h"
43#include "RooErrorHandler.h"
44#include "RooArgSet.h"
45#include "RooMsgService.h"
46#include "RooRealVar.h"
47#include "RooArgList.h"
48#include "RooWorkspace.h"
49
50#include "TH1.h"
51
52#include <array>
53#include <sstream>
54#include <stdexcept>
55#include <iostream>
56
58
59
60////////////////////////////////////////////////////////////////////////////////
61
63 : _normIntMgr(this)
64{
65 _dataSet.removeSelfFromDir(); // files must not delete _dataSet.
66}
67
68
69////////////////////////////////////////////////////////////////////////////////
70/// Create a function which returns binewise-values
71/// This class contains N RooAbsReals's, one for each
72/// bin from the given RooRealVar.
73///
74/// The value of the function in the ith bin is
75/// given by:
76///
77/// F(i) = gamma_i * nominal(i)
78///
79/// Where the nominal values are simply fixed
80/// numbers (default = 1.0 for all i)
81ParamHistFunc::ParamHistFunc(const char* name, const char* title,
82 const RooArgList& vars, const RooArgList& paramSet) :
83 RooAbsReal(name, title),
84 _normIntMgr(this),
85 _dataVars("!dataVars","data Vars", this),
86 _paramSet("!paramSet","bin parameters", this),
87 _dataSet( (std::string(name)+"_dataSet").c_str(), "", vars)
88{
89
90 // Create the dataset that stores the binning info:
91
92 // _dataSet = RooDataSet("
93
94 _dataSet.removeSelfFromDir(); // files must not delete _dataSet.
95
96 // Set the binning
97 // //_binning = var.getBinning().clone() ;
98
99 // Create the set of parameters
100 // controlling the height of each bin
101
102 // Get the number of bins
103 _numBins = GetNumBins( vars );
104
105 // Add the parameters (with checking)
107 addParamSet( paramSet );
108}
109
110
111////////////////////////////////////////////////////////////////////////////////
112/// Create a function which returns bin-wise values.
113/// This class allows to multiply bin contents of histograms
114/// with the values of a set of RooAbsReal.
115///
116/// The value of the function in the ith bin is
117/// given by:
118/// \f[
119/// F(i) = \gamma_{i} * \mathrm{nominal}(i)
120/// \f]
121///
122/// Where the nominal values are taken from the histogram,
123/// and the \f$ \gamma_{i} \f$ can be set from the outside.
124ParamHistFunc::ParamHistFunc(const char* name, const char* title,
125 const RooArgList& vars, const RooArgList& paramSet,
126 const TH1* Hist ) :
127 RooAbsReal(name, title),
128 _normIntMgr(this),
129 // _dataVar("!dataVar","data Var", this, (RooRealVar&) var),
130 _dataVars("!dataVars","data Vars", this),
131 _paramSet("!paramSet","bin parameters", this),
132 _dataSet( (std::string(name)+"_dataSet").c_str(), "", vars, Hist)
133{
134
135 _dataSet.removeSelfFromDir(); // files must not delete _dataSet.
136
137 // Get the number of bins
138 _numBins = GetNumBins( vars );
139
140 // Add the parameters (with checking)
142 addParamSet( paramSet );
143}
144
145
147
148 // A helper method to get the number of bins
149
150 if( vars.empty() ) return 0;
151
152 Int_t numBins = 1;
153
154 for (auto comp : vars) {
155 if (!dynamic_cast<RooRealVar*>(comp)) {
156 auto errorMsg = std::string("ParamHistFunc::GetNumBins") + vars.GetName() + ") ERROR: component "
157 + comp->GetName() + " in vars list is not of type RooRealVar";
158 oocoutE(nullptr, InputArguments) << errorMsg << std::endl;
159 throw std::runtime_error(errorMsg);
160 }
161 auto var = static_cast<RooRealVar*>(comp);
162
163 Int_t varNumBins = var->numBins();
164 numBins *= varNumBins;
165 }
166
167 return numBins;
168}
169
170
171////////////////////////////////////////////////////////////////////////////////
172
174 RooAbsReal(other, name),
175 _normIntMgr(other._normIntMgr, this),
176 _dataVars("!dataVars", this, other._dataVars ),
177 _paramSet("!paramSet", this, other._paramSet),
178 _numBins( other._numBins ),
179 _dataSet( other._dataSet )
180{
181 _dataSet.removeSelfFromDir(); // files must not delete _dataSet.
182
183 // Copy constructor
184 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
185}
186
187
188////////////////////////////////////////////////////////////////////////////////
189/// Get the parameter associated with the index.
190/// The index follows RooDataHist indexing conventions.
191/// It uses the binMap to convert the RooDataSet style index
192/// into the TH1 style index (which is how they are stored
193/// internally in the '_paramSet' vector).
195
196 auto const& n = _numBinsPerDim;
197
198 // check if _numBins needs to be filled
199 if(n.x == 0) {
201 }
202
203 // Unravel the index to 3D coordinates. We can't use the index directly,
204 // because in the parameter set the dimensions are ordered in reverse order
205 // compared to the RooDataHist (for historical reasons).
206 const int i = index / n.yz;
207 const int tmp = index % n.yz;
208 const int j = tmp / n.z;
209 const int k = tmp % n.z;
210
211 const int idx = i + j * n.x + k * n.xy;
212 if (idx >= _numBins) {
213 throw std::runtime_error("invalid index");
214 }
215 return static_cast<RooAbsReal &>(_paramSet[idx]);
216}
217
218
219////////////////////////////////////////////////////////////////////////////////
220
223 return getParameter( index );
224}
225
226
228 RooAbsReal& var = getParameter( index );
229 var.setAttribute( "Constant", varConst );
230}
231
232
233void ParamHistFunc::setConstant( bool constant ) {
234 for( int i=0; i < numBins(); ++i) {
235 setParamConst(i, constant);
236 }
237}
238
239
240////////////////////////////////////////////////////////////////////////////////
241
243 int num_hist_bins = shape->GetNbinsX()*shape->GetNbinsY()*shape->GetNbinsZ();
244
245 if( num_hist_bins != numBins() ) {
246 std::cout << "Error - ParamHistFunc: cannot set Shape of ParamHistFunc: " << GetName()
247 << " using histogram: " << shape->GetName()
248 << ". Bins don't match" << std::endl;
249 throw std::runtime_error("setShape");
250 }
251
252
253 Int_t TH1BinNumber = 0;
254 for( Int_t i = 0; i < numBins(); ++i) {
255
256 TH1BinNumber++;
257
258 while( shape->IsBinUnderflow(TH1BinNumber) || shape->IsBinOverflow(TH1BinNumber) ){
259 TH1BinNumber++;
260 }
261
262 RooRealVar* param = dynamic_cast<RooRealVar*>(&_paramSet[i]);
263 if(!param) {
264 std::cout << "Error - ParamHisFunc: cannot set Shape of ParamHistFunc: " << GetName()
265 << " - param is not RooRealVar" << std::endl;
266 throw std::runtime_error("setShape");
267 }
268 param->setVal( shape->GetBinContent(TH1BinNumber) );
269 }
270
271}
272
273
274////////////////////////////////////////////////////////////////////////////////
275/// Create the list of RooRealVar
276/// parameters which represent the
277/// height of the histogram bins.
278/// The list 'vars' represents the
279/// observables (corresponding to histogram bins)
280/// that these newly created parameters will
281/// be mapped to. (ie, we create one parameter
282/// per observable in vars and per bin in each observable)
283
284/// Store them in a list using:
285/// _paramSet.add( createParamSet() );
286/// This list is stored in the "TH1" index order
288 const RooArgList& vars) {
289
290
291 // Get the number of bins
292 // in the nominal histogram
293
294 RooArgList paramSet;
295
296 Int_t numVars = vars.size();
297 Int_t numBins = GetNumBins( vars );
298
299 if( numVars == 0 ) {
300 std::cout << "Warning - ParamHistFunc::createParamSet() :"
301 << " No Variables provided. Not making constraint terms."
302 << std::endl;
303 return paramSet;
304 }
305
306 else if( numVars == 1 ) {
307
308 // For each bin, create a RooRealVar
309 for( Int_t i = 0; i < numBins; ++i) {
310
311 std::stringstream VarNameStream;
312 VarNameStream << Prefix << "_bin_" << i;
313 std::string VarName = VarNameStream.str();
314
315 RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 );
316 // "Hard-Code" a minimum of 0.0
317 gamma.setMin( 0.0 );
318 gamma.setConstant( false );
319
320 w.import( gamma, RooFit::RecycleConflictNodes() );
321 RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName );
322
323 paramSet.add( *gamma_wspace );
324
325 }
326 }
327
328 else if( numVars == 2 ) {
329
330 // Create a vector of indices
331 // all starting at 0
332 std::vector< Int_t > Indices(numVars, 0);
333
334 RooRealVar* varx = static_cast<RooRealVar*>(vars.at(0));
335 RooRealVar* vary = static_cast<RooRealVar*>(vars.at(1));
336
337 // For each bin, create a RooRealVar
338 for( Int_t j = 0; j < vary->numBins(); ++j) {
339 for( Int_t i = 0; i < varx->numBins(); ++i) {
340
341 // Ordering is important:
342 // To match TH1, list goes over x bins
343 // first, then y
344
345 std::stringstream VarNameStream;
346 VarNameStream << Prefix << "_bin_" << i << "_" << j;
347 std::string VarName = VarNameStream.str();
348
349 RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 );
350 // "Hard-Code" a minimum of 0.0
351 gamma.setMin( 0.0 );
352 gamma.setConstant( false );
353
354 w.import( gamma, RooFit::RecycleConflictNodes() );
355 RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName );
356
357 paramSet.add( *gamma_wspace );
358
359 }
360 }
361 }
362
363 else if( numVars == 3 ) {
364
365 // Create a vector of indices
366 // all starting at 0
367 std::vector< Int_t > Indices(numVars, 0);
368
369 RooRealVar* varx = static_cast<RooRealVar*>(vars.at(0));
370 RooRealVar* vary = static_cast<RooRealVar*>(vars.at(1));
371 RooRealVar* varz = static_cast<RooRealVar*>(vars.at(2));
372
373 // For each bin, create a RooRealVar
374 for( Int_t k = 0; k < varz->numBins(); ++k) {
375 for( Int_t j = 0; j < vary->numBins(); ++j) {
376 for( Int_t i = 0; i < varx->numBins(); ++i) {
377
378 // Ordering is important:
379 // To match TH1, list goes over x bins
380 // first, then y, then z
381
382 std::stringstream VarNameStream;
383 VarNameStream << Prefix << "_bin_" << i << "_" << j << "_" << k;
384 std::string VarName = VarNameStream.str();
385
386 RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 );
387 // "Hard-Code" a minimum of 0.0
388 gamma.setMin( 0.0 );
389 gamma.setConstant( false );
390
391 w.import( gamma, RooFit::RecycleConflictNodes() );
392 RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName );
393
394 paramSet.add( *gamma_wspace );
395
396 }
397 }
398 }
399 }
400
401 else {
402 std::cout << " Error: ParamHistFunc doesn't support dimensions > 3D " << std::endl;
403 }
404
405 return paramSet;
406
407}
408
409
410////////////////////////////////////////////////////////////////////////////////
411/// Create the list of RooRealVar parameters which scale the
412/// height of histogram bins.
413/// The list `vars` represents the observables (corresponding to histogram bins)
414/// that these newly created parameters will
415/// be mapped to. *I.e.*, we create one parameter
416/// per observable in `vars` and per bin in each observable.
417///
418/// The new parameters are initialised to 1 with an uncertainty of +/- 1.,
419/// their range is set to the function arguments.
420///
421/// Store the parameters in a list using:
422/// ```
423/// _paramSet.add( createParamSet() );
424/// ```
425/// This list is stored in the "TH1" index order.
427 const RooArgList& vars,
428 double gamma_min, double gamma_max) {
429
430
431
433
434 for (auto comp : params) {
435 auto var = static_cast<RooRealVar*>(comp);
436
437 var->setMin( gamma_min );
438 var->setMax( gamma_max );
439 }
440
441 return params;
442
443}
444
445
446////////////////////////////////////////////////////////////////////////////////
447/// Create the list of RooRealVar
448/// parameters which represent the
449/// height of the histogram bins.
450/// Store them in a list
452 double gamma_min, double gamma_max) {
453
454 // Get the number of bins
455 // in the nominal histogram
456
457 RooArgList paramSet;
458
459 if( gamma_max <= gamma_min ) {
460
461 std::cout << "Warning: gamma_min <= gamma_max: Using default values (0, 10)" << std::endl;
462
463 gamma_min = 0.0;
464 gamma_max = 10.0;
465
466 }
467
468 double gamma_nominal = 1.0;
469
470 if( gamma_nominal < gamma_min ) {
471 gamma_nominal = gamma_min;
472 }
473
474 if( gamma_nominal > gamma_max ) {
475 gamma_nominal = gamma_max;
476 }
477
478 // For each bin, create a RooRealVar
479 for( Int_t i = 0; i < numBins; ++i) {
480
481 std::stringstream VarNameStream;
482 VarNameStream << Prefix << "_bin_" << i;
483 std::string VarName = VarNameStream.str();
484
485 auto gamma = std::make_unique<RooRealVar>(VarName.c_str(), VarName.c_str(),
486 gamma_nominal, gamma_min, gamma_max);
487 gamma->setConstant( false );
488 paramSet.addOwned(std::move(gamma));
489
490 }
491
492 return paramSet;
493
494}
495
496
498 int numVars = vars.size();
499
500 if (numVars > 3 || numVars < 1) {
501 std::cout << "ParamHistFunc() - Only works for 1-3 variables (1d-3d)" << std::endl;
502 throw -1;
503 }
504
505 int numBinsX = numVars >= 1 ? static_cast<RooRealVar const&>(*vars[0]).numBins() : 1;
506 int numBinsY = numVars >= 2 ? static_cast<RooRealVar const&>(*vars[1]).numBins() : 1;
507 int numBinsZ = numVars >= 3 ? static_cast<RooRealVar const&>(*vars[2]).numBins() : 1;
508
509 return {numBinsX, numBinsY, numBinsZ};
510}
511
512
513////////////////////////////////////////////////////////////////////////////////
514/// Get the index of the gamma parameter associated with the current bin.
515/// e.g. `RooRealVar& currentParam = getParameter( getCurrentBin() );`
517 // We promise that our coordinates and the data hist coordinates have same layout.
518 return _dataSet.getIndex(_dataVars, /*fast=*/true);
519}
520
521
522////////////////////////////////////////////////////////////////////////////////
523
525 // return 0 for success
526 // return 1 for failure
527
528 // Check that the supplied list has
529 // the right number of arguments:
530
531 Int_t numVarBins = GetNumBins(_dataVars);
532 Int_t numElements = params.size();
533
534 if( numVarBins != numElements ) {
535 std::cout << "ParamHistFunc::addParamSet - ERROR - "
536 << "Supplied list of parameters " << params.GetName()
537 << " has " << numElements << " elements but the ParamHistFunc"
538 << GetName() << " has " << numVarBins << " bins."
539 << std::endl;
540 return 1;
541
542 }
543
544 // Check that the elements
545 // are actually RooAbsreal's
546 // If so, add them to the
547 // list of params
548
550
551 return 0;
552}
553
554
555////////////////////////////////////////////////////////////////////////////////
556/// Find the bin corresponding to the current value of the observable, and evaluate
557/// the associated parameter.
559{
560 return getParameter().getVal();
561}
562
563////////////////////////////////////////////////////////////////////////////////
564/// Find all bins corresponding to the values of the observables in `evalData`, and evaluate
565/// the associated parameters.
566/// \param[in,out] evalData Input/output data for evaluating the ParamHistFunc.
567/// \param[in] normSet Normalisation set passed on to objects that are serving values to us.
569{
570 std::span<double> output = ctx.output();
571 std::size_t size = output.size();
572
573 auto const& n = _numBinsPerDim;
574 // check if _numBins needs to be filled
575 if(n.x == 0) {
577 }
578
579 // Different from the evaluate() funnction that first retrieves the indices
580 // corresponding to the RooDataHist and then transforms them, we can use the
581 // right bin multiplicators to begin with.
582 std::array<int, 3> idxMult{{1, n.x, n.xy}};
583
584 // As a working buffer for the bin indices, we use the tail of the output
585 // buffer. We can't use the same starting pointer, otherwise we would
586 // overwrite the later bin indices as we fill the output.
587 auto indexBuffer = reinterpret_cast<int*>(output.data() + size) - size;
588 std::fill(indexBuffer, indexBuffer + size, 0); // output buffer for bin indices needs to be zero-initialized
589
590 // Use the vectorized RooAbsBinning::binNumbers() to update the total bin
591 // index for each dimension, using the `coef` parameter to multiply with the
592 // right index multiplication factor for each dimension.
593 for (std::size_t iVar = 0; iVar < _dataVars.size(); ++iVar) {
594 _dataSet.getBinnings()[iVar]->binNumbers(ctx.at(&_dataVars[iVar]).data(), indexBuffer, size, idxMult[iVar]);
595 }
596
597 // Finally, look up the parameters and get their values to fill the output buffer
598 for (std::size_t i = 0; i < size; ++i) {
599 output[i] = static_cast<RooAbsReal const&>(_paramSet[indexBuffer[i]]).getVal();
600 }
601}
602
603////////////////////////////////////////////////////////////////////////////////
604/// Advertise that all integrals can be handled internally.
605
607 const RooArgSet* normSet,
608 const char* /*rangeName*/) const
609{
610 // Handle trivial no-integration scenario
611 if (allVars.empty()) return 0 ;
612 if (_forceNumInt) return 0 ;
613
614
615 // Select subset of allVars that are actual dependents
616 analVars.add(allVars) ;
617
618 // Check if this configuration was created before
619 Int_t sterileIdx(-1) ;
620 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObj(normSet,&analVars,&sterileIdx,(const char*)nullptr)) ;
621 if (cache) {
622 return _normIntMgr.lastIndex()+1 ;
623 }
624
625 // Create new cache element
626 cache = new CacheElem ;
627
628 // Store cache element
629 Int_t code = _normIntMgr.setObj(normSet,&analVars,cache,nullptr) ;
630
631 return code+1 ;
632
633}
634
635
636////////////////////////////////////////////////////////////////////////////////
637/// Implement analytical integrations by doing appropriate weighting from component integrals
638/// functions to integrators of components
639
640double ParamHistFunc::analyticalIntegralWN(Int_t /*code*/, const RooArgSet* /*normSet2*/,
641 const char* /*rangeName*/) const
642{
643 double value(0) ;
644
645 // Simply loop over bins,
646 // get the height, and
647 // multiply by the bind width
648 auto binVolumes = _dataSet.binVolumes(0, _dataSet.numEntries());
649
650 for (unsigned int i=0; i < _paramSet.size(); ++i) {
651 const auto& param = static_cast<const RooAbsReal&>(_paramSet[i]);
652
653 // Get the gamma's value
654 const double paramVal = param.getVal();
655
656 // Finally, get the subtotal
657 value += paramVal * binVolumes[i];
658 }
659
660 return value;
661
662}
663
664
665
666////////////////////////////////////////////////////////////////////////////////
667/// Return sampling hint for making curves of (projections) of this function
668/// as the recursive division strategy of RooCurve cannot deal efficiently
669/// with the vertical lines that occur in a non-interpolated histogram
670
671std::list<double>* ParamHistFunc::plotSamplingHint(RooAbsRealLValue& obs, double xlo,
672 double xhi) const
673{
674 // copied and edited from RooHistFunc
675 RooAbsLValue* lvarg = &obs;
676
677 // Retrieve position of all bin boundaries
678 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
679 double* boundaries = binning->array() ;
680
681 std::list<double>* hint = new std::list<double> ;
682
683 // Widen range slightly
684 xlo = xlo - 0.01*(xhi-xlo) ;
685 xhi = xhi + 0.01*(xhi-xlo) ;
686
687 double delta = (xhi-xlo)*1e-8 ;
688
689 // Construct array with pairs of points positioned epsilon to the left and
690 // right of the bin boundaries
691 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
692 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
693 hint->push_back(boundaries[i]-delta) ;
694 hint->push_back(boundaries[i]+delta) ;
695 }
696 }
697 return hint ;
698}
699
700
701////////////////////////////////////////////////////////////////////////////////
702/// Return sampling hint for making curves of (projections) of this function
703/// as the recursive division strategy of RooCurve cannot deal efficiently
704/// with the vertical lines that occur in a non-interpolated histogram
705
706std::list<double> *ParamHistFunc::binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
707{
708 // copied and edited from RooHistFunc
709 RooAbsLValue *lvarg = &obs;
710
711 // look for variable in the DataHist, and if found, return the binning
712 std::string varName = dynamic_cast<TObject *>(lvarg)->GetName();
713 RooArgSet const &vars = *_dataSet.get(); // guaranteed to be in the same order as the binnings vector
714 auto &binnings = _dataSet.getBinnings();
715 for (size_t i = 0; i < vars.size(); i++) {
716 if (varName == vars[i]->GetName()) {
717 // found the variable, return its binning
718 double *boundaries = binnings.at(i)->array();
719 std::list<double> *hint = new std::list<double>;
720 for (int j = 0; j < binnings.at(i)->numBoundaries(); j++) {
721 if (boundaries[j] >= xlo && boundaries[j] <= xhi) {
722 hint->push_back(boundaries[j]);
723 }
724 }
725 return hint;
726 }
727 }
728 // variable not found, return null
729 return nullptr;
730}
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define oocoutE(o, a)
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
A class which maps the current values of a RooRealVar (or a set of RooRealVars) to one of a number of...
static NumBins getNumBinsPerDim(RooArgSet const &vars)
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
RooDataHist _dataSet
static Int_t GetNumBins(const RooArgSet &vars)
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
void setConstant(bool constant)
Int_t getCurrentBin() const
Get the index of the gamma parameter associated with the current bin.
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise that all integrals can be handled internally.
RooObjCacheManager _normIntMgr
! The integration cache manager
Int_t numBins() const
double evaluate() const override
Find the bin corresponding to the current value of the observable, and evaluate the associated parame...
void doEval(RooFit::EvalContext &) const override
Find all bins corresponding to the values of the observables in evalData, and evaluate the associated...
RooAbsReal & getParameter() const
Int_t addParamSet(const RooArgList &params)
NumBins _numBinsPerDim
static RooArgList createParamSet(RooWorkspace &w, const std::string &, const RooArgList &Vars)
Create the list of RooRealVar parameters which represent the height of the histogram bins.
void setParamConst(Int_t, bool=true)
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
RooListProxy _paramSet
interpolation parameters
void setShape(TH1 *shape)
RooListProxy _dataVars
The RooRealVars.
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
Abstract base class for RooRealVar binning definitions.
virtual Int_t numBoundaries() const =0
virtual double * array() const =0
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
bool addTyped(const RooAbsCollection &list, bool silent=false)
Adds elements of a given RooAbsCollection to the container if they match the specified type.
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract base class for objects that are lvalues, i.e.
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const =0
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Int_t numBins(const char *rangeName=nullptr) const override
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
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:538
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
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=nullptr)
Setter function without integration set.
Int_t lastIndex() const
Return index of slot used in last get or set operation.
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=nullptr, const TNamed *isetRangeName=nullptr)
Getter function without integration set.
std::vector< std::unique_ptr< const RooAbsBinning > > const & getBinnings() const
Int_t getIndex(const RooAbsCollection &coord, bool fast=false) const
Calculate bin number of the given coordinates.
void removeSelfFromDir()
std::span< const double > binVolumes(std::size_t first, std::size_t len) const
Retrieve all bin volumes. Bins are indexed according to getIndex().
Definition RooDataHist.h:95
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
void setMin(const char *name, double value)
Set minimum of name range to given value.
Persistable container for RooFit projects.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
virtual Int_t GetNbinsY() const
Definition TH1.h:299
virtual Int_t GetNbinsZ() const
Definition TH1.h:300
virtual Int_t GetNbinsX() const
Definition TH1.h:298
Bool_t IsBinUnderflow(Int_t bin, Int_t axis=0) const
Return true if the bin is underflow.
Definition TH1.cxx:5243
Bool_t IsBinOverflow(Int_t bin, Int_t axis=0) const
Return true if the bin is overflow.
Definition TH1.cxx:5211
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5090
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
RooCmdArg RecycleConflictNodes(bool flag=true)
const Int_t n
Definition legend1.C:16
static int Prefix[4096]
Definition gifdecode.c:12
static void output()