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$ \theta_i \f$:
7* \f[
8* A = \mathrm{nominal} + \sum_i I_i(\theta_i;\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i).
9* \f]
10* for additive interpolation modes (interpCodes 0, 2, 3, and 4), or:
11* \f[
12* A = \mathrm{nominal}\prod_i I_i(\theta_i;\mathrm{low}_i/\mathrm{nominal}, 1, \mathrm{high}_i/\mathrm{nominal}).
13* \f]
14* for multiplicative interpolation modes (interpCodes 1, 5, and 6). The interpCodes determine the function \f$ I_i \f$ (see table below).
15*
16* Note that a PiecewiseInterpolation with \f$ \mathrm{nominal}=1 \f$, N variations, and a multiplicative interpolation mode is equivalent to N
17* PiecewiseInterpolations each with a single variation and the same interpolation code, all inside a RooProduct.
18*
19* If an \f$ \theta_i \f$ is zero, the distribution is identical to the nominal distribution, at
20* \f$ \pm 1 \f$ it is identical to the up/down distribution for that specific \f$ i \f$.
21*
22* PiecewiseInterpolation will behave identically (except for differences in the interpCode assignments) to a FlexibleInterpVar if both its nominal, and high and low variation sets
23* are all RooRealVar.
24*
25* The class supports several interpolation methods, which can be selected for each parameter separately
26* using setInterpCode(). The default interpolation code is 0. The table below provides details of the interpCodes:
27
28| **interpCode** | **Name** | **Description** |
29|----------------|----------|-----------------|
30| 0 (default) | Additive Piecewise Linear | \f$ I_0(\theta;x_{-},x_0,x_{+}) = \theta(x_{+} - x_0) \f$ for \f$ \theta>=0 \f$, otherwise \f$ \theta(x_0 - x_{-}) \f$. Not recommended except if using a symmetric variation, because of discontinuities in derivatives. |
31| 1 | Multiplicative Piecewise Exponential | \f$ I_1(\theta;x_{-},x_0,x_{+}) = (x_{+}/x_0)^{\theta} \f$ for \f$ \theta>=0 \f$, otherwise \f$ (x_{-}/x_0)^{-\theta} \f$. |
32| 2 | Additive Quadratic Interp. + Linear Extrap. | Deprecated by interpCode 4. |
33| 4 | Additive Poly Interp. + Linear Extrap. | \f$ I_4(\theta;x_{-},x_0,x_{+}) = I_0(\theta;x_{-},x_0,x_{+}) \f$ if \f$ |\theta|>=1 \f$, otherwise \f$ \theta(\frac{x_{+}-x_{-}}{2}+\theta\frac{x_{+}+x_{-}-2x_{0}}{16}(15+\theta^2(3\alpha^2-10))) \f$ (6th-order polynomial through origin for with matching 0th,1st,2nd derivatives at boundary). |
34| 5 | Multiplicative Poly Interp. + Exponential Extrap. | \f$ I_5(\theta;x_{-},x_0,x_{+}) = I_1(\theta;x_{-},x_0,x_{+}) \f$ if \f$ |\theta|>=1 \f$, otherwise 6th-order polynomial for \f$ |\theta_i|<1 \f$ with matching 0th,1st,2nd derivatives at boundary. Recommended for normalization factors. In FlexibleInterpVar this is interpCode=4. |
35| 6 | Multiplicative Poly Interp. + Linear Extrap. | \f$ I_6(\theta;x_{-},x_0,x_{+}) = 1+I_4(\theta;x_{-},x_0,x_{+}). \f$ Recommended for normalization factors that must not have roots (i.e. be equal to 0) outside of \f$ |\theta_i|<1 \f$. |
36
37*/
38
40
42
43#include "Riostream.h"
44#include "TBuffer.h"
45
46#include "RooAbsReal.h"
47#include "RooAbsPdf.h"
48#include "RooErrorHandler.h"
49#include "RooArgSet.h"
50#include "RooRealVar.h"
51#include "RooMsgService.h"
52#include "RooNumIntConfig.h"
53#include "RooTrace.h"
54#include "RooDataHist.h"
55#include "RooHistFunc.h"
56
57#include <exception>
58#include <cmath>
59#include <algorithm>
60
61using std::endl, std::cout;
62
64
65////////////////////////////////////////////////////////////////////////////////
66
68{
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Construct a new interpolation. The value of the function will be
74/// \f[
75/// A = \sum_i \mathrm{Interpolate}(\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i).
76/// \f]
77/// \param name Name of the object.
78/// \param title Title (for e.g. plotting)
79/// \param nominal Nominal value of the function.
80/// \param lowSet Set of down variations.
81/// \param highSet Set of up variations.
82/// \param paramSet Parameters that control the interpolation.
83PiecewiseInterpolation::PiecewiseInterpolation(const char *name, const char *title, const RooAbsReal &nominal,
84 const RooArgList &lowSet, const RooArgList &highSet,
85 const RooArgList &paramSet)
86 : RooAbsReal(name, title),
87 _normIntMgr(this),
88 _nominal("!nominal", "nominal value", this, (RooAbsReal &)nominal),
89 _lowSet("!lowSet", "low-side variation", this),
90 _highSet("!highSet", "high-side variation", this),
91 _paramSet("!paramSet", "high-side variation", this),
92 _positiveDefinite(false)
93
94{
95 // KC: check both sizes
96 if (lowSet.size() != highSet.size()) {
97 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: input lists should be of equal length" << endl ;
99 }
100
101 for (auto *comp : lowSet) {
102 if (!dynamic_cast<RooAbsReal*>(comp)) {
103 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
104 << " in first list is not of type RooAbsReal" << endl ;
106 }
107 _lowSet.add(*comp) ;
108 }
109
110
111 for (auto *comp : highSet) {
112 if (!dynamic_cast<RooAbsReal*>(comp)) {
113 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
114 << " in first list is not of type RooAbsReal" << endl ;
116 }
117 _highSet.add(*comp) ;
118 }
119
120
121 for (auto *comp : paramSet) {
122 if (!dynamic_cast<RooAbsReal*>(comp)) {
123 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
124 << " in first list is not of type RooAbsReal" << endl ;
126 }
127 _paramSet.add(*comp) ;
128 _interpCode.push_back(0); // default code: linear interpolation
129 }
130
131
132 // Choose special integrator by default
133 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
135}
136
137
138
139////////////////////////////////////////////////////////////////////////////////
140/// Copy constructor
141
143 RooAbsReal(other, name),
144 _normIntMgr(other._normIntMgr, this),
145 _nominal("!nominal",this,other._nominal),
146 _lowSet("!lowSet",this,other._lowSet),
147 _highSet("!highSet",this,other._highSet),
148 _paramSet("!paramSet",this,other._paramSet),
149 _positiveDefinite(other._positiveDefinite),
150 _interpCode(other._interpCode)
151{
152 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
154}
155
156
157
158////////////////////////////////////////////////////////////////////////////////
159/// Destructor
160
162{
164}
165
166
167
168
169////////////////////////////////////////////////////////////////////////////////
170/// Calculate and return current value of self
171
173{
174 ///////////////////
175 double nominal = _nominal;
176 double sum(nominal) ;
177
178 for (unsigned int i=0; i < _paramSet.size(); ++i) {
179 auto param = static_cast<RooAbsReal*>(_paramSet.at(i));
180 auto low = static_cast<RooAbsReal*>(_lowSet.at(i));
181 auto high = static_cast<RooAbsReal*>(_highSet.at(i));
183 sum += flexibleInterpSingle(_interpCode[i], low->getVal(), high->getVal(), 1.0, nominal, param->getVal(), sum);
184 }
185
186 if(_positiveDefinite && (sum<0)){
187 sum = 0;
188 // cout <<"sum < 0 forcing positive definite"<<endl;
189 // int code = 1;
190 // RooArgSet* myset = new RooArgSet();
191 // cout << "integral = " << analyticalIntegralWN(code, myset) << endl;
192 } else if(sum<0){
193 cxcoutD(Tracing) <<"PiecewiseInterpolation::evaluate - sum < 0, not forcing positive definite"<<endl;
194 }
195 return sum;
196
197}
198
199namespace {
200
201inline double broadcast(std::span<const double> const &s, std::size_t i)
202{
203 return s.size() > 1 ? s[i] : s[0];
204}
205
206} // namespace
207
208////////////////////////////////////////////////////////////////////////////////
209/// Interpolate between input distributions for all values of the observable in `evalData`.
210/// \param[in,out] evalData Struct holding spans pointing to input data. The results of this function will be stored here.
211/// \param[in] normSet Arguments to normalise over.
213{
214 std::span<double> sum = ctx.output();
215
216 auto nominal = ctx.at(_nominal);
217
218 for (std::size_t j = 0; j < sum.size(); ++j) {
219 sum[j] = broadcast(nominal, j);
220 }
221
222 for (unsigned int i = 0; i < _paramSet.size(); ++i) {
223 auto param = ctx.at(_paramSet.at(i));
224 auto low = ctx.at(_lowSet.at(i));
225 auto high = ctx.at(_highSet.at(i));
226
227 for (std::size_t j = 0; j < sum.size(); ++j) {
229 sum[j] += flexibleInterpSingle(_interpCode[i], broadcast(low, j), broadcast(high, j), 1.0, broadcast(nominal, j),
230 broadcast(param, j), sum[j]);
231 }
232 }
233
234 if (_positiveDefinite) {
235 for (std::size_t j = 0; j < sum.size(); ++j) {
236 if (sum[j] < 0.)
237 sum[j] = 0.;
238 }
239 }
240}
241
242////////////////////////////////////////////////////////////////////////////////
243
245{
246 if(allVars.size()==1){
247 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
248 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
249 int nbins = (static_cast<RooRealVar*>(allVars.first()))->numBins();
250 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
251 return true;
252 }else{
253 cout << "Currently BinIntegrator only knows how to deal with 1-d "<<endl;
254 return false;
255 }
256 return false;
257}
258
259////////////////////////////////////////////////////////////////////////////////
260/// Advertise that all integrals can be handled internally.
261
263 const RooArgSet* normSet, const char* /*rangeName*/) const
264{
265 /*
266 cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " <<endl;
267 cout << "all vars = "<<endl;
268 allVars.Print("v");
269 cout << "anal vars = "<<endl;
270 analVars.Print("v");
271 cout << "normset vars = "<<endl;
272 if(normSet2)
273 normSet2->Print("v");
274 */
275
276
277 // Handle trivial no-integration scenario
278 if (allVars.empty()) return 0 ;
279 if (_forceNumInt) return 0 ;
280
281
282 // Force using numeric integration
283 // use special numeric integrator
284 return 0;
285
286
287 // KC: check if interCode=0 for all
288 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it) {
289 if (!_interpCode.empty() && _interpCode[it - _paramSet.begin()] != 0) {
290 // can't factorize integral
291 cout << "can't factorize integral" << endl;
292 return 0;
293 }
294 }
295
296 // Select subset of allVars that are actual dependents
297 analVars.add(allVars) ;
298
299 // Check if this configuration was created before
300 Int_t sterileIdx(-1) ;
301 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObj(normSet,&analVars,&sterileIdx)) ;
302 if (cache) {
303 return _normIntMgr.lastIndex()+1 ;
304 }
305
306 // Create new cache element
307 cache = new CacheElem ;
308
309 // Make list of function projection and normalization integrals
310 RooAbsReal *func ;
311
312 // do variations
313 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it)
314 {
315 auto i = it - _paramSet.begin();
316 func = static_cast<RooAbsReal *>(_lowSet.at(i));
317 cache->_lowIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
318
319 func = static_cast<RooAbsReal *>(_highSet.at(i));
320 cache->_highIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
321 }
322
323 // Store cache element
324 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,nullptr) ;
325
326 return code+1 ;
327}
328
329
330
331
332////////////////////////////////////////////////////////////////////////////////
333/// Implement analytical integrations by doing appropriate weighting from component integrals
334/// functions to integrators of components
335
336double PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
337{
338 /*
339 cout <<"Enter analytic Integral"<<endl;
340 printDirty(true);
341 // _nominal.arg().setDirtyInhibit(true) ;
342 _nominal.arg().setShapeDirty() ;
343 RooAbsReal* temp ;
344 RooFIter lowIter(_lowSet.fwdIterator()) ;
345 while((temp=(RooAbsReal*)lowIter.next())) {
346 // temp->setDirtyInhibit(true) ;
347 temp->setShapeDirty() ;
348 }
349 RooFIter highIter(_highSet.fwdIterator()) ;
350 while((temp=(RooAbsReal*)highIter.next())) {
351 // temp->setDirtyInhibit(true) ;
352 temp->setShapeDirty() ;
353 }
354 */
355
356 /*
357 RooAbsArg::setDirtyInhibit(true);
358 printDirty(true);
359 cout <<"done setting dirty inhibit = true"<<endl;
360
361 // old integral, only works for linear and not positive definite
362 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
363
364
365 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
366 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
367 cout <<"iset = "<<endl;
368 iset->Print("v");
369
370 double sum = 0;
371 RooArgSet* vars = getVariables();
372 vars->remove(_paramSet);
373 _paramSet.Print("v");
374 vars->Print("v");
375 if(vars->size()==1){
376 RooRealVar* obs = (RooRealVar*) vars->first();
377 for(int i=0; i<obs->numBins(); ++i){
378 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
379 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
380 cout << "obs = " << obs->getVal() << " sum = " << sum << endl;
381 }
382 } else{
383 cout <<"only know how to deal with 1 observable right now"<<endl;
384 }
385 */
386
387 /*
388 _nominal.arg().setDirtyInhibit(false) ;
389 RooFIter lowIter2(_lowSet.fwdIterator()) ;
390 while((temp=(RooAbsReal*)lowIter2.next())) {
391 temp->setDirtyInhibit(false) ;
392 }
393 RooFIter highIter2(_highSet.fwdIterator()) ;
394 while((temp=(RooAbsReal*)highIter2.next())) {
395 temp->setDirtyInhibit(false) ;
396 }
397 */
398
399 /*
400 RooAbsArg::setDirtyInhibit(false);
401 printDirty(true);
402 cout <<"done"<<endl;
403 cout << "sum = " <<sum<<endl;
404 //return sum;
405 */
406
407 // old integral, only works for linear and not positive definite
408 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObjByIndex(code-1)) ;
409 if( cache==nullptr ) {
410 std::cout << "Error: Cache Element is nullptr" << std::endl;
411 throw std::exception();
412 }
413
414 // old integral, only works for linear and not positive definite
415
416 RooAbsReal *low;
417 RooAbsReal *high;
418 double value(0);
419 double nominal(0);
420
421 // get nominal
422 int i=0;
423 for (auto funcInt : static_range_cast<RooAbsReal*>(cache->_funcIntList)) {
424 value += funcInt->getVal() ;
425 nominal = value;
426 i++;
427 }
428 if(i==0 || i>1) { cout << "problem, wrong number of nominal functions"<<endl; }
429
430 // now get low/high variations
431 // KC: old interp code with new iterator
432
433 i = 0;
434 for (auto const *param : static_range_cast<RooAbsReal *>(_paramSet)) {
435 low = static_cast<RooAbsReal *>(cache->_lowIntList.at(i));
436 high = static_cast<RooAbsReal *>(cache->_highIntList.at(i));
437
438 if(param->getVal() > 0) {
439 value += param->getVal()*(high->getVal() - nominal);
440 } else {
441 value += param->getVal()*(nominal - low->getVal());
442 }
443 ++i;
444 }
445
446 /* // MB : old bit of interpolation code
447 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
448 low = (RooAbsReal*)lowIntIter->Next() ;
449 high = (RooAbsReal*)highIntIter->Next() ;
450
451 if(param->getVal()>0) {
452 value += param->getVal()*(high->getVal() - nominal );
453 } else {
454 value += param->getVal()*(nominal - low->getVal());
455 }
456 ++i;
457 }
458 */
459
460 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
461 while( (param=(RooAbsReal*)paramIter.next()) ) {
462 low = (RooAbsReal*)lowIntIter.next() ;
463 high = (RooAbsReal*)highIntIter.next() ;
464
465 if(_interpCode.empty() || _interpCode.at(i)==0){
466 // piece-wise linear
467 if(param->getVal()>0)
468 value += param->getVal()*(high->getVal() - nominal );
469 else
470 value += param->getVal()*(nominal - low->getVal());
471 } else if(_interpCode.at(i)==1){
472 // piece-wise log
473 if(param->getVal()>=0)
474 value *= pow(high->getVal()/nominal, +param->getVal());
475 else
476 value *= pow(low->getVal()/nominal, -param->getVal());
477 } else if(_interpCode.at(i)==2){
478 // parabolic with linear
479 double a = 0.5*(high->getVal()+low->getVal())-nominal;
480 double b = 0.5*(high->getVal()-low->getVal());
481 double c = 0;
482 if(param->getVal()>1 ){
483 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
484 } else if(param->getVal()<-1 ) {
485 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
486 } else {
487 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
488 }
489 } else if(_interpCode.at(i)==3){
490 //parabolic version of log-normal
491 double a = 0.5*(high->getVal()+low->getVal())-nominal;
492 double b = 0.5*(high->getVal()-low->getVal());
493 double c = 0;
494 if(param->getVal()>1 ){
495 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
496 } else if(param->getVal()<-1 ) {
497 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
498 } else {
499 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
500 }
501
502 } else {
503 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
504 << " with unknown interpolation code" << endl ;
505 }
506 ++i;
507 }
508 */
509
510 // cout << "value = " << value <<endl;
511 return value;
512}
513
514void PiecewiseInterpolation::setInterpCode(RooAbsReal &param, int code, bool /*silent*/)
515{
516 int index = _paramSet.index(&param);
517 if (index < 0) {
518 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName() << " is not in list"
519 << std::endl;
520 return;
521 }
523}
524
526{
527 for (std::size_t i = 0; i < _interpCode.size(); ++i) {
528 setInterpCodeForParam(i, code);
529 }
530}
531
533{
534 RooAbsArg const &param = _paramSet[iParam];
535 if (code < 0 || code > 6) {
536 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
537 << " with unknown interpolation code " << code << ", keeping current code "
538 << _interpCode[iParam] << std::endl;
539 return;
540 }
541 if (code == 3) {
542 // In the past, code 3 was equivalent to code 2, which confused users.
543 // Now, we just say that code 3 doesn't exist and default to code 2 in
544 // that case for backwards compatible behavior.
545 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
546 << " with unknown interpolation code " << code << ", defaulting to code 2" << std::endl;
547 code = 2;
548 }
549 _interpCode.at(iParam) = code;
551}
552
553////////////////////////////////////////////////////////////////////////////////
554
556 for(unsigned int i=0; i<_interpCode.size(); ++i){
557 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
558 }
559}
560
561
562////////////////////////////////////////////////////////////////////////////////
563/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
564
565std::list<double>* PiecewiseInterpolation::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
566{
567 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
568}
569
570
571////////////////////////////////////////////////////////////////////////////////
572/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
573
575{
576 return _nominal.arg().isBinnedDistribution(obs) ;
577}
578
579
580
581////////////////////////////////////////////////////////////////////////////////
582
583std::list<double>* PiecewiseInterpolation::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
584{
585 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
586}
587
588////////////////////////////////////////////////////////////////////////////////
589/// Stream an object of class PiecewiseInterpolation.
590
592{
593 if (R__b.IsReading()) {
595 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
596 if (_interpCode.empty()) _interpCode.resize(_paramSet.size());
597 } else {
599 }
600}
601
602
603/*
604////////////////////////////////////////////////////////////////////////////////
605/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
606/// product operator construction
607
608void PiecewiseInterpolation::printMetaArgs(ostream& os) const
609{
610 _lowIter->Reset() ;
611 if (_highIter) {
612 _highIter->Reset() ;
613 }
614
615 bool first(true) ;
616
617 RooAbsArg* arg1, *arg2 ;
618 if (_highSet.size()!=0) {
619
620 while((arg1=(RooAbsArg*)_lowIter->Next())) {
621 if (!first) {
622 os << " + " ;
623 } else {
624 first = false ;
625 }
626 arg2=(RooAbsArg*)_highIter->Next() ;
627 os << arg1->GetName() << " * " << arg2->GetName() ;
628 }
629
630 } else {
631
632 while((arg1=(RooAbsArg*)_lowIter->Next())) {
633 if (!first) {
634 os << " + " ;
635 } else {
636 first = false ;
637 }
638 os << arg1->GetName() ;
639 }
640
641 }
642
643 os << " " ;
644}
645
646*/
#define coutI(a)
#define cxcoutD(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
#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
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.
static TClass * Class()
~PiecewiseInterpolation() override
Destructor.
void setInterpCodeForParam(int iParam, int code)
void setInterpCode(RooAbsReal &param, int code, bool silent=true)
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...
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
friend void RooRefArray::Streamer(TBuffer &)
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:431
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
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.
const_iterator end() const
Storage_t::size_type size() const
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...
const_iterator begin() const
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.
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
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:348
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...
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 ...
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.
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.
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
double flexibleInterpSingle(unsigned int code, double low, double high, double boundary, double nominal, double paramVal, double res)
Definition MathFuncs.h:213
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345