Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooLegendre.h
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#ifndef ROO_LEGENDRE
15#define ROO_LEGENDRE
16
17#include "RooAbsReal.h"
18#include "RooRealProxy.h"
19
20class RooLegendre : public RooAbsReal {
21public:
22 RooLegendre() ;
23 // an (associated) Legendre polynomial, P_l^m(x)
24 // note: P_l(x) == P_l^0(x)
25 RooLegendre(const char *name, const char *title, RooAbsReal& ctheta, int l, int m=0);
26 // product of two associated Legendre polynomials, P_l1^m1(ctheta) * P_l2^m2(ctheta)
27 RooLegendre(const char *name, const char *title, RooAbsReal& ctheta, int l1, int m1, int l2, int m2);
28
29 RooLegendre(const RooLegendre& other, const char *name = nullptr);
30 TObject* clone(const char* newname) const override { return new RooLegendre(*this, newname); }
31
32 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
33 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
34
35 Int_t getMaxVal( const RooArgSet& vars) const override;
36 double maxVal( Int_t code) const override;
37
38protected: // allow RooSpHarmonic access...
40 int _l1,_m1;
41 int _l2,_m2;
42
43 double evaluate() const override;
44 void doEval(RooFit::EvalContext &) const override;
45
46 ClassDefOverride(RooLegendre,1) // Legendre polynomial
47};
48
49#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
char name[80]
Definition TGX11.cxx:110
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Compute the associated Legendre polynomials using ROOT::Math::assoc_legendre().
Definition RooLegendre.h:20
double evaluate() const override
Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0.
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise capability to determine maximum value of function for given set of observables.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
this was verified to match mathematica for l1 in [0,2], m1 in [0,l1], l2 in [l1,4],...
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooRealProxy _ctheta
Definition RooLegendre.h:39
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
TObject * clone(const char *newname) const override
Definition RooLegendre.h:30
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
Mother of all ROOT objects.
Definition TObject.h:41
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4