RooSpHarmonic.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitModels                                                     *
00004  *    File: $Id: RooSpHarmonic.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 //   Implementation of the so-called real spherical harmonics, using the orthonormal normalization,
00019 // which are related to spherical harmonics as:
00020 //
00021 //  Y_{l0} = Y_l^0   (m=0)
00022 //  Y_{lm} = \frac{1}{\sqrt{2}}  \left( Y_l^m     + (-1)^m     Y_l^{-m}   \right) (m>0)
00023 //  Y_{lm} = \frac{1}{i\sqrt{2}} \left( Y_l^{|m|} - (-1)^{|m|} Y_l^{-|m|} \right) (m<0)
00024 //
00025 // which implies:
00026 //
00027 // Y_{l0}(\cos\theta,\phi) =          N_{l0}   P_l^0    (\cos\theta)              (m=0)
00028 // Y_{lm}(\cos\theta,\phi) = \sqrt{2} N_{lm}   P_l^m    (\cos\theta) cos(|m|\phi) (m>0)
00029 // Y_{lm}(\cos\theta,\phi) = \sqrt{2} N_{l|m|} P_l^{|m|}(\cos\theta) sin(|m|\phi) (m<0)
00030 //
00031 // where
00032 //  N_{lm} = \sqrt{ \frac{2l+1}{4\pi} \frac{ (l-m)! }{ (l+m)! } } 
00033 //
00034 // Note that the normalization corresponds to the orthonormal case,
00035 // and thus we have Int d\cos\theta d\phi Y_{lm} Y_{l'm'} = \delta_{ll'} \delta{mm'}
00036 //
00037 // Note that in addition, this code can also represent the product of two
00038 // (real) spherical harmonics -- it actually uses the fact that Y_{00} = \sqrt{\frac{1}{4\pi}}
00039 // in order to represent a single spherical harmonics by multiplying it
00040 // by \sqrt{4\pi} Y_00, as this makes it trivial to compute the analytical
00041 // integrals, using the orthogonality properties of Y_l^m...
00042 //
00043 // END_HTML
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   // TODO: check that phi.max - phi.min = 2 pi... ctheta.max = +1, and ctheta.min = -1
00120   // we don't support indefinite integrals... maybe one day, when there is a use for it.....
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 }

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