00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #include "RooFit.h"
00024
00025 #include "Riostream.h"
00026 #include "Riostream.h"
00027 #include <math.h>
00028
00029 #include "RooChebychev.h"
00030 #include "RooAbsReal.h"
00031 #include "RooRealVar.h"
00032 #include "RooArgList.h"
00033
00034 ClassImp(RooChebychev)
00035 ;
00036
00037
00038
00039 RooChebychev::RooChebychev()
00040 {
00041 }
00042
00043
00044
00045 RooChebychev::RooChebychev(const char* name, const char* title,
00046 RooAbsReal& x, const RooArgList& coefList):
00047 RooAbsPdf(name, title),
00048 _x("x", "Dependent", this, x),
00049 _coefList("coefficients","List of coefficients",this)
00050 {
00051
00052 TIterator* coefIter = coefList.createIterator() ;
00053 RooAbsArg* coef ;
00054 while((coef = (RooAbsArg*)coefIter->Next())) {
00055 if (!dynamic_cast<RooAbsReal*>(coef)) {
00056 cout << "RooChebychev::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName()
00057 << " is not of type RooAbsReal" << endl ;
00058 assert(0) ;
00059 }
00060 _coefList.add(*coef) ;
00061 }
00062 delete coefIter ;
00063 }
00064
00065
00066
00067
00068 RooChebychev::RooChebychev(const RooChebychev& other, const char* name) :
00069 RooAbsPdf(other, name),
00070 _x("x", this, other._x),
00071 _coefList("coefList",this,other._coefList)
00072 {
00073 }
00074
00075 inline static double p0(double ,double a) { return a; }
00076 inline static double p1(double t,double a,double b) { return a*t+b; }
00077 inline static double p2(double t,double a,double b,double c) { return p1(t,p1(t,a,b),c); }
00078 inline static double p3(double t,double a,double b,double c,double d) { return p2(t,p1(t,a,b),c,d); }
00079 inline static double p4(double t,double a,double b,double c,double d,double e) { return p3(t,p1(t,a,b),c,d,e); }
00080
00081
00082
00083 Double_t RooChebychev::evaluate() const
00084 {
00085
00086 Double_t xmin = _x.min(); Double_t xmax = _x.max();
00087 Double_t x(-1+2*(_x-xmin)/(xmax-xmin));
00088 Double_t x2(x*x);
00089 Double_t sum(0) ;
00090 switch (_coefList.getSize()) {
00091 case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x*p3(x2,64,-112,56,-7);
00092 case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*p3(x2,32,-48,18,-1);
00093 case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x*p2(x2,16,-20,5);
00094 case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*p2(x2,8,-8,1);
00095 case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x*p1(x2,4,-3);
00096 case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*p1(x2,2,-1);
00097 case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x;
00098 case 0: sum+=1;
00099 }
00100 return sum;
00101 }
00102
00103
00104
00105 Int_t RooChebychev::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
00106 {
00107
00108 if (rangeName && strlen(rangeName)) {
00109 return 0 ;
00110 }
00111
00112 if (matchArgs(allVars, analVars, _x)) return 1;
00113 return 0;
00114 }
00115
00116
00117
00118 Double_t RooChebychev::analyticalIntegral(Int_t code, const char* rangeName) const
00119 {
00120 assert(code==1) ;
00121
00122 Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
00123 Double_t norm(0) ;
00124 switch(_coefList.getSize()) {
00125
00126 case 7: case 6: norm+=((RooAbsReal&)_coefList[5]).getVal()*(-1 + 18./3. - 48./5. + 32./7.);
00127
00128 case 5: case 4: norm+=((RooAbsReal&)_coefList[3]).getVal()*( 1 - 8./3. + 8./5.);
00129
00130 case 3: case 2: norm+=((RooAbsReal&)_coefList[1]).getVal()*(-1 + 2./3.);
00131 case 1: case 0: norm+= 1;
00132 }
00133 norm *= xmax-xmin;
00134 return norm;
00135 }