Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooPolynomial.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$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/** \class RooPolynomial
18 \ingroup Roofit
19
20RooPolynomial implements a polynomial p.d.f of the form
21\f[ f(x) = \mathcal{N} \cdot \sum_{i} a_{i} * x^i \f]
22By default, the coefficient \f$ a_0 \f$ is chosen to be 1, as polynomial
23probability density functions have one degree of freedom
24less than polynomial functions due to the normalisation condition. \f$ \mathcal{N} \f$
25is a normalisation constant that is automatically calculated when the polynomial is used
26in computations.
27
28The sum can be truncated at the low end. See the main constructor
29RooPolynomial::RooPolynomial(const char*, const char*, RooAbsReal&, const RooArgList&, Int_t)
30**/
31
32#include "RooPolynomial.h"
33#include "RooArgList.h"
34#include "RooMsgService.h"
35#include "RooPolyVar.h"
36
37#include "TError.h"
38#include <vector>
39
41
42////////////////////////////////////////////////////////////////////////////////
43/// Create a polynomial in the variable `x`.
44/// \param[in] name Name of the PDF
45/// \param[in] title Title for plotting the PDF
46/// \param[in] x The variable of the polynomial
47/// \param[in] coefList The coefficients \f$ a_i \f$
48/// \param[in] lowestOrder [optional] Truncate the sum such that it skips the lower orders:
49/// \f[
50/// 1. + \sum_{i=0}^{\mathrm{coefList.size()}} a_{i} * x^{(i + \mathrm{lowestOrder})}
51/// \f]
52///
53/// This means that
54/// \code{.cpp}
55/// RooPolynomial pol("pol", "pol", x, RooArgList(a, b), lowestOrder = 2)
56/// \endcode
57/// computes
58/// \f[
59/// \mathrm{pol}(x) = 1 * x^0 + (0 * x^{\ldots}) + a * x^2 + b * x^3.
60/// \f]
61
62
63RooPolynomial::RooPolynomial(const char* name, const char* title,
64 RooAbsReal& x, const RooArgList& coefList, Int_t lowestOrder) :
65 RooAbsPdf(name, title),
66 _x("x", "Dependent", this, x),
67 _coefList("coefList","List of coefficients",this),
68 _lowestOrder(lowestOrder)
69{
70 // Check lowest order
71 if (_lowestOrder<0) {
72 coutE(InputArguments) << "RooPolynomial::ctor(" << GetName()
73 << ") WARNING: lowestOrder must be >=0, setting value to 0" << std::endl;
74 _lowestOrder=0 ;
75 }
76
77 for (auto *coef : coefList) {
78 if (!dynamic_cast<RooAbsReal*>(coef)) {
79 coutE(InputArguments) << "RooPolynomial::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName()
80 << " is not of type RooAbsReal" << std::endl;
81 R__ASSERT(0) ;
82 }
83 _coefList.add(*coef) ;
84 }
85}
86
87////////////////////////////////////////////////////////////////////////////////
88
89RooPolynomial::RooPolynomial(const char* name, const char* title,
90 RooAbsReal& x) :
91 RooAbsPdf(name, title),
92 _x("x", "Dependent", this, x),
93 _coefList("coefList","List of coefficients",this),
94 _lowestOrder(1)
95{ }
96
97////////////////////////////////////////////////////////////////////////////////
98/// Copy constructor
99
101 RooAbsPdf(other, name),
102 _x("x", this, other._x),
103 _coefList("coefList",this,other._coefList),
104 _lowestOrder(other._lowestOrder)
105{ }
106
107////////////////////////////////////////////////////////////////////////////////
108
110{
111 // Calculate and return value of polynomial
112
113 const unsigned sz = _coefList.getSize();
114 const int lowestOrder = _lowestOrder;
115 if (!sz) return lowestOrder ? 1. : 0.;
116 _wksp.clear();
117 _wksp.reserve(sz);
118 {
119 const RooArgSet* nset = _coefList.nset();
120 for (auto *c : static_range_cast<RooAbsReal *>(_coefList)) {
121 _wksp.push_back(c->getVal(nset));
122 }
123 }
124 const double x = _x;
125 double retVal = _wksp[sz - 1];
126 for (unsigned i = sz - 1; i--; ) retVal = _wksp[i] + x * retVal;
127 return retVal * std::pow(x, lowestOrder) + (lowestOrder ? 1.0 : 0.0);
128}
129
130/// Compute multiple values of Polynomial.
131void RooPolynomial::computeBatch(cudaStream_t *stream, double *output, size_t nEvents,
132 RooFit::Detail::DataMap const &dataMap) const
133{
134 return RooPolyVar::computeBatchImpl(stream, output, nEvents, dataMap, _x.arg(), _coefList, _lowestOrder);
135}
136
137////////////////////////////////////////////////////////////////////////////////
138/// Advertise to RooFit that this function can be analytically integrated.
139Int_t RooPolynomial::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
140{
141 if (matchArgs(allVars, analVars, _x)) return 1;
142 return 0;
143}
144
145////////////////////////////////////////////////////////////////////////////////
146/// Do the analytical integral according to the code that was returned by getAnalyticalIntegral().
147double RooPolynomial::analyticalIntegral(Int_t code, const char* rangeName) const
148{
149 R__ASSERT(code==1) ;
150
151 const double xmin = _x.min(rangeName), xmax = _x.max(rangeName);
152 const int lowestOrder = _lowestOrder;
153 const unsigned sz = _coefList.getSize();
154 if (!sz) return xmax - xmin;
155 _wksp.clear();
156 _wksp.reserve(sz);
157 {
158 const RooArgSet* nset = _coefList.nset();
159 unsigned i = 1 + lowestOrder;
160 for (auto *c : static_range_cast<RooAbsReal *>(_coefList)) {
161 _wksp.push_back(c->getVal(nset) / double(i));
162 ++i;
163 }
164 }
165 double min = _wksp[sz - 1], max = _wksp[sz - 1];
166 for (unsigned i = sz - 1; i--; )
167 min = _wksp[i] + xmin * min, max = _wksp[i] + xmax * max;
168 return max * std::pow(xmax, 1 + lowestOrder) - min * std::pow(xmin, 1 + lowestOrder) +
169 (lowestOrder ? (xmax - xmin) : 0.);
170}
#define c(i)
Definition RSha256.hxx:101
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:117
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Int_t getSize() const
Return the number of elements in the collection.
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
static void computeBatchImpl(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &, RooAbsReal const &x, RooArgList const &coefs, int lowestOrder)
RooPolynomial implements a polynomial p.d.f of the form.
double evaluate() const override
do not persist
RooRealProxy _x
std::vector< double > _wksp
RooAbsReal const & x() const
Get the x variable.
RooListProxy _coefList
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Do the analytical integral according to the code that was returned by getAnalyticalIntegral().
int lowestOrder() const
Return the order for the first coefficient in the list.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Advertise to RooFit that this function can be analytically integrated.
RooArgList const & coefList() const
Get the coefficient list.
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Polynomial.
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
const T & arg() const
Return reference to object held in proxy.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
static void output()