Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
61
62////////////////////////////////////////////////////////////////////////////////
63
68
69////////////////////////////////////////////////////////////////////////////////
70/// Construct a new interpolation. The value of the function will be
71/// \f[
72/// A = \sum_i \mathrm{Interpolate}(\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i).
73/// \f]
74/// \param name Name of the object.
75/// \param title Title (for e.g. plotting)
76/// \param nominal Nominal value of the function.
77/// \param lowSet Set of down variations.
78/// \param highSet Set of up variations.
79/// \param paramSet Parameters that control the interpolation.
80PiecewiseInterpolation::PiecewiseInterpolation(const char *name, const char *title, const RooAbsReal &nominal,
81 const RooArgList &lowSet, const RooArgList &highSet,
82 const RooArgList &paramSet)
83 : RooAbsReal(name, title),
84 _normIntMgr(this),
85 _nominal("!nominal", "nominal value", this, (RooAbsReal &)nominal),
86 _lowSet("!lowSet", "low-side variation", this),
87 _highSet("!highSet", "high-side variation", this),
88 _paramSet("!paramSet", "high-side variation", this),
89 _positiveDefinite(false)
90
91{
92 // KC: check both sizes
93 if (lowSet.size() != highSet.size()) {
94 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: input lists should be of equal length" << std::endl ;
96 }
97
98 for (auto *comp : lowSet) {
99 if (!dynamic_cast<RooAbsReal*>(comp)) {
100 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
101 << " in first list is not of type RooAbsReal" << std::endl ;
103 }
104 _lowSet.add(*comp) ;
105 }
106
107
108 for (auto *comp : highSet) {
109 if (!dynamic_cast<RooAbsReal*>(comp)) {
110 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
111 << " in first list is not of type RooAbsReal" << std::endl ;
113 }
114 _highSet.add(*comp) ;
115 }
116
117
118 for (auto *comp : paramSet) {
119 if (!dynamic_cast<RooAbsReal*>(comp)) {
120 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
121 << " in first list is not of type RooAbsReal" << std::endl ;
123 }
124 _paramSet.add(*comp) ;
125 _interpCode.push_back(0); // default code: linear interpolation
126 }
127
128
129 // Choose special integrator by default
130 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
132}
133
134
135
136////////////////////////////////////////////////////////////////////////////////
137/// Copy constructor
138
141 _normIntMgr(other._normIntMgr, this),
142 _nominal("!nominal",this,other._nominal),
143 _lowSet("!lowSet",this,other._lowSet),
144 _highSet("!highSet",this,other._highSet),
145 _paramSet("!paramSet",this,other._paramSet),
146 _positiveDefinite(other._positiveDefinite),
147 _interpCode(other._interpCode)
148{
149 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
151}
152
153
154
155////////////////////////////////////////////////////////////////////////////////
156/// Destructor
157
162
163
164
165
166////////////////////////////////////////////////////////////////////////////////
167/// Calculate and return current value of self
168
170{
171 ///////////////////
172 double nominal = _nominal;
173 double sum(nominal) ;
174
175 for (unsigned int i=0; i < _paramSet.size(); ++i) {
176 auto param = static_cast<RooAbsReal*>(_paramSet.at(i));
177 auto low = static_cast<RooAbsReal*>(_lowSet.at(i));
178 auto high = static_cast<RooAbsReal*>(_highSet.at(i));
180 sum += flexibleInterpSingle(_interpCode[i], low->getVal(), high->getVal(), 1.0, nominal, param->getVal(), sum);
181 }
182
183 if(_positiveDefinite && (sum<0)){
184 sum = 0;
185 // std::cout <<"sum < 0 forcing positive definite"<< std::endl;
186 // int code = 1;
187 // RooArgSet* myset = new RooArgSet();
188 // std::cout << "integral = " << analyticalIntegralWN(code, myset) << std::endl;
189 } else if(sum<0){
190 cxcoutD(Tracing) <<"PiecewiseInterpolation::evaluate - sum < 0, not forcing positive definite"<< std::endl;
191 }
192 return sum;
193
194}
195
196namespace {
197
198inline double broadcast(std::span<const double> const &s, std::size_t i)
199{
200 return s.size() > 1 ? s[i] : s[0];
201}
202
203} // namespace
204
205////////////////////////////////////////////////////////////////////////////////
206/// Interpolate between input distributions for all values of the observable in `evalData`.
207/// \param[in,out] ctx Struct holding spans pointing to input data. The results of this function will be stored here.
209{
210 std::span<double> sum = ctx.output();
211
212 auto nominal = ctx.at(_nominal);
213
214 for (std::size_t j = 0; j < sum.size(); ++j) {
215 sum[j] = broadcast(nominal, j);
216 }
217
218 for (unsigned int i = 0; i < _paramSet.size(); ++i) {
219 auto param = ctx.at(_paramSet.at(i));
220 auto low = ctx.at(_lowSet.at(i));
221 auto high = ctx.at(_highSet.at(i));
222
223 for (std::size_t j = 0; j < sum.size(); ++j) {
225 sum[j] += flexibleInterpSingle(_interpCode[i], broadcast(low, j), broadcast(high, j), 1.0, broadcast(nominal, j),
226 broadcast(param, j), sum[j]);
227 }
228 }
229
230 if (_positiveDefinite) {
231 for (std::size_t j = 0; j < sum.size(); ++j) {
232 if (sum[j] < 0.)
233 sum[j] = 0.;
234 }
235 }
236}
237
238////////////////////////////////////////////////////////////////////////////////
239
241{
242 if(allVars.size()==1){
243 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
244 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
245 int nbins = (static_cast<RooRealVar*>(allVars.first()))->numBins();
246 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
247 return true;
248 }else{
249 std::cout << "Currently BinIntegrator only knows how to deal with 1-d "<< std::endl;
250 return false;
251 }
252 return false;
253}
254
255////////////////////////////////////////////////////////////////////////////////
256/// Advertise that all integrals can be handled internally.
257
259 const RooArgSet* normSet, const char* /*rangeName*/) const
260{
261 /*
262 std::cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " << std::endl;
263 std::cout << "all vars = "<< std::endl;
264 allVars.Print("v");
265 std::cout << "anal vars = "<< std::endl;
266 analVars.Print("v");
267 std::cout << "normset vars = "<< std::endl;
268 if(normSet2)
269 normSet2->Print("v");
270 */
271
272
273 // Handle trivial no-integration scenario
274 if (allVars.empty()) return 0 ;
275 if (_forceNumInt) return 0 ;
276
277
278 // Force using numeric integration
279 // use special numeric integrator
280 return 0;
281
282
283 // KC: check if interCode=0 for all
284 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it) {
285 if (!_interpCode.empty() && _interpCode[it - _paramSet.begin()] != 0) {
286 // can't factorize integral
287 std::cout << "can't factorize integral" << std::endl;
288 return 0;
289 }
290 }
291
292 // Select subset of allVars that are actual dependents
293 analVars.add(allVars) ;
294
295 // Check if this configuration was created before
296 Int_t sterileIdx(-1) ;
297 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObj(normSet,&analVars,&sterileIdx)) ;
298 if (cache) {
299 return _normIntMgr.lastIndex()+1 ;
300 }
301
302 // Create new cache element
303 cache = new CacheElem ;
304
305 // Make list of function projection and normalization integrals
306 RooAbsReal *func ;
307
308 // do variations
309 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it)
310 {
311 auto i = it - _paramSet.begin();
312 func = static_cast<RooAbsReal *>(_lowSet.at(i));
313 cache->_lowIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
314
315 func = static_cast<RooAbsReal *>(_highSet.at(i));
316 cache->_highIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
317 }
318
319 // Store cache element
320 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,nullptr) ;
321
322 return code+1 ;
323}
324
325
326
327
328////////////////////////////////////////////////////////////////////////////////
329/// Implement analytical integrations by doing appropriate weighting from component integrals
330/// functions to integrators of components
331
332double PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
333{
334 /*
335 std::cout <<"Enter analytic Integral"<< std::endl;
336 printDirty(true);
337 // _nominal.arg().setDirtyInhibit(true) ;
338 _nominal.arg().setShapeDirty() ;
339 RooAbsReal* temp ;
340 RooFIter lowIter(_lowSet.fwdIterator()) ;
341 while((temp=(RooAbsReal*)lowIter.next())) {
342 // temp->setDirtyInhibit(true) ;
343 temp->setShapeDirty() ;
344 }
345 RooFIter highIter(_highSet.fwdIterator()) ;
346 while((temp=(RooAbsReal*)highIter.next())) {
347 // temp->setDirtyInhibit(true) ;
348 temp->setShapeDirty() ;
349 }
350 */
351
352 /*
353 RooAbsArg::setDirtyInhibit(true);
354 printDirty(true);
355 std::cout <<"done setting dirty inhibit = true"<< std::endl;
356
357 // old integral, only works for linear and not positive definite
358 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
359
360
361 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
362 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
363 std::cout <<"iset = "<< std::endl;
364 iset->Print("v");
365
366 double sum = 0;
367 RooArgSet* vars = getVariables();
368 vars->remove(_paramSet);
369 _paramSet.Print("v");
370 vars->Print("v");
371 if(vars->size()==1){
372 RooRealVar* obs = (RooRealVar*) vars->first();
373 for(int i=0; i<obs->numBins(); ++i){
374 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
375 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
376 std::cout << "obs = " << obs->getVal() << " sum = " << sum << std::endl;
377 }
378 } else{
379 std::cout <<"only know how to deal with 1 observable right now"<< std::endl;
380 }
381 */
382
383 /*
384 _nominal.arg().setDirtyInhibit(false) ;
385 RooFIter lowIter2(_lowSet.fwdIterator()) ;
386 while((temp=(RooAbsReal*)lowIter2.next())) {
387 temp->setDirtyInhibit(false) ;
388 }
389 RooFIter highIter2(_highSet.fwdIterator()) ;
390 while((temp=(RooAbsReal*)highIter2.next())) {
391 temp->setDirtyInhibit(false) ;
392 }
393 */
394
395 /*
396 RooAbsArg::setDirtyInhibit(false);
397 printDirty(true);
398 std::cout <<"done"<< std::endl;
399 std::cout << "sum = " <<sum<< std::endl;
400 //return sum;
401 */
402
403 // old integral, only works for linear and not positive definite
404 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObjByIndex(code-1)) ;
405 if( cache==nullptr ) {
406 std::cout << "Error: Cache Element is nullptr" << std::endl;
407 throw std::exception();
408 }
409
410 // old integral, only works for linear and not positive definite
411
412 RooAbsReal *low;
413 RooAbsReal *high;
414 double value(0);
415 double nominal(0);
416
417 // get nominal
418 int i=0;
420 value += funcInt->getVal() ;
421 nominal = value;
422 i++;
423 }
424 if(i==0 || i>1) { std::cout << "problem, wrong number of nominal functions"<< std::endl; }
425
426 // now get low/high variations
427 // KC: old interp code with new iterator
428
429 i = 0;
430 for (auto const *param : static_range_cast<RooAbsReal *>(_paramSet)) {
431 low = static_cast<RooAbsReal *>(cache->_lowIntList.at(i));
432 high = static_cast<RooAbsReal *>(cache->_highIntList.at(i));
433
434 if(param->getVal() > 0) {
435 value += param->getVal()*(high->getVal() - nominal);
436 } else {
437 value += param->getVal()*(nominal - low->getVal());
438 }
439 ++i;
440 }
441
442 /* // MB : old bit of interpolation code
443 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
444 low = (RooAbsReal*)lowIntIter->Next() ;
445 high = (RooAbsReal*)highIntIter->Next() ;
446
447 if(param->getVal()>0) {
448 value += param->getVal()*(high->getVal() - nominal );
449 } else {
450 value += param->getVal()*(nominal - low->getVal());
451 }
452 ++i;
453 }
454 */
455
456 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
457 while( (param=(RooAbsReal*)paramIter.next()) ) {
458 low = (RooAbsReal*)lowIntIter.next() ;
459 high = (RooAbsReal*)highIntIter.next() ;
460
461 if(_interpCode.empty() || _interpCode.at(i)==0){
462 // piece-wise linear
463 if(param->getVal()>0)
464 value += param->getVal()*(high->getVal() - nominal );
465 else
466 value += param->getVal()*(nominal - low->getVal());
467 } else if(_interpCode.at(i)==1){
468 // piece-wise log
469 if(param->getVal()>=0)
470 value *= pow(high->getVal()/nominal, +param->getVal());
471 else
472 value *= pow(low->getVal()/nominal, -param->getVal());
473 } else if(_interpCode.at(i)==2){
474 // parabolic with linear
475 double a = 0.5*(high->getVal()+low->getVal())-nominal;
476 double b = 0.5*(high->getVal()-low->getVal());
477 double c = 0;
478 if(param->getVal()>1 ){
479 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
480 } else if(param->getVal()<-1 ) {
481 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
482 } else {
483 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
484 }
485 } else if(_interpCode.at(i)==3){
486 //parabolic version of log-normal
487 double a = 0.5*(high->getVal()+low->getVal())-nominal;
488 double b = 0.5*(high->getVal()-low->getVal());
489 double c = 0;
490 if(param->getVal()>1 ){
491 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
492 } else if(param->getVal()<-1 ) {
493 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
494 } else {
495 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
496 }
497
498 } else {
499 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
500 << " with unknown interpolation code" << std::endl ;
501 }
502 ++i;
503 }
504 */
505
506 // std::cout << "value = " << value << std::endl;
507 return value;
508}
509
510void PiecewiseInterpolation::setInterpCode(RooAbsReal &param, int code, bool /*silent*/)
511{
512 int index = _paramSet.index(&param);
513 if (index < 0) {
514 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName() << " is not in list"
515 << std::endl;
516 return;
517 }
519}
520
522{
523 for (std::size_t i = 0; i < _interpCode.size(); ++i) {
524 setInterpCodeForParam(i, code);
525 }
526}
527
529{
530 RooAbsArg const &param = _paramSet[iParam];
531 if (code < 0 || code > 6) {
532 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
533 << " with unknown interpolation code " << code << ", keeping current code "
534 << _interpCode[iParam] << std::endl;
535 return;
536 }
537 if (code == 3) {
538 // In the past, code 3 was equivalent to code 2, which confused users.
539 // Now, we just say that code 3 doesn't exist and default to code 2 in
540 // that case for backwards compatible behavior.
541 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
542 << " with unknown interpolation code " << code << ", defaulting to code 2" << std::endl;
543 code = 2;
544 }
545 _interpCode.at(iParam) = code;
547}
548
549////////////////////////////////////////////////////////////////////////////////
550
552 for(unsigned int i=0; i<_interpCode.size(); ++i){
553 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) << std::endl;
554 }
555}
556
557
558////////////////////////////////////////////////////////////////////////////////
559/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
560
561std::list<double>* PiecewiseInterpolation::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
562{
563 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
564}
565
566
567////////////////////////////////////////////////////////////////////////////////
568/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
569
571{
572 return _nominal.arg().isBinnedDistribution(obs) ;
573}
574
575
576
577////////////////////////////////////////////////////////////////////////////////
578
579std::list<double>* PiecewiseInterpolation::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
580{
581 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
582}
583
584////////////////////////////////////////////////////////////////////////////////
585/// Stream an object of class PiecewiseInterpolation.
586
588{
589 if (R__b.IsReading()) {
590 R__b.ReadClassBuffer(PiecewiseInterpolation::Class(),this);
591 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
592 if (_interpCode.empty()) _interpCode.resize(_paramSet.size());
593 } else {
594 R__b.WriteClassBuffer(PiecewiseInterpolation::Class(),this);
595 }
596}
597
598
599/*
600////////////////////////////////////////////////////////////////////////////////
601/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
602/// product operator construction
603
604void PiecewiseInterpolation::printMetaArgs(ostream& os) const
605{
606 _lowIter->Reset() ;
607 if (_highIter) {
608 _highIter->Reset() ;
609 }
610
611 bool first(true) ;
612
613 RooAbsArg* arg1, *arg2 ;
614 if (_highSet.size()!=0) {
615
616 while((arg1=(RooAbsArg*)_lowIter->Next())) {
617 if (!first) {
618 os << " + " ;
619 } else {
620 first = false ;
621 }
622 arg2=(RooAbsArg*)_highIter->Next() ;
623 os << arg1->GetName() << " * " << arg2->GetName() ;
624 }
625
626 } else {
627
628 while((arg1=(RooAbsArg*)_lowIter->Next())) {
629 if (!first) {
630 os << " + " ;
631 } else {
632 first = false ;
633 }
634 os << arg1->GetName() ;
635 }
636
637 }
638
639 os << " " ;
640}
641
642*/
#define coutI(a)
#define cxcoutD(a)
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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.
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.
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
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.
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
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
double flexibleInterpSingle(unsigned int code, double low, double high, double boundary, double nominal, double paramVal, double res)
Definition MathFuncs.h:232
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345