ROOT  6.06/09
Reference Guide
RooPolyVar.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofitcore:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // Class RooPolyVar is a RooAbsReal implementing a polynomial in terms
21 // of a list of RooAbsReal coefficients
22 // <pre>
23 // f(x) = sum_i a_i * x
24 // </pre>
25 // Class RooPolyvar implements analytical integrals of all polynomials
26 // it can define.
27 // END_HTML
28 //
29 
30 #include <cmath>
31 
32 #include "RooPolyVar.h"
33 #include "RooArgList.h"
34 #include "RooMsgService.h"
35 //#include "Riostream.h"
36 
37 #include "TError.h"
38 
39 using namespace std;
40 
42 ;
43 
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 /// Default constructor
47 
48 RooPolyVar::RooPolyVar() : _lowestOrder(0)
49 { }
50 
51 
52 ////////////////////////////////////////////////////////////////////////////////
53 /// Construct polynomial in x with coefficients in coefList. If
54 /// lowestOrder is not zero, then the first element in coefList is
55 /// interpreted as as the 'lowestOrder' coefficients and all
56 /// subsequent coeffient elements are shifted by a similar amount.
57 RooPolyVar::RooPolyVar(const char* name, const char* title,
58  RooAbsReal& x, const RooArgList& coefList, Int_t lowestOrder) :
59  RooAbsReal(name, title),
60  _x("x", "Dependent", this, x),
61  _coefList("coefList","List of coefficients",this),
62  _lowestOrder(lowestOrder)
63 {
64  // Check lowest order
65  if (_lowestOrder<0) {
66  coutE(InputArguments) << "RooPolyVar::ctor(" << GetName()
67  << ") WARNING: lowestOrder must be >=0, setting value to 0" << endl ;
68  _lowestOrder=0 ;
69  }
70 
71  RooFIter coefIter = coefList.fwdIterator() ;
72  RooAbsArg* coef ;
73  while((coef = (RooAbsArg*)coefIter.next())) {
74  if (!dynamic_cast<RooAbsReal*>(coef)) {
75  coutE(InputArguments) << "RooPolyVar::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName()
76  << " is not of type RooAbsReal" << endl ;
77  R__ASSERT(0) ;
78  }
79  _coefList.add(*coef) ;
80  }
81 }
82 
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 /// Constructor of flat polynomial function
86 
87 RooPolyVar::RooPolyVar(const char* name, const char* title,
88  RooAbsReal& x) :
89  RooAbsReal(name, title),
90  _x("x", "Dependent", this, x),
91  _coefList("coefList","List of coefficients",this),
92  _lowestOrder(1)
93 { }
94 
95 
96 
97 ////////////////////////////////////////////////////////////////////////////////
98 /// Copy constructor
99 
100 RooPolyVar::RooPolyVar(const RooPolyVar& other, const char* name) :
101  RooAbsReal(other, name),
102  _x("x", this, other._x),
103  _coefList("coefList",this,other._coefList),
104  _lowestOrder(other._lowestOrder)
105 { }
106 
107 
108 
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// Destructor
112 
114 { }
115 
116 
117 
118 
119 ////////////////////////////////////////////////////////////////////////////////
120 /// Calculate and return value of polynomial
121 
123 {
124  const unsigned sz = _coefList.getSize();
125  const int lowestOrder = _lowestOrder;
126  if (!sz) return lowestOrder ? 1. : 0.;
127  _wksp.clear();
128  _wksp.reserve(sz);
129  {
130  const RooArgSet* nset = _coefList.nset();
132  RooAbsReal* c;
133  while ((c = (RooAbsReal*) it.next())) _wksp.push_back(c->getVal(nset));
134  }
135  const Double_t x = _x;
136  Double_t retVal = _wksp[sz - 1];
137  for (unsigned i = sz - 1; i--; ) retVal = _wksp[i] + x * retVal;
138  return retVal * std::pow(x, lowestOrder);
139 }
140 
141 
142 
143 ////////////////////////////////////////////////////////////////////////////////
144 /// Advertise that we can internally integrate over x
145 
146 Int_t RooPolyVar::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
147 {
148  if (matchArgs(allVars, analVars, _x)) return 1;
149  return 0;
150 }
151 
152 
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 /// Calculate and return analytical integral over x
156 
157 Double_t RooPolyVar::analyticalIntegral(Int_t code, const char* rangeName) const
158 {
159  R__ASSERT(code==1) ;
160 
161  const Double_t xmin = _x.min(rangeName), xmax = _x.max(rangeName);
162  const int lowestOrder = _lowestOrder;
163  const unsigned sz = _coefList.getSize();
164  if (!sz) return xmax - xmin;
165  _wksp.clear();
166  _wksp.reserve(sz);
167  {
168  const RooArgSet* nset = _coefList.nset();
170  unsigned i = 1 + lowestOrder;
171  RooAbsReal* c;
172  while ((c = (RooAbsReal*) it.next())) {
173  _wksp.push_back(c->getVal(nset) / Double_t(i));
174  ++i;
175  }
176  }
177  Double_t min = _wksp[sz - 1], max = _wksp[sz - 1];
178  for (unsigned i = sz - 1; i--; )
179  min = _wksp[i] + xmin * min, max = _wksp[i] + xmax * max;
180  return max * std::pow(xmax, 1 + lowestOrder) - min * std::pow(xmin, 1 + lowestOrder);
181 }
const RooArgSet * nset() const
Definition: RooAbsProxy.h:47
#define coutE(a)
Definition: RooMsgService.h:35
float xmin
Definition: THbookFile.cxx:93
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
virtual ~RooPolyVar()
Destructor.
Definition: RooPolyVar.cxx:113
RooFIter fwdIterator() const
Double_t evaluate() const
do not persist
Definition: RooPolyVar.cxx:122
#define R__ASSERT(e)
Definition: TError.h:98
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Advertise that we can internally integrate over x.
Definition: RooPolyVar.cxx:146
Double_t x[n]
Definition: legend1.C:17
double pow(double, double)
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Calculate and return analytical integral over x.
Definition: RooPolyVar.cxx:157
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
float xmax
Definition: THbookFile.cxx:93
RooAbsArg * next()
Int_t _lowestOrder
Definition: RooPolyVar.h:47
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
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
#define name(a, b)
Definition: linkTestLib0.cpp:5
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
RooListProxy _coefList
Definition: RooPolyVar.h:46
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
std::vector< Double_t > _wksp
Definition: RooPolyVar.h:49
RooRealProxy _x
Definition: RooPolyVar.h:45
RooPolyVar()
Default constructor.
Definition: RooPolyVar.cxx:48
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
ClassImp(RooPolyVar)