Logo ROOT  
Reference Guide
RooSpHarmonic.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_SPHARMONIC
15#define ROO_SPHARMONIC
16
17#include "RooLegendre.h"
18#include "RooRealProxy.h"
19
20class RooSpHarmonic : public RooLegendre {
21public:
23 RooSpHarmonic(const char *name, const char *title, RooAbsReal& ctheta, RooAbsReal& phi, int l, int m);
24 RooSpHarmonic(const char *name, const char *title, RooAbsReal& ctheta, RooAbsReal& phi, int l1, int m1, int l2, int m2);
25
26 RooSpHarmonic(const RooSpHarmonic& other, const char* name = 0);
27 virtual TObject* clone(const char* newname) const { return new RooSpHarmonic(*this, newname); }
28 inline virtual ~RooSpHarmonic() { }
29
30 virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
31 virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
32
33 virtual Int_t getMaxVal( const RooArgSet& vars) const;
34 virtual Double_t maxVal( Int_t code) const;
35
36private:
38 double _n;
40
41 Double_t evaluate() const;
42
43 ClassDef(RooSpHarmonic,1) // SpHarmonic polynomial
44};
45
46#endif
double Double_t
Definition: RtypesCore.h:57
#define ClassDef(name, id)
Definition: Rtypes.h:322
char name[80]
Definition: TGX11.cxx:109
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Compute the associated Legendre polynomials using ROOT::Math::assoc_legendre().
Definition: RooLegendre.h:20
Implementation of the so-called real spherical harmonics, using the orthonormal normalization,...
Definition: RooSpHarmonic.h:20
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
TODO: check that phi.max - phi.min = 2 pi... ctheta.max = +1, and ctheta.min = -1 we don't support in...
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 TObject * clone(const char *newname) const
Definition: RooSpHarmonic.h:27
RooRealProxy _phi
Definition: RooSpHarmonic.h:37
virtual ~RooSpHarmonic()
Definition: RooSpHarmonic.h:28
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
Note: P_0^0 = 1, so P_l^m = P_l^m P_0^0.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Advertise capability to determine maximum value of function for given set of observables.
Mother of all ROOT objects.
Definition: TObject.h:37
static constexpr double m2
auto * m
Definition: textangle.C:8
auto * l
Definition: textangle.C:4