ROOT  6.06/09
Reference Guide
RooBernstein.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * Kyle Cranmer
7  * *
8  *****************************************************************************/
9 
10 //////////////////////////////////////////////////////////////////////////////
11 //
12 // BEGIN_HTML
13 // Bernstein basis polynomials are positive-definite in the range [0,1].
14 // In this implementation, we extend [0,1] to be the range of the parameter.
15 // There are n+1 Bernstein basis polynomials of degree n.
16 // Thus, by providing N coefficients that are positive-definite, there
17 // is a natural way to have well bahaved polynomail PDFs.
18 // For any n, the n+1 basis polynomials 'form a partition of unity', eg.
19 // they sum to one for all values of x. See
20 // http://www.idav.ucdavis.edu/education/CAGDNotes/Bernstein-Polynomials.pdf
21 // END_HTML
22 //
23 
24 #include "RooFit.h"
25 
26 #include "Riostream.h"
27 #include "Riostream.h"
28 #include <math.h>
29 #include "TMath.h"
30 #include "RooBernstein.h"
31 #include "RooAbsReal.h"
32 #include "RooRealVar.h"
33 #include "RooArgList.h"
34 
35 #include "TError.h"
36 
37 using namespace std;
38 
40 ;
41 
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 
46 {
47 }
48 
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 /// Constructor
52 
53 RooBernstein::RooBernstein(const char* name, const char* title,
54  RooAbsReal& x, const RooArgList& coefList):
55  RooAbsPdf(name, title),
56  _x("x", "Dependent", this, x),
57  _coefList("coefficients","List of coefficients",this)
58 {
59  TIterator* coefIter = coefList.createIterator() ;
60  RooAbsArg* coef ;
61  while((coef = (RooAbsArg*)coefIter->Next())) {
62  if (!dynamic_cast<RooAbsReal*>(coef)) {
63  cout << "RooBernstein::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName()
64  << " is not of type RooAbsReal" << endl ;
65  R__ASSERT(0) ;
66  }
67  _coefList.add(*coef) ;
68  }
69  delete coefIter ;
70 }
71 
72 
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 
76 RooBernstein::RooBernstein(const RooBernstein& other, const char* name) :
77  RooAbsPdf(other, name),
78  _x("x", this, other._x),
79  _coefList("coefList",this,other._coefList)
80 {
81 }
82 
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 
87 {
88  Double_t xmin = _x.min();
89  Double_t x = (_x - xmin) / (_x.max() - xmin); // rescale to [0,1]
90  Int_t degree = _coefList.getSize() - 1; // n+1 polys of degree n
92 
93  if(degree == 0) {
94 
95  return ((RooAbsReal *)iter.next())->getVal();
96 
97  } else if(degree == 1) {
98 
99  Double_t a0 = ((RooAbsReal *)iter.next())->getVal(); // c0
100  Double_t a1 = ((RooAbsReal *)iter.next())->getVal() - a0; // c1 - c0
101  return a1 * x + a0;
102 
103  } else if(degree == 2) {
104 
105  Double_t a0 = ((RooAbsReal *)iter.next())->getVal(); // c0
106  Double_t a1 = 2 * (((RooAbsReal *)iter.next())->getVal() - a0); // 2 * (c1 - c0)
107  Double_t a2 = ((RooAbsReal *)iter.next())->getVal() - a1 - a0; // c0 - 2 * c1 + c2
108  return (a2 * x + a1) * x + a0;
109 
110  } else if(degree > 2) {
111 
112  Double_t t = x;
113  Double_t s = 1 - x;
114 
115  Double_t result = ((RooAbsReal *)iter.next())->getVal() * s;
116  for(Int_t i = 1; i < degree; i++) {
117  result = (result + t * TMath::Binomial(degree, i) * ((RooAbsReal *)iter.next())->getVal()) * s;
118  t *= x;
119  }
120  result += t * ((RooAbsReal *)iter.next())->getVal();
121 
122  return result;
123  }
124 
125  // in case list of arguments passed is empty
126  return TMath::SignalingNaN();
127 }
128 
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 /// No analytical calculation available (yet) of integrals over subranges
132 
133 Int_t RooBernstein::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
134 {
135  if (rangeName && strlen(rangeName)) {
136  return 0 ;
137  }
138 
139  if (matchArgs(allVars, analVars, _x)) return 1;
140  return 0;
141 }
142 
143 
144 ////////////////////////////////////////////////////////////////////////////////
145 
146 Double_t RooBernstein::analyticalIntegral(Int_t code, const char* rangeName) const
147 {
148  R__ASSERT(code==1) ;
149  Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
150  Int_t degree= _coefList.getSize()-1; // n+1 polys of degree n
151  Double_t norm(0) ;
152 
154  Double_t temp=0;
155  for (int i=0; i<=degree; ++i){
156  // for each of the i Bernstein basis polynomials
157  // represent it in the 'power basis' (the naive polynomial basis)
158  // where the integral is straight forward.
159  temp = 0;
160  for (int j=i; j<=degree; ++j){ // power basis≈ß
161  temp += pow(-1.,j-i) * TMath::Binomial(degree, j) * TMath::Binomial(j,i) / (j+1);
162  }
163  temp *= ((RooAbsReal*)iter.next())->getVal(); // include coeff
164  norm += temp; // add this basis's contribution to total
165  }
166 
167  norm *= xmax-xmin;
168  return norm;
169 }
RooRealProxy _x
Definition: RooBernstein.h:39
float xmin
Definition: THbookFile.cxx:93
ClassImp(RooBernstein)
RooFIter fwdIterator() const
RooListProxy _coefList
Definition: RooBernstein.h:40
#define R__ASSERT(e)
Definition: TError.h:98
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Double_t evaluate() const
Iterator abstract base class.
Definition: TIterator.h:32
Double_t x[n]
Definition: legend1.C:17
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
double pow(double, double)
TIterator * createIterator(Bool_t dir=kIterForward) const
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
No analytical calculation available (yet) of integrals over subranges.
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
Double_t Binomial(Int_t n, Int_t k)
Calculate the binomial coefficient n over k.
Definition: TMath.cxx:2072
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
float xmax
Definition: THbookFile.cxx:93
Double_t SignalingNaN()
Definition: TMath.h:642
RooAbsArg * next()
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
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
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
virtual TObject * Next()=0
double result[121]
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Int_t getSize() const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57