Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RooStepFunction.cxx
Go to the documentation of this file.
1
2/*****************************************************************************
3 * Project: RooFit *
4 * Package: RooFitBabar *
5 * @(#)root/roofit:$Id$
6 * Author: *
7 * Tristan du Pree, Nikhef, Amsterdam, tdupree@nikhef.nl *
8 * Wouter Verkerke, Nikhef, Amsterdam, verkerke@nikhef.nl
9 * *
10 * Copyright (c) 2009, NIKHEF. 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 RooStepFunction
18 \ingroup Roofit
19
20The Step Function is a binned function whose parameters
21are the heights of each bin.
22
23This function may be used to describe oddly shaped distributions. A RooStepFunction
24has free parameters. In particular, any statistical uncertainty
25used to model this efficiency may be understood with these free parameters.
26
27Note that in contrast to RooParametricStepFunction, a RooStepFunction is NOT a PDF,
28but a not-normalized function (RooAbsReal)
29**/
30
31#include <RooStepFunction.h>
32
33#include <RooArgList.h>
34#include <RooCurve.h>
36#include <RooMath.h>
37#include <RooMsgService.h>
38
39
40namespace {
41
42// Get the values for each element in a RooListProxy, using the correct normalization set.
43std::span<double const> values(RooListProxy const &listProxy, std::vector<double> &valueCache)
44{
45 valueCache.resize(listProxy.size());
46 for (std::size_t i = 0; i < listProxy.size(); ++i) {
47 valueCache[i] = static_cast<RooAbsReal &>(listProxy[i]).getVal(listProxy.nset());
48 }
49 return valueCache;
50}
51
52} // namespace
53
54////////////////////////////////////////////////////////////////////////////////
55/// Constructor
56
57RooStepFunction::RooStepFunction(const char* name, const char* title,
58 RooAbsReal& x, const RooArgList& coefList, const RooArgList& boundaryList, bool interpolate) :
59 RooAbsReal(name, title),
60 _x("x", "Dependent", this, x),
61 _coefList("coefList","List of coefficients",this),
62 _boundaryList("boundaryList","List of boundaries",this),
63 _interpolate(interpolate)
64{
65 _coefList.addTyped<RooAbsReal>(coefList);
67
68 if (_boundaryList.size()!=_coefList.size()+1) {
69 coutE(InputArguments) << "RooStepFunction::ctor(" << GetName() << ") ERROR: Number of boundaries must be number of coefficients plus 1" << std::endl ;
70 throw std::invalid_argument("RooStepFunction::ctor() ERROR: Number of boundaries must be number of coefficients plus 1") ;
71 }
72
73}
74
75////////////////////////////////////////////////////////////////////////////////
76/// Copy constructor
77
80 _x("x", this, other._x),
81 _coefList("coefList",this,other._coefList),
82 _boundaryList("boundaryList",this,other._boundaryList),
83 _interpolate(other._interpolate)
84{
85}
86
87
88////////////////////////////////////////////////////////////////////////////////
89/// Transfer contents to std::vector for use below
90
92{
93 std::span<const double> b = values(_boundaryList, _boundaryCache);
94 int nb = b.size();
95
96 // Return zero if outside any boundaries
97 if ((_x<b[0]) || (_x>b[nb-1])) return 0 ;
98
99 if (!_interpolate) {
100
101 // No interpolation -- Return values bin-by-bin
102 for (int i=0;i<nb-1;i++){
103 if (_x>b[i]&&_x<=b[i+1]) {
104 return (static_cast<RooAbsReal*>(_coefList.at(i)))->getVal() ;
105 }
106 }
107 return 0 ;
108
109 }
110
111 std::vector<double> c(_coefList.size()+3) ;
112
113 // Interpolation
114
115 // Make array of (b[0],bin centers,b[last])
116 c[0] = b[0] ; c[nb] = b[nb-1] ;
117 for (int i=0 ; i<nb-1 ; i++) {
118 c[i+1] = (b[i]+b[i+1])/2 ;
119 }
120
121 // Make array of (0,coefficient values,0)
122 int nc(0) ;
123 std::vector<double> y(_coefList.size()+3) ;
124 y[nc++] = 0 ;
125 for(auto * coef : static_range_cast<RooAbsReal*>(_coefList)) {
126 y[nc++] = coef->getVal() ;
127 }
128 y[nc++] = 0 ;
129
130 for (int i=0;i<nc-1;i++){
131 if (_x>c[i]&&_x<=c[i+1]) {
132 double xx[2] ; xx[0]=c[i] ; xx[1]=c[i+1] ;
133 double yy[2] ; yy[0]=y[i] ; yy[1]=y[i+1] ;
134 return RooMath::interpolate(xx,yy,2,_x) ;
135 }
136 }
137 return 0;
138}
139
140
141std::list<double> *RooStepFunction::plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const
142{
143 if (obs.namePtr() != _x->namePtr()) {
144 return nullptr;
145 }
146
147 // Use the helper function from RooCurve to make sure to get sampling hints
148 // that work with the RooFitPlotting.
150}
151
152int RooStepFunction::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
153{
154 // We only support analytical integration if we integrate over x and there is
155 // no interpolation.
156 if (!_interpolate && matchArgs(allVars, analVars, _x))
157 return 1;
158 return 0;
159}
160
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define coutE(a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
const TNamed * namePtr() const
De-duplicated pointer to this object's name.
Definition RooAbsArg.h:504
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 base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
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
static double interpolate(double yArr[], Int_t nOrder, double x)
Definition RooMath.cxx:78
The Step Function is a binned function whose parameters are the heights of each bin.
std::vector< double > _boundaryCache
RooListProxy _coefList
int getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooListProxy _boundaryList
double evaluate() const override
Transfer contents to std::vector for use below.
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...
std::vector< double > _coefCache
double analyticalIntegral(int code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy _x
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.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
double stepFunctionIntegral(double xmin, double xmax, std::size_t nBins, double const *boundaries, double const *coefs)
Definition MathFuncs.h:776