Logo ROOT   6.16/01
Reference Guide
RooChebychev.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * GR, Gerhard Raven, UC San Diego, Gerhard.Raven@slac.stanford.edu
7 * *
8 * Copyright (c) 2000-2005, Regents of the University of California *
9 * and Stanford University. All rights reserved. *
10 * *
11 * Redistribution and use in source and binary forms, *
12 * with or without modification, are permitted according to the terms *
13 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14 *****************************************************************************/
15
16/** \class RooChebychev
17 \ingroup Roofit
18
19Chebychev polynomial p.d.f. of the first kind
20**/
21
22#include <cmath>
23#include <iostream>
24
25#include "RooFit.h"
26
27#include "Riostream.h"
28
29#include "RooChebychev.h"
30#include "RooAbsReal.h"
31#include "RooRealVar.h"
32#include "RooArgList.h"
33#include "RooNameReg.h"
34
35#include "TError.h"
36
37#if defined(__my_func__)
38#undef __my_func__
39#endif
40#if defined(WIN32)
41#define __my_func__ __FUNCTION__
42#else
43#define __my_func__ __func__
44#endif
45
47
48////////////////////////////////////////////////////////////////////////////////
49
50RooChebychev::RooChebychev() : _refRangeName(0)
51{
52}
53
54////////////////////////////////////////////////////////////////////////////////
55/// Constructor
56
57RooChebychev::RooChebychev(const char* name, const char* title,
58 RooAbsReal& x, const RooArgList& coefList):
59 RooAbsPdf(name, title),
60 _x("x", "Dependent", this, x),
61 _coefList("coefficients","List of coefficients",this),
62 _refRangeName(0)
63{
64 TIterator* coefIter = coefList.createIterator() ;
65 RooAbsArg* coef ;
66 while((coef = (RooAbsArg*)coefIter->Next())) {
67 if (!dynamic_cast<RooAbsReal*>(coef)) {
68 std::cerr << "RooChebychev::ctor(" << GetName() <<
69 ") ERROR: coefficient " << coef->GetName() <<
70 " is not of type RooAbsReal" << std::endl ;
71 R__ASSERT(0) ;
72 }
73 _coefList.add(*coef) ;
74 }
75
76 delete coefIter ;
77}
78
79////////////////////////////////////////////////////////////////////////////////
80
81RooChebychev::RooChebychev(const RooChebychev& other, const char* name) :
82 RooAbsPdf(other, name),
83 _x("x", this, other._x),
84 _coefList("coefList",this,other._coefList),
85 _refRangeName(other._refRangeName)
86{
87}
88
89//inline static double p0(double ,double a) { return a; }
90inline static double p1(double t,double a,double b) { return a*t+b; }
91inline static double p2(double t,double a,double b,double c) { return p1(t,p1(t,a,b),c); }
92inline static double p3(double t,double a,double b,double c,double d) { return p2(t,p1(t,a,b),c,d); }
93//inline static double p4(double t,double a,double b,double c,double d,double e) { return p3(t,p1(t,a,b),c,d,e); }
94
95////////////////////////////////////////////////////////////////////////////////
96
97void RooChebychev::selectNormalizationRange(const char* rangeName, Bool_t force)
98{
99 if (rangeName && (force || !_refRangeName)) {
101 }
102 if (!rangeName) {
103 _refRangeName = 0 ;
104 }
105}
106
107////////////////////////////////////////////////////////////////////////////////
108
110{
112 Double_t x(-1+2*(_x-xmin)/(xmax-xmin));
113 Double_t x2(x*x);
114 Double_t sum(0) ;
115 switch (_coefList.getSize()) {
116 case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x*p3(x2,64,-112,56,-7);
117 case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*p3(x2,32,-48,18,-1);
118 case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x*p2(x2,16,-20,5);
119 case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*p2(x2,8,-8,1);
120 case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x*p1(x2,4,-3);
121 case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*p1(x2,2,-1);
122 case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x;
123 case 0: sum+=1; break;
124 default: std::cerr << "In " << __my_func__ << " (" << __FILE__ << ", line " <<
125 __LINE__ << "): Higher order Chebychev polynomials currently "
126 "unimplemented." << std::endl;
127 R__ASSERT(false);
128 }
129 return sum;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
134Int_t RooChebychev::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /* rangeName */) const
135{
136 if (matchArgs(allVars, analVars, _x)) return 1;
137 return 0;
138}
139
140////////////////////////////////////////////////////////////////////////////////
141
142Double_t RooChebychev::analyticalIntegral(Int_t code, const char* rangeName) const
143{
144 R__ASSERT(1 == code);
145
146 // the full range of the function is mapped to the normalised [-1, 1] range
147
148 const Double_t xminfull(_x.min(_refRangeName?_refRangeName->GetName():0)) ;
149 const Double_t xmaxfull(_x.max(_refRangeName?_refRangeName->GetName():0)) ;
150
151 const Double_t fullRange = xmaxfull - xminfull;
152
153 // define limits of the integration range on a normalised scale
154 Double_t minScaled = -1., maxScaled = +1.;
155
156 minScaled = -1. + 2. * (_x.min(rangeName) - xminfull) / fullRange;
157 maxScaled = +1. - 2. * (xmaxfull - _x.max(rangeName)) / fullRange;
158
159 // return half of the range since the normalised range runs from -1 to 1
160 // which has a range of two
161 double val = 0.5 * fullRange * (evalAnaInt(maxScaled) - evalAnaInt(minScaled));
162 //std::cout << " integral = " << val << std::endl;
163 return val;
164}
165
166////////////////////////////////////////////////////////////////////////////////
167
169{
170 const Double_t x2 = x * x;
171 Double_t sum = 0.;
172 switch (_coefList.getSize()) {
173 case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x2*p3(x2,8.,-112./6.,14.,-7./2.);
174 case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*x*p3(x2,32./7.,-48./5.,6.,-1.);
175 case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x2*p2(x2,16./6.,-5.,2.5);
176 case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*x*p2(x2,8./5.,-8./3.,1.);
177 case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x2*p1(x2,1.,-3./2.);
178 case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*x*p1(x2,2./3.,-1.);
179 case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x2*.5;
180 case 0: sum+=x; break;
181
182 default: std::cerr << "In " << __my_func__ << " (" << __FILE__ << ", line " <<
183 __LINE__ << "): Higher order Chebychev polynomials currently "
184 "unimplemented." << std::endl;
185 R__ASSERT(false);
186 }
187 return sum;
188}
189
190#undef __my_func__
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
#define c(i)
Definition: RSha256.hxx:101
static double p3(double t, double a, double b, double c, double d)
static double p1(double t, double a, double b)
#define __my_func__
static double p2(double t, double a, double b, double c)
static const double x2[5]
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:363
#define R__ASSERT(e)
Definition: TError.h:96
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
Int_t getSize() const
TIterator * createIterator(Bool_t dir=kIterForward) const
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
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:28
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
Double_t evalAnaInt(const Double_t x) const
virtual void selectNormalizationRange(const char *rangeName=0, Bool_t force=kFALSE)
Interface function to force use of a given normalization range to interpret function value.
RooRealProxy _x
Definition: RooChebychev.h:43
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
TNamed * _refRangeName
Definition: RooChebychev.h:45
RooListProxy _coefList
Definition: RooChebychev.h:44
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
const TNamed * constPtr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:91
static RooNameReg & instance()
Return reference to singleton instance.
Definition: RooNameReg.cxx:64
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Double_t x[n]
Definition: legend1.C:17
auto * a
Definition: textangle.C:12
static long int sum(long int i)
Definition: Factory.cxx:2258