ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 \file RooChebychev.cxx
18 \class RooChebychev
19 \ingroup Roofit
20 
21 Chebychev polynomial p.d.f. of the first kind
22 **/
23 
24 #include <cmath>
25 #include <iostream>
26 
27 #include "RooFit.h"
28 
29 #include "Riostream.h"
30 
31 #include "RooChebychev.h"
32 #include "RooAbsReal.h"
33 #include "RooRealVar.h"
34 #include "RooArgList.h"
35 #include "RooNameReg.h"
36 
37 #include "TError.h"
38 
39 #if defined(__my_func__)
40 #undef __my_func__
41 #endif
42 #if defined(WIN32)
43 #define __my_func__ __FUNCTION__
44 #else
45 #define __my_func__ __func__
46 #endif
47 
49 ;
50 
51 ////////////////////////////////////////////////////////////////////////////////
52 
53 RooChebychev::RooChebychev() : _refRangeName(0)
54 {
55 }
56 
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 /// Constructor
60 
61 RooChebychev::RooChebychev(const char* name, const char* title,
62  RooAbsReal& x, const RooArgList& coefList):
63  RooAbsPdf(name, title),
64  _x("x", "Dependent", this, x),
65  _coefList("coefficients","List of coefficients",this),
66  _refRangeName(0)
67 {
68  TIterator* coefIter = coefList.createIterator() ;
69  RooAbsArg* coef ;
70  while((coef = (RooAbsArg*)coefIter->Next())) {
71  if (!dynamic_cast<RooAbsReal*>(coef)) {
72  std::cerr << "RooChebychev::ctor(" << GetName() <<
73  ") ERROR: coefficient " << coef->GetName() <<
74  " is not of type RooAbsReal" << std::endl ;
75  R__ASSERT(0) ;
76  }
77  _coefList.add(*coef) ;
78  }
79 
80  delete coefIter ;
81 }
82 
83 
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 
87 RooChebychev::RooChebychev(const RooChebychev& other, const char* name) :
88  RooAbsPdf(other, name),
89  _x("x", this, other._x),
90  _coefList("coefList",this,other._coefList),
91  _refRangeName(other._refRangeName)
92 {
93 }
94 
95 //inline static double p0(double ,double a) { return a; }
96 inline static double p1(double t,double a,double b) { return a*t+b; }
97 inline static double p2(double t,double a,double b,double c) { return p1(t,p1(t,a,b),c); }
98 inline static double p3(double t,double a,double b,double c,double d) { return p2(t,p1(t,a,b),c,d); }
99 //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); }
100 
101 
102 ////////////////////////////////////////////////////////////////////////////////
103 
104 void RooChebychev::selectNormalizationRange(const char* rangeName, Bool_t force)
105 {
106  if (rangeName && (force || !_refRangeName)) {
108  }
109  if (!rangeName) {
110  _refRangeName = 0 ;
111  }
112 }
113 
114 
115 ////////////////////////////////////////////////////////////////////////////////
116 
118 {
120  Double_t x(-1+2*(_x-xmin)/(xmax-xmin));
121  Double_t x2(x*x);
122  Double_t sum(0) ;
123  switch (_coefList.getSize()) {
124  case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x*p3(x2,64,-112,56,-7);
125  case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*p3(x2,32,-48,18,-1);
126  case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x*p2(x2,16,-20,5);
127  case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*p2(x2,8,-8,1);
128  case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x*p1(x2,4,-3);
129  case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*p1(x2,2,-1);
130  case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x;
131  case 0: sum+=1; break;
132  default: std::cerr << "In " << __my_func__ << " (" << __FILE__ << ", line " <<
133  __LINE__ << "): Higher order Chebychev polynomials currently "
134  "unimplemented." << std::endl;
135  R__ASSERT(false);
136  }
137  return sum;
138 }
139 
140 
141 ////////////////////////////////////////////////////////////////////////////////
142 
143 Int_t RooChebychev::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /* rangeName */) const
144 {
145  if (matchArgs(allVars, analVars, _x)) return 1;
146  return 0;
147 }
148 
149 
150 ////////////////////////////////////////////////////////////////////////////////
151 
152 Double_t RooChebychev::analyticalIntegral(Int_t code, const char* rangeName) const
153 {
154  R__ASSERT(1 == code);
155 
156  // the full range of the function is mapped to the normalised [-1, 1] range
157 
158  const Double_t xminfull(_x.min(_refRangeName?_refRangeName->GetName():0)) ;
159  const Double_t xmaxfull(_x.max(_refRangeName?_refRangeName->GetName():0)) ;
160 
161  const Double_t fullRange = xmaxfull - xminfull;
162 
163  // define limits of the integration range on a normalised scale
164  Double_t minScaled = -1., maxScaled = +1.;
165 
166  minScaled = -1. + 2. * (_x.min(rangeName) - xminfull) / fullRange;
167  maxScaled = +1. - 2. * (xmaxfull - _x.max(rangeName)) / fullRange;
168 
169  // return half of the range since the normalised range runs from -1 to 1
170  // which has a range of two
171  double val = 0.5 * fullRange * (evalAnaInt(maxScaled) - evalAnaInt(minScaled));
172  //std::cout << " integral = " << val << std::endl;
173  return val;
174 }
175 
177 {
178  const Double_t x2 = x * x;
179  Double_t sum = 0.;
180  switch (_coefList.getSize()) {
181  case 7: sum+=((RooAbsReal&)_coefList[6]).getVal()*x2*p3(x2,8.,-112./6.,14.,-7./2.);
182  case 6: sum+=((RooAbsReal&)_coefList[5]).getVal()*x*p3(x2,32./7.,-48./5.,6.,-1.);
183  case 5: sum+=((RooAbsReal&)_coefList[4]).getVal()*x2*p2(x2,16./6.,-5.,2.5);
184  case 4: sum+=((RooAbsReal&)_coefList[3]).getVal()*x*p2(x2,8./5.,-8./3.,1.);
185  case 3: sum+=((RooAbsReal&)_coefList[2]).getVal()*x2*p1(x2,1.,-3./2.);
186  case 2: sum+=((RooAbsReal&)_coefList[1]).getVal()*x*p1(x2,2./3.,-1.);
187  case 1: sum+=((RooAbsReal&)_coefList[0]).getVal()*x2*.5;
188  case 0: sum+=x; break;
189 
190  default: std::cerr << "In " << __my_func__ << " (" << __FILE__ << ", line " <<
191  __LINE__ << "): Higher order Chebychev polynomials currently "
192  "unimplemented." << std::endl;
193  R__ASSERT(false);
194  }
195  return sum;
196 }
197 
198 #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)
return c
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
int d
Definition: tornado.py:11
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:64
TThread * t[5]
Definition: threadsh1.C:13
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
#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:91
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Int_t getSize() const
Chebychev polynomial p.d.f.
Definition: RooChebychev.h:25
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