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
200{
201 std::size_t n = _interpCode.size();
202
203 std::string resName = "total_" + ctx.getTmpVarName();
204 for (std::size_t i = 0; i < n; ++i) {
205 if (_interpCode[i] != _interpCode[0]) {
206 coutE(InputArguments) << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
207 "different interpolation codes for the same class object "
208 << endl;
209 }
210 }
211
212 // The PiecewiseInterpolation class is used in the context of HistFactory
213 // models, where is is always used the same way: all RooAbsReals in _lowSet,
214 // _histSet, and also nominal are 1D RooHistFuncs with with same structure.
215 //
216 // Therefore, we can make a big optimization: we get the bin index only once
217 // here in the generated code for PiecewiseInterpolation. Then, we also
218 // rearrange the histogram data in such a way that we can always pass the
219 // same arrays to the free function that implements the interpolation, just
220 // with a dynamic offset calculated from the bin index.
221 RooDataHist const &nomHist = dynamic_cast<RooHistFunc const &>(*_nominal).dataHist();
222 int nBins = nomHist.numEntries();
223 std::vector<double> valsNominal;
224 std::vector<double> valsLow;
225 std::vector<double> valsHigh;
226 for (int i = 0; i < nBins; ++i) {
227 valsNominal.push_back(nomHist.weight(i));
228 }
229 for (int i = 0; i < nBins; ++i) {
230 for (std::size_t iParam = 0; iParam < n; ++iParam) {
231 valsLow.push_back(dynamic_cast<RooHistFunc const &>(_lowSet[iParam]).dataHist().weight(i));
232 valsHigh.push_back(dynamic_cast<RooHistFunc const &>(_highSet[iParam]).dataHist().weight(i));
233 }
234 }
235 std::string idxName = ctx.getTmpVarName();
236 std::string valsNominalStr = ctx.buildArg(valsNominal);
237 std::string valsLowStr = ctx.buildArg(valsLow);
238 std::string valsHighStr = ctx.buildArg(valsHigh);
239 std::string nStr = std::to_string(n);
240 std::string code;
241
242 std::string lowName = ctx.getTmpVarName();
243 std::string highName = ctx.getTmpVarName();
244 std::string nominalName = ctx.getTmpVarName();
245 code += "unsigned int " + idxName + " = " + nomHist.calculateTreeIndexForCodeSquash(this, ctx, dynamic_cast<RooHistFunc const &>(*_nominal).variables()) + ";\n";
246 code += "double const* " + lowName + " = " + valsLowStr + " + " + nStr + " * " + idxName + ";\n";
247 code += "double const* " + highName + " = " + valsHighStr + " + " + nStr + " * " + idxName + ";\n";
248 code += "double " + nominalName + " = *(" + valsNominalStr + " + " + idxName + ");\n";
249
250 std::string funcCall = ctx.buildCall("RooFit::Detail::MathFuncs::flexibleInterp", _interpCode[0], _paramSet, n,
251 lowName, highName, 1.0, nominalName, 0.0);
252 code += "double " + resName + " = " + funcCall + ";\n";
253
255 code += resName + " = " + resName + " < 0 ? 0 : " + resName + ";\n";
256
257 ctx.addToCodeBody(this, code);
258 ctx.addResult(this, resName);
259}
260
261namespace {
262
263inline double broadcast(std::span<const double> const &s, std::size_t i)
264{
265 return s.size() > 1 ? s[i] : s[0];
266}
267
268} // namespace
269
270////////////////////////////////////////////////////////////////////////////////
271/// Interpolate between input distributions for all values of the observable in `evalData`.
272/// \param[in,out] evalData Struct holding spans pointing to input data. The results of this function will be stored here.
273/// \param[in] normSet Arguments to normalise over.
275{
276 std::span<double> sum = ctx.output();
277
278 auto nominal = ctx.at(_nominal);
279
280 for (std::size_t j = 0; j < sum.size(); ++j) {
281 sum[j] = broadcast(nominal, j);
282 }
283
284 for (unsigned int i = 0; i < _paramSet.size(); ++i) {
285 auto param = ctx.at(_paramSet.at(i));
286 auto low = ctx.at(_lowSet.at(i));
287 auto high = ctx.at(_highSet.at(i));
288
289 for (std::size_t j = 0; j < sum.size(); ++j) {
291 sum[j] += flexibleInterpSingle(_interpCode[i], broadcast(low, j), broadcast(high, j), 1.0, broadcast(nominal, j),
292 broadcast(param, j), sum[j]);
293 }
294 }
295
296 if (_positiveDefinite) {
297 for (std::size_t j = 0; j < sum.size(); ++j) {
298 if (sum[j] < 0.)
299 sum[j] = 0.;
300 }
301 }
302}
303
304////////////////////////////////////////////////////////////////////////////////
305
307{
308 if(allVars.size()==1){
309 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
310 temp->specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
311 int nbins = (static_cast<RooRealVar*>(allVars.first()))->numBins();
312 temp->specialIntegratorConfig(true)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
313 return true;
314 }else{
315 cout << "Currently BinIntegrator only knows how to deal with 1-d "<<endl;
316 return false;
317 }
318 return false;
319}
320
321////////////////////////////////////////////////////////////////////////////////
322/// Advertise that all integrals can be handled internally.
323
325 const RooArgSet* normSet, const char* /*rangeName*/) const
326{
327 /*
328 cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " <<endl;
329 cout << "all vars = "<<endl;
330 allVars.Print("v");
331 cout << "anal vars = "<<endl;
332 analVars.Print("v");
333 cout << "normset vars = "<<endl;
334 if(normSet2)
335 normSet2->Print("v");
336 */
337
338
339 // Handle trivial no-integration scenario
340 if (allVars.empty()) return 0 ;
341 if (_forceNumInt) return 0 ;
342
343
344 // Force using numeric integration
345 // use special numeric integrator
346 return 0;
347
348
349 // KC: check if interCode=0 for all
350 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it) {
351 if (!_interpCode.empty() && _interpCode[it - _paramSet.begin()] != 0) {
352 // can't factorize integral
353 cout << "can't factorize integral" << endl;
354 return 0;
355 }
356 }
357
358 // Select subset of allVars that are actual dependents
359 analVars.add(allVars) ;
360
361 // Check if this configuration was created before
362 Int_t sterileIdx(-1) ;
363 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObj(normSet,&analVars,&sterileIdx)) ;
364 if (cache) {
365 return _normIntMgr.lastIndex()+1 ;
366 }
367
368 // Create new cache element
369 cache = new CacheElem ;
370
371 // Make list of function projection and normalization integrals
372 RooAbsReal *func ;
373
374 // do variations
375 for (auto it = _paramSet.begin(); it != _paramSet.end(); ++it)
376 {
377 auto i = it - _paramSet.begin();
378 func = static_cast<RooAbsReal *>(_lowSet.at(i));
379 cache->_lowIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
380
381 func = static_cast<RooAbsReal *>(_highSet.at(i));
382 cache->_highIntList.addOwned(std::unique_ptr<RooAbsReal>{func->createIntegral(analVars)});
383 }
384
385 // Store cache element
386 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,nullptr) ;
387
388 return code+1 ;
389}
390
391
392
393
394////////////////////////////////////////////////////////////////////////////////
395/// Implement analytical integrations by doing appropriate weighting from component integrals
396/// functions to integrators of components
397
398double PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
399{
400 /*
401 cout <<"Enter analytic Integral"<<endl;
402 printDirty(true);
403 // _nominal.arg().setDirtyInhibit(true) ;
404 _nominal.arg().setShapeDirty() ;
405 RooAbsReal* temp ;
406 RooFIter lowIter(_lowSet.fwdIterator()) ;
407 while((temp=(RooAbsReal*)lowIter.next())) {
408 // temp->setDirtyInhibit(true) ;
409 temp->setShapeDirty() ;
410 }
411 RooFIter highIter(_highSet.fwdIterator()) ;
412 while((temp=(RooAbsReal*)highIter.next())) {
413 // temp->setDirtyInhibit(true) ;
414 temp->setShapeDirty() ;
415 }
416 */
417
418 /*
419 RooAbsArg::setDirtyInhibit(true);
420 printDirty(true);
421 cout <<"done setting dirty inhibit = true"<<endl;
422
423 // old integral, only works for linear and not positive definite
424 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
425
426
427 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
428 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
429 cout <<"iset = "<<endl;
430 iset->Print("v");
431
432 double sum = 0;
433 RooArgSet* vars = getVariables();
434 vars->remove(_paramSet);
435 _paramSet.Print("v");
436 vars->Print("v");
437 if(vars->size()==1){
438 RooRealVar* obs = (RooRealVar*) vars->first();
439 for(int i=0; i<obs->numBins(); ++i){
440 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
441 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
442 cout << "obs = " << obs->getVal() << " sum = " << sum << endl;
443 }
444 } else{
445 cout <<"only know how to deal with 1 observable right now"<<endl;
446 }
447 */
448
449 /*
450 _nominal.arg().setDirtyInhibit(false) ;
451 RooFIter lowIter2(_lowSet.fwdIterator()) ;
452 while((temp=(RooAbsReal*)lowIter2.next())) {
453 temp->setDirtyInhibit(false) ;
454 }
455 RooFIter highIter2(_highSet.fwdIterator()) ;
456 while((temp=(RooAbsReal*)highIter2.next())) {
457 temp->setDirtyInhibit(false) ;
458 }
459 */
460
461 /*
462 RooAbsArg::setDirtyInhibit(false);
463 printDirty(true);
464 cout <<"done"<<endl;
465 cout << "sum = " <<sum<<endl;
466 //return sum;
467 */
468
469 // old integral, only works for linear and not positive definite
470 CacheElem* cache = static_cast<CacheElem*>(_normIntMgr.getObjByIndex(code-1)) ;
471 if( cache==nullptr ) {
472 std::cout << "Error: Cache Element is nullptr" << std::endl;
473 throw std::exception();
474 }
475
476 // old integral, only works for linear and not positive definite
477
478 RooAbsReal *low;
479 RooAbsReal *high;
480 double value(0);
481 double nominal(0);
482
483 // get nominal
484 int i=0;
485 for (auto funcInt : static_range_cast<RooAbsReal*>(cache->_funcIntList)) {
486 value += funcInt->getVal() ;
487 nominal = value;
488 i++;
489 }
490 if(i==0 || i>1) { cout << "problem, wrong number of nominal functions"<<endl; }
491
492 // now get low/high variations
493 // KC: old interp code with new iterator
494
495 i = 0;
496 for (auto const *param : static_range_cast<RooAbsReal *>(_paramSet)) {
497 low = static_cast<RooAbsReal *>(cache->_lowIntList.at(i));
498 high = static_cast<RooAbsReal *>(cache->_highIntList.at(i));
499
500 if(param->getVal() > 0) {
501 value += param->getVal()*(high->getVal() - nominal);
502 } else {
503 value += param->getVal()*(nominal - low->getVal());
504 }
505 ++i;
506 }
507
508 /* // MB : old bit of interpolation code
509 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
510 low = (RooAbsReal*)lowIntIter->Next() ;
511 high = (RooAbsReal*)highIntIter->Next() ;
512
513 if(param->getVal()>0) {
514 value += param->getVal()*(high->getVal() - nominal );
515 } else {
516 value += param->getVal()*(nominal - low->getVal());
517 }
518 ++i;
519 }
520 */
521
522 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
523 while( (param=(RooAbsReal*)paramIter.next()) ) {
524 low = (RooAbsReal*)lowIntIter.next() ;
525 high = (RooAbsReal*)highIntIter.next() ;
526
527 if(_interpCode.empty() || _interpCode.at(i)==0){
528 // piece-wise linear
529 if(param->getVal()>0)
530 value += param->getVal()*(high->getVal() - nominal );
531 else
532 value += param->getVal()*(nominal - low->getVal());
533 } else if(_interpCode.at(i)==1){
534 // piece-wise log
535 if(param->getVal()>=0)
536 value *= pow(high->getVal()/nominal, +param->getVal());
537 else
538 value *= pow(low->getVal()/nominal, -param->getVal());
539 } else if(_interpCode.at(i)==2){
540 // parabolic with linear
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 } else if(_interpCode.at(i)==3){
552 //parabolic version of log-normal
553 double a = 0.5*(high->getVal()+low->getVal())-nominal;
554 double b = 0.5*(high->getVal()-low->getVal());
555 double c = 0;
556 if(param->getVal()>1 ){
557 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
558 } else if(param->getVal()<-1 ) {
559 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
560 } else {
561 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
562 }
563
564 } else {
565 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
566 << " with unknown interpolation code" << endl ;
567 }
568 ++i;
569 }
570 */
571
572 // cout << "value = " << value <<endl;
573 return value;
574}
575
576void PiecewiseInterpolation::setInterpCode(RooAbsReal &param, int code, bool /*silent*/)
577{
578 int index = _paramSet.index(&param);
579 if (index < 0) {
580 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName() << " is not in list"
581 << std::endl;
582 return;
583 }
585}
586
588{
589 for (std::size_t i = 0; i < _interpCode.size(); ++i) {
590 setInterpCodeForParam(i, code);
591 }
592}
593
595{
596 RooAbsArg const &param = _paramSet[iParam];
597 if (code < 0 || code > 6) {
598 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
599 << " with unknown interpolation code " << code << ", keeping current code "
600 << _interpCode[iParam] << std::endl;
601 return;
602 }
603 if (code == 3) {
604 // In the past, code 3 was equivalent to code 2, which confused users.
605 // Now, we just say that code 3 doesn't exist and default to code 2 in
606 // that case for backwards compatible behavior.
607 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
608 << " with unknown interpolation code " << code << ", defaulting to code 2" << std::endl;
609 code = 2;
610 }
611 _interpCode.at(iParam) = code;
613}
614
615////////////////////////////////////////////////////////////////////////////////
616
618 for(unsigned int i=0; i<_interpCode.size(); ++i){
619 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
620 }
621}
622
623
624////////////////////////////////////////////////////////////////////////////////
625/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
626
627std::list<double>* PiecewiseInterpolation::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
628{
629 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
630}
631
632
633////////////////////////////////////////////////////////////////////////////////
634/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
635
637{
638 return _nominal.arg().isBinnedDistribution(obs) ;
639}
640
641
642
643////////////////////////////////////////////////////////////////////////////////
644
645std::list<double>* PiecewiseInterpolation::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
646{
647 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
648}
649
650////////////////////////////////////////////////////////////////////////////////
651/// Stream an object of class PiecewiseInterpolation.
652
654{
655 if (R__b.IsReading()) {
657 specialIntegratorConfig(true)->method1D().setLabel("RooBinIntegrator") ;
658 if (_interpCode.empty()) _interpCode.resize(_paramSet.size());
659 } else {
661 }
662}
663
664
665/*
666////////////////////////////////////////////////////////////////////////////////
667/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
668/// product operator construction
669
670void PiecewiseInterpolation::printMetaArgs(ostream& os) const
671{
672 _lowIter->Reset() ;
673 if (_highIter) {
674 _highIter->Reset() ;
675 }
676
677 bool first(true) ;
678
679 RooAbsArg* arg1, *arg2 ;
680 if (_highSet.size()!=0) {
681
682 while((arg1=(RooAbsArg*)_lowIter->Next())) {
683 if (!first) {
684 os << " + " ;
685 } else {
686 first = false ;
687 }
688 arg2=(RooAbsArg*)_highIter->Next() ;
689 os << arg1->GetName() << " * " << arg2->GetName() ;
690 }
691
692 } else {
693
694 while((arg1=(RooAbsArg*)_lowIter->Next())) {
695 if (!first) {
696 os << " + " ;
697 } else {
698 first = false ;
699 }
700 os << arg1->GetName() ;
701 }
702
703 }
704
705 os << " " ;
706}
707
708*/
#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)
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...
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
friend void RooRefArray::Streamer(TBuffer &)
void setValueDirty()
Mark the element dirty. This forces a re-evaluation when a value is requested.
Definition RooAbsArg.h:462
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
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
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:539
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
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...
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
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
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:213
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345