00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #include "RooFit.h"
00047 #include "Riostream.h"
00048 #include <math.h>
00049
00050 #include "RooSpHarmonic.h"
00051 #include "Math/SpecFunc.h"
00052 #include "TMath.h"
00053
00054 ClassImp(RooSpHarmonic)
00055 ;
00056
00057
00058
00059 namespace {
00060 inline double N(int l, int m=0) {
00061 double n = sqrt( double(2*l+1)/(4*TMath::Pi())*TMath::Factorial(l-m)/TMath::Factorial(l+m) );
00062 return m==0 ? n : TMath::Sqrt2() * n;
00063 }
00064 }
00065
00066
00067 RooSpHarmonic::RooSpHarmonic()
00068 {
00069 }
00070
00071
00072 RooSpHarmonic::RooSpHarmonic(const char* name, const char* title, RooAbsReal& ctheta, RooAbsReal& phi, int l, int m)
00073 : RooLegendre(name, title,ctheta,l,m<0?-m:m)
00074 , _phi("phi", "phi", this, phi)
00075 , _n( double(4)/(2*sqrt(TMath::Pi())) )
00076 , _sgn1( m==0 ? 0 : m<0 ? -1 : +1 )
00077 , _sgn2( 0 )
00078 {
00079 }
00080
00081
00082 RooSpHarmonic::RooSpHarmonic(const char* name, const char* title, RooAbsReal& ctheta, RooAbsReal& phi, int l1, int m1, int l2, int m2)
00083 : RooLegendre(name, title,ctheta,l1, m1<0?-m1:m1,l2,m2<0?-m2:m2)
00084 , _phi("phi", "phi", this, phi)
00085 , _n(1)
00086 , _sgn1( m1==0 ? 0 : m1<0 ? -1 : +1 )
00087 , _sgn2( m2==0 ? 0 : m2<0 ? -1 : +1 )
00088 {
00089 }
00090
00091
00092 RooSpHarmonic::RooSpHarmonic(const RooSpHarmonic& other, const char* name)
00093 : RooLegendre(other, name)
00094 , _phi("phi", this,other._phi)
00095 , _n(other._n)
00096 , _sgn1( other._sgn1 )
00097 , _sgn2( other._sgn2 )
00098 {
00099 }
00100
00101
00102 Double_t RooSpHarmonic::evaluate() const
00103 {
00104 double n = _n*N(_l1,_m1)*N(_l2,_m2)*RooLegendre::evaluate();
00105 if (_sgn1!=0) n *= (_sgn1<0 ? sin(_m1*_phi) : cos(_m1*_phi) );
00106 if (_sgn2!=0) n *= (_sgn2<0 ? sin(_m2*_phi) : cos(_m2*_phi) );
00107 return n;
00108 }
00109
00110 namespace {
00111 bool fullRange(const RooRealProxy& x ,const char* range) {
00112 return ( x.min(range) == x.min() && x.max(range) == x.max() ) ;
00113 }
00114 }
00115
00116
00117 Int_t RooSpHarmonic::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
00118 {
00119
00120
00121 bool noRange = ( rangeName == 0 || strlen(rangeName)==0 );
00122 bool phiOK = ( noRange || fullRange(_phi,rangeName) );
00123 bool cthetaOK = ( noRange || fullRange(_ctheta,rangeName) );
00124 if (cthetaOK && phiOK && matchArgs(allVars, analVars, _ctheta,_phi)) return 3;
00125 if ( phiOK && matchArgs(allVars, analVars, _phi)) return 2;
00126 return RooLegendre::getAnalyticalIntegral(allVars,analVars,rangeName);
00127 }
00128
00129
00130 Double_t RooSpHarmonic::analyticalIntegral(Int_t code, const char* range) const
00131 {
00132 if (code==3) {
00133 return (_l1==_l2 && _sgn1*_m1==_sgn2*_m2 ) ? _n : 0 ;
00134 } else if (code == 2) {
00135 if (_m1!=0 || _m2!=0) return 0;
00136 return _n*N(_l1)*N(_l2)*2*TMath::Pi()*RooLegendre::evaluate();
00137 } else {
00138 double n = _n*N(_l1,_m1)*N(_l2,_m2)*RooLegendre::analyticalIntegral(code,range);
00139 if (_sgn1!=0) n *= (_sgn1<0 ? sin(_m1*_phi) : cos(_m1*_phi) );
00140 if (_sgn2!=0) n *= (_sgn2<0 ? sin(_m2*_phi) : cos(_m2*_phi) );
00141 return n;
00142 }
00143 }
00144
00145 Int_t RooSpHarmonic::getMaxVal( const RooArgSet& vars) const {
00146 return RooLegendre::getMaxVal(vars);
00147 }
00148
00149 Double_t RooSpHarmonic::maxVal( Int_t code) const {
00150 double n = _n*N(_l1,_m1)*N(_l2,_m2);
00151 return n*RooLegendre::maxVal(code);
00152 }