00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include "RooFit.h"
00022 #include "Riostream.h"
00023 #include <math.h>
00024 #include <string>
00025
00026 #include "RooLegendre.h"
00027 #include "RooAbsReal.h"
00028 #include "Math/SpecFunc.h"
00029 #include "TMath.h"
00030
00031 ClassImp(RooLegendre)
00032 ;
00033
00034
00035
00036 namespace {
00037 inline double a(int p, int l, int m) {
00038 double r = TMath::Factorial(l+m)/TMath::Factorial(m+p)/TMath::Factorial(p)/TMath::Factorial(l-m-2*p);
00039 r /= pow(2.,m+2*p);
00040 return p%2==0 ? r : -r ;
00041 }
00042 }
00043
00044
00045 RooLegendre::RooLegendre()
00046 {
00047 }
00048
00049
00050 RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l, int m)
00051 : RooAbsReal(name, title)
00052 , _ctheta("ctheta", "ctheta", this, ctheta)
00053 , _l1(l),_m1(m),_l2(0),_m2(0)
00054 {
00055
00056
00057
00058
00059
00060 }
00061
00062
00063 RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l1, int m1, int l2, int m2)
00064 : RooAbsReal(name, title)
00065 , _ctheta("ctheta", "ctheta", this, ctheta)
00066 , _l1(l1),_m1(m1),_l2(l2),_m2(m2)
00067 {
00068 }
00069
00070
00071 RooLegendre::RooLegendre(const RooLegendre& other, const char* name)
00072 : RooAbsReal(other, name)
00073 , _ctheta("ctheta", this, other._ctheta)
00074 , _l1(other._l1), _m1(other._m1)
00075 , _l2(other._l2), _m2(other._m2)
00076 {
00077 }
00078
00079
00080 Double_t RooLegendre::evaluate() const
00081 {
00082
00083
00084 #ifdef R__HAS_MATHMORE
00085 double r = 1;
00086 if (_l1!=0||_m1!=0) r *= ROOT::Math::assoc_legendre(_l1,_m1,_ctheta);
00087 if (_l2!=0||_m2!=0) r *= ROOT::Math::assoc_legendre(_l2,_m2,_ctheta);
00088 if ((_m1+_m2)%2==1) r = -r;
00089 return r;
00090 #else
00091 throw std::string("RooLegendre: ERROR: This class require installation of the MathMore library") ;
00092 return 0 ;
00093 #endif
00094 }
00095
00096
00097 namespace {
00098 bool fullRange(const RooRealProxy& x ,const char* range)
00099 { return range==0 || strlen(range)==0
00100 || ( x.min(range) == x.min() && x.max(range) == x.max() ) ;
00101 }
00102 }
00103 Int_t RooLegendre::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
00104 {
00105
00106 if (fullRange(_ctheta,rangeName) && matchArgs(allVars, analVars, _ctheta)) return 1;
00107 return 0;
00108 }
00109
00110
00111 Double_t RooLegendre::analyticalIntegral(Int_t code, const char* ) const
00112 {
00113
00114
00115 assert(code==1) ;
00116 if ( _m1==_m2 ) return ( _l1 == _l2) ? TMath::Factorial(_l1+_m2)/TMath::Factorial(_l1-_m1)*double(2)/(2*_l1+1) : 0.;
00117 if ( (_l1+_l2-_m1-_m2)%2 != 0 ) return 0;
00118
00119
00120
00121
00122
00123
00124
00125
00126 double r=0;
00127 for (int p1=0; 2*p1 <= _l1-_m1 ;++p1) {
00128 double a1 = a(p1,_l1,_m1);
00129 for (int p2=0; 2*p2 <= _l2-_m2 ; ++p2) {
00130 double a2 = a(p2,_l2,_m2);
00131 r+= a1*a2*TMath::Gamma( double(_l1+_l2-_m1-_m2-2*p1-2*p2+1)/2 )*TMath::Gamma( double(_m1+_m2+2*p1+2*p2+2)/2 );
00132 }
00133 }
00134 r /= TMath::Gamma( double(_l1+_l2+3)/2 );
00135
00136 if ((_m1+_m2)%2==1) r = -r;
00137 return r;
00138 }
00139
00140 Int_t RooLegendre::getMaxVal( const RooArgSet& ) const {
00141 if (_m1==0&&_m2==0) return 1;
00142
00143 if (_l1<3&&_l2<3) return 1;
00144 return 0;
00145 }
00146
00147 namespace {
00148 inline double maxSingle(int i, int j) {
00149 assert(j<=i);
00150
00151 if (j==0) return 1;
00152 assert(i<3);
00153
00154 if (i<2) return 1;
00155
00156 static const double m2[3] = { 3,3 };
00157 return m2[j-1];
00158 }
00159 }
00160 Double_t RooLegendre::maxVal( Int_t ) const {
00161 return maxSingle(_l1,_m1)*maxSingle(_l2,_m2);
00162 }