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#include "RooDataHist.h"
37#include "RooHistFunc.h"
38
39#include <exception>
40#include <cmath>
41#include <algorithm>
42
43using std::endl, std::cout;
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.
65PiecewiseInterpolation::PiecewiseInterpolation(const char *name, const char *title, const RooAbsReal &nominal,
66 const RooArgList &lowSet, const RooArgList &highSet,
67 const RooArgList &paramSet)
68 : RooAbsReal(name, title),
69 _normIntMgr(this),
70 _nominal("!nominal", "nominal value", this, (RooAbsReal &)nominal),
71 _lowSet("!lowSet", "low-side variation", this),
72 _highSet("!highSet", "high-side variation", this),
73 _paramSet("!paramSet", "high-side variation", this),
74 _positiveDefinite(false)
75
76{
77 // KC: check both sizes
78 if (lowSet.size() != highSet.size()) {
79 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: input lists should be of equal length" << endl ;
81 }
82
83 for (auto *comp : lowSet) {
84 if (!dynamic_cast<RooAbsReal*>(comp)) {
85 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
86 << " in first list is not of type RooAbsReal" << endl ;
88 }
89 _lowSet.add(*comp) ;
90 }
91
92
93 for (auto *comp : highSet) {
94 if (!dynamic_cast<RooAbsReal*>(comp)) {
95 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
96 << " in first list is not of type RooAbsReal" << endl ;
98 }
99 _highSet.add(*comp) ;
100 }
101
102
103 for (auto *comp : paramSet) {
104 if (!dynamic_cast<RooAbsReal*>(comp)) {
105 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
106 << " in first list is not of type RooAbsReal" << endl ;
108 }
109 _paramSet.add(*comp) ;
110 _interpCode.push_back(0); // default code: linear interpolation
111 }
112
113
114 // Choose special integrator by default
115 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
117}
118
119
120
121////////////////////////////////////////////////////////////////////////////////
122/// Copy constructor
123
125 RooAbsReal(other, name),
126 _normIntMgr(other._normIntMgr, this),
127 _nominal("!nominal",this,other._nominal),
128 _lowSet("!lowSet",this,other._lowSet),
129 _highSet("!highSet",this,other._highSet),
130 _paramSet("!paramSet",this,other._paramSet),
131 _positiveDefinite(other._positiveDefinite),
132 _interpCode(other._interpCode)
133{
134 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
136}
137
138
139
140////////////////////////////////////////////////////////////////////////////////
141/// Destructor
142
144{
146}
147
148
149
150
151////////////////////////////////////////////////////////////////////////////////
152/// Calculate and return current value of self
153
155{
156 ///////////////////
157 double nominal = _nominal;
158 double sum(nominal) ;
159
160 for (unsigned int i=0; i < _paramSet.size(); ++i) {
161 auto param = static_cast<RooAbsReal*>(_paramSet.at(i));
162 auto low = static_cast<RooAbsReal*>(_lowSet.at(i));
163 auto high = static_cast<RooAbsReal*>(_highSet.at(i));
164 Int_t icode = _interpCode[i] ;
165
166 if(icode < 0 || icode > 5) {
167 coutE(InputArguments) << "PiecewiseInterpolation::evaluate ERROR: " << param->GetName()
168 << " with unknown interpolation code" << icode << endl ;
169 }
171 sum += flexibleInterpSingle(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 std::size_t 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
204 // The PiecewiseInterpolation class is used in the context of HistFactory
205 // models, where is is always used the same way: all RooAbsReals in _lowSet,
206 // _histSet, and also nominal are 1D RooHistFuncs with with same structure.
207 //
208 // Therefore, we can make a big optimization: we get the bin index only once
209 // here in the generated code for PiecewiseInterpolation. Then, we also
210 // rearrange the histogram data in such a way that we can always pass the
211 // same arrays to the free function that implements the interpolation, just
212 // with a dynamic offset calculated from the bin index.
213 RooDataHist const &nomHist = dynamic_cast<RooHistFunc const &>(*_nominal).dataHist();
214 int nBins = nomHist.numEntries();
215 std::vector<double> valsNominal;
216 std::vector<double> valsLow;
217 std::vector<double> valsHigh;
218 for (int i = 0; i < nBins; ++i) {
219 valsNominal.push_back(nomHist.weight(i));
220 }
221 for (int i = 0; i < nBins; ++i) {
222 for (std::size_t iParam = 0; iParam < n; ++iParam) {
223 valsLow.push_back(dynamic_cast<RooHistFunc const &>(_lowSet[iParam]).dataHist().weight(i));
224 valsHigh.push_back(dynamic_cast<RooHistFunc const &>(_highSet[iParam]).dataHist().weight(i));
225 }
226 }
227 std::string idxName =
228 nomHist.calculateTreeIndexForCodeSquash(this, ctx, dynamic_cast<RooHistFunc const &>(*_nominal).variables());
229 std::string valsNominalStr = ctx.buildArg(valsNominal);
230 std::string valsLowStr = ctx.buildArg(valsLow);
231 std::string valsHighStr = ctx.buildArg(valsHigh);
232 std::string nStr = std::to_string(n);
233 std::string code;
234
235 std::string lowName = ctx.getTmpVarName();
236 std::string highName = ctx.getTmpVarName();
237 std::string nominalName = ctx.getTmpVarName();
238 code += "double * " + lowName + " = " + valsLowStr + " + " + nStr + " * " + idxName + ";\n";
239 code += "double * " + highName + " = " + valsHighStr + " + " + nStr + " * " + idxName + ";\n";
240 code += "double " + nominalName + " = *(" + valsNominalStr + " + " + idxName + ");\n";
241
242 std::string funcCall = ctx.buildCall("RooFit::Detail::MathFuncs::flexibleInterp", _interpCode[0], _paramSet, n,
243 lowName, highName, 1.0, nominalName, 0.0);
244 code += "double " + resName + " = " + funcCall + ";\n";
245
247 code += resName + " = " + resName + " < 0 ? 0 : " + resName + ";\n";
248
249 ctx.addToCodeBody(this, code);
250 ctx.addResult(this, resName);
251}
252
253////////////////////////////////////////////////////////////////////////////////
254/// Interpolate between input distributions for all values of the observable in `evalData`.
255/// \param[in,out] evalData Struct holding spans pointing to input data. The results of this function will be stored here.
256/// \param[in] normSet Arguments to normalise over.
258{
259 std::span<double> sum = ctx.output();
260
261 auto nominal = ctx.at(_nominal);
262 for(unsigned int j=0; j < nominal.size(); ++j) {
263 sum[j] = nominal[j];
264 }
265
266 for (unsigned int i=0; i < _paramSet.size(); ++i) {
267 const double param = static_cast<RooAbsReal*>(_paramSet.at(i))->getVal();
268 auto low = ctx.at(_lowSet.at(i));
269 auto high = ctx.at(_highSet.at(i));
270 const int icode = _interpCode[i];
271
272 if (icode < 0 || icode > 5) {
273 coutE(InputArguments) << "PiecewiseInterpolation::doEval(): " << _paramSet[i].GetName()
274 << " with unknown interpolation code" << icode << std::endl;
275 throw std::invalid_argument("PiecewiseInterpolation::doEval() got invalid interpolation code " + std::to_string(icode));
276 }
277
278 for (unsigned int j=0; j < nominal.size(); ++j) {
280 sum[j] += flexibleInterpSingle(icode, low[j], high[j], 1.0, nominal[j], param, sum[j]);
281 }
282 }
283
284 if (_positiveDefinite) {
285 for(unsigned int j=0; j < nominal.size(); ++j) {
286 if (sum[j] < 0.)
287 sum[j] = 0.;
288 }
289 }
290}
291
292////////////////////////////////////////////////////////////////////////////////
293
295{
296 if(allVars.size()==1){
297 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
298 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
299 int nbins = (static_cast<RooRealVar*>(allVars.first()))->numBins();
300 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
301 return true;
302 }else{
303 cout << "Currently BinIntegrator only knows how to deal with 1-d "<<endl;
304 return false;
305 }
306 return false;
307}
308
309////////////////////////////////////////////////////////////////////////////////
310/// Advertise that all integrals can be handled internally.
311
313 const RooArgSet* normSet, const char* /*rangeName*/) const
314{
315 /*
316 cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " <<endl;
317 cout << "all vars = "<<endl;
318 allVars.Print("v");
319 cout << "anal vars = "<<endl;
320 analVars.Print("v");
321 cout << "normset vars = "<<endl;
322 if(normSet2)
323 normSet2->Print("v");
324 */
325
326
327 // Handle trivial no-integration scenario
328 if (allVars.empty()) return 0 ;
329 if (_forceNumInt) return 0 ;
330
331
332 // Force using numeric integration
333 // use special numeric integrator
334 return 0;
335
336
337 // KC: check if interCode=0 for all
338 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it) {
339 if (!_interpCode.empty() && _interpCode[it - _paramSet.begin()] != 0) {
340 // can't factorize integral
341 cout << "can't factorize integral" << endl;
342 return 0;
343 }
344 }
345
346 // Select subset of allVars that are actual dependents
347 analVars.add(allVars) ;
348
349 // Check if this configuration was created before
350 Int_t sterileIdx(-1) ;
351 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObj(normSet,&analVars,&sterileIdx)) ;
352 if (cache) {
353 return _normIntMgr.lastIndex()+1 ;
354 }
355
356 // Create new cache element
357 cache = new CacheElem ;
358
359 // Make list of function projection and normalization integrals
360 RooAbsReal *func ;
361
362 // do variations
363 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it)
364 {
365 auto i = it - _paramSet.begin();
366 func = static_cast<RooAbsReal *>(_lowSet.at(i));
367 cache->_lowIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
368
369 func = static_cast<RooAbsReal *>(_highSet.at(i));
370 cache->_highIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
371 }
372
373 // Store cache element
374 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,nullptr) ;
375
376 return code+1 ;
377}
378
379
380
381
382////////////////////////////////////////////////////////////////////////////////
383/// Implement analytical integrations by doing appropriate weighting from component integrals
384/// functions to integrators of components
385
386double PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
387{
388 /*
389 cout <<"Enter analytic Integral"<<endl;
390 printDirty(true);
391 // _nominal.arg().setDirtyInhibit(true) ;
392 _nominal.arg().setShapeDirty() ;
393 RooAbsReal* temp ;
394 RooFIter lowIter(_lowSet.fwdIterator()) ;
395 while((temp=(RooAbsReal*)lowIter.next())) {
396 // temp->setDirtyInhibit(true) ;
397 temp->setShapeDirty() ;
398 }
399 RooFIter highIter(_highSet.fwdIterator()) ;
400 while((temp=(RooAbsReal*)highIter.next())) {
401 // temp->setDirtyInhibit(true) ;
402 temp->setShapeDirty() ;
403 }
404 */
405
406 /*
407 RooAbsArg::setDirtyInhibit(true);
408 printDirty(true);
409 cout <<"done setting dirty inhibit = true"<<endl;
410
411 // old integral, only works for linear and not positive definite
412 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
413
414
415 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
416 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
417 cout <<"iset = "<<endl;
418 iset->Print("v");
419
420 double sum = 0;
421 RooArgSet* vars = getVariables();
422 vars->remove(_paramSet);
423 _paramSet.Print("v");
424 vars->Print("v");
425 if(vars->size()==1){
426 RooRealVar* obs = (RooRealVar*) vars->first();
427 for(int i=0; i<obs->numBins(); ++i){
428 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
429 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
430 cout << "obs = " << obs->getVal() << " sum = " << sum << endl;
431 }
432 } else{
433 cout <<"only know how to deal with 1 observable right now"<<endl;
434 }
435 */
436
437 /*
438 _nominal.arg().setDirtyInhibit(false) ;
439 RooFIter lowIter2(_lowSet.fwdIterator()) ;
440 while((temp=(RooAbsReal*)lowIter2.next())) {
441 temp->setDirtyInhibit(false) ;
442 }
443 RooFIter highIter2(_highSet.fwdIterator()) ;
444 while((temp=(RooAbsReal*)highIter2.next())) {
445 temp->setDirtyInhibit(false) ;
446 }
447 */
448
449 /*
450 RooAbsArg::setDirtyInhibit(false);
451 printDirty(true);
452 cout <<"done"<<endl;
453 cout << "sum = " <<sum<<endl;
454 //return sum;
455 */
456
457 // old integral, only works for linear and not positive definite
458 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObjByIndex(code-1)) ;
459 if( cache==nullptr ) {
460 std::cout << "Error: Cache Element is nullptr" << std::endl;
461 throw std::exception();
462 }
463
464 // old integral, only works for linear and not positive definite
465
466 RooAbsReal *low;
467 RooAbsReal *high;
468 double value(0);
469 double nominal(0);
470
471 // get nominal
472 int i=0;
473 for (auto funcInt : static_range_cast<RooAbsReal*>(cache->_funcIntList)) {
474 value += funcInt->getVal() ;
475 nominal = value;
476 i++;
477 }
478 if(i==0 || i>1) { cout << "problem, wrong number of nominal functions"<<endl; }
479
480 // now get low/high variations
481 // KC: old interp code with new iterator
482
483 i = 0;
484 for (auto const *param : static_range_cast<RooAbsReal *>(_paramSet)) {
485 low = static_cast<RooAbsReal *>(cache->_lowIntList.at(i));
486 high = static_cast<RooAbsReal *>(cache->_highIntList.at(i));
487
488 if(param->getVal() > 0) {
489 value += param->getVal()*(high->getVal() - nominal);
490 } else {
491 value += param->getVal()*(nominal - low->getVal());
492 }
493 ++i;
494 }
495
496 /* // MB : old bit of interpolation code
497 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
498 low = (RooAbsReal*)lowIntIter->Next() ;
499 high = (RooAbsReal*)highIntIter->Next() ;
500
501 if(param->getVal()>0) {
502 value += param->getVal()*(high->getVal() - nominal );
503 } else {
504 value += param->getVal()*(nominal - low->getVal());
505 }
506 ++i;
507 }
508 */
509
510 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
511 while( (param=(RooAbsReal*)paramIter.next()) ) {
512 low = (RooAbsReal*)lowIntIter.next() ;
513 high = (RooAbsReal*)highIntIter.next() ;
514
515 if(_interpCode.empty() || _interpCode.at(i)==0){
516 // piece-wise linear
517 if(param->getVal()>0)
518 value += param->getVal()*(high->getVal() - nominal );
519 else
520 value += param->getVal()*(nominal - low->getVal());
521 } else if(_interpCode.at(i)==1){
522 // piece-wise log
523 if(param->getVal()>=0)
524 value *= pow(high->getVal()/nominal, +param->getVal());
525 else
526 value *= pow(low->getVal()/nominal, -param->getVal());
527 } else if(_interpCode.at(i)==2){
528 // parabolic with linear
529 double a = 0.5*(high->getVal()+low->getVal())-nominal;
530 double b = 0.5*(high->getVal()-low->getVal());
531 double c = 0;
532 if(param->getVal()>1 ){
533 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
534 } else if(param->getVal()<-1 ) {
535 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
536 } else {
537 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
538 }
539 } else if(_interpCode.at(i)==3){
540 //parabolic version of log-normal
541 double a = 0.5*(high->getVal()+low->getVal())-nominal;
542 double b = 0.5*(high->getVal()-low->getVal());
543 double c = 0;
544 if(param->getVal()>1 ){
545 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
546 } else if(param->getVal()<-1 ) {
547 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
548 } else {
549 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
550 }
551
552 } else {
553 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
554 << " with unknown interpolation code" << endl ;
555 }
556 ++i;
557 }
558 */
559
560 // cout << "value = " << value <<endl;
561 return value;
562}
563
564
565////////////////////////////////////////////////////////////////////////////////
566
567void PiecewiseInterpolation::setInterpCode(RooAbsReal& param, int code, bool silent){
568 int index = _paramSet.index(&param);
569 if(index<0){
570 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
571 << " is not in list" << endl ;
572 } else {
573 if(!silent){
574 coutW(InputArguments) << "PiecewiseInterpolation::setInterpCode : " << param.GetName()
575 << " is now " << code << endl ;
576 }
577 _interpCode.at(index) = code;
578 }
579}
580
581
582////////////////////////////////////////////////////////////////////////////////
583
585 for(unsigned int i=0; i<_interpCode.size(); ++i){
586 _interpCode.at(i) = code;
587 }
588}
589
590
591////////////////////////////////////////////////////////////////////////////////
592
594 for(unsigned int i=0; i<_interpCode.size(); ++i){
595 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
596 }
597}
598
599
600////////////////////////////////////////////////////////////////////////////////
601/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
602
603std::list<double>* PiecewiseInterpolation::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
604{
605 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
606}
607
608
609////////////////////////////////////////////////////////////////////////////////
610/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
611
613{
614 return _nominal.arg().isBinnedDistribution(obs) ;
615}
616
617
618
619////////////////////////////////////////////////////////////////////////////////
620
621std::list<double>* PiecewiseInterpolation::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
622{
623 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
624}
625
626////////////////////////////////////////////////////////////////////////////////
627/// Stream an object of class PiecewiseInterpolation.
628
630{
631 if (R__b.IsReading()) {
633 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
634 if (_interpCode.empty()) _interpCode.resize(_paramSet.size());
635 } else {
637 }
638}
639
640
641/*
642////////////////////////////////////////////////////////////////////////////////
643/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
644/// product operator construction
645
646void PiecewiseInterpolation::printMetaArgs(ostream& os) const
647{
648 _lowIter->Reset() ;
649 if (_highIter) {
650 _highIter->Reset() ;
651 }
652
653 bool first(true) ;
654
655 RooAbsArg* arg1, *arg2 ;
656 if (_highSet.size()!=0) {
657
658 while((arg1=(RooAbsArg*)_lowIter->Next())) {
659 if (!first) {
660 os << " + " ;
661 } else {
662 first = false ;
663 }
664 arg2=(RooAbsArg*)_highIter->Next() ;
665 os << arg1->GetName() << " * " << arg2->GetName() ;
666 }
667
668 } else {
669
670 while((arg1=(RooAbsArg*)_lowIter->Next())) {
671 if (!first) {
672 os << " + " ;
673 } else {
674 first = false ;
675 }
676 os << arg1->GetName() ;
677 }
678
679 }
680
681 os << " " ;
682}
683
684*/
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.")
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
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...
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
std::string calculateTreeIndexForCodeSquash(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
double weight(std::size_t i) const
Return weight of i-th bin.
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::string buildArg(RooAbsCollection const &x)
Function to save a RooListProxy as an array in the squashed code.
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
A real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:31
RooDataHist & dataHist()
Return RooDataHist that is represented.
Definition RooHistFunc.h:45
RooArgSet const & variables() const
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 flexibleInterpSingle(unsigned int code, double low, double high, double boundary, double nominal, double paramVal, double res)
Definition MathFuncs.h:184
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345