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)
106 addVarSet( vars );
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)
141 addVarSet( vars );
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 return static_cast<RooAbsReal&>(_paramSet[i + j * n.x + k * n.xy]);
212}
213
214
215////////////////////////////////////////////////////////////////////////////////
216
219 return getParameter( index );
220}
221
222
224 RooAbsReal& var = getParameter( index );
225 var.setAttribute( "Constant", varConst );
226}
227
228
229void ParamHistFunc::setConstant( bool constant ) {
230 for( int i=0; i < numBins(); ++i) {
231 setParamConst(i, constant);
232 }
233}
234
235
236////////////////////////////////////////////////////////////////////////////////
237
239 int num_hist_bins = shape->GetNbinsX()*shape->GetNbinsY()*shape->GetNbinsZ();
240
241 if( num_hist_bins != numBins() ) {
242 std::cout << "Error - ParamHistFunc: cannot set Shape of ParamHistFunc: " << GetName()
243 << " using histogram: " << shape->GetName()
244 << ". Bins don't match" << std::endl;
245 throw std::runtime_error("setShape");
246 }
247
248
249 Int_t TH1BinNumber = 0;
250 for( Int_t i = 0; i < numBins(); ++i) {
251
252 TH1BinNumber++;
253
254 while( shape->IsBinUnderflow(TH1BinNumber) || shape->IsBinOverflow(TH1BinNumber) ){
255 TH1BinNumber++;
256 }
257
258 RooRealVar* param = dynamic_cast<RooRealVar*>(&_paramSet[i]);
259 if(!param) {
260 std::cout << "Error - ParamHisFunc: cannot set Shape of ParamHistFunc: " << GetName()
261 << " - param is not RooRealVar" << std::endl;
262 throw std::runtime_error("setShape");
263 }
264 param->setVal( shape->GetBinContent(TH1BinNumber) );
265 }
266
267}
268
269
270////////////////////////////////////////////////////////////////////////////////
271/// Create the list of RooRealVar
272/// parameters which represent the
273/// height of the histogram bins.
274/// The list 'vars' represents the
275/// observables (corresponding to histogram bins)
276/// that these newly created parameters will
277/// be mapped to. (ie, we create one parameter
278/// per observable in vars and per bin in each observable)
279
280/// Store them in a list using:
281/// _paramSet.add( createParamSet() );
282/// This list is stored in the "TH1" index order
284 const RooArgList& vars) {
285
286
287 // Get the number of bins
288 // in the nominal histogram
289
290 RooArgList paramSet;
291
292 Int_t numVars = vars.getSize();
293 Int_t numBins = GetNumBins( vars );
294
295 if( numVars == 0 ) {
296 std::cout << "Warning - ParamHistFunc::createParamSet() :"
297 << " No Variables provided. Not making constraint terms."
298 << std::endl;
299 return paramSet;
300 }
301
302 else if( numVars == 1 ) {
303
304 // For each bin, create a RooRealVar
305 for( Int_t i = 0; i < numBins; ++i) {
306
307 std::stringstream VarNameStream;
308 VarNameStream << Prefix << "_bin_" << i;
309 std::string VarName = VarNameStream.str();
310
311 RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 );
312 // "Hard-Code" a minimum of 0.0
313 gamma.setMin( 0.0 );
314 gamma.setConstant( false );
315
316 w.import( gamma, RooFit::RecycleConflictNodes() );
317 RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName );
318
319 paramSet.add( *gamma_wspace );
320
321 }
322 }
323
324 else if( numVars == 2 ) {
325
326 // Create a vector of indices
327 // all starting at 0
328 std::vector< Int_t > Indices(numVars, 0);
329
330 RooRealVar* varx = (RooRealVar*) vars.at(0);
331 RooRealVar* vary = (RooRealVar*) vars.at(1);
332
333 // For each bin, create a RooRealVar
334 for( Int_t j = 0; j < vary->numBins(); ++j) {
335 for( Int_t i = 0; i < varx->numBins(); ++i) {
336
337 // Ordering is important:
338 // To match TH1, list goes over x bins
339 // first, then y
340
341 std::stringstream VarNameStream;
342 VarNameStream << Prefix << "_bin_" << i << "_" << j;
343 std::string VarName = VarNameStream.str();
344
345 RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 );
346 // "Hard-Code" a minimum of 0.0
347 gamma.setMin( 0.0 );
348 gamma.setConstant( false );
349
350 w.import( gamma, RooFit::RecycleConflictNodes() );
351 RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName );
352
353 paramSet.add( *gamma_wspace );
354
355 }
356 }
357 }
358
359 else if( numVars == 3 ) {
360
361 // Create a vector of indices
362 // all starting at 0
363 std::vector< Int_t > Indices(numVars, 0);
364
365 RooRealVar* varx = (RooRealVar*) vars.at(0);
366 RooRealVar* vary = (RooRealVar*) vars.at(1);
367 RooRealVar* varz = (RooRealVar*) vars.at(2);
368
369 // For each bin, create a RooRealVar
370 for( Int_t k = 0; k < varz->numBins(); ++k) {
371 for( Int_t j = 0; j < vary->numBins(); ++j) {
372 for( Int_t i = 0; i < varx->numBins(); ++i) {
373
374 // Ordering is important:
375 // To match TH1, list goes over x bins
376 // first, then y, then z
377
378 std::stringstream VarNameStream;
379 VarNameStream << Prefix << "_bin_" << i << "_" << j << "_" << k;
380 std::string VarName = VarNameStream.str();
381
382 RooRealVar gamma( VarName.c_str(), VarName.c_str(), 1.0 );
383 // "Hard-Code" a minimum of 0.0
384 gamma.setMin( 0.0 );
385 gamma.setConstant( false );
386
387 w.import( gamma, RooFit::RecycleConflictNodes() );
388 RooRealVar* gamma_wspace = (RooRealVar*) w.var( VarName );
389
390 paramSet.add( *gamma_wspace );
391
392 }
393 }
394 }
395 }
396
397 else {
398 std::cout << " Error: ParamHistFunc doesn't support dimensions > 3D " << std::endl;
399 }
400
401 return paramSet;
402
403}
404
405
406////////////////////////////////////////////////////////////////////////////////
407/// Create the list of RooRealVar parameters which scale the
408/// height of histogram bins.
409/// The list `vars` represents the observables (corresponding to histogram bins)
410/// that these newly created parameters will
411/// be mapped to. *I.e.*, we create one parameter
412/// per observable in `vars` and per bin in each observable.
413///
414/// The new parameters are initialised to 1 with an uncertainty of +/- 1.,
415/// their range is set to the function arguments.
416///
417/// Store the parameters in a list using:
418/// ```
419/// _paramSet.add( createParamSet() );
420/// ```
421/// This list is stored in the "TH1" index order.
423 const RooArgList& vars,
424 double gamma_min, double gamma_max) {
425
426
427
429
430 for (auto comp : params) {
431 auto var = static_cast<RooRealVar*>(comp);
432
433 var->setMin( gamma_min );
434 var->setMax( gamma_max );
435 }
436
437 return params;
438
439}
440
441
442////////////////////////////////////////////////////////////////////////////////
443/// Create the list of RooRealVar
444/// parameters which represent the
445/// height of the histogram bins.
446/// Store them in a list
448 double gamma_min, double gamma_max) {
449
450 // Get the number of bins
451 // in the nominal histogram
452
453 RooArgList paramSet;
454
455 if( gamma_max <= gamma_min ) {
456
457 std::cout << "Warning: gamma_min <= gamma_max: Using default values (0, 10)" << std::endl;
458
459 gamma_min = 0.0;
460 gamma_max = 10.0;
461
462 }
463
464 double gamma_nominal = 1.0;
465
466 if( gamma_nominal < gamma_min ) {
467 gamma_nominal = gamma_min;
468 }
469
470 if( gamma_nominal > gamma_max ) {
471 gamma_nominal = gamma_max;
472 }
473
474 // For each bin, create a RooRealVar
475 for( Int_t i = 0; i < numBins; ++i) {
476
477 std::stringstream VarNameStream;
478 VarNameStream << Prefix << "_bin_" << i;
479 std::string VarName = VarNameStream.str();
480
481 RooRealVar* gamma = new RooRealVar( VarName.c_str(), VarName.c_str(),
482 gamma_nominal, gamma_min, gamma_max );
483 gamma->setConstant( false );
484 paramSet.add( *gamma );
485
486 }
487
488 return paramSet;
489
490}
491
492
494 int numVars = vars.size();
495
496 if (numVars > 3 || numVars < 1) {
497 std::cout << "ParamHistFunc() - Only works for 1-3 variables (1d-3d)" << std::endl;
498 throw -1;
499 }
500
501 int numBinsX = numVars >= 1 ? static_cast<RooRealVar const&>(*vars[0]).numBins() : 1;
502 int numBinsY = numVars >= 2 ? static_cast<RooRealVar const&>(*vars[1]).numBins() : 1;
503 int numBinsZ = numVars >= 3 ? static_cast<RooRealVar const&>(*vars[2]).numBins() : 1;
504
505 return {numBinsX, numBinsY, numBinsZ};
506}
507
508
509////////////////////////////////////////////////////////////////////////////////
510/// Get the index of the gamma parameter associated with the current bin.
511/// e.g. `RooRealVar& currentParam = getParameter( getCurrentBin() );`
513 // We promise that our coordinates and the data hist coordinates have same layout.
514 return _dataSet.getIndex(_dataVars, /*fast=*/true);
515}
516
517
518////////////////////////////////////////////////////////////////////////////////
519/// return 0 for success
520/// return 1 for failure
521/// Check that the elements
522/// are actually RooRealVar's
523/// If so, add them to the
524/// list of vars
526 for(auto const& comp : vars) {
527 if (!dynamic_cast<RooRealVar*>(comp)) {
528 auto errorMsg = std::string("ParamHistFunc::(") + GetName() + ") ERROR: component "
529 + comp->GetName() + " in variables list is not of type RooRealVar";
530 coutE(InputArguments) << errorMsg << std::endl;
531 throw std::runtime_error(errorMsg);
532 }
533 _dataVars.add( *comp );
534 }
535 return 0;
536}
537
538
539////////////////////////////////////////////////////////////////////////////////
540
542 // return 0 for success
543 // return 1 for failure
544
545 // Check that the supplied list has
546 // the right number of arguments:
547
548 Int_t numVarBins = GetNumBins(_dataVars);
549 Int_t numElements = params.getSize();
550
551 if( numVarBins != numElements ) {
552 std::cout << "ParamHistFunc::addParamSet - ERROR - "
553 << "Supplied list of parameters " << params.GetName()
554 << " has " << numElements << " elements but the ParamHistFunc"
555 << GetName() << " has " << numVarBins << " bins."
556 << std::endl;
557 return 1;
558
559 }
560
561 // Check that the elements
562 // are actually RooAbsreal's
563 // If so, add them to the
564 // list of params
565
566 for (const auto comp : params) {
567 if (!dynamic_cast<const RooAbsReal*>(comp)) {
568 auto errorMsg = std::string("ParamHistFunc::(") + GetName() + ") ERROR: component "
569 + comp->GetName() + " in parameter list is not of type RooAbsReal.";
570 coutE(InputArguments) << errorMsg << std::endl;
571 throw std::runtime_error(errorMsg);
572 }
573
574 _paramSet.add( *comp );
575 }
576
577 return 0;
578}
579
580
581////////////////////////////////////////////////////////////////////////////////
582/// Find the bin corresponding to the current value of the observable, and evaluate
583/// the associated parameter.
585{
586 return getParameter().getVal();
587}
588
589
590////////////////////////////////////////////////////////////////////////////////
591/// Find all bins corresponding to the values of the observables in `evalData`, and evaluate
592/// the associated parameters.
593/// \param[in,out] evalData Input/output data for evaluating the ParamHistFunc.
594/// \param[in] normSet Normalisation set passed on to objects that are serving values to us.
595void ParamHistFunc::computeBatch(cudaStream_t*, double* output, size_t size, RooFit::Detail::DataMap const& dataMap) const {
596
597 auto const& n = _numBinsPerDim;
598 // check if _numBins needs to be filled
599 if(n.x == 0) {
601 }
602
603 // Different from the evaluate() funnction that first retrieves the indices
604 // corresponding to the RooDataHist and then transforms them, we can use the
605 // right bin multiplicators to begin with.
606 std::array<int, 3> idxMult{{1, n.x, n.xy}};
607
608 // As a working buffer for the bin indices, we use the tail of the output
609 // buffer. We can't use the same starting pointer, otherwise we would
610 // overwrite the later bin indices as we fill the output.
611 auto indexBuffer = reinterpret_cast<int*>(output + size) - size;
612 std::fill(indexBuffer, indexBuffer + size, 0); // output buffer for bin indices needs to be zero-initialized
613
614 // Use the vectorized RooAbsBinning::binNumbers() to update the total bin
615 // index for each dimension, using the `coef` parameter to multiply with the
616 // right index multiplication factor for each dimension.
617 for (std::size_t iVar = 0; iVar < _dataVars.size(); ++iVar) {
618 _dataSet.getBinnings()[iVar]->binNumbers(dataMap.at(&_dataVars[iVar]).data(), indexBuffer, size, idxMult[iVar]);
619 }
620
621 // Finally, look up the parameters and get their values to fill the output buffer
622 for (std::size_t i = 0; i < size; ++i) {
623 output[i] = static_cast<RooAbsReal const&>(_paramSet[indexBuffer[i]]).getVal();
624 }
625}
626
627////////////////////////////////////////////////////////////////////////////////
628/// Advertise that all integrals can be handled internally.
629
631 const RooArgSet* normSet,
632 const char* /*rangeName*/) const
633{
634 // Handle trivial no-integration scenario
635 if (allVars.empty()) return 0 ;
636 if (_forceNumInt) return 0 ;
637
638
639 // Select subset of allVars that are actual dependents
640 analVars.add(allVars) ;
641
642 // Check if this configuration was created before
643 Int_t sterileIdx(-1) ;
644 CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx,(const char*)0) ;
645 if (cache) {
646 return _normIntMgr.lastIndex()+1 ;
647 }
648
649 // Create new cache element
650 cache = new CacheElem ;
651
652 // Store cache element
653 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,0) ;
654
655 return code+1 ;
656
657}
658
659
660////////////////////////////////////////////////////////////////////////////////
661/// Implement analytical integrations by doing appropriate weighting from component integrals
662/// functions to integrators of components
663
664double ParamHistFunc::analyticalIntegralWN(Int_t /*code*/, const RooArgSet* /*normSet2*/,
665 const char* /*rangeName*/) const
666{
667 double value(0) ;
668
669 // Simply loop over bins,
670 // get the height, and
671 // multiply by the bind width
672 auto binVolumes = _dataSet.binVolumes(0, _dataSet.numEntries());
673
674 for (unsigned int i=0; i < _paramSet.size(); ++i) {
675 const auto& param = static_cast<const RooAbsReal&>(_paramSet[i]);
676
677 // Get the gamma's value
678 const double paramVal = param.getVal();
679
680 // Finally, get the subtotal
681 value += paramVal * binVolumes[i];
682 }
683
684 return value;
685
686}
687
688
689
690////////////////////////////////////////////////////////////////////////////////
691/// Return sampling hint for making curves of (projections) of this function
692/// as the recursive division strategy of RooCurve cannot deal efficiently
693/// with the vertical lines that occur in a non-interpolated histogram
694
695std::list<double>* ParamHistFunc::plotSamplingHint(RooAbsRealLValue& obs, double xlo,
696 double xhi) const
697{
698 // copied and edited from RooHistFunc
699 RooAbsLValue* lvarg = &obs;
700
701 // Retrieve position of all bin boundaries
702 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
703 double* boundaries = binning->array() ;
704
705 std::list<double>* hint = new std::list<double> ;
706
707 // Widen range slighty
708 xlo = xlo - 0.01*(xhi-xlo) ;
709 xhi = xhi + 0.01*(xhi-xlo) ;
710
711 double delta = (xhi-xlo)*1e-8 ;
712
713 // Construct array with pairs of points positioned epsilon to the left and
714 // right of the bin boundaries
715 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
716 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
717 hint->push_back(boundaries[i]-delta) ;
718 hint->push_back(boundaries[i]+delta) ;
719 }
720 }
721 return hint ;
722}
723
724
725////////////////////////////////////////////////////////////////////////////////
726/// Return sampling hint for making curves of (projections) of this function
727/// as the recursive division strategy of RooCurve cannot deal efficiently
728/// with the vertical lines that occur in a non-interpolated histogram
729
730std::list<double>* ParamHistFunc::binBoundaries(RooAbsRealLValue& obs, double xlo,
731 double xhi) const
732{
733 // copied and edited from RooHistFunc
734 RooAbsLValue* lvarg = &obs;
735
736 // Retrieve position of all bin boundaries
737 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
738 double* boundaries = binning->array() ;
739
740 std::list<double>* hint = new std::list<double> ;
741
742 // Construct array with pairs of points positioned epsilon to the left and
743 // right of the bin boundaries
744 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
745 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
746 hint->push_back(boundaries[i]) ;
747 }
748 }
749
750 return hint ;
751}
#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 coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
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...
Int_t addVarSet(const RooArgList &vars)
return 0 for success return 1 for failure Check that the elements are actually RooRealVar's If so,...
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 computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const override
Find all bins corresponding to the values of the observables in evalData, and evaluate the associated...
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
virtual Int_t numBoundaries() const =0
virtual double * array() const =0
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
Int_t getSize() const
Return the number of elements in the collection.
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
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
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
Int_t numBins(const char *rangeName=nullptr) const override
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:512
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
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.
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
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.
RooSpan< 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:89
void removeSelfFromDir()
RooSpan< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:21
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
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.
The RooWorkspace is a persistable container for RooFit projects.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual Int_t GetNbinsY() const
Definition TH1.h:296
virtual Int_t GetNbinsZ() const
Definition TH1.h:297
virtual Int_t GetNbinsX() const
Definition TH1.h:295
Bool_t IsBinUnderflow(Int_t bin, Int_t axis=0) const
Return true if the bin is underflow.
Definition TH1.cxx:5178
Bool_t IsBinOverflow(Int_t bin, Int_t axis=0) const
Return true if the bin is overflow.
Definition TH1.cxx:5146
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5025
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
RooCmdArg RecycleConflictNodes(bool flag=true)
const Int_t n
Definition legend1.C:16
static int Prefix[4096]
Definition gifdecode.c:12
static void output()