Logo ROOT  
Reference Guide
 
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>
35#include <RooMsgService.h>
36#include <RooMath.h>
37#include <RooRealVar.h>
38
40
41////////////////////////////////////////////////////////////////////////////////
42/// Constructor
43
44RooStepFunction::RooStepFunction(const char* name, const char* title,
45 RooAbsReal& x, const RooArgList& coefList, const RooArgList& boundaryList, bool interpolate) :
46 RooAbsReal(name, title),
47 _x("x", "Dependent", this, x),
48 _coefList("coefList","List of coefficients",this),
49 _boundaryList("boundaryList","List of boundaries",this),
50 _interpolate(interpolate)
51{
52 _coefList.addTyped<RooAbsReal>(coefList);
53 _boundaryList.addTyped<RooAbsReal>(boundaryList);
54
55 if (_boundaryList.size()!=_coefList.size()+1) {
56 coutE(InputArguments) << "RooStepFunction::ctor(" << GetName() << ") ERROR: Number of boundaries must be number of coefficients plus 1" << std::endl ;
57 throw std::invalid_argument("RooStepFunction::ctor() ERROR: Number of boundaries must be number of coefficients plus 1") ;
58 }
59
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Copy constructor
64
66 RooAbsReal(other, name),
67 _x("x", this, other._x),
68 _coefList("coefList",this,other._coefList),
69 _boundaryList("boundaryList",this,other._boundaryList),
70 _interpolate(other._interpolate)
71{
72}
73
74
75////////////////////////////////////////////////////////////////////////////////
76/// Transfer contents to std::vector for use below
77
79{
80 std::vector<double> b(_boundaryList.size()) ;
81 std::vector<double> c(_coefList.size()+3) ;
82 Int_t nb(0) ;
83 for (auto * boundary : static_range_cast<RooAbsReal*>(_boundaryList)) {
84 b[nb++] = boundary->getVal() ;
85 }
86
87 // Return zero if outside any boundaries
88 if ((_x<b[0]) || (_x>b[nb-1])) return 0 ;
89
90 if (!_interpolate) {
91
92 // No interpolation -- Return values bin-by-bin
93 for (Int_t i=0;i<nb-1;i++){
94 if (_x>b[i]&&_x<=b[i+1]) {
95 return (static_cast<RooAbsReal*>(_coefList.at(i)))->getVal() ;
96 }
97 }
98 return 0 ;
99
100 } else {
101
102 // Interpolation
103
104 // Make array of (b[0],bin centers,b[last])
105 c[0] = b[0] ; c[nb] = b[nb-1] ;
106 for (Int_t i=0 ; i<nb-1 ; i++) {
107 c[i+1] = (b[i]+b[i+1])/2 ;
108 }
109
110 // Make array of (0,coefficient values,0)
111 Int_t nc(0) ;
112 std::vector<double> y(_coefList.size()+3) ;
113 y[nc++] = 0 ;
114 for(auto * coef : static_range_cast<RooAbsReal*>(_coefList)) {
115 y[nc++] = coef->getVal() ;
116 }
117 y[nc++] = 0 ;
118
119 for (Int_t i=0;i<nc-1;i++){
120 if (_x>c[i]&&_x<=c[i+1]) {
121 double xx[2] ; xx[0]=c[i] ; xx[1]=c[i+1] ;
122 double yy[2] ; yy[0]=y[i] ; yy[1]=y[i+1] ;
123 return RooMath::interpolate(xx,yy,2,_x) ;
124 }
125 }
126 return 0;
127 }
128}
129
130
131std::list<double> *RooStepFunction::plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const
132{
133 if (obs.namePtr() != _x->namePtr()) {
134 return nullptr;
135 }
136
137 // Retrieve position of all bin boundaries
138 std::vector<double> boundaries;
140 for (auto *boundary : static_range_cast<RooAbsReal *>(_boundaryList)) {
141 boundaries.push_back(boundary->getVal());
142 }
143
144 // Use the helper function from RooCurve to make sure to get sampling hints
145 // that work with the RooFitPlotting.
147}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:382
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
void reserve(Storage_t::size_type count)
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
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
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:893
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.
const RooArgList & boundaries()
RooListProxy _coefList
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...
RooRealProxy _x
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17