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
23#include "RooFit.h"
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 "RunContext.h"
37
38#include <exception>
39#include <math.h>
40#include <algorithm>
41
42using namespace std;
43
45;
46
47
48////////////////////////////////////////////////////////////////////////////////
49
51{
54}
55
56
57
58////////////////////////////////////////////////////////////////////////////////
59/// Construct a new interpolation. The value of the function will be
60/// \f[
61/// A = \sum_i \mathrm{Interpolate}(\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i).
62/// \f]
63/// \param name Name of the object.
64/// \param title Title (for e.g. plotting)
65/// \param nominal Nominal value of the function.
66/// \param lowSet Set of down variations.
67/// \param highSet Set of up variations.
68/// \param paramSet Parameters that control the interpolation.
69/// \param takeOwnership If true, the PiecewiseInterpolation object will take ownership of the arguments in the low, high and parameter sets.
70PiecewiseInterpolation::PiecewiseInterpolation(const char* name, const char* title, const RooAbsReal& nominal,
71 const RooArgList& lowSet,
72 const RooArgList& highSet,
73 const RooArgList& paramSet,
74 Bool_t takeOwnership) :
75 RooAbsReal(name, title),
76 _normIntMgr(this),
77 _nominal("!nominal","nominal value", this, (RooAbsReal&)nominal),
78 _lowSet("!lowSet","low-side variation",this),
79 _highSet("!highSet","high-side variation",this),
80 _paramSet("!paramSet","high-side variation",this),
81 _positiveDefinite(false)
82
83{
84 // KC: check both sizes
85 if (lowSet.getSize() != highSet.getSize()) {
86 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: input lists should be of equal length" << endl ;
88 }
89
90 RooFIter inputIter1 = lowSet.fwdIterator() ;
91 RooAbsArg* comp ;
92 while((comp = inputIter1.next())) {
93 if (!dynamic_cast<RooAbsReal*>(comp)) {
94 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
95 << " in first list is not of type RooAbsReal" << endl ;
97 }
98 _lowSet.add(*comp) ;
99 if (takeOwnership) {
100 _ownedList.addOwned(*comp) ;
101 }
102 }
103
104
105 RooFIter inputIter2 = highSet.fwdIterator() ;
106 while((comp = inputIter2.next())) {
107 if (!dynamic_cast<RooAbsReal*>(comp)) {
108 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
109 << " in first list is not of type RooAbsReal" << endl ;
111 }
112 _highSet.add(*comp) ;
113 if (takeOwnership) {
114 _ownedList.addOwned(*comp) ;
115 }
116 }
117
118
119 RooFIter inputIter3 = paramSet.fwdIterator() ;
120 while((comp = inputIter3.next())) {
121 if (!dynamic_cast<RooAbsReal*>(comp)) {
122 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
123 << " in first list is not of type RooAbsReal" << endl ;
125 }
126 _paramSet.add(*comp) ;
127 if (takeOwnership) {
128 _ownedList.addOwned(*comp) ;
129 }
130 _interpCode.push_back(0); // default code: linear interpolation
131 }
132
133
134 // Choose special integrator by default
135 specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
137}
138
139
140
141////////////////////////////////////////////////////////////////////////////////
142/// Copy constructor
143
145 RooAbsReal(other, name),
146 _normIntMgr(other._normIntMgr, this),
147 _nominal("!nominal",this,other._nominal),
148 _lowSet("!lowSet",this,other._lowSet),
149 _highSet("!highSet",this,other._highSet),
150 _paramSet("!paramSet",this,other._paramSet),
151 _positiveDefinite(other._positiveDefinite),
152 _interpCode(other._interpCode)
153{
154 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
156}
157
158
159
160////////////////////////////////////////////////////////////////////////////////
161/// Destructor
162
164{
166}
167
168
169
170
171////////////////////////////////////////////////////////////////////////////////
172/// Calculate and return current value of self
173
175{
176 ///////////////////
177 Double_t nominal = _nominal;
178 Double_t sum(nominal) ;
179
180 for (unsigned int i=0; i < _paramSet.size(); ++i) {
181 auto param = static_cast<RooAbsReal*>(_paramSet.at(i));
182 auto low = static_cast<RooAbsReal*>(_lowSet.at(i));
183 auto high = static_cast<RooAbsReal*>(_highSet.at(i));
184 Int_t icode = _interpCode[i] ;
185
186 switch(icode) {
187 case 0: {
188 // piece-wise linear
189 if(param->getVal()>0)
190 sum += param->getVal()*(high->getVal() - nominal );
191 else
192 sum += param->getVal()*(nominal - low->getVal());
193 break ;
194 }
195 case 1: {
196 // pice-wise log
197 if(param->getVal()>=0)
198 sum *= pow(high->getVal()/nominal, +param->getVal());
199 else
200 sum *= pow(low->getVal()/nominal, -param->getVal());
201 break ;
202 }
203 case 2: {
204 // parabolic with linear
205 double a = 0.5*(high->getVal()+low->getVal())-nominal;
206 double b = 0.5*(high->getVal()-low->getVal());
207 double c = 0;
208 if(param->getVal()>1 ){
209 sum += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
210 } else if(param->getVal()<-1 ) {
211 sum += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
212 } else {
213 sum += a*pow(param->getVal(),2) + b*param->getVal()+c;
214 }
215 break ;
216 }
217 case 3: {
218 //parabolic version of log-normal
219 double a = 0.5*(high->getVal()+low->getVal())-nominal;
220 double b = 0.5*(high->getVal()-low->getVal());
221 double c = 0;
222 if(param->getVal()>1 ){
223 sum += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
224 } else if(param->getVal()<-1 ) {
225 sum += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
226 } else {
227 sum += a*pow(param->getVal(),2) + b*param->getVal()+c;
228 }
229 break ;
230 }
231 case 4: {
232
233 // WVE ****************************************************************
234 // WVE *** THIS CODE IS CRITICAL TO HISTFACTORY FIT CPU PERFORMANCE ***
235 // WVE *** Do not modify unless you know what you are doing... ***
236 // WVE ****************************************************************
237
238 double x = param->getVal();
239 if (x>1) {
240 sum += x*(high->getVal() - nominal );
241 } else if (x<-1) {
242 sum += x*(nominal - low->getVal());
243 } else {
244 double eps_plus = high->getVal() - nominal;
245 double eps_minus = nominal - low->getVal();
246 double S = 0.5 * (eps_plus + eps_minus);
247 double A = 0.0625 * (eps_plus - eps_minus);
248
249 //fcns+der+2nd_der are eq at bd
250
251 double val = nominal + x * (S + x * A * ( 15 + x * x * (-10 + x * x * 3 ) ) );
252
253
254 if (val < 0) val = 0;
255 sum += val-nominal;
256 }
257 break ;
258
259 // WVE ****************************************************************
260 }
261 case 5: {
262
263 double x0 = 1.0;//boundary;
264 double x = param->getVal();
265
266 if (x > x0 || x < -x0)
267 {
268 if(x>0)
269 sum += x*(high->getVal() - nominal );
270 else
271 sum += x*(nominal - low->getVal());
272 }
273 else if (nominal != 0)
274 {
275 double eps_plus = high->getVal() - nominal;
276 double eps_minus = nominal - low->getVal();
277 double S = (eps_plus + eps_minus)/2;
278 double A = (eps_plus - eps_minus)/2;
279
280 //fcns+der are eq at bd
281 double a = S;
282 double b = 3*A/(2*x0);
283 //double c = 0;
284 double d = -A/(2*x0*x0*x0);
285
286 double val = nominal + a*x + b*pow(x, 2) + 0/*c*pow(x, 3)*/ + d*pow(x, 4);
287 if (val < 0) val = 0;
288
289 //cout << "Using interp code 5, val = " << val << endl;
290
291 sum += val-nominal;
292 }
293 break ;
294 }
295 default: {
296 coutE(InputArguments) << "PiecewiseInterpolation::evaluate ERROR: " << param->GetName()
297 << " with unknown interpolation code" << icode << endl ;
298 break ;
299 }
300 }
301 }
302
303 if(_positiveDefinite && (sum<0)){
304 sum = 0;
305 // cout <<"sum < 0 forcing positive definite"<<endl;
306 // int code = 1;
307 // RooArgSet* myset = new RooArgSet();
308 // cout << "integral = " << analyticalIntegralWN(code, myset) << endl;
309 } else if(sum<0){
310 cxcoutD(Tracing) <<"PiecewiseInterpolation::evaluate - sum < 0, not forcing positive definite"<<endl;
311 }
312 return sum;
313
314}
315
316
317////////////////////////////////////////////////////////////////////////////////
318/// Interpolate between input distributions for all values of the observable in `evalData`.
319/// \param[in/out] evalData Struct holding spans pointing to input data. The results of this function will be stored here.
320/// \param[in] normSet Arguments to normalise over.
321void PiecewiseInterpolation::computeBatch(cudaStream_t*, double* sum, size_t /*size*/, RooFit::Detail::DataMap const& dataMap) const {
322 auto nominal = dataMap.at(_nominal);
323 for(unsigned int j=0; j < nominal.size(); ++j) {
324 sum[j] = nominal[j];
325 }
326
327 for (unsigned int i=0; i < _paramSet.size(); ++i) {
328 const double param = static_cast<RooAbsReal*>(_paramSet.at(i))->getVal();
329 auto low = dataMap.at(_lowSet.at(i));
330 auto high = dataMap.at(_highSet.at(i));
331 const int icode = _interpCode[i];
332
333 switch(icode) {
334 case 0: {
335 // piece-wise linear
336 for (unsigned int j=0; j < nominal.size(); ++j) {
337 if(param >0)
338 sum[j] += param * (high[j] - nominal[j]);
339 else
340 sum[j] += param * (nominal[j] - low[j] );
341 }
342 break;
343 }
344 case 1: {
345 // pice-wise log
346 for (unsigned int j=0; j < nominal.size(); ++j) {
347 if(param >=0)
348 sum[j] *= pow(high[j]/ nominal[j], +param);
349 else
350 sum[j] *= pow(low[j] / nominal[j], -param);
351 }
352 break;
353 }
354 case 2:
355 // parabolic with linear
356 for (unsigned int j=0; j < nominal.size(); ++j) {
357 const double a = 0.5*(high[j]+low[j])-nominal[j];
358 const double b = 0.5*(high[j]-low[j]);
359 const double c = 0;
360 if (param > 1.) {
361 sum[j] += (2*a+b)*(param -1)+high[j]-nominal[j];
362 } else if (param < -1.) {
363 sum[j] += -1*(2*a-b)*(param +1)+low[j]-nominal[j];
364 } else {
365 sum[j] += a*pow(param ,2) + b*param +c;
366 }
367 }
368 break;
369 case 3: {
370 //parabolic version of log-normal
371 for (unsigned int j=0; j < nominal.size(); ++j) {
372 const double a = 0.5*(high[j]+low[j])-nominal[j];
373 const double b = 0.5*(high[j]-low[j]);
374 const double c = 0;
375 if (param > 1.) {
376 sum[j] += (2*a+b)*(param -1)+high[j]-nominal[j];
377 } else if (param < -1.) {
378 sum[j] += -1*(2*a-b)*(param +1)+low[j]-nominal[j];
379 } else {
380 sum[j] += a*pow(param ,2) + b*param +c;
381 }
382 }
383 break;
384 }
385 case 4:
386 for (unsigned int j=0; j < nominal.size(); ++j) {
387 const double x = param;
388 if (x > 1.) {
389 sum[j] += x * (high[j] - nominal[j]);
390 } else if (x < -1.) {
391 sum[j] += x * (nominal[j] - low[j]);
392 } else {
393 const double eps_plus = high[j] - nominal[j];
394 const double eps_minus = nominal[j] - low[j];
395 const double S = 0.5 * (eps_plus + eps_minus);
396 const double A = 0.0625 * (eps_plus - eps_minus);
397
398 double val = nominal[j] + x * (S + x * A * ( 15. + x * x * (-10. + x * x * 3. ) ) );
399
400 if (val < 0.) val = 0.;
401 sum[j] += val - nominal[j];
402 }
403 }
404 break;
405 case 5:
406 for (unsigned int j=0; j < nominal.size(); ++j) {
407 if (param > 1. || param < -1.) {
408 if(param>0)
409 sum[j] += param * (high[j] - nominal[j]);
410 else
411 sum[j] += param * (nominal[j] - low[j] );
412 } else if (nominal[j] != 0) {
413 const double eps_plus = high[j] - nominal[j];
414 const double eps_minus = nominal[j] - low[j];
415 const double S = (eps_plus + eps_minus)/2;
416 const double A = (eps_plus - eps_minus)/2;
417
418 //fcns+der are eq at bd
419 const double a = S;
420 const double b = 3*A/(2*1.);
421 //double c = 0;
422 const double d = -A/(2*1.*1.*1.);
423
424 double val = nominal[j] + a * param + b * pow(param, 2) + d * pow(param, 4);
425 if (val < 0.) val = 0.;
426
427 sum[j] += val - nominal[j];
428 }
429 }
430 break;
431 default:
432 coutE(InputArguments) << "PiecewiseInterpolation::evaluateSpan(): " << _paramSet[i].GetName()
433 << " with unknown interpolation code" << icode << std::endl;
434 throw std::invalid_argument("PiecewiseInterpolation::evaluateSpan() got invalid interpolation code " + std::to_string(icode));
435 break;
436 }
437 }
438
439 if (_positiveDefinite) {
440 for(unsigned int j=0; j < nominal.size(); ++j) {
441 if (sum[j] < 0.)
442 sum[j] = 0.;
443 }
444 }
445}
446
447////////////////////////////////////////////////////////////////////////////////
448
450{
451 if(allVars.getSize()==1){
452 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
453 temp->specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
454 int nbins = ((RooRealVar*) allVars.first())->numBins();
455 temp->specialIntegratorConfig(kTRUE)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
456 return true;
457 }else{
458 cout << "Currently BinIntegrator only knows how to deal with 1-d "<<endl;
459 return false;
460 }
461 return false;
462}
463
464////////////////////////////////////////////////////////////////////////////////
465/// Advertise that all integrals can be handled internally.
466
468 const RooArgSet* normSet, const char* /*rangeName*/) const
469{
470 /*
471 cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " <<endl;
472 cout << "all vars = "<<endl;
473 allVars.Print("v");
474 cout << "anal vars = "<<endl;
475 analVars.Print("v");
476 cout << "normset vars = "<<endl;
477 if(normSet2)
478 normSet2->Print("v");
479 */
480
481
482 // Handle trivial no-integration scenario
483 if (allVars.getSize()==0) return 0 ;
484 if (_forceNumInt) return 0 ;
485
486
487 // Force using numeric integration
488 // use special numeric integrator
489 return 0;
490
491
492 // KC: check if interCode=0 for all
493 RooFIter paramIterExtra(_paramSet.fwdIterator()) ;
494 int i=0;
495 while( paramIterExtra.next() ) {
496 if(!_interpCode.empty() && _interpCode[i]!=0){
497 // can't factorize integral
498 cout <<"can't factorize integral"<<endl;
499 return 0;
500 }
501 ++i;
502 }
503
504 // Select subset of allVars that are actual dependents
505 analVars.add(allVars) ;
506 // RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
507 // RooArgSet* normSet = getObservables();
508 // RooArgSet* normSet = 0;
509
510
511 // Check if this configuration was created before
512 Int_t sterileIdx(-1) ;
513 CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx) ;
514 if (cache) {
515 return _normIntMgr.lastIndex()+1 ;
516 }
517
518 // Create new cache element
519 cache = new CacheElem ;
520
521 // Make list of function projection and normalization integrals
522 RooAbsReal *func ;
523 // const RooArgSet* nset = _paramList.nset() ;
524
525 // do nominal
526 func = (RooAbsReal*)(&_nominal.arg()) ;
527 RooAbsReal* funcInt = func->createIntegral(analVars) ;
528 cache->_funcIntList.addOwned(*funcInt) ;
529
530 // do variations
531 RooFIter lowIter(_lowSet.fwdIterator()) ;
532 RooFIter highIter(_highSet.fwdIterator()) ;
533 RooFIter paramIter(_paramSet.fwdIterator()) ;
534
535 // int i=0;
536 i=0;
537 while(paramIter.next() ) {
538 func = (RooAbsReal*)lowIter.next() ;
539 funcInt = func->createIntegral(analVars) ;
540 cache->_lowIntList.addOwned(*funcInt) ;
541
542 func = (RooAbsReal*)highIter.next() ;
543 funcInt = func->createIntegral(analVars) ;
544 cache->_highIntList.addOwned(*funcInt) ;
545 ++i;
546 }
547
548 // Store cache element
549 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,0) ;
550
551 return code+1 ;
552}
553
554
555
556
557////////////////////////////////////////////////////////////////////////////////
558/// Implement analytical integrations by doing appropriate weighting from component integrals
559/// functions to integrators of components
560
561Double_t PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
562{
563 /*
564 cout <<"Enter analytic Integral"<<endl;
565 printDirty(true);
566 // _nominal.arg().setDirtyInhibit(kTRUE) ;
567 _nominal.arg().setShapeDirty() ;
568 RooAbsReal* temp ;
569 RooFIter lowIter(_lowSet.fwdIterator()) ;
570 while((temp=(RooAbsReal*)lowIter.next())) {
571 // temp->setDirtyInhibit(kTRUE) ;
572 temp->setShapeDirty() ;
573 }
574 RooFIter highIter(_highSet.fwdIterator()) ;
575 while((temp=(RooAbsReal*)highIter.next())) {
576 // temp->setDirtyInhibit(kTRUE) ;
577 temp->setShapeDirty() ;
578 }
579 */
580
581 /*
582 RooAbsArg::setDirtyInhibit(kTRUE);
583 printDirty(true);
584 cout <<"done setting dirty inhibit = true"<<endl;
585
586 // old integral, only works for linear and not positive definite
587 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
588
589
590 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
591 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
592 cout <<"iset = "<<endl;
593 iset->Print("v");
594
595 double sum = 0;
596 RooArgSet* vars = getVariables();
597 vars->remove(_paramSet);
598 _paramSet.Print("v");
599 vars->Print("v");
600 if(vars->getSize()==1){
601 RooRealVar* obs = (RooRealVar*) vars->first();
602 for(int i=0; i<obs->numBins(); ++i){
603 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
604 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
605 cout << "obs = " << obs->getVal() << " sum = " << sum << endl;
606 }
607 } else{
608 cout <<"only know how to deal with 1 observable right now"<<endl;
609 }
610 */
611
612 /*
613 _nominal.arg().setDirtyInhibit(kFALSE) ;
614 RooFIter lowIter2(_lowSet.fwdIterator()) ;
615 while((temp=(RooAbsReal*)lowIter2.next())) {
616 temp->setDirtyInhibit(kFALSE) ;
617 }
618 RooFIter highIter2(_highSet.fwdIterator()) ;
619 while((temp=(RooAbsReal*)highIter2.next())) {
620 temp->setDirtyInhibit(kFALSE) ;
621 }
622 */
623
624 /*
625 RooAbsArg::setDirtyInhibit(kFALSE);
626 printDirty(true);
627 cout <<"done"<<endl;
628 cout << "sum = " <<sum<<endl;
629 //return sum;
630 */
631
632 // old integral, only works for linear and not positive definite
633 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
634 if( cache==NULL ) {
635 std::cout << "Error: Cache Element is NULL" << std::endl;
636 throw std::exception();
637 }
638
639 // old integral, only works for linear and not positive definite
640 RooFIter funcIntIter = cache->_funcIntList.fwdIterator() ;
641 RooFIter lowIntIter = cache->_lowIntList.fwdIterator() ;
642 RooFIter highIntIter = cache->_highIntList.fwdIterator() ;
643 RooAbsReal *funcInt(0), *low(0), *high(0), *param(0) ;
644 Double_t value(0) ;
645 Double_t nominal(0);
646
647 // get nominal
648 int i=0;
649 while(( funcInt = (RooAbsReal*)funcIntIter.next())) {
650 value += funcInt->getVal() ;
651 nominal = value;
652 i++;
653 }
654 if(i==0 || i>1) { cout << "problem, wrong number of nominal functions"<<endl; }
655
656 // now get low/high variations
657 i = 0;
658 RooFIter paramIter(_paramSet.fwdIterator()) ;
659
660 // KC: old interp code with new iterator
661 while( (param=(RooAbsReal*)paramIter.next()) ) {
662 low = (RooAbsReal*)lowIntIter.next() ;
663 high = (RooAbsReal*)highIntIter.next() ;
664
665 if(param->getVal()>0) {
666 value += param->getVal()*(high->getVal() - nominal );
667 } else {
668 value += param->getVal()*(nominal - low->getVal());
669 }
670 ++i;
671 }
672
673 /* // MB : old bit of interpolation code
674 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
675 low = (RooAbsReal*)lowIntIter->Next() ;
676 high = (RooAbsReal*)highIntIter->Next() ;
677
678 if(param->getVal()>0) {
679 value += param->getVal()*(high->getVal() - nominal );
680 } else {
681 value += param->getVal()*(nominal - low->getVal());
682 }
683 ++i;
684 }
685 */
686
687 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
688 while( (param=(RooAbsReal*)paramIter.next()) ) {
689 low = (RooAbsReal*)lowIntIter.next() ;
690 high = (RooAbsReal*)highIntIter.next() ;
691
692 if(_interpCode.empty() || _interpCode.at(i)==0){
693 // piece-wise linear
694 if(param->getVal()>0)
695 value += param->getVal()*(high->getVal() - nominal );
696 else
697 value += param->getVal()*(nominal - low->getVal());
698 } else if(_interpCode.at(i)==1){
699 // pice-wise log
700 if(param->getVal()>=0)
701 value *= pow(high->getVal()/nominal, +param->getVal());
702 else
703 value *= pow(low->getVal()/nominal, -param->getVal());
704 } else if(_interpCode.at(i)==2){
705 // parabolic with linear
706 double a = 0.5*(high->getVal()+low->getVal())-nominal;
707 double b = 0.5*(high->getVal()-low->getVal());
708 double c = 0;
709 if(param->getVal()>1 ){
710 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
711 } else if(param->getVal()<-1 ) {
712 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
713 } else {
714 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
715 }
716 } else if(_interpCode.at(i)==3){
717 //parabolic version of log-normal
718 double a = 0.5*(high->getVal()+low->getVal())-nominal;
719 double b = 0.5*(high->getVal()-low->getVal());
720 double c = 0;
721 if(param->getVal()>1 ){
722 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
723 } else if(param->getVal()<-1 ) {
724 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
725 } else {
726 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
727 }
728
729 } else {
730 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
731 << " with unknown interpolation code" << endl ;
732 }
733 ++i;
734 }
735 */
736
737 // cout << "value = " << value <<endl;
738 return value;
739}
740
741
742////////////////////////////////////////////////////////////////////////////////
743
744void PiecewiseInterpolation::setInterpCode(RooAbsReal& param, int code, bool silent){
745 int index = _paramSet.index(&param);
746 if(index<0){
747 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
748 << " is not in list" << endl ;
749 } else {
750 if(!silent){
751 coutW(InputArguments) << "PiecewiseInterpolation::setInterpCode : " << param.GetName()
752 << " is now " << code << endl ;
753 }
754 _interpCode.at(index) = code;
755 }
756}
757
758
759////////////////////////////////////////////////////////////////////////////////
760
762 for(unsigned int i=0; i<_interpCode.size(); ++i){
763 _interpCode.at(i) = code;
764 }
765}
766
767
768////////////////////////////////////////////////////////////////////////////////
769
771 for(unsigned int i=0; i<_interpCode.size(); ++i){
772 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
773 }
774}
775
776
777////////////////////////////////////////////////////////////////////////////////
778/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
779
781{
782 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
783}
784
785
786////////////////////////////////////////////////////////////////////////////////
787/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
788
790{
791 return _nominal.arg().isBinnedDistribution(obs) ;
792}
793
794
795
796////////////////////////////////////////////////////////////////////////////////
797
799{
800 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
801}
802
803////////////////////////////////////////////////////////////////////////////////
804/// Stream an object of class PiecewiseInterpolation.
805
807{
808 if (R__b.IsReading()) {
809 R__b.ReadClassBuffer(PiecewiseInterpolation::Class(),this);
810 specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
811 if (_interpCode.empty()) _interpCode.resize(_paramSet.getSize());
812 } else {
813 R__b.WriteClassBuffer(PiecewiseInterpolation::Class(),this);
814 }
815}
816
817
818/*
819////////////////////////////////////////////////////////////////////////////////
820/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
821/// product operator construction
822
823void PiecewiseInterpolation::printMetaArgs(ostream& os) const
824{
825 _lowIter->Reset() ;
826 if (_highIter) {
827 _highIter->Reset() ;
828 }
829
830 Bool_t first(kTRUE) ;
831
832 RooAbsArg* arg1, *arg2 ;
833 if (_highSet.getSize()!=0) {
834
835 while((arg1=(RooAbsArg*)_lowIter->Next())) {
836 if (!first) {
837 os << " + " ;
838 } else {
839 first = kFALSE ;
840 }
841 arg2=(RooAbsArg*)_highIter->Next() ;
842 os << arg1->GetName() << " * " << arg2->GetName() ;
843 }
844
845 } else {
846
847 while((arg1=(RooAbsArg*)_lowIter->Next())) {
848 if (!first) {
849 os << " + " ;
850 } else {
851 first = kFALSE ;
852 }
853 os << arg1->GetName() ;
854 }
855
856 }
857
858 os << " " ;
859}
860
861*/
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#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
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
The PiecewiseInterpolation is a class that can morph distributions into each other,...
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
RooListProxy _lowSet
Low-side variation.
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
RooListProxy _highSet
High-side variation.
void setInterpCode(RooAbsReal &param, int code, bool silent=false)
Bool_t _positiveDefinite
protect against negative and 0 bins.
Bool_t setBinIntegrator(RooArgSet &allVars)
RooObjCacheManager _normIntMgr
! The integration cache manager
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Advertise that all integrals can be handled internally.
double evaluate() const
Calculate and return current value of self.
void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const
Interpolate between input distributions for all values of the observable in evalData.
RooListProxy _paramSet
interpolation parameters
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
virtual ~PiecewiseInterpolation()
Destructor.
RooArgList _ownedList
List of owned components.
RooRealProxy _nominal
The nominal value.
virtual Bool_t isBinnedDistribution(const RooArgSet &obs) const
WVE note: assumes nominal and alternates have identical structure, must add explicit check.
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
friend void RooRefArray::Streamer(TBuffer &)
RooAbsCacheElement is the abstract base class for objects to be stored in RooAbsCache cache manager o...
Int_t getSize() const
Int_t index(const RooAbsArg *arg) const
Returns index of given arg, or -1 if arg is not in the collection.
RooFIter fwdIterator() const
One-time forward iterator.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Add an argument and transfer the ownership to the collection.
Storage_t::size_type size() const
Bool_t setRealValue(const char *name, Double_t newVal=0, Bool_t verbose=kFALSE)
Set value of a RooAbsRealLValye stored in set with given name to newVal No error messages are printed...
RooAbsArg * first() const
const char * GetName() const
Returns name of object.
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
Bool_t _forceNumInt
Definition RooAbsReal.h:487
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Retrieve bin boundaries if this distribution is binned in obs.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
virtual Bool_t isBinnedDistribution(const RooArgSet &) const
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooAbsReal.h:345
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:35
T * getObj(const RooArgSet *nset, Int_t *sterileIndex=0, const TNamed *isetRangeName=0)
T * getObjByIndex(Int_t index) const
Retrieve payload object by slot index.
Int_t lastIndex() const
Int_t setObj(const RooArgSet *nset, T *obj, const TNamed *isetRangeName=0)
virtual Bool_t setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
static void softAbort()
A one-time forward iterator working on RooLinkedList or RooAbsCollection.
RooAbsArg * next()
Return next element or nullptr if at end.
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
RooCategory & method1D()
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
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=0)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345