RooLegendre.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  *    File: $Id: RooLegendre.cxx 37380 2010-12-07 21:26:55Z wouter $
00005  * Authors:                                                                  *
00006  *   GR, Gerhard Raven,   Nikhef & VU, Gerhard.Raven@nikhef.nl
00007  *                                                                           *
00008  * Copyright (c) 2010, Nikhef & VU. All rights reserved.
00009  *                                                                           *
00010  * Redistribution and use in source and binary forms,                        *
00011  * with or without modification, are permitted according to the terms        *
00012  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00013  *****************************************************************************/
00014 
00015 //////////////////////////////////////////////////////////////////////////////
00016 //
00017 // BEGIN_HTML
00018 // END_HTML
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   //TODO: for now, we assume that ctheta has a range [-1,1]
00056   // should map the ctheta range onto this interval, and adjust integrals...
00057 
00058   //TODO: we assume m>=0
00059   //      should map m<0 back to m>=0...
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   // TODO: check that 0<=m_i<=l_i; on the other hand, assoc_legendre already does that ;-)
00083   // Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0
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   // don't support indefinite integrals...
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   // this was verified to match mathematica for 
00114   // l1 in [0,2], m1 in [0,l1], l2 in [l1,4], m2 in [0,l2]
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; // these combinations are odd under x -> -x
00118 
00119   // from B.R. Wong, "On the overlap integral of associated Legendre Polynomials" 1998 J. Phys. A: Math. Gen. 31 1101
00120   // TODO: update to the result of 
00121   //       H. A. Mavromatis
00122   //       "A single-sum expression for the overlap integral of two associated Legendre polynomials"
00123   //       1999 J. Phys. A: Math. Gen. 32 2601
00124   //       http://iopscience.iop.org/0305-4470/32/13/011/pdf/0305-4470_32_13_011.pdf
00125   //       For that we need Wigner 3-j, which Lorenzo has added for Root 5.28... (note: check Condon-Shortly convention in this paper!)
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& /*vars*/) const {
00141     if (_m1==0&&_m2==0) return 1;
00142     // does anyone know the analytical expression for the  max values in case m!=0??
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         //   x0 : 1 (ordinary Legendre)
00151         if (j==0) return 1;
00152         assert(i<3);
00153         //   11: 1
00154         if (i<2) return 1;
00155         //   21: 3   22: 3
00156         static const double m2[3] = { 3,3 };
00157         return m2[j-1];
00158     }
00159 }
00160 Double_t RooLegendre::maxVal( Int_t /*code*/) const {
00161     return maxSingle(_l1,_m1)*maxSingle(_l2,_m2);
00162 }

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