#include <RooBinnedPdf.hh>
Public Member Functions | |
RooBinnedPdf (const char *name, const char *title, RooAbsReal &x, const RooArgList &coefList, const TArrayD &limits) | |
RooBinnedPdf (const RooBinnedPdf &other, const char *name=0) | |
virtual TObject * | clone (const char *newname) const |
virtual | ~RooBinnedPdf () |
Int_t | getAnalyticalIntegral (RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const |
Double_t | analyticalIntegral (Int_t code, const char *rangeName=0) const |
Int_t | getnBins () |
Double_t * | getLimits () |
Protected Member Functions | |
Double_t | evaluate () const |
Protected Attributes | |
RooRealProxy | _x |
RooListProxy | _coefList |
TArrayD | _limits |
Int_t | _nBins |
Private Member Functions | |
Double_t | localEval (const Double_t) const |
Definition at line 30 of file RooBinnedPdf.hh.
RooBinnedPdf::RooBinnedPdf | ( | const char * | name, | |
const char * | title, | |||
RooAbsReal & | x, | |||
const RooArgList & | coefList, | |||
const TArrayD & | limits | |||
) |
Definition at line 65 of file RooBinnedPdf.cc.
References _coefList, _limits, and _nBins.
Referenced by clone().
00066 : 00067 RooAbsPdf(name, title), 00068 _x("x", "Dependent", this, x), 00069 _coefList("coefList","List of coefficients",this), 00070 _nBins(limits.GetSize()-1) 00071 { 00072 // Check lowest order 00073 if (_nBins<0) { 00074 cout << "RooBinnedPdf::ctor(" << GetName() 00075 << ") WARNING: nBins must be >=0, setting value to 0" << endl ; 00076 _nBins=0 ; 00077 } 00078 00079 TIterator* coefIter = coefList.createIterator() ; 00080 RooAbsArg* coef ; 00081 while(coef = (RooAbsArg*)coefIter->Next()) { 00082 if (!dynamic_cast<RooAbsReal*>(coef)) { 00083 cout << "RooBinnedPdf::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName() 00084 << " is not of type RooAbsReal" << endl ; 00085 assert(0) ; 00086 } 00087 _coefList.add(*coef) ; 00088 } 00089 delete coefIter ; 00090 00091 // Bin limits 00092 limits.Copy(_limits); 00093 }
RooBinnedPdf::RooBinnedPdf | ( | const RooBinnedPdf & | other, | |
const char * | name = 0 | |||
) |
RooBinnedPdf::~RooBinnedPdf | ( | ) | [virtual] |
Double_t RooBinnedPdf::analyticalIntegral | ( | Int_t | code, | |
const char * | rangeName = 0 | |||
) | const |
Definition at line 121 of file RooBinnedPdf.cc.
References _limits, _nBins, _x, and localEval().
00122 { 00123 assert(code==1) ; 00124 00125 // Old method: 00126 // Double_t sum(1.0) ; 00127 // return sum; 00128 00129 Double_t integral(0.0); Double_t width(0.0); 00130 Double_t average(0.0); 00131 Double_t min(_x.min(rangeName)); Double_t max(_x.max(rangeName)); 00132 Int_t first(0); Int_t last(_nBins-1); 00133 00134 if (min >= _limits[0] && max <= _limits[_nBins]){ 00135 00136 // Get first and last bin 00137 int i=0 ; 00138 for (i=0; i<_nBins && min>=_limits[i];++i){;} 00139 first = i-1; if(first < 0) first = 0; 00140 for (i=0; i<_nBins && max>=_limits[i];++i){;} 00141 last = i-1; if(last > _nBins-1) last = _nBins-1; 00142 00143 if(first != last) { 00144 width = _limits[first+1] - _limits[first]; 00145 average = _limits[first] + (width/2.0); 00146 width = _limits[first+1] - min; 00147 integral += (width * localEval(average)); 00148 // integral += (width * static_cast<RooRealVar*>(_coefList.at(first))->getVal()); 00149 00150 width = _limits[last+1] - _limits[last]; 00151 average = _limits[last] + (width/2.0); 00152 width = max - _limits[last]; 00153 integral += (width * localEval(average)); 00154 // integral += (width * static_cast<RooRealVar*>(_coefList.at(last))->getVal()); 00155 00156 for(i=first+1; i < last; i++) { 00157 width = _limits[i+1] - _limits[i]; 00158 average = _limits[i] + (width/2.0); 00159 integral += (width * localEval(average)); 00160 // integral += (width * static_cast<RooRealVar*>(_coefList.at(i))->getVal()); 00161 } 00162 } 00163 else { 00164 width = _limits[first+1] - _limits[first]; 00165 average = _limits[first] + (width/2.0); 00166 width = max - min; 00167 integral += (width * localEval(average)); 00168 // integral += (width * static_cast<RooRealVar*>(_coefList.at(first))->getVal()); 00169 } 00170 } 00171 return integral; 00172 }
virtual TObject* RooBinnedPdf::clone | ( | const char * | newname | ) | const [inline, virtual] |
Definition at line 41 of file RooBinnedPdf.hh.
References RooBinnedPdf().
00041 { return new RooBinnedPdf(*this, newname); }
Double_t RooBinnedPdf::evaluate | ( | ) | const [protected] |
Int_t RooBinnedPdf::getAnalyticalIntegral | ( | RooArgSet & | allVars, | |
RooArgSet & | analVars, | |||
const char * | rangeName = 0 | |||
) | const |
Definition at line 113 of file RooBinnedPdf.cc.
References _x.
00114 { 00115 if (matchArgs(allVars, analVars, _x)) return 1; 00116 return 0; 00117 }
Double_t * RooBinnedPdf::getLimits | ( | ) |
Definition at line 212 of file RooBinnedPdf.cc.
References _limits.
00212 { 00213 Double_t* limoutput = _limits.GetArray(); 00214 return limoutput; 00215 }
Int_t RooBinnedPdf::getnBins | ( | ) |
Definition at line 208 of file RooBinnedPdf.cc.
References _nBins.
00208 { 00209 return _nBins; 00210 }
Double_t RooBinnedPdf::localEval | ( | const | Double_t | ) | const [private] |
Definition at line 181 of file RooBinnedPdf.cc.
References _coefList, _limits, and _nBins.
Referenced by analyticalIntegral(), and evaluate().
00182 { 00183 Double_t value(0); 00184 if (xval >= _limits[0] && xval < _limits[_nBins]){ 00185 double sum(0),binval(0); 00186 int i=0 ; 00187 for (i=0; i<_nBins-1 && xval>=_limits[i];++i) { 00188 binval = (1-sum)*static_cast<RooRealVar*>(_coefList.at(i))->getVal() ; 00189 //binval = static_cast<RooRealVar*>(_coefList.at(i))->getVal() ; 00190 sum += binval ; 00191 } 00192 if( xval>=_limits[_nBins-1] ) { // the last bin 00193 binval = 1-sum ; 00194 i = _nBins ; 00195 } 00196 double binwidth = _limits[i] - _limits[i-1]; 00197 value = binval/binwidth ; 00198 if (value<=0){ 00199 cout << "RooBinnedPdf: sum of values gt 1.0 -- beware!!" 00200 << value << " " << binval << " " << sum << " " << i << " " << xval << endl; 00201 value = 0.000001; 00202 } 00203 } 00204 return value; 00205 }
RooListProxy RooBinnedPdf::_coefList [protected] |
TArrayD RooBinnedPdf::_limits [protected] |
Definition at line 53 of file RooBinnedPdf.hh.
Referenced by analyticalIntegral(), getLimits(), localEval(), and RooBinnedPdf().
Int_t RooBinnedPdf::_nBins [protected] |
Definition at line 54 of file RooBinnedPdf.hh.
Referenced by analyticalIntegral(), getnBins(), localEval(), and RooBinnedPdf().
RooRealProxy RooBinnedPdf::_x [protected] |
Definition at line 51 of file RooBinnedPdf.hh.
Referenced by analyticalIntegral(), evaluate(), and getAnalyticalIntegral().