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
37#include <exception>
38#include <math.h>
39
40using namespace std;
41
43;
44
45
46////////////////////////////////////////////////////////////////////////////////
47
49{
52}
53
54
55
56////////////////////////////////////////////////////////////////////////////////
57/// Construct a new interpolation. The value of the function will be
58/// \f[
59/// A = \sum_i \mathrm{Interpolate}(\mathrm{low}_i, \mathrm{nominal}, \mathrm{high}_i).
60/// \f]
61/// \param name Name of the object.
62/// \param title Title (for e.g. plotting)
63/// \param nominal Nominal value of the function.
64/// \param lowSet Set of down variations.
65/// \param highSet Set of up variations.
66/// \param paramSet Parameters that control the interpolation.
67/// \param takeOwnership If true, the PiecewiseInterpolation object will take ownership of the arguments in the low, high and parameter sets.
68PiecewiseInterpolation::PiecewiseInterpolation(const char* name, const char* title, const RooAbsReal& nominal,
69 const RooArgList& lowSet,
70 const RooArgList& highSet,
71 const RooArgList& paramSet,
72 Bool_t takeOwnership) :
73 RooAbsReal(name, title),
74 _nominal("!nominal","nominal value", this, (RooAbsReal&)nominal),
75 _lowSet("!lowSet","low-side variation",this),
76 _highSet("!highSet","high-side variation",this),
77 _paramSet("!paramSet","high-side variation",this),
78 _positiveDefinite(false)
79
80{
81 // KC: check both sizes
82 if (lowSet.getSize() != highSet.getSize()) {
83 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: input lists should be of equal length" << endl ;
85 }
86
87 RooFIter inputIter1 = lowSet.fwdIterator() ;
88 RooAbsArg* comp ;
89 while((comp = inputIter1.next())) {
90 if (!dynamic_cast<RooAbsReal*>(comp)) {
91 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
92 << " in first list is not of type RooAbsReal" << endl ;
94 }
95 _lowSet.add(*comp) ;
96 if (takeOwnership) {
97 _ownedList.addOwned(*comp) ;
98 }
99 }
100
101
102 RooFIter inputIter2 = highSet.fwdIterator() ;
103 while((comp = inputIter2.next())) {
104 if (!dynamic_cast<RooAbsReal*>(comp)) {
105 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
106 << " in first list is not of type RooAbsReal" << endl ;
108 }
109 _highSet.add(*comp) ;
110 if (takeOwnership) {
111 _ownedList.addOwned(*comp) ;
112 }
113 }
114
115
116 RooFIter inputIter3 = paramSet.fwdIterator() ;
117 while((comp = inputIter3.next())) {
118 if (!dynamic_cast<RooAbsReal*>(comp)) {
119 coutE(InputArguments) << "PiecewiseInterpolation::ctor(" << GetName() << ") ERROR: component " << comp->GetName()
120 << " in first list is not of type RooAbsReal" << endl ;
122 }
123 _paramSet.add(*comp) ;
124 if (takeOwnership) {
125 _ownedList.addOwned(*comp) ;
126 }
127 _interpCode.push_back(0); // default code: linear interpolation
128 }
129
130
131 // Choose special integrator by default
132 specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
134}
135
136
137
138////////////////////////////////////////////////////////////////////////////////
139/// Copy constructor
140
142 RooAbsReal(other, name),
143 _nominal("!nominal",this,other._nominal),
144 _lowSet("!lowSet",this,other._lowSet),
145 _highSet("!highSet",this,other._highSet),
146 _paramSet("!paramSet",this,other._paramSet),
147 _positiveDefinite(other._positiveDefinite),
148 _interpCode(other._interpCode)
149{
150 // Member _ownedList is intentionally not copy-constructed -- ownership is not transferred
152}
153
154
155
156////////////////////////////////////////////////////////////////////////////////
157/// Destructor
158
160{
162}
163
164
165
166
167////////////////////////////////////////////////////////////////////////////////
168/// Calculate and return current value of self
169
171{
172 ///////////////////
173 Double_t nominal = _nominal;
174 Double_t sum(nominal) ;
175
176 for (unsigned int i=0; i < _paramSet.size(); ++i) {
177 auto param = static_cast<RooAbsReal*>(_paramSet.at(i));
178 auto low = static_cast<RooAbsReal*>(_lowSet.at(i));
179 auto high = static_cast<RooAbsReal*>(_highSet.at(i));
180 Int_t icode = _interpCode[i] ;
181
182 switch(icode) {
183 case 0: {
184 // piece-wise linear
185 if(param->getVal()>0)
186 sum += param->getVal()*(high->getVal() - nominal );
187 else
188 sum += param->getVal()*(nominal - low->getVal());
189 break ;
190 }
191 case 1: {
192 // pice-wise log
193 if(param->getVal()>=0)
194 sum *= pow(high->getVal()/nominal, +param->getVal());
195 else
196 sum *= pow(low->getVal()/nominal, -param->getVal());
197 break ;
198 }
199 case 2: {
200 // parabolic with linear
201 double a = 0.5*(high->getVal()+low->getVal())-nominal;
202 double b = 0.5*(high->getVal()-low->getVal());
203 double c = 0;
204 if(param->getVal()>1 ){
205 sum += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
206 } else if(param->getVal()<-1 ) {
207 sum += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
208 } else {
209 sum += a*pow(param->getVal(),2) + b*param->getVal()+c;
210 }
211 break ;
212 }
213 case 3: {
214 //parabolic version of log-normal
215 double a = 0.5*(high->getVal()+low->getVal())-nominal;
216 double b = 0.5*(high->getVal()-low->getVal());
217 double c = 0;
218 if(param->getVal()>1 ){
219 sum += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
220 } else if(param->getVal()<-1 ) {
221 sum += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
222 } else {
223 sum += a*pow(param->getVal(),2) + b*param->getVal()+c;
224 }
225 break ;
226 }
227 case 4: {
228
229 // WVE ****************************************************************
230 // WVE *** THIS CODE IS CRITICAL TO HISTFACTORY FIT CPU PERFORMANCE ***
231 // WVE *** Do not modify unless you know what you are doing... ***
232 // WVE ****************************************************************
233
234 double x = param->getVal();
235 if (x>1) {
236 sum += x*(high->getVal() - nominal );
237 } else if (x<-1) {
238 sum += x*(nominal - low->getVal());
239 } else {
240 double eps_plus = high->getVal() - nominal;
241 double eps_minus = nominal - low->getVal();
242 double S = 0.5 * (eps_plus + eps_minus);
243 double A = 0.0625 * (eps_plus - eps_minus);
244
245 //fcns+der+2nd_der are eq at bd
246
247 double val = nominal + x * (S + x * A * ( 15 + x * x * (-10 + x * x * 3 ) ) );
248
249
250 if (val < 0) val = 0;
251 sum += val-nominal;
252 }
253 break ;
254
255 // WVE ****************************************************************
256 }
257 case 5: {
258
259 double x0 = 1.0;//boundary;
260 double x = param->getVal();
261
262 if (x > x0 || x < -x0)
263 {
264 if(x>0)
265 sum += x*(high->getVal() - nominal );
266 else
267 sum += x*(nominal - low->getVal());
268 }
269 else if (nominal != 0)
270 {
271 double eps_plus = high->getVal() - nominal;
272 double eps_minus = nominal - low->getVal();
273 double S = (eps_plus + eps_minus)/2;
274 double A = (eps_plus - eps_minus)/2;
275
276 //fcns+der are eq at bd
277 double a = S;
278 double b = 3*A/(2*x0);
279 //double c = 0;
280 double d = -A/(2*x0*x0*x0);
281
282 double val = nominal + a*x + b*pow(x, 2) + 0/*c*pow(x, 3)*/ + d*pow(x, 4);
283 if (val < 0) val = 0;
284
285 //cout << "Using interp code 5, val = " << val << endl;
286
287 sum += val-nominal;
288 }
289 break ;
290 }
291 default: {
292 coutE(InputArguments) << "PiecewiseInterpolation::evaluate ERROR: " << param->GetName()
293 << " with unknown interpolation code" << icode << endl ;
294 break ;
295 }
296 }
297 }
298
299 if(_positiveDefinite && (sum<0)){
300 sum = 1e-6;
301 sum = 0;
302 // cout <<"sum < 0 forcing positive definite"<<endl;
303 // int code = 1;
304 // RooArgSet* myset = new RooArgSet();
305 // cout << "integral = " << analyticalIntegralWN(code, myset) << endl;
306 } else if(sum<0){
307 cxcoutD(Tracing) <<"PiecewiseInterpolation::evaluate - sum < 0, not forcing positive definite"<<endl;
308 }
309 return sum;
310
311}
312
313////////////////////////////////////////////////////////////////////////////////
314
316{
317 if(allVars.getSize()==1){
318 RooAbsReal* temp = const_cast<PiecewiseInterpolation*>(this);
319 temp->specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
320 int nbins = ((RooRealVar*) allVars.first())->numBins();
321 temp->specialIntegratorConfig(kTRUE)->getConfigSection("RooBinIntegrator").setRealValue("numBins",nbins);
322 return true;
323 }else{
324 cout << "Currently BinIntegrator only knows how to deal with 1-d "<<endl;
325 return false;
326 }
327 return false;
328}
329
330////////////////////////////////////////////////////////////////////////////////
331/// Advertise that all integrals can be handled internally.
332
334 const RooArgSet* normSet, const char* /*rangeName*/) const
335{
336 /*
337 cout << "---------------------------\nin PiecewiseInterpolation get analytic integral " <<endl;
338 cout << "all vars = "<<endl;
339 allVars.Print("v");
340 cout << "anal vars = "<<endl;
341 analVars.Print("v");
342 cout << "normset vars = "<<endl;
343 if(normSet2)
344 normSet2->Print("v");
345 */
346
347
348 // Handle trivial no-integration scenario
349 if (allVars.getSize()==0) return 0 ;
350 if (_forceNumInt) return 0 ;
351
352
353 // Force using numeric integration
354 // use special numeric integrator
355 return 0;
356
357
358 // KC: check if interCode=0 for all
359 RooFIter paramIterExtra(_paramSet.fwdIterator()) ;
360 int i=0;
361 while( paramIterExtra.next() ) {
362 if(!_interpCode.empty() && _interpCode[i]!=0){
363 // can't factorize integral
364 cout <<"can't factorize integral"<<endl;
365 return 0;
366 }
367 ++i;
368 }
369
370 // Select subset of allVars that are actual dependents
371 analVars.add(allVars) ;
372 // RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
373 // RooArgSet* normSet = getObservables();
374 // RooArgSet* normSet = 0;
375
376
377 // Check if this configuration was created before
378 Int_t sterileIdx(-1) ;
379 CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx) ;
380 if (cache) {
381 return _normIntMgr.lastIndex()+1 ;
382 }
383
384 // Create new cache element
385 cache = new CacheElem ;
386
387 // Make list of function projection and normalization integrals
388 RooAbsReal *func ;
389 // const RooArgSet* nset = _paramList.nset() ;
390
391 // do nominal
392 func = (RooAbsReal*)(&_nominal.arg()) ;
393 RooAbsReal* funcInt = func->createIntegral(analVars) ;
394 cache->_funcIntList.addOwned(*funcInt) ;
395
396 // do variations
397 RooFIter lowIter(_lowSet.fwdIterator()) ;
398 RooFIter highIter(_highSet.fwdIterator()) ;
399 RooFIter paramIter(_paramSet.fwdIterator()) ;
400
401 // int i=0;
402 i=0;
403 while(paramIter.next() ) {
404 func = (RooAbsReal*)lowIter.next() ;
405 funcInt = func->createIntegral(analVars) ;
406 cache->_lowIntList.addOwned(*funcInt) ;
407
408 func = (RooAbsReal*)highIter.next() ;
409 funcInt = func->createIntegral(analVars) ;
410 cache->_highIntList.addOwned(*funcInt) ;
411 ++i;
412 }
413
414 // Store cache element
415 Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,0) ;
416
417 return code+1 ;
418}
419
420
421
422
423////////////////////////////////////////////////////////////////////////////////
424/// Implement analytical integrations by doing appropriate weighting from component integrals
425/// functions to integrators of components
426
427Double_t PiecewiseInterpolation::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
428{
429 /*
430 cout <<"Enter analytic Integral"<<endl;
431 printDirty(true);
432 // _nominal.arg().setDirtyInhibit(kTRUE) ;
433 _nominal.arg().setShapeDirty() ;
434 RooAbsReal* temp ;
435 RooFIter lowIter(_lowSet.fwdIterator()) ;
436 while((temp=(RooAbsReal*)lowIter.next())) {
437 // temp->setDirtyInhibit(kTRUE) ;
438 temp->setShapeDirty() ;
439 }
440 RooFIter highIter(_highSet.fwdIterator()) ;
441 while((temp=(RooAbsReal*)highIter.next())) {
442 // temp->setDirtyInhibit(kTRUE) ;
443 temp->setShapeDirty() ;
444 }
445 */
446
447 /*
448 RooAbsArg::setDirtyInhibit(kTRUE);
449 printDirty(true);
450 cout <<"done setting dirty inhibit = true"<<endl;
451
452 // old integral, only works for linear and not positive definite
453 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
454
455
456 std::unique_ptr<RooArgSet> vars2( getParameters(RooArgSet()) );
457 std::unique_ptr<RooArgSet> iset( _normIntMgr.nameSet2ByIndex(code-1)->select(*vars2) );
458 cout <<"iset = "<<endl;
459 iset->Print("v");
460
461 double sum = 0;
462 RooArgSet* vars = getVariables();
463 vars->remove(_paramSet);
464 _paramSet.Print("v");
465 vars->Print("v");
466 if(vars->getSize()==1){
467 RooRealVar* obs = (RooRealVar*) vars->first();
468 for(int i=0; i<obs->numBins(); ++i){
469 obs->setVal( obs->getMin() + (.5+i)*(obs->getMax()-obs->getMin())/obs->numBins());
470 sum+=evaluate()*(obs->getMax()-obs->getMin())/obs->numBins();
471 cout << "obs = " << obs->getVal() << " sum = " << sum << endl;
472 }
473 } else{
474 cout <<"only know how to deal with 1 observable right now"<<endl;
475 }
476 */
477
478 /*
479 _nominal.arg().setDirtyInhibit(kFALSE) ;
480 RooFIter lowIter2(_lowSet.fwdIterator()) ;
481 while((temp=(RooAbsReal*)lowIter2.next())) {
482 temp->setDirtyInhibit(kFALSE) ;
483 }
484 RooFIter highIter2(_highSet.fwdIterator()) ;
485 while((temp=(RooAbsReal*)highIter2.next())) {
486 temp->setDirtyInhibit(kFALSE) ;
487 }
488 */
489
490 /*
491 RooAbsArg::setDirtyInhibit(kFALSE);
492 printDirty(true);
493 cout <<"done"<<endl;
494 cout << "sum = " <<sum<<endl;
495 //return sum;
496 */
497
498 // old integral, only works for linear and not positive definite
499 CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
500 if( cache==NULL ) {
501 std::cout << "Error: Cache Element is NULL" << std::endl;
502 throw std::exception();
503 }
504
505 // old integral, only works for linear and not positive definite
506 RooFIter funcIntIter = cache->_funcIntList.fwdIterator() ;
507 RooFIter lowIntIter = cache->_lowIntList.fwdIterator() ;
508 RooFIter highIntIter = cache->_highIntList.fwdIterator() ;
509 RooAbsReal *funcInt(0), *low(0), *high(0), *param(0) ;
510 Double_t value(0) ;
511 Double_t nominal(0);
512
513 // get nominal
514 int i=0;
515 while(( funcInt = (RooAbsReal*)funcIntIter.next())) {
516 value += funcInt->getVal() ;
517 nominal = value;
518 i++;
519 }
520 if(i==0 || i>1) { cout << "problem, wrong number of nominal functions"<<endl; }
521
522 // now get low/high variations
523 i = 0;
524 RooFIter paramIter(_paramSet.fwdIterator()) ;
525
526 // KC: old interp code with new iterator
527 while( (param=(RooAbsReal*)paramIter.next()) ) {
528 low = (RooAbsReal*)lowIntIter.next() ;
529 high = (RooAbsReal*)highIntIter.next() ;
530
531 if(param->getVal()>0) {
532 value += param->getVal()*(high->getVal() - nominal );
533 } else {
534 value += param->getVal()*(nominal - low->getVal());
535 }
536 ++i;
537 }
538
539 /* // MB : old bit of interpolation code
540 while( (param=(RooAbsReal*)_paramIter->Next()) ) {
541 low = (RooAbsReal*)lowIntIter->Next() ;
542 high = (RooAbsReal*)highIntIter->Next() ;
543
544 if(param->getVal()>0) {
545 value += param->getVal()*(high->getVal() - nominal );
546 } else {
547 value += param->getVal()*(nominal - low->getVal());
548 }
549 ++i;
550 }
551 */
552
553 /* KC: the code below is wrong. Can't pull out a constant change to a non-linear shape deformation.
554 while( (param=(RooAbsReal*)paramIter.next()) ) {
555 low = (RooAbsReal*)lowIntIter.next() ;
556 high = (RooAbsReal*)highIntIter.next() ;
557
558 if(_interpCode.empty() || _interpCode.at(i)==0){
559 // piece-wise linear
560 if(param->getVal()>0)
561 value += param->getVal()*(high->getVal() - nominal );
562 else
563 value += param->getVal()*(nominal - low->getVal());
564 } else if(_interpCode.at(i)==1){
565 // pice-wise log
566 if(param->getVal()>=0)
567 value *= pow(high->getVal()/nominal, +param->getVal());
568 else
569 value *= pow(low->getVal()/nominal, -param->getVal());
570 } else if(_interpCode.at(i)==2){
571 // parabolic with linear
572 double a = 0.5*(high->getVal()+low->getVal())-nominal;
573 double b = 0.5*(high->getVal()-low->getVal());
574 double c = 0;
575 if(param->getVal()>1 ){
576 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
577 } else if(param->getVal()<-1 ) {
578 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
579 } else {
580 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
581 }
582 } else if(_interpCode.at(i)==3){
583 //parabolic version of log-normal
584 double a = 0.5*(high->getVal()+low->getVal())-nominal;
585 double b = 0.5*(high->getVal()-low->getVal());
586 double c = 0;
587 if(param->getVal()>1 ){
588 value += (2*a+b)*(param->getVal()-1)+high->getVal()-nominal;
589 } else if(param->getVal()<-1 ) {
590 value += -1*(2*a-b)*(param->getVal()+1)+low->getVal()-nominal;
591 } else {
592 value += a*pow(param->getVal(),2) + b*param->getVal()+c;
593 }
594
595 } else {
596 coutE(InputArguments) << "PiecewiseInterpolation::analyticalIntegralWN ERROR: " << param->GetName()
597 << " with unknown interpolation code" << endl ;
598 }
599 ++i;
600 }
601 */
602
603 // cout << "value = " << value <<endl;
604 return value;
605}
606
607
608////////////////////////////////////////////////////////////////////////////////
609
611 int index = _paramSet.index(&param);
612 if(index<0){
613 coutE(InputArguments) << "PiecewiseInterpolation::setInterpCode ERROR: " << param.GetName()
614 << " is not in list" << endl ;
615 } else {
616 coutW(InputArguments) << "PiecewiseInterpolation::setInterpCode : " << param.GetName()
617 << " is now " << code << endl ;
618 _interpCode.at(index) = code;
619 }
620}
621
622
623////////////////////////////////////////////////////////////////////////////////
624
626 for(unsigned int i=0; i<_interpCode.size(); ++i){
627 _interpCode.at(i) = code;
628 }
629}
630
631
632////////////////////////////////////////////////////////////////////////////////
633
635 for(unsigned int i=0; i<_interpCode.size(); ++i){
636 coutI(InputArguments) <<"interp code for " << _paramSet.at(i)->GetName() << " = " << _interpCode.at(i) <<endl;
637 }
638}
639
640
641////////////////////////////////////////////////////////////////////////////////
642/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
643
645{
646 return _nominal.arg().binBoundaries(obs,xlo,xhi) ;
647}
648
649
650////////////////////////////////////////////////////////////////////////////////
651/// WVE note: assumes nominal and alternates have identical structure, must add explicit check
652
654{
655 return _nominal.arg().isBinnedDistribution(obs) ;
656}
657
658
659
660////////////////////////////////////////////////////////////////////////////////
661
663{
664 return _nominal.arg().plotSamplingHint(obs,xlo,xhi) ;
665}
666
667////////////////////////////////////////////////////////////////////////////////
668/// Stream an object of class PiecewiseInterpolation.
669
671{
672 if (R__b.IsReading()) {
673 R__b.ReadClassBuffer(PiecewiseInterpolation::Class(),this);
674 specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
675 if (_interpCode.empty()) _interpCode.resize(_paramSet.getSize());
676 } else {
677 R__b.WriteClassBuffer(PiecewiseInterpolation::Class(),this);
678 }
679}
680
681
682/*
683////////////////////////////////////////////////////////////////////////////////
684/// Customized printing of arguments of a PiecewiseInterpolation to more intuitively reflect the contents of the
685/// product operator construction
686
687void PiecewiseInterpolation::printMetaArgs(ostream& os) const
688{
689 _lowIter->Reset() ;
690 if (_highIter) {
691 _highIter->Reset() ;
692 }
693
694 Bool_t first(kTRUE) ;
695
696 RooAbsArg* arg1, *arg2 ;
697 if (_highSet.getSize()!=0) {
698
699 while((arg1=(RooAbsArg*)_lowIter->Next())) {
700 if (!first) {
701 os << " + " ;
702 } else {
703 first = kFALSE ;
704 }
705 arg2=(RooAbsArg*)_highIter->Next() ;
706 os << arg1->GetName() << " * " << arg2->GetName() ;
707 }
708
709 } else {
710
711 while((arg1=(RooAbsArg*)_lowIter->Next())) {
712 if (!first) {
713 os << " + " ;
714 } else {
715 first = kFALSE ;
716 }
717 os << arg1->GetName() ;
718 }
719
720 }
721
722 os << " " ;
723}
724
725*/
#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 e(i)
Definition RSha256.hxx:103
#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:91
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
double pow(double, double)
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.
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...
Bool_t setBinIntegrator(RooArgSet &allVars)
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Advertise that all integrals can be handled internally.
Double_t evaluate() const
Calculate and return current value of self.
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...
void setInterpCode(RooAbsReal &param, int code)
virtual ~PiecewiseInterpolation()
Destructor.
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:72
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
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
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.
RooFIter fwdIterator() const
One-time forward iterator.
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
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:61
Bool_t _forceNumInt
Definition RooAbsReal.h:478
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:91
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 ...
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:341
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:70
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
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.
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