Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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/** \class RooLegendre
16 \ingroup Roofit
17
18 Compute the associated Legendre polynomials using ROOT::Math::assoc_legendre().
19
20 Since the Legendre polynomials have a value range of [-1, 1], these cannot be implemented as a PDF.
21 They can be used in sums, though, for example using a RooRealSumFunc of RooLegendre plus an offset.
22**/
23
24#include "RooLegendre.h"
25#include "RooBatchCompute.h"
26#include "RooAbsReal.h"
27
28#include "Math/SpecFunc.h"
29#include "TMath.h"
30
31#include <cmath>
32#include <string>
33#include <algorithm>
34
35using namespace std;
36
38
39////////////////////////////////////////////////////////////////////////////////
40
41namespace {
42 inline double a(int p, int l, int m) {
44 r /= pow(2.,m+2*p);
45 return p%2==0 ? r : -r ;
46 }
47
48 void throwIfNoMathMore() {
49#ifndef R__HAS_MATHMORE
50 throw std::runtime_error("RooLegendre needs functions from MathMore. It is not available in this root build.");
51#endif
52 }
53
54 void checkCoeffs(int m1, int l1, int m2, int l2) {
55 if (m1 < 0 || m2 < 0) {
56 throw std::invalid_argument("RooLegendre: m coefficients need to be >= 0.");
57 }
58 if (l1 < m1 || l2 < m2) {
59 throw std::invalid_argument("RooLegendre: m coefficients need to be smaller than corresponding l.");
60 }
61 }
62}
63
64////////////////////////////////////////////////////////////////////////////////
65
67 _l1(1),_m1(1),_l2(0),_m2(0)
68{
69 throwIfNoMathMore();
70}
71
72////////////////////////////////////////////////////////////////////////////////
73///TODO: for now, we assume that ctheta has a range [-1,1]
74/// should map the ctheta range onto this interval, and adjust integrals...
75
76RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l, int m)
77 : RooAbsReal(name, title)
78 , _ctheta("ctheta", "ctheta", this, ctheta)
79 , _l1(l),_m1(m),_l2(0),_m2(0)
80{
81 checkCoeffs(_m1, _l1, _m2, _l2);
82
83 throwIfNoMathMore();
84}
85
86////////////////////////////////////////////////////////////////////////////////
87
88RooLegendre::RooLegendre(const char* name, const char* title, RooAbsReal& ctheta, int l1, int m1, int l2, int m2)
89 : RooAbsReal(name, title)
90 , _ctheta("ctheta", "ctheta", this, ctheta)
91 , _l1(l1),_m1(m1),_l2(l2),_m2(m2)
92{
93 checkCoeffs(_m1, _l1, _m2, _l2);
94
95 throwIfNoMathMore();
96}
97
98////////////////////////////////////////////////////////////////////////////////
99
100RooLegendre::RooLegendre(const RooLegendre& other, const char* name)
101 : RooAbsReal(other, name)
102 , _ctheta("ctheta", this, other._ctheta)
103 , _l1(other._l1), _m1(other._m1)
104 , _l2(other._l2), _m2(other._m2)
105{
106}
107
108////////////////////////////////////////////////////////////////////////////////
109/// Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0
110
112{
113#ifdef R__HAS_MATHMORE
114 double r = 1;
115 double ctheta = std::max(-1., std::min((double)_ctheta, +1.));
116 if (_l1!=0||_m1!=0) r *= ROOT::Math::assoc_legendre(_l1,_m1,ctheta);
117 if (_l2!=0||_m2!=0) r *= ROOT::Math::assoc_legendre(_l2,_m2,ctheta);
118 if ((_m1+_m2)%2==1) r = -r;
119 return r;
120#else
121 throwIfNoMathMore();
122 return 0.;
123#endif
124}
125
126
127////////////////////////////////////////////////////////////////////////////////
128
129namespace {
130//Author: Emmanouil Michalainas, CERN 26 August 2019
131
132void compute( size_t batchSize, const int l1, const int m1, const int l2, const int m2,
133 double * __restrict output,
134 double const * __restrict TH)
135{
136#ifdef R__HAS_MATHMORE
137 double legendre1=1.0, legendreMinus1=1.0;
138 if (l1+m1 > 0) {
139 legendre1 = ROOT::Math::internal::legendre(l1,m1,1.0);
140 legendreMinus1 = ROOT::Math::internal::legendre(l1,m1,-1.0);
141 }
142 if (l2+m2 > 0) {
143 legendre1 *= ROOT::Math::internal::legendre(l2,m2,1.0);
144 legendreMinus1 *= ROOT::Math::internal::legendre(l2,m2,-1.0);
145 }
146
147 for (size_t i=0; i<batchSize; i++) {
148 if (TH[i] <= -1.0) {
149 output[i] = legendreMinus1;
150 } else if (TH[i] >= 1.0) {
151 output[i] = legendre1;
152 }
153 else {
154 output[i] = 1.0;
155 if (l1+m1 > 0) {
156 output[i] *= ROOT::Math::internal::legendre(l1,m1,TH[i]);
157 }
158 if (l2+m2 > 0) {
159 output[i] *= ROOT::Math::internal::legendre(l2,m2,TH[i]);
160 }
161 }
162 }
163
164#else
165 (void) batchSize, (void) l1, (void)m1, (void)l2, (void)m2, (void)output, (void)TH;
166 throwIfNoMathMore();
167#endif
168}
169};
170
172 RooSpan<const double> cthetaData = _ctheta->getValues(evalData, normSet);
173 size_t batchSize = cthetaData.size();
174 auto output = evalData.makeBatch(this, batchSize);
175 compute(batchSize, _l1, _m1, _l2, _m2, output.data(), cthetaData.data());
176 return output;
177}
178
179
180////////////////////////////////////////////////////////////////////////////////
181
182namespace {
183 Bool_t fullRange(const RooRealProxy& x ,const char* range)
184 {
185 return range == 0 || strlen(range) == 0
186 ? std::fabs(x.min() + 1.) < 1.e-8 && std::fabs(x.max() - 1.) < 1.e-8
187 : std::fabs(x.min(range) + 1.) < 1.e-8 && std::fabs(x.max(range) - 1.) < 1.e-8;
188 }
189}
190
191////////////////////////////////////////////////////////////////////////////////
192
193Int_t RooLegendre::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
194{
195 // don't support indefinite integrals...
196 if (fullRange(_ctheta,rangeName) && matchArgs(allVars, analVars, _ctheta)) return 1;
197 return 0;
198}
199
200////////////////////////////////////////////////////////////////////////////////
201/// this was verified to match mathematica for
202/// l1 in [0,2], m1 in [0,l1], l2 in [l1,4], m2 in [0,l2]
203
205{
206 R__ASSERT(code==1) ;
207 if ( _m1==_m2 ) return ( _l1 == _l2) ? TMath::Factorial(_l1+_m2)/TMath::Factorial(_l1-_m1)*double(2)/(2*_l1+1) : 0.;
208 if ( (_l1+_l2-_m1-_m2)%2 != 0 ) return 0; // these combinations are odd under x -> -x
209
210 // from B.R. Wong, "On the overlap integral of associated Legendre Polynomials" 1998 J. Phys. A: Math. Gen. 31 1101
211 // TODO: update to the result of
212 // H. A. Mavromatis
213 // "A single-sum expression for the overlap integral of two associated Legendre polynomials"
214 // 1999 J. Phys. A: Math. Gen. 32 2601
215 // http://iopscience.iop.org/0305-4470/32/13/011/pdf/0305-4470_32_13_011.pdf
216 // For that we need Wigner 3-j, which Lorenzo has added for Root 5.28... (note: check Condon-Shortly convention in this paper!)
217 double r=0;
218 for (int p1=0; 2*p1 <= _l1-_m1 ;++p1) {
219 double a1 = a(p1,_l1,_m1);
220 for (int p2=0; 2*p2 <= _l2-_m2 ; ++p2) {
221 double a2 = a(p2,_l2,_m2);
222 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 );
223 }
224 }
225 r /= TMath::Gamma( double(_l1+_l2+3)/2 );
226
227 if ((_m1+_m2)%2==1) r = -r;
228 return r;
229}
230
231////////////////////////////////////////////////////////////////////////////////
232
233Int_t RooLegendre::getMaxVal( const RooArgSet& /*vars*/) const {
234 if (_m1==0&&_m2==0) return 1;
235 // does anyone know the analytical expression for the max values in case m!=0??
236 if (_l1<3&&_l2<3) return 1;
237 return 0;
238}
239
240namespace {
241 inline double maxSingle(int i, int j) {
242 R__ASSERT(j<=i);
243 // x0 : 1 (ordinary Legendre)
244 if (j==0) return 1;
245 R__ASSERT(i<3);
246 // 11: 1
247 if (i<2) return 1;
248 // 21: 3 22: 3
249 static const double m2[3] = { 3,3 };
250 return m2[j-1];
251 }
252}
254 return maxSingle(_l1,_m1)*maxSingle(_l2,_m2);
255}
double
ROOT::R::TRInterface & r
Definition Object.C:4
#define a(i)
Definition RSha256.hxx:99
#define ClassImp(name)
Definition Rtypes.h:364
#define R__ASSERT(e)
Definition TError.h:120
char name[80]
Definition TGX11.cxx:110
double pow(double, double)
typedef void((*Func_t)())
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this trans...
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Compute the associated Legendre polynomials using ROOT::Math::assoc_legendre().
Definition RooLegendre.h:20
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 evaluate() const
Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0.
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],...
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
RooRealProxy _ctheta
Definition RooLegendre.h:40
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Evaluate this object for a batch/span of data points.
A simple container to hold a batch of data values.
Definition RooSpan.h:34
constexpr std::span< T >::pointer data() const
Definition RooSpan.h:106
constexpr std::span< T >::index_type size() const noexcept
Definition RooSpan.h:121
double assoc_legendre(unsigned l, unsigned m, double x)
Computes the associated Legendre polynomials.
Double_t x[n]
Definition legend1.C:17
double legendre(unsigned l, unsigned m, double x)
Double_t Factorial(Int_t i)
Compute factorial(n).
Definition TMath.cxx:247
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition TMath.cxx:348
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31
RooSpan< double > makeBatch(const RooAbsReal *owner, std::size_t size)
Create a writable batch.
auto * m
Definition textangle.C:8
auto * l
Definition textangle.C:4
static void output(int code)
Definition gifencode.c:226