Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
21The coefficient that goes with \f$ T_0(x)=1 \f$ (i.e. the constant polynomial) is
22implicitly assumed to be 1, and the list of coefficients supplied by callers
23starts with the coefficient that goes with \f$ T_1(x)=x \f$ (i.e. the linear term).
24**/
25
26#include "RooChebychev.h"
27#include "RooRealVar.h"
28#include "RooArgList.h"
29#include "RooNameReg.h"
30#include "RooBatchCompute.h"
31
34
35#include <cmath>
36
38
39////////////////////////////////////////////////////////////////////////////////
40
42
43////////////////////////////////////////////////////////////////////////////////
44/// Constructor
45
46RooChebychev::RooChebychev(const char* name, const char* title,
47 RooAbsReal& x, const RooArgList& coefList):
48 RooAbsPdf(name, title),
49 _x("x", "Dependent", this, x),
50 _coefList("coefficients","List of coefficients",this)
51{
52 _coefList.addTyped<RooAbsReal>(coefList);
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
57RooChebychev::RooChebychev(const RooChebychev& other, const char* name) :
58 RooAbsPdf(other, name),
59 _x("x", this, other._x),
60 _coefList("coefList",this,other._coefList),
61 _refRangeName(other._refRangeName)
62{
63}
64
65////////////////////////////////////////////////////////////////////////////////
66
67void RooChebychev::selectNormalizationRange(const char* rangeName, bool force)
68{
69 if (rangeName && (force || !_refRangeName)) {
70 _refRangeName = const_cast<TNamed*>(RooNameReg::instance().constPtr(rangeName));
71 }
72 if (!rangeName) {
73 _refRangeName = nullptr ;
74 }
75}
76
77////////////////////////////////////////////////////////////////////////////////
78
80{
81 // first bring the range of the variable _x to the normalised range [-1, 1]
82 // calculate sum_k c_k T_k(x) where x is given in the normalised range,
83 // c_0 = 1, and the higher coefficients are given in _coefList
84 double xmax = _x.max(_refRangeName ? _refRangeName->GetName() : nullptr);
85 double xmin = _x.min(_refRangeName ? _refRangeName->GetName() : nullptr);
86
87 std::vector<double> coeffs;
88 for (auto it : _coefList)
89 coeffs.push_back(static_cast<const RooAbsReal &>(*it).getVal());
91}
92
94{
95 // first bring the range of the variable _x to the normalised range [-1, 1]
96 // calculate sum_k c_k T_k(x) where x is given in the normalised range,
97 // c_0 = 1, and the higher coefficients are given in _coefList
98 double xmax = _x.max(_refRangeName ? _refRangeName->GetName() : nullptr);
99 double xmin = _x.min(_refRangeName ? _refRangeName->GetName() : nullptr);
100
101 ctx.addResult(this,
102 ctx.buildCall("RooFit::Detail::EvaluateFuncs::chebychevEvaluate", _coefList, _coefList.size(), _x, xmin, xmax));
103}
104
105////////////////////////////////////////////////////////////////////////////////
106/// Compute multiple values of Chebychev.
107void RooChebychev::computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &dataMap) const
108{
109 std::vector<double> extraArgs;
110 extraArgs.reserve(_coefList.size() + 2);
111 for (auto *coef : _coefList) {
112 extraArgs.push_back(static_cast<const RooAbsReal *>(coef)->getVal());
113 }
114 extraArgs.push_back(_x.min(_refRangeName ? _refRangeName->GetName() : nullptr));
115 extraArgs.push_back(_x.max(_refRangeName ? _refRangeName->GetName() : nullptr));
116 RooBatchCompute::compute(dataMap.config(this), RooBatchCompute::Chebychev, output, nEvents, {dataMap.at(_x)},
117 extraArgs);
118}
119
120////////////////////////////////////////////////////////////////////////////////
121
122
123Int_t RooChebychev::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /* rangeName */) const
124{
125 if (matchArgs(allVars, analVars, _x)) return 1;
126 return 0;
127}
128
129////////////////////////////////////////////////////////////////////////////////
130
131double RooChebychev::analyticalIntegral(Int_t code, const char* rangeName) const
132{
133 assert(1 == code); (void)code;
134
135 double xmax = _x.max(_refRangeName ? _refRangeName->GetName() : nullptr);
136 double xmaxFull = _x.max(rangeName);
137 double xmin = _x.min(_refRangeName ? _refRangeName->GetName() : nullptr);
138 double xminFull = _x.min(rangeName);
139 unsigned int sz = _coefList.size();
140
141 std::vector<double> coeffs;
142 for (auto it : _coefList)
143 coeffs.push_back(static_cast<const RooAbsReal &>(*it).getVal());
144
145 return RooFit::Detail::AnalyticalIntegrals::chebychevIntegral(coeffs.data(), sz, xmin, xmax, xminFull, xmaxFull);
146}
147
148std::string RooChebychev::buildCallToAnalyticIntegral(Int_t /* code */, const char *rangeName,
150{
151 double xmax = _x.max(_refRangeName ? _refRangeName->GetName() : nullptr);
152 double xmaxFull = _x.max(rangeName);
153 double xmin = _x.min(_refRangeName ? _refRangeName->GetName() : nullptr);
154 double xminFull = _x.min(rangeName);
155 unsigned int sz = _coefList.size();
156
157 return ctx.buildCall("RooFit::Detail::AnalyticalIntegrals::chebychevIntegral", _coefList, sz, xmin, xmax, xminFull, xmaxFull);
158}
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Storage_t::size_type size() const
bool addTyped(const RooAbsCollection &list, bool silent=false)
Adds elements of a given RooAbsCollection to the container if they match the specified type.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
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
Chebychev polynomial p.d.f.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy _x
std::string buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
void selectNormalizationRange(const char *rangeName=nullptr, bool force=false) override
Interface function to force use of a given normalization range to interpret function value.
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
TNamed * _refRangeName
void computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Chebychev.
RooListProxy _coefList
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
A class to maintain the context for squashing of RooFit models into code.
std::string buildCall(std::string const &funcname, Args_t const &...args)
Build the code to call the function with name funcname, passing some arguments.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition DataMap.cxx:73
const TNamed * constPtr(const char *stringPtr)
Return a unique TNamed pointer for given C++ string.
static RooNameReg & instance()
Return reference to singleton instance.
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, RestrictArr output, size_t size, VarSpan vars, ArgSpan extraArgs={})
double chebychevIntegral(double const *coeffs, unsigned int nCoeffs, double xMin, double xMax, double xMinFull, double xMaxFull)
double chebychevEvaluate(double *coeffs, unsigned int nCoeffs, double x_in, double xMin, double xMax)
static void output()