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