ROOT  6.06/09
Reference Guide
RooSpHarmonic.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * File: $Id$
5  * Authors: *
6  * GR, Gerhard Raven, Nikhef & VU, Gerhard.Raven@nikhef.nl
7  * *
8  * Copyright (c) 2010, Nikhef & VU. All rights reserved.
9  * *
10  * Redistribution and use in source and binary forms, *
11  * with or without modification, are permitted according to the terms *
12  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
13  *****************************************************************************/
14 
15 //////////////////////////////////////////////////////////////////////////////
16 //
17 // BEGIN_HTML
18 // Implementation of the so-called real spherical harmonics, using the orthonormal normalization,
19 // which are related to spherical harmonics as:
20 //
21 // Y_{l0} = Y_l^0 (m=0)
22 // Y_{lm} = \frac{1}{\sqrt{2}} \left( Y_l^m + (-1)^m Y_l^{-m} \right) (m>0)
23 // Y_{lm} = \frac{1}{i\sqrt{2}} \left( Y_l^{|m|} - (-1)^{|m|} Y_l^{-|m|} \right) (m<0)
24 //
25 // which implies:
26 //
27 // Y_{l0}(\cos\theta,\phi) = N_{l0} P_l^0 (\cos\theta) (m=0)
28 // Y_{lm}(\cos\theta,\phi) = \sqrt{2} N_{lm} P_l^m (\cos\theta) cos(|m|\phi) (m>0)
29 // Y_{lm}(\cos\theta,\phi) = \sqrt{2} N_{l|m|} P_l^{|m|}(\cos\theta) sin(|m|\phi) (m<0)
30 //
31 // where
32 // N_{lm} = \sqrt{ \frac{2l+1}{4\pi} \frac{ (l-m)! }{ (l+m)! } }
33 //
34 // Note that the normalization corresponds to the orthonormal case,
35 // and thus we have Int d\cos\theta d\phi Y_{lm} Y_{l'm'} = \delta_{ll'} \delta{mm'}
36 //
37 // Note that in addition, this code can also represent the product of two
38 // (real) spherical harmonics -- it actually uses the fact that Y_{00} = \sqrt{\frac{1}{4\pi}}
39 // in order to represent a single spherical harmonics by multiplying it
40 // by \sqrt{4\pi} Y_00, as this makes it trivial to compute the analytical
41 // integrals, using the orthogonality properties of Y_l^m...
42 //
43 // END_HTML
44 //
45 
46 #include "RooFit.h"
47 #include "Riostream.h"
48 #include <math.h>
49 
50 #include "RooSpHarmonic.h"
51 #include "Math/SpecFunc.h"
52 #include "TMath.h"
53 
54 using namespace std;
55 
57 ;
58 
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 
62 namespace {
63  inline double N(int l, int m=0) {
64  double n = sqrt( double(2*l+1)/(4*TMath::Pi())*TMath::Factorial(l-m)/TMath::Factorial(l+m) );
65  return m==0 ? n : TMath::Sqrt2() * n;
66  }
67 }
68 
69 ////////////////////////////////////////////////////////////////////////////////
70 
72  _n(0),
73  _sgn1(0),
74  _sgn2(0)
75 {
76 }
77 
78 ////////////////////////////////////////////////////////////////////////////////
79 
80 RooSpHarmonic::RooSpHarmonic(const char* name, const char* title, RooAbsReal& ctheta, RooAbsReal& phi, int l, int m)
81  : RooLegendre(name, title,ctheta,l,m<0?-m:m)
82  , _phi("phi", "phi", this, phi)
83  , _n( 2*sqrt(TMath::Pi()))
84  , _sgn1( m==0 ? 0 : m<0 ? -1 : +1 )
85  , _sgn2( 0 )
86 {
87 }
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 
91 RooSpHarmonic::RooSpHarmonic(const char* name, const char* title, RooAbsReal& ctheta, RooAbsReal& phi, int l1, int m1, int l2, int m2)
92  : RooLegendre(name, title,ctheta,l1, m1<0?-m1:m1,l2,m2<0?-m2:m2)
93  , _phi("phi", "phi", this, phi)
94  , _n(1)
95  , _sgn1( m1==0 ? 0 : m1<0 ? -1 : +1 )
96  , _sgn2( m2==0 ? 0 : m2<0 ? -1 : +1 )
97 {
98 }
99 
100 ////////////////////////////////////////////////////////////////////////////////
101 
103  : RooLegendre(other, name)
104  , _phi("phi", this,other._phi)
105  , _n(other._n)
106  , _sgn1( other._sgn1 )
107  , _sgn2( other._sgn2 )
108 {
109 }
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 
114 {
115  double n = _n*N(_l1,_m1)*N(_l2,_m2)*RooLegendre::evaluate();
116  if (_sgn1!=0) n *= (_sgn1<0 ? sin(_m1*_phi) : cos(_m1*_phi) );
117  if (_sgn2!=0) n *= (_sgn2<0 ? sin(_m2*_phi) : cos(_m2*_phi) );
118  return n;
119 }
120 
121 //_____________________________________________________________________________
122 namespace {
123  Bool_t fullRange(const RooRealProxy& x, const char* range, Bool_t phi)
124  {
125  if (phi) {
126  return range == 0 || strlen(range) == 0
127  ? std::fabs(x.max() - x.min() - TMath::TwoPi()) < 1.e-8
128  : std::fabs(x.max(range) - x.min(range) - TMath::TwoPi()) < 1.e-8;
129  }
130 
131  return range == 0 || strlen(range) == 0
132  ? std::fabs(x.min() + 1.) < 1.e-8 && std::fabs(x.max() - 1.) < 1.e-8
133  : std::fabs(x.min(range) + 1.) < 1.e-8 && std::fabs(x.max(range) - 1.) < 1.e-8;
134  }
135 }
136 
137 ////////////////////////////////////////////////////////////////////////////////
138 /// TODO: check that phi.max - phi.min = 2 pi... ctheta.max = +1, and ctheta.min = -1
139 /// we don't support indefinite integrals... maybe one day, when there is a use for it.....
140 
141 Int_t RooSpHarmonic::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
142 {
143  // we don't support indefinite integrals... maybe one day, when there is a use for it.....
144  Bool_t cthetaOK = fullRange(_ctheta, rangeName, kFALSE);
145  Bool_t phiOK = fullRange(_phi, rangeName, kTRUE );
146  if (cthetaOK && phiOK && matchArgs(allVars, analVars, _ctheta, _phi)) return 3;
147  if ( phiOK && matchArgs(allVars, analVars, _phi)) return 2;
148  return RooLegendre::getAnalyticalIntegral(allVars, analVars, rangeName);
149 }
150 
151 ////////////////////////////////////////////////////////////////////////////////
152 
153 Double_t RooSpHarmonic::analyticalIntegral(Int_t code, const char* range) const
154 {
155  if (code==3) {
156  return (_l1==_l2 && _sgn1*_m1==_sgn2*_m2 ) ? _n : 0 ;
157  } else if (code == 2) {
158  if ( _sgn1*_m1 != _sgn2*_m2) return 0;
159  return ( _m1==0 ? 2 : 1 ) * TMath::Pi()*_n*N(_l1,_m1)*N(_l2,_m2)*RooLegendre::evaluate();
160  } else {
161  double n = _n*N(_l1,_m1)*N(_l2,_m2)*RooLegendre::analyticalIntegral(code,range);
162  if (_sgn1!=0) n *= (_sgn1<0 ? sin(_m1*_phi) : cos(_m1*_phi) );
163  if (_sgn2!=0) n *= (_sgn2<0 ? sin(_m2*_phi) : cos(_m2*_phi) );
164  return n;
165  }
166 }
167 
169  return RooLegendre::getMaxVal(vars);
170 }
171 
173  double n = _n*N(_l1,_m1)*N(_l2,_m2);
174  return n*RooLegendre::maxVal(code);
175 }
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Double_t Sqrt2()
Definition: TMath.h:51
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
#define N
RooRealProxy _ctheta
Definition: RooLegendre.h:40
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
this was verified to match mathematica for l1 in [0,2], m1 in [0,l1], l2 in [l1,4], m2 in [0,l2]
STL namespace.
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
double cos(double)
double sqrt(double)
RooRealProxy _phi
Definition: RooSpHarmonic.h:37
Double_t TwoPi()
Definition: TMath.h:45
double sin(double)
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
Double_t evaluate() const
TODO: check that 0<=m_i<=l_i; on the other hand, assoc_legendre already does that ;-) Note: P_0^0 = 1...
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
this was verified to match mathematica for l1 in [0,2], m1 in [0,l1], l2 in [l1,4], m2 in [0,l2]
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
double Pi()
Mathematical constants.
Definition: Math.h:68
TMarker * m
Definition: textangle.C:8
Double_t Pi()
Definition: TMath.h:44
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
ClassImp(RooSpHarmonic)
#define name(a, b)
Definition: linkTestLib0.cpp:5
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
Double_t evaluate() const
TODO: check that 0<=m_i<=l_i; on the other hand, assoc_legendre already does that ;-) Note: P_0^0 = 1...
Definition: RooLegendre.cxx:93
Double_t Factorial(Int_t i)
Compute factorial(n).
Definition: TMath.cxx:250
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
TODO: check that phi.max - phi.min = 2 pi...
const Int_t n
Definition: legend1.C:16
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57