Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RooParametricStepFunction.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitBabar *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * Aaron Roodman, Stanford Linear Accelerator Center, Stanford University *
7 * *
8 * Copyright (c) 2004, Stanford University. All rights reserved. *
9 *
10 * Redistribution and use in source and binary forms, *
11 * with or without modification, are permitted according to the terms *
12 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
13 *****************************************************************************/
14
15/** \class RooParametricStepFunction
16 \ingroup Roofit
17
18The Parametric Step Function PDF is a binned distribution whose parameters
19are the heights of each bin. This PDF was first used in BaBar's B0->pi0pi0
20paper BABAR Collaboration (B. Aubert et al.) Phys.Rev.Lett.91:241801,2003.
21
22This PDF may be used to describe oddly shaped distributions. It differs
23from a RooKeysPdf or a RooHistPdf in that a RooParametricStepFunction
24has free parameters. In particular, any statistical uncertainty in
25sample used to model this PDF may be understood with these free parameters;
26this is not possible with non-parametric PDFs.
27
28The RooParametricStepFunction has Nbins-1 free parameters. Note that
29the limits of the dependent variable must match the low and hi bin limits.
30
31Here is an example showing how to use the RooParametricStepFunction to fit toy
32data generated from a uniform distribution:
33
34~~~ {.cpp}
35// Define some constant parameters
36const int nBins = 10;
37const double xMin = 0.0;
38const double xMax = 10.0;
39const double binWidth = (xMax - xMin) / nBins;
40const std::size_t nEvents = 10000;
41
42// Fill the bin boundaries
43std::vector<double> binBoundaries(nBins + 1);
44
45for(int i = 0; i <= nBins; ++i) {
46 binBoundaries[i] = i * binWidth;
47}
48
49// The RooParametricStepFunction needs a TArrayD
50TArrayD binBoundariesTArr{int(binBoundaries.size()), binBoundaries.data()};
51
52RooRealVar x{"x", "x", xMin, xMax};
53
54// There are nBins-1 free coefficient parameters, whose sum must be <= 1.0.
55// We all set them to 0.1, such that the resulting step function pdf is
56// initially uniform.
57RooArgList coefList;
58for(int i = 0; i < nBins - 1; ++i) {
59 auto name = std::string("coef_") + std::to_string(i);
60 coefList.addOwned(std::make_unique<RooRealVar>(name.c_str(), name.c_str(), 0.1, 0.0, 1.0));
61}
62
63// Create the RooParametricStepFunction pdf
64RooParametricStepFunction pdf{"pdf", "pdf", x, coefList, binBoundariesTArr, int(nBins)};
65
66// Generate some toy data, following our uniform step function pdf
67std::unique_ptr<RooAbsData> data{pdf.generate(x, nEvents)};
68
69// Fit the step function to the toy data
70pdf.fitTo(*data);
71
72// Now we plot the data and the pdf, the latter should not be uniform
73// anymore because the coefficients were fit to the toy data
74RooPlot *xframe = x.frame(RooFit::Title("Fitting uniform toy data with a step function p.d.f."));
75data->plotOn(xframe);
76pdf.plotOn(xframe);
77xframe->Draw();
78~~~
79*/
80
82
83#include <RooCurve.h>
84#include <RooRealVar.h>
85
86
87////////////////////////////////////////////////////////////////////////////////
88/// Constructor
89
91 RooAbsReal& x, const RooArgList& coefList, TArrayD const& limits, Int_t nBins) :
92 RooAbsPdf(name, title),
93 _x("x", "Dependent", this, x),
94 _coefList("coefList","List of coefficients",this),
95 _nBins(nBins)
96{
97
98 // Check lowest order
99 if (_nBins<0) {
100 std::cout << "RooParametricStepFunction::ctor(" << GetName()
101 << ") WARNING: nBins must be >=0, setting value to 0" << std::endl ;
102 _nBins=0 ;
103 }
104
105 _coefList.addTyped<RooAbsReal>(coefList);
106
107 // Bin limits
108 limits.Copy(_limits);
109
110}
111
112////////////////////////////////////////////////////////////////////////////////
113/// Copy constructor
114
116 RooAbsPdf(other, name),
117 _x("x", this, other._x),
118 _coefList("coefList",this,other._coefList),
119 _nBins(other._nBins)
120{
121 (other._limits).Copy(_limits);
122}
123
124////////////////////////////////////////////////////////////////////////////////
125
126Int_t RooParametricStepFunction::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
127{
128 if (matchArgs(allVars, analVars, _x)) return 1;
129 return 0;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
134double RooParametricStepFunction::analyticalIntegral(Int_t /*code*/, const char* rangeName) const
135{
136 // Case without range is trivial: p.d.f is by construction normalized
137 if (!rangeName) {
138 return 1.0 ;
139 }
140
141 // Case with ranges, calculate integral explicitly
142 double xmin = _x.min(rangeName) ;
143 double xmax = _x.max(rangeName) ;
144
145 double sum=0 ;
146 Int_t i ;
147 for (i=1 ; i<=_nBins ; i++) {
148 double binVal = (i<_nBins) ? (static_cast<RooRealVar*>(_coefList.at(i-1))->getVal()) : lastBinValue() ;
149 if (_limits[i-1]>=xmin && _limits[i]<=xmax) {
150 // Bin fully in the integration domain
151 sum += (_limits[i]-_limits[i-1])*binVal ;
152 } else if (_limits[i-1]<xmin && _limits[i]>xmax) {
153 // Domain is fully contained in this bin
154 sum += (xmax-xmin)*binVal ;
155 // Exit here, this is the last bin to be processed by construction
156 return sum ;
157 } else if (_limits[i-1]<xmin && _limits[i]<=xmax && _limits[i]>xmin) {
158 // Lower domain boundary is in bin
159 sum += (_limits[i]-xmin)*binVal ;
160 } else if (_limits[i-1]>=xmin && _limits[i]>xmax && _limits[i-1]<xmax) {
161 sum += (xmax-_limits[i-1])*binVal ;
162 // Upper domain boundary is in bin
163 // Exit here, this is the last bin to be processed by construction
164 return sum ;
165 }
166 }
167
168 return sum;
169}
170
171////////////////////////////////////////////////////////////////////////////////
172
174{
175 double sum(0.);
176 double binSize(0.);
177 for (Int_t j=1;j<_nBins;j++){
178 RooRealVar* tmp = static_cast<RooRealVar*>(_coefList.at(j-1));
179 binSize = _limits[j] - _limits[j-1];
180 sum = sum + tmp->getVal()*binSize;
181 }
182 binSize = _limits[_nBins] - _limits[_nBins-1];
183 return (1.0 - sum)/binSize;
184}
185
186////////////////////////////////////////////////////////////////////////////////
187
189{
190 double value(0.);
191 if (_x >= _limits[0] && _x < _limits[_nBins]){
192
193 for (Int_t i=1;i<=_nBins;i++){
194 if (_x < _limits[i]){
195 // in Bin i-1 (starting with Bin 0)
196 if (i<_nBins) {
197 // not in last Bin
198 value = static_cast<RooRealVar*>(_coefList.at(i-1))->getVal();
199 break;
200 } else {
201 value = lastBinValue();
202 if (value<=0.0){
203 value = 0.000001;
204 // std::cout << "RooParametricStepFunction: sum of values gt 1.0 -- beware!!" << std::endl;
205 }
206 break;
207 }
208 }
209 }
210
211 }
212 return value;
213
214}
215
216
217std::list<double>* RooParametricStepFunction::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
218{
219 if(obs.namePtr() != _x->namePtr()) {
220 return nullptr;
221 }
222
223 // Retrieve position of all bin boundaries
224 std::span<const double> boundaries{_limits.GetArray(), static_cast<std::size_t>(_limits.GetSize())};
225
226 // Use the helper function from RooCurve to make sure to get sampling hints
227 // that work with the RooFitPlotting.
228 return RooCurve::plotSamplingHintForBinBoundaries(boundaries, xlo, xhi);
229}
int Int_t
Signed integer 4 bytes (int).
Definition RtypesCore.h:59
char name[80]
Definition TGX11.cxx:148
float xmin
float xmax
const TNamed * namePtr() const
De-duplicated pointer to this object's name.
Definition RooAbsArg.h:502
RooAbsPdf()
Default constructor.
friend class RooAbsReal
Definition RooAbsPdf.h:342
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
bool matchArgs(const RooArgSet &allDeps, RooArgSet &analDeps, const RooArgProxy &a, const Proxies &... proxies) const
Definition RooAbsReal.h:425
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:24
static std::list< double > * plotSamplingHintForBinBoundaries(std::span< const double > boundaries, double xlo, double xhi)
Returns sampling hints for a histogram with given boundaries.
Definition RooCurve.cxx:897
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Array of doubles (64 bits per element).
Definition TArrayD.h:27
void Copy(TArrayD &array) const
Definition TArrayD.h:42
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Double_t x[n]
Definition legend1.C:17
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2338