ROOT  6.06/09
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 //////////////////////////////////////////////////////////////////////////////
17 //
18 // BEGIN_HTML
19 // Chebychev polynomial p.d.f. of the first kind
20 // END_HTML
21 //
22 
23 #include <cmath>
24 #include <iostream>
25 
26 #include "RooFit.h"
27 
28 #include "Riostream.h"
29 
30 #include "RooChebychev.h"
31 #include "RooAbsReal.h"
32 #include "RooRealVar.h"
33 #include "RooArgList.h"
34 #include "RooNameReg.h"
35 
36 #include "TError.h"
37 
38 #if defined(__my_func__)
39 #undef __my_func__
40 #endif
41 #if defined(WIN32)
42 #define __my_func__ __FUNCTION__
43 #else
44 #define __my_func__ __func__
45 #endif
46 
48 ;
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 
52 RooChebychev::RooChebychev() : _refRangeName(0)
53 {
54 }
55 
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 /// Constructor
59 
60 RooChebychev::RooChebychev(const char* name, const char* title,
61  RooAbsReal& x, const RooArgList& coefList):
62  RooAbsPdf(name, title),
63  _x("x", "Dependent", this, x),
64  _coefList("coefficients","List of coefficients",this),
65  _refRangeName(0)
66 {
67  TIterator* coefIter = coefList.createIterator() ;
68  RooAbsArg* coef ;
69  while((coef = (RooAbsArg*)coefIter->Next())) {
70  if (!dynamic_cast<RooAbsReal*>(coef)) {
71  std::cerr << "RooChebychev::ctor(" << GetName() <<
72  ") ERROR: coefficient " << coef->GetName() <<
73  " is not of type RooAbsReal" << std::endl ;
74  R__ASSERT(0) ;
75  }
76  _coefList.add(*coef) ;
77  }
78 
79  delete coefIter ;
80 }
81 
82 
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 
86 RooChebychev::RooChebychev(const RooChebychev& other, const char* name) :
87  RooAbsPdf(other, name),
88  _x("x", this, other._x),
89  _coefList("coefList",this,other._coefList),
90  _refRangeName(other._refRangeName)
91 {
92 }
93 
94 //inline static double p0(double ,double a) { return a; }
95 inline static double p1(double t,double a,double b) { return a*t+b; }
96 inline static double p2(double t,double a,double b,double c) { return p1(t,p1(t,a,b),c); }
97 inline static double p3(double t,double a,double b,double c,double d) { return p2(t,p1(t,a,b),c,d); }
98 //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); }
99 
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 
103 void RooChebychev::selectNormalizationRange(const char* rangeName, Bool_t force)
104 {
105  if (rangeName && (force || !_refRangeName)) {
107  }
108  if (!rangeName) {
109  _refRangeName = 0 ;
110  }
111 }
112 
113 
114 ////////////////////////////////////////////////////////////////////////////////
115 
117 {
119  Double_t x(-1+2*(_x-xmin)/(xmax-xmin));
120  Double_t x2(x*x);
121  Double_t sum(0) ;
122  switch (_coefList.getSize()) {
123  case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x*p3(x2,64,-112,56,-7);
124  case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*p3(x2,32,-48,18,-1);
125  case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x*p2(x2,16,-20,5);
126  case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*p2(x2,8,-8,1);
127  case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x*p1(x2,4,-3);
128  case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*p1(x2,2,-1);
129  case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x;
130  case 0: sum+=1; break;
131  default: std::cerr << "In " << __my_func__ << " (" << __FILE__ << ", line " <<
132  __LINE__ << "): Higher order Chebychev polynomials currently "
133  "unimplemented." << std::endl;
134  R__ASSERT(false);
135  }
136  return sum;
137 }
138 
139 
140 ////////////////////////////////////////////////////////////////////////////////
141 
142 Int_t RooChebychev::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /* rangeName */) const
143 {
144  if (matchArgs(allVars, analVars, _x)) return 1;
145  return 0;
146 }
147 
148 
149 ////////////////////////////////////////////////////////////////////////////////
150 
151 Double_t RooChebychev::analyticalIntegral(Int_t code, const char* rangeName) const
152 {
153  R__ASSERT(1 == code);
154 
155  // the full range of the function is mapped to the normalised [-1, 1] range
156 
157  const Double_t xminfull(_x.min(_refRangeName?_refRangeName->GetName():0)) ;
158  const Double_t xmaxfull(_x.max(_refRangeName?_refRangeName->GetName():0)) ;
159 
160  const Double_t fullRange = xmaxfull - xminfull;
161 
162  // define limits of the integration range on a normalised scale
163  Double_t minScaled = -1., maxScaled = +1.;
164 
165  minScaled = -1. + 2. * (_x.min(rangeName) - xminfull) / fullRange;
166  maxScaled = +1. - 2. * (xmaxfull - _x.max(rangeName)) / fullRange;
167 
168  // return half of the range since the normalised range runs from -1 to 1
169  // which has a range of two
170  double val = 0.5 * fullRange * (evalAnaInt(maxScaled) - evalAnaInt(minScaled));
171  //std::cout << " integral = " << val << std::endl;
172  return val;
173 }
174 
176 {
177  const Double_t x2 = x * x;
178  Double_t sum = 0.;
179  switch (_coefList.getSize()) {
180  case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x2*p3(x2,8.,-112./6.,14.,-7./2.);
181  case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*x*p3(x2,32./7.,-48./5.,6.,-1.);
182  case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x2*p2(x2,16./6.,-5.,2.5);
183  case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*x*p2(x2,8./5.,-8./3.,1.);
184  case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x2*p1(x2,1.,-3./2.);
185  case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*x*p1(x2,2./3.,-1.);
186  case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x2*.5;
187  case 0: sum+=x; break;
188 
189  default: std::cerr << "In " << __my_func__ << " (" << __FILE__ << ", line " <<
190  __LINE__ << "): Higher order Chebychev polynomials currently "
191  "unimplemented." << std::endl;
192  R__ASSERT(false);
193  }
194  return sum;
195 }
196 
197 #undef __my_func__
float xmin
Definition: THbookFile.cxx:93
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
static double p3(double t, double a, double b, double c, double d)
Double_t evalAnaInt(const Double_t x) const
ClassImp(RooChebychev)
#define R__ASSERT(e)
Definition: TError.h:98
RooRealProxy _x
Definition: RooChebychev.h:43
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
Iterator abstract base class.
Definition: TIterator.h:32
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
static double p2(double t, double a, double b, double c)
TIterator * createIterator(Bool_t dir=kIterForward) const
TNamed * _refRangeName
Definition: RooChebychev.h:45
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
RooListProxy _coefList
Definition: RooChebychev.h:44
static RooNameReg & instance()
Return reference to singleton instance.
Definition: RooNameReg.cxx:63
#define __my_func__
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
float xmax
Definition: THbookFile.cxx:93
static double p1(double t, double a, double b)
Double_t evaluate() const
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
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
virtual TObject * Next()=0
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...
const TNamed * constPtr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
Definition: RooNameReg.cxx:90
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Int_t getSize() const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57