RooBernstein.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  * @(#)root/roofit:$Id: RooBernstein.cxx 28259 2009-04-16 16:21:16Z wouter $
00005  * Authors:                                                                  *
00006  *   Kyle Cranmer
00007  *                                                                           *
00008  *****************************************************************************/
00009 
00010 //////////////////////////////////////////////////////////////////////////////
00011 //
00012 // BEGIN_HTML
00013 // Bernstein basis polynomials are positive-definite in the range [0,1].
00014 // In this implementation, we extend [0,1] to be the range of the parameter.
00015 // There are n+1 Bernstein basis polynomials of degree n.
00016 // Thus, by providing N coefficients that are positive-definite, there 
00017 // is a natural way to have well bahaved polynomail PDFs.
00018 // For any n, the n+1 basis polynomials 'form a partition of unity', eg.
00019 //  they sum to one for all values of x. See
00020 // http://www.idav.ucdavis.edu/education/CAGDNotes/Bernstein-Polynomials.pdf
00021 // END_HTML
00022 //
00023 
00024 #include "RooFit.h"
00025 
00026 #include "Riostream.h"
00027 #include "Riostream.h"
00028 #include <math.h>
00029 #include "TMath.h"
00030 #include "RooBernstein.h"
00031 #include "RooAbsReal.h"
00032 #include "RooRealVar.h"
00033 #include "RooArgList.h"
00034 
00035 ClassImp(RooBernstein)
00036 ;
00037 
00038 
00039 //_____________________________________________________________________________
00040 RooBernstein::RooBernstein()
00041 {
00042 }
00043 
00044 
00045 //_____________________________________________________________________________
00046 RooBernstein::RooBernstein(const char* name, const char* title, 
00047                            RooAbsReal& x, const RooArgList& coefList): 
00048   RooAbsPdf(name, title),
00049   _x("x", "Dependent", this, x),
00050   _coefList("coefficients","List of coefficients",this)
00051 {
00052   // Constructor
00053   TIterator* coefIter = coefList.createIterator() ;
00054   RooAbsArg* coef ;
00055   while((coef = (RooAbsArg*)coefIter->Next())) {
00056     if (!dynamic_cast<RooAbsReal*>(coef)) {
00057       cout << "RooBernstein::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName() 
00058            << " is not of type RooAbsReal" << endl ;
00059       assert(0) ;
00060     }
00061     _coefList.add(*coef) ;
00062   }
00063   delete coefIter ;
00064 }
00065 
00066 
00067 
00068 //_____________________________________________________________________________
00069 RooBernstein::RooBernstein(const RooBernstein& other, const char* name) :
00070   RooAbsPdf(other, name), 
00071   _x("x", this, other._x), 
00072   _coefList("coefList",this,other._coefList)
00073 {
00074 }
00075 
00076 
00077 //_____________________________________________________________________________
00078 Double_t RooBernstein::evaluate() const 
00079 {
00080 
00081   Double_t xmin = _x.min(); Double_t xmax = _x.max();
00082   Int_t degree= _coefList.getSize()-1; // n+1 polys of degree n
00083 
00084   Double_t temp=0, tempx=0;
00085   for (int i=0; i<=degree; ++i){
00086     tempx = (_x-xmin)/(xmax-xmin); // rescale to [0,1]
00087     temp += ((RooAbsReal&) _coefList[i]).getVal() *
00088       TMath::Binomial(degree, i) * pow(tempx,i) * pow(1-tempx,degree-i);
00089   }
00090   return temp;
00091 
00092 }
00093 
00094 
00095 //_____________________________________________________________________________
00096 Int_t RooBernstein::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const 
00097 {
00098   // No analytical calculation available (yet) of integrals over subranges
00099   if (rangeName && strlen(rangeName)) {
00100     return 0 ;
00101   }
00102 
00103   if (matchArgs(allVars, analVars, _x)) return 1;
00104   return 0;
00105 }
00106 
00107 
00108 //_____________________________________________________________________________
00109 Double_t RooBernstein::analyticalIntegral(Int_t code, const char* rangeName) const 
00110 {
00111   assert(code==1) ;
00112   Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
00113   Int_t degree= _coefList.getSize()-1; // n+1 polys of degree n
00114   Double_t norm(0) ;
00115 
00116   Double_t temp=0;
00117   for (int i=0; i<=degree; ++i){
00118     // for each of the i Bernstein basis polynomials
00119     // represent it in the 'power basis' (the naive polynomial basis)
00120     // where the integral is straight forward.
00121     temp = 0;
00122     for (int j=i; j<=degree; ++j){ // power basisŧ
00123       temp += pow(-1.,j-i) * TMath::Binomial(degree, j) * TMath::Binomial(j,i) / (j+1);
00124     }
00125     temp *= ((RooAbsReal&) _coefList[i]).getVal(); // include coeff
00126     norm += temp; // add this basis's contribution to total
00127   }
00128 
00129   norm *= xmax-xmin;
00130   return norm;
00131 }

Generated on Tue Jul 5 14:54:33 2011 for ROOT_528-00b_version by  doxygen 1.5.1