ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
RooLegendre.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 \file RooLegendre.cxx
17 \class RooLegendre
18 \ingroup Roofit
19 
20 **/
21 
22 #include "RooFit.h"
23 #include "Riostream.h"
24 #include <math.h>
25 #include <string>
26 #include <algorithm>
27 
28 #include "RooLegendre.h"
29 #include "RooAbsReal.h"
30 #include "Math/SpecFunc.h"
31 #include "TMath.h"
32 
33 #include "TError.h"
34 
35 using namespace std;
36 
38 ;
39 
40 
41 ////////////////////////////////////////////////////////////////////////////////
42 
43 namespace {
44  inline double a(int p, int l, int m) {
46  r /= pow(2.,m+2*p);
47  return p%2==0 ? r : -r ;
48  }
49 }
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 
54  _l1(1),_m1(1),_l2(0),_m2(0)
55 {
56 }
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 ///TODO: for now, we assume that ctheta has a range [-1,1]
60 /// should map the ctheta range onto this interval, and adjust integrals...
61 
62 RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l, int m)
63  : RooAbsReal(name, title)
64  , _ctheta("ctheta", "ctheta", this, ctheta)
65  , _l1(l),_m1(m),_l2(0),_m2(0)
66 {
67  //TODO: we assume m>=0
68  // should map m<0 back to m>=0...
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 
73 RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l1, int m1, int l2, int m2)
74  : RooAbsReal(name, title)
75  , _ctheta("ctheta", "ctheta", this, ctheta)
76  , _l1(l1),_m1(m1),_l2(l2),_m2(m2)
77 {
78 }
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 
82 RooLegendre::RooLegendre(const RooLegendre& other, const char* name)
83  : RooAbsReal(other, name)
84  , _ctheta("ctheta", this, other._ctheta)
85  , _l1(other._l1), _m1(other._m1)
86  , _l2(other._l2), _m2(other._m2)
87 {
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// TODO: check that 0<=m_i<=l_i; on the other hand, assoc_legendre already does that ;-)
92 /// Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0
93 
95 {
96 #ifdef R__HAS_MATHMORE
97  double r = 1;
98  double ctheta = std::max(-1., std::min((double)_ctheta, +1.));
99  if (_l1!=0||_m1!=0) r *= ROOT::Math::assoc_legendre(_l1,_m1,ctheta);
100  if (_l2!=0||_m2!=0) r *= ROOT::Math::assoc_legendre(_l2,_m2,ctheta);
101  if ((_m1+_m2)%2==1) r = -r;
102  return r;
103 #else
104  throw std::string("RooLegendre: ERROR: This class require installation of the MathMore library") ;
105  return 0 ;
106 #endif
107 }
108 
109 ////////////////////////////////////////////////////////////////////////////////
110 
111 namespace {
112  Bool_t fullRange(const RooRealProxy& x ,const char* range)
113  {
114  return range == 0 || strlen(range) == 0
115  ? std::fabs(x.min() + 1.) < 1.e-8 && std::fabs(x.max() - 1.) < 1.e-8
116  : std::fabs(x.min(range) + 1.) < 1.e-8 && std::fabs(x.max(range) - 1.) < 1.e-8;
117  }
118 }
119 
120 //_____________________________________________________________________________
121 Int_t RooLegendre::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
122 {
123  // don't support indefinite integrals...
124  if (fullRange(_ctheta,rangeName) && matchArgs(allVars, analVars, _ctheta)) return 1;
125  return 0;
126 }
127 
128 ////////////////////////////////////////////////////////////////////////////////
129 /// this was verified to match mathematica for
130 /// l1 in [0,2], m1 in [0,l1], l2 in [l1,4], m2 in [0,l2]
131 
133 {
134  R__ASSERT(code==1) ;
135  if ( _m1==_m2 ) return ( _l1 == _l2) ? TMath::Factorial(_l1+_m2)/TMath::Factorial(_l1-_m1)*double(2)/(2*_l1+1) : 0.;
136  if ( (_l1+_l2-_m1-_m2)%2 != 0 ) return 0; // these combinations are odd under x -> -x
137 
138  // from B.R. Wong, "On the overlap integral of associated Legendre Polynomials" 1998 J. Phys. A: Math. Gen. 31 1101
139  // TODO: update to the result of
140  // H. A. Mavromatis
141  // "A single-sum expression for the overlap integral of two associated Legendre polynomials"
142  // 1999 J. Phys. A: Math. Gen. 32 2601
143  // http://iopscience.iop.org/0305-4470/32/13/011/pdf/0305-4470_32_13_011.pdf
144  // For that we need Wigner 3-j, which Lorenzo has added for Root 5.28... (note: check Condon-Shortly convention in this paper!)
145  double r=0;
146  for (int p1=0; 2*p1 <= _l1-_m1 ;++p1) {
147  double a1 = a(p1,_l1,_m1);
148  for (int p2=0; 2*p2 <= _l2-_m2 ; ++p2) {
149  double a2 = a(p2,_l2,_m2);
150  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 );
151  }
152  }
153  r /= TMath::Gamma( double(_l1+_l2+3)/2 );
154 
155  if ((_m1+_m2)%2==1) r = -r;
156  return r;
157 }
158 
159 Int_t RooLegendre::getMaxVal( const RooArgSet& /*vars*/) const {
160  if (_m1==0&&_m2==0) return 1;
161  // does anyone know the analytical expression for the max values in case m!=0??
162  if (_l1<3&&_l2<3) return 1;
163  return 0;
164 }
165 
166 namespace {
167  inline double maxSingle(int i, int j) {
168  R__ASSERT(j<=i);
169  // x0 : 1 (ordinary Legendre)
170  if (j==0) return 1;
171  R__ASSERT(i<3);
172  // 11: 1
173  if (i<2) return 1;
174  // 21: 3 22: 3
175  static const double m2[3] = { 3,3 };
176  return m2[j-1];
177  }
178 }
180  return maxSingle(_l1,_m1)*maxSingle(_l2,_m2);
181 }
ClassImp(RooLegendre)
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
#define R__ASSERT(e)
Definition: TError.h:98
RooRealProxy _ctheta
Definition: RooLegendre.h:40
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
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]
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition: TMath.cxx:352
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
Double_t x[n]
Definition: legend1.C:17
TLine l1(2.5, 4.5, 15.5, 4.5)
static double p2(double t, double a, double b, double c)
double pow(double, double)
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
ROOT::R::TRInterface & r
Definition: Object.C:4
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
TMarker * m
Definition: textangle.C:8
double assoc_legendre(unsigned l, unsigned m, double x)
Computes the associated Legendre polynomials.
TLine * l
Definition: textangle.C:4
static double p1(double t, double a, double b)
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
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
#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:94
Double_t Factorial(Int_t i)
Compute factorial(n).
Definition: TMath.cxx:250
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57