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