Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
PiecewiseInterpolation.cxx
Go to the documentation of this file.
1/** \class PiecewiseInterpolation
2 * \ingroup HistFactory
3 * The PiecewiseInterpolation is a class that can morph distributions into each other, which
4 * is useful to estimate systematic uncertainties. Given a nominal distribution and one or
5 * more altered or distorted ones, it computes a new shape depending on the value of the nuisance
6 * parameters \f$ \alpha_i \f$:
7 * \f[
8 * A = \sum_i \mathrm{Interpolate}(\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i, \alpha_i).
9 * \f]
10 * If an \f$ \alpha_i \f$ is zero, the distribution is identical to the nominal distribution, at
11 * \f$ \pm 1 \f$ it is identical to the up/down distribution for that specific \f$ i \f$.
12 *
13 * The class supports several interpolation methods, which can be selected for each parameter separately
14 * using setInterpCode(). The default interpolation code is 4. This performs
15 * - \f$ |\alpha | > 1 \f$: Linear extrapolation.
16 * - \f$ |\alpha | < 1 \f$: Polynomial interpolation. A sixth-order polynomial is used. Its coefficients
17 * are chosen such that function, first, and second derivative at \f$ \alpha \pm 1 \f$ match the values
18 * that the extrapolation procedure uses.
19 */
20
22
24
25#include "Riostream.h"
26#include "TBuffer.h"
27
28#include "RooAbsReal.h"
29#include "RooAbsPdf.h"
30#include "RooErrorHandler.h"
31#include "RooArgSet.h"
32#include "RooRealVar.h"
33#include "RooMsgService.h"
34#include "RooNumIntConfig.h"
35#include "RooTrace.h"
36
37#include <exception>
38#include <cmath>
39#include <algorithm>
40
41using std::endl, std::cout;
42
44
46
47////////////////////////////////////////////////////////////////////////////////
48
50{
52}
53
54////////////////////////////////////////////////////////////////////////////////
55/// Construct a new interpolation. The value of the function will be
56/// \f[
57/// A = \sum_i \mathrm{Interpolate}(\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i).
58/// \f]
59/// \param name Name of the object.
60/// \param title Title (for e.g. plotting)
61/// \param nominal Nominal value of the function.
62/// \param lowSet Set of down variations.
63/// \param highSet Set of up variations.
64/// \param paramSet Parameters that control the interpolation.
65/// \param takeOwnership If true, the PiecewiseInterpolation object will take ownership of the arguments in the low, high and parameter sets.
66PiecewiseInterpolation::PiecewiseInterpolation(const char *name, const char *title, const RooAbsReal &nominal,
67 const RooArgList &lowSet, const RooArgList &highSet,
68 const RooArgList &paramSet)
69 : RooAbsReal(name, title),
70 _normIntMgr(this),
71 _nominal("!nominal", "nominal value", this, (RooAbsReal &)nominal),
72 _lowSet("!lowSet", "low-side variation", this),
73 _highSet("!highSet", "high-side variation", this),
74 _paramSet("!paramSet", "high-side variation", this),
75 _positiveDefinite(false)
76
77{
78 // KC: check both sizes
79 if (lowSet.size() != highSet.size()) {
80 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: input lists should be of equal length" << endl ;
82 }
83
84 for (auto *comp : lowSet) {
85 if (!dynamic_cast<RooAbsReal*>(comp)) {
86 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
87 << " in first list is not of type RooAbsReal" << endl ;
89 }
90 _lowSet.add(*comp) ;
91 }
92
93
94 for (auto *comp : highSet) {
95 if (!dynamic_cast<RooAbsReal*>(comp)) {
96 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
97 << " in first list is not of type RooAbsReal" << endl ;
99 }
100 _highSet.add(*comp) ;
101 }
102
103
104 for (auto *comp : paramSet) {
105 if (!dynamic_cast<RooAbsReal*>(comp)) {
106 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
107 << " in first list is not of type RooAbsReal" << endl ;
109 }
110 _paramSet.add(*comp) ;
111 _interpCode.push_back(0); // default code: linear interpolation
112 }
113
114
115 // Choose special integrator by default
116 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
118}
119
120
121
122////////////////////////////////////////////////////////////////////////////////
123/// Copy constructor
124
126 RooAbsReal(other, name),
127 _normIntMgr(other._normIntMgr, this),
128 _nominal("!nominal",this,other._nominal),
129 _lowSet("!lowSet",this,other._lowSet),
130 _highSet("!highSet",this,other._highSet),
131 _paramSet("!paramSet",this,other._paramSet),
132 _positiveDefinite(other._positiveDefinite),
133 _interpCode(other._interpCode)
134{
135 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
137}
138
139
140
141////////////////////////////////////////////////////////////////////////////////
142/// Destructor
143
145{
147}
148
149
150
151
152////////////////////////////////////////////////////////////////////////////////
153/// Calculate and return current value of self
154
156{
157 ///////////////////
158 double nominal = _nominal;
159 double sum(nominal) ;
160
161 for (unsigned int i=0; i < _paramSet.size(); ++i) {
162 auto param = static_cast<RooAbsReal*>(_paramSet.at(i));
163 auto low = static_cast<RooAbsReal*>(_lowSet.at(i));
164 auto high = static_cast<RooAbsReal*>(_highSet.at(i));
165 Int_t icode = _interpCode[i] ;
166
167 if(icode < 0 || icode > 5) {
168 coutE(InputArguments) << "PiecewiseInterpolation::evaluate ERROR: " << param->GetName()
169 << " with unknown interpolation code" << icode << endl ;
170 }
171 sum += flexibleInterp(icode, low->getVal(), high->getVal(), 1.0, nominal, param->getVal(), sum);
172 }
173
174 if(_positiveDefinite && (sum<0)){
175 sum = 0;
176 // cout <<"sum < 0 forcing positive definite"<<endl;
177 // int code = 1;
178 // RooArgSet* myset = new RooArgSet();
179 // cout << "integral = " << analyticalIntegralWN(code, myset) << endl;
180 } else if(sum<0){
181 cxcoutD(Tracing) <<"PiecewiseInterpolation::evaluate - sum < 0, not forcing positive definite"<<endl;
182 }
183 return sum;
184
185}
186
188{
189 unsigned int n = _interpCode.size();
190
191 std::string resName = "total_" + ctx.getTmpVarName();
192 for (std::size_t i = 0; i < n; ++i) {
193 if (_interpCode[i] < 0 || _interpCode[i] > 5) {
194 coutE(InputArguments) << "PiecewiseInterpolation::evaluate ERROR: " << _paramSet[i].GetName()
195 << " with unknown interpolation code" << _interpCode[i] << endl;
196 }
197 if (_interpCode[i] != _interpCode[0]) {
198 coutE(InputArguments) << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
199 "different interpolation codes for the same class object "
200 << endl;
201 }
202 }
203 std::string funcCall;
204 funcCall = ctx.buildCall("RooFit::Detail::EvaluateFuncs::piecewiseInterpolationEvaluate", _interpCode[0], _lowSet,
206 std::string code = "double " + resName + " = " + funcCall + ";\n";
207
209 code += resName + " = " + resName + " < 0 ? 0 : " + resName + ";\n";
210
211 ctx.addToCodeBody(this, code);
212 ctx.addResult(this, resName);
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Interpolate between input distributions for all values of the observable in `evalData`.
217/// \param[in,out] evalData Struct holding spans pointing to input data. The results of this function will be stored here.
218/// \param[in] normSet Arguments to normalise over.
220{
221 std::span<double> sum = ctx.output();
222
223 auto nominal = ctx.at(_nominal);
224 for(unsigned int j=0; j < nominal.size(); ++j) {
225 sum[j] = nominal[j];
226 }
227
228 for (unsigned int i=0; i < _paramSet.size(); ++i) {
229 const double param = static_cast<RooAbsReal*>(_paramSet.at(i))->getVal();
230 auto low = ctx.at(_lowSet.at(i));
231 auto high = ctx.at(_highSet.at(i));
232 const int icode = _interpCode[i];
233
234 if (icode < 0 || icode > 5) {
235 coutE(InputArguments) << "PiecewiseInterpolation::doEval(): " << _paramSet[i].GetName()
236 << " with unknown interpolation code" << icode << std::endl;
237 throw std::invalid_argument("PiecewiseInterpolation::doEval() got invalid interpolation code " + std::to_string(icode));
238 }
239
240 for (unsigned int j=0; j < nominal.size(); ++j) {
241 sum[j] += flexibleInterp(icode, low[j], high[j], 1.0, nominal[j], param, sum[j]);
242 }
243 }
244
245 if (_positiveDefinite) {
246 for(unsigned int j=0; j < nominal.size(); ++j) {
247 if (sum[j] < 0.)
248 sum[j] = 0.;
249 }
250 }
251}
252
253////////////////////////////////////////////////////////////////////////////////
254
256{
257 if(allVars.size()==1){
258 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
259 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
260 int nbins = (static_cast<RooRealVar*>(allVars.first()))->numBins();
261 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
262 return true;
263 }else{
264 cout << "Currently BinIntegrator only knows how to deal with 1-d "<<endl;
265 return false;
266 }
267 return false;
268}
269
270////////////////////////////////////////////////////////////////////////////////
271/// Advertise that all integrals can be handled internally.
272
274 const RooArgSet* normSet, const char* /*rangeName*/) const
275{
276 /*
277 cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " <<endl;
278 cout << "all vars = "<<endl;
279 allVars.Print("v");
280 cout << "anal vars = "<<endl;
281 analVars.Print("v");
282 cout << "normset vars = "<<endl;
283 if(normSet2)
284 normSet2->Print("v");
285 */
286
287
288 // Handle trivial no-integration scenario
289 if (allVars.empty()) return 0 ;
290 if (_forceNumInt) return 0 ;
291
292
293 // Force using numeric integration
294 // use special numeric integrator
295 return 0;
296
297
298 // KC: check if interCode=0 for all
299 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it) {
300 if (!_interpCode.empty() && _interpCode[it - _paramSet.begin()] != 0) {
301 // can't factorize integral
302 cout << "can't factorize integral" << endl;
303 return 0;
304 }
305 }
306
307 // Select subset of allVars that are actual dependents
308 analVars.add(allVars) ;
309
310 // Check if this configuration was created before
311 Int_t sterileIdx(-1) ;
312 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObj(normSet,&analVars,&sterileIdx)) ;
313 if (cache) {
314 return _normIntMgr.lastIndex()+1 ;
315 }
316
317 // Create new cache element
318 cache = new CacheElem ;
319
320 // Make list of function projection and normalization integrals
321 RooAbsReal *func ;
322
323 // do variations
324 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it)
325 {
326 auto i = it - _paramSet.begin();
327 func = static_cast<RooAbsReal *>(_lowSet.at(i));
328 cache->_lowIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
329
330 func = static_cast<RooAbsReal *>(_highSet.at(i));
331 cache->_highIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
332 }
333
334 // Store cache element
335 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,nullptr) ;
336
337 return code+1 ;
338}
339
340
341
342
343////////////////////////////////////////////////////////////////////////////////
344/// Implement analytical integrations by doing appropriate weighting from component integrals
345/// functions to integrators of components
346
347double PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
348{
349 /*
350 cout <<"Enter analytic Integral"<<endl;
351 printDirty(true);
352 // _nominal.arg().setDirtyInhibit(true) ;
353 _nominal.arg().setShapeDirty() ;
354 RooAbsReal* temp ;
355 RooFIter lowIter(_lowSet.fwdIterator()) ;
356 while((temp=(RooAbsReal*)lowIter.next())) {
357 // temp->setDirtyInhibit(true) ;
358 temp->setShapeDirty() ;
359 }
360 RooFIter highIter(_highSet.fwdIterator()) ;
361 while((temp=(RooAbsReal*)highIter.next())) {
362 // temp->setDirtyInhibit(true) ;
363 temp->setShapeDirty() ;
364 }
365 */
366
367 /*
368 RooAbsArg::setDirtyInhibit(true);
369 printDirty(true);
370 cout <<"done setting dirty inhibit = true"<<endl;
371
372 // old integral, only works for linear and not positive definite
373 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
374
375
376 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
377 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
378 cout <<"iset = "<<endl;
379 iset->Print("v");
380
381 double sum = 0;
382 RooArgSet* vars = getVariables();
383 vars->remove(_paramSet);
384 _paramSet.Print("v");
385 vars->Print("v");
386 if(vars->size()==1){
387 RooRealVar* obs = (RooRealVar*) vars->first();
388 for(int i=0; i<obs->numBins(); ++i){
389 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
390 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
391 cout << "obs = " << obs->getVal() << " sum = " << sum << endl;
392 }
393 } else{
394 cout <<"only know how to deal with 1 observable right now"<<endl;
395 }
396 */
397
398 /*
399 _nominal.arg().setDirtyInhibit(false) ;
400 RooFIter lowIter2(_lowSet.fwdIterator()) ;
401 while((temp=(RooAbsReal*)lowIter2.next())) {
402 temp->setDirtyInhibit(false) ;
403 }
404 RooFIter highIter2(_highSet.fwdIterator()) ;
405 while((temp=(RooAbsReal*)highIter2.next())) {
406 temp->setDirtyInhibit(false) ;
407 }
408 */
409
410 /*
411 RooAbsArg::setDirtyInhibit(false);
412 printDirty(true);
413 cout <<"done"<<endl;
414 cout << "sum = " <<sum<<endl;
415 //return sum;
416 */
417
418 // old integral, only works for linear and not positive definite
419 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObjByIndex(code-1)) ;
420 if( cache==nullptr ) {
421 std::cout << "Error: Cache Element is nullptr" << std::endl;
422 throw std::exception();
423 }
424
425 // old integral, only works for linear and not positive definite
426
427 RooAbsReal *low;
428 RooAbsReal *high;
429 double value(0);
430 double nominal(0);
431
432 // get nominal
433 int i=0;
434 for (auto funcInt : static_range_cast<RooAbsReal*>(cache->_funcIntList)) {
435 value += funcInt->getVal() ;
436 nominal = value;
437 i++;
438 }
439 if(i==0 || i>1) { cout << "problem, wrong number of nominal functions"<<endl; }
440
441 // now get low/high variations
442 // KC: old interp code with new iterator
443
444 i = 0;
445 for (auto const *param : static_range_cast<RooAbsReal *>(_paramSet)) {
446 low = static_cast<RooAbsReal *>(cache->_lowIntList.at(i));
447 high = static_cast<RooAbsReal *>(cache->_highIntList.at(i));
448
449 if(param->getVal() > 0) {
450 value += param->getVal()*(high->getVal() - nominal);
451 } else {
452 value += param->getVal()*(nominal - low->getVal());
453 }
454 ++i;
455 }
456
457 /* // MB : old bit of interpolation code
458 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
459 low = (RooAbsReal*)lowIntIter->Next() ;
460 high = (RooAbsReal*)highIntIter->Next() ;
461
462 if(param->getVal()>0) {
463 value += param->getVal()*(high->getVal() - nominal );
464 } else {
465 value += param->getVal()*(nominal - low->getVal());
466 }
467 ++i;
468 }
469 */
470
471 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
472 while( (param=(RooAbsReal*)paramIter.next()) ) {
473 low = (RooAbsReal*)lowIntIter.next() ;
474 high = (RooAbsReal*)highIntIter.next() ;
475
476 if(_interpCode.empty() || _interpCode.at(i)==0){
477 // piece-wise linear
478 if(param->getVal()>0)
479 value += param->getVal()*(high->getVal() - nominal );
480 else
481 value += param->getVal()*(nominal - low->getVal());
482 } else if(_interpCode.at(i)==1){
483 // piece-wise log
484 if(param->getVal()>=0)
485 value *= pow(high->getVal()/nominal, +param->getVal());
486 else
487 value *= pow(low->getVal()/nominal, -param->getVal());
488 } else if(_interpCode.at(i)==2){
489 // parabolic with linear
490 double a = 0.5*(high->getVal()+low->getVal())-nominal;
491 double b = 0.5*(high->getVal()-low->getVal());
492 double c = 0;
493 if(param->getVal()>1 ){
494 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
495 } else if(param->getVal()<-1 ) {
496 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
497 } else {
498 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
499 }
500 } else if(_interpCode.at(i)==3){
501 //parabolic version of log-normal
502 double a = 0.5*(high->getVal()+low->getVal())-nominal;
503 double b = 0.5*(high->getVal()-low->getVal());
504 double c = 0;
505 if(param->getVal()>1 ){
506 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
507 } else if(param->getVal()<-1 ) {
508 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
509 } else {
510 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
511 }
512
513 } else {
514 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
515 << " with unknown interpolation code" << endl ;
516 }
517 ++i;
518 }
519 */
520
521 // cout << "value = " << value <<endl;
522 return value;
523}
524
525
526////////////////////////////////////////////////////////////////////////////////
527
528void PiecewiseInterpolation::setInterpCode(RooAbsReal& param, int code, bool silent){
529 int index = _paramSet.index(&param);
530 if(index<0){
531 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
532 << " is not in list" << endl ;
533 } else {
534 if(!silent){
535 coutW(InputArguments) << "PiecewiseInterpolation::setInterpCode : " << param.GetName()
536 << " is now " << code << endl ;
537 }
538 _interpCode.at(index) = code;
539 }
540}
541
542
543////////////////////////////////////////////////////////////////////////////////
544
546 for(unsigned int i=0; i<_interpCode.size(); ++i){
547 _interpCode.at(i) = code;
548 }
549}
550
551
552////////////////////////////////////////////////////////////////////////////////
553
555 for(unsigned int i=0; i<_interpCode.size(); ++i){
556 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
557 }
558}
559
560
561////////////////////////////////////////////////////////////////////////////////
562/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
563
564std::list<double>* PiecewiseInterpolation::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
565{
566 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
567}
568
569
570////////////////////////////////////////////////////////////////////////////////
571/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
572
574{
575 return _nominal.arg().isBinnedDistribution(obs) ;
576}
577
578
579
580////////////////////////////////////////////////////////////////////////////////
581
582std::list<double>* PiecewiseInterpolation::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
583{
584 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
585}
586
587////////////////////////////////////////////////////////////////////////////////
588/// Stream an object of class PiecewiseInterpolation.
589
591{
592 if (R__b.IsReading()) {
594 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
595 if (_interpCode.empty()) _interpCode.resize(_paramSet.size());
596 } else {
598 }
599}
600
601
602/*
603////////////////////////////////////////////////////////////////////////////////
604/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
605/// product operator construction
606
607void PiecewiseInterpolation::printMetaArgs(ostream& os) const
608{
609 _lowIter->Reset() ;
610 if (_highIter) {
611 _highIter->Reset() ;
612 }
613
614 bool first(true) ;
615
616 RooAbsArg* arg1, *arg2 ;
617 if (_highSet.size()!=0) {
618
619 while((arg1=(RooAbsArg*)_lowIter->Next())) {
620 if (!first) {
621 os << " + " ;
622 } else {
623 first = false ;
624 }
625 arg2=(RooAbsArg*)_highIter->Next() ;
626 os << arg1->GetName() << " * " << arg2->GetName() ;
627 }
628
629 } else {
630
631 while((arg1=(RooAbsArg*)_lowIter->Next())) {
632 if (!first) {
633 os << " + " ;
634 } else {
635 first = false ;
636 }
637 os << arg1->GetName() ;
638 }
639
640 }
641
642 os << " " ;
643}
644
645*/
RooSetProxy _paramSet
Parameters of the test statistic (=parameters of the input function)
#define coutI(a)
#define cxcoutD(a)
#define coutW(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
#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
The PiecewiseInterpolation is a class that can morph distributions into each other,...
bool _positiveDefinite
protect against negative and 0 bins.
RooListProxy _lowSet
Low-side variation.
RooListProxy _highSet
High-side variation.
bool isBinnedDistribution(const RooArgSet &obs) const override
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
void setInterpCode(RooAbsReal &param, int code, bool silent=false)
static TClass * Class()
~PiecewiseInterpolation() override
Destructor.
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
RooObjCacheManager _normIntMgr
! The integration cache manager
bool setBinIntegrator(RooArgSet &allVars)
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise that all integrals can be handled internally.
RooListProxy _paramSet
interpolation parameters
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
RooRealProxy _nominal
The nominal value.
double evaluate() const override
Calculate and return current value of self.
void doEval(RooFit::EvalContext &) const override
Interpolate between input distributions for all values of the observable in evalData.
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...
friend void RooRefArray::Streamer(TBuffer &)
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
Storage_t::size_type size() const
TIterator begin()
RooAbsArg * first() const
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
bool setRealValue(const char *name, double newVal=0.0, bool verbose=false)
Set value of a RooAbsRealLValue stored in set with given name to newVal No error messages are printed...
TIterator end() and range-based for loops.")
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
virtual std::list< double > * binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
RooFit::OwningPtr< RooAbsReal > createIntegral(const RooArgSet &iset, 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
Create an object that represents the integral of the function over one or more observables listed in ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
bool _forceNumInt
Force numerical integration if flag set.
Definition RooAbsReal.h:546
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
virtual bool isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:353
virtual std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
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.
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
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 setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
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...
static void softAbort()
Soft abort function that interrupts macro execution but doesn't kill ROOT.
A class to maintain the context for squashing of RooFit models into code.
std::string buildCall(std::string const &funcname, Args_t const &...args)
Build the code to call the function with name funcname, passing some arguments.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
void addToCodeBody(RooAbsArg const *klass, std::string const &in)
Adds the input string to the squashed code body.
std::string getTmpVarName() const
Get a unique variable name to be used in the generated code.
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
RooCategory & method1D()
Variable that can be changed from the outside.
Definition RooRealVar.h:37
const T & arg() const
Return reference to object held in proxy.
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const Int_t n
Definition legend1.C:16
double flexibleInterp(unsigned int code, double low, double high, double boundary, double nominal, double paramVal, double res)
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345