MINT2
RooCubicSplinePdf.cpp
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id: RooBSpline.cxx 45780 2012-08-31 15:45:27Z moneta $
5  * Authors: *
6  * Gerhard Raven
7  * *
8  *****************************************************************************/
9 
11 //
12 // BEGIN_HTML
13 // BSpline basis polynomials are positive-definite in the range [0,1].
14 // In this implementation, we extend [0,1] to be the range of the parameter.
15 // There are n+1 BSpline basis polynomials of degree n.
16 // Thus, by providing N coefficients that are positive-definite, there
17 // is a natural way to have well bahaved polynomail PDFs.
18 // For any n, the n+1 basis polynomials 'form a partition of unity', eg.
19 // they sum to one for all values of x. See
20 // END_HTML
21 //
22 // STD & STL
23 #include <cmath>
24 #include <complex>
25 #include <sstream>
26 
27 // ROOT
28 #include "TH1.h"
29 #include "TGraph.h"
30 #include "TGraphErrors.h"
31 
32 // RooFit
33 #include "RooMsgService.h"
34 #include "RooMath.h"
35 #include "RooAbsReal.h"
36 #include "RooRealVar.h"
37 #include "RooConstVar.h"
38 #include "RooArgList.h"
39 
40 // P2VV
41 #include "Mint/RooCubicSplinePdf.h"
43 
44 using namespace std;
45 
46 //ClassImp(RooCubicSplinePdf);
47 
48 //_____________________________________________________________________________
49 void RooCubicSplinePdf::init(const char* name,
50  const vector<double>& heights,
51  const vector<double>& errors,
52  double smooth, bool constCoeffs) {
53  vector<double> values(heights);
54  if ( smooth > 0 ) {
55  assert(int(errors.size()) == _aux.size());
56  _aux.smooth( values, errors, smooth );
57  }
58  _aux.computeCoefficients( values );
59  for (unsigned int i=0;i<values.size();++i) {
60  if (constCoeffs) {
61  _coefList.add( RooFit::RooConst( values[i] ) );
62  } else {
63  stringstream name_str;
64  name_str << name << "_smoothed_bin_" << i;
65  string n(name_str.str());
66  RooRealVar* coeff = new RooRealVar(n.c_str(), n.c_str(), values[i], 0.0, 1.0);
67  if (i == 0 || i == values.size() - 1) coeff->setConstant(true);
68  _coefList.add(*coeff);
69  addOwnedComponents( *coeff );
70  }
71  }
72 }
73 
74 //_____________________________________________________________________________
76  : _aux(0,0)
77 {
78 }
79 
80 //_____________________________________________________________________________
81 RooCubicSplinePdf::RooCubicSplinePdf(const char* name, const char* title, RooRealVar& x,
82  const vector<double>& _knots,
83  const vector<double>& values,
84  const vector<double>& errors,
85  double smooth, bool constCoeffs) :
86  RooAbsPdf(name, title),
87  _x("x", "Dependent", this, x),
88  _coefList("coefficients","List of coefficients",this),
89  _aux(_knots.begin(),_knots.end())
90 {
91  init(name, values, errors, smooth, constCoeffs);
92 }
93 
94 //_____________________________________________________________________________
95 RooCubicSplinePdf::RooCubicSplinePdf(const char* name, const char* title,
96  RooRealVar& x, const TGraph* graph, bool constCoeffs) :
97  RooAbsPdf(name, title),
98  _x("x", "Dependent", this, x),
99  _coefList("coefficients","List of coefficients",this),
100  _aux(0,0)
101 {
102  int nPoints = graph->GetN();
103  vector<double> centres, values;
104  for (int i=0;i<nPoints ;++i) {
105  double xx,yy;
106  graph->GetPoint(i,xx,yy);
107  centres.push_back(xx);
108  values.push_back(yy);
109  }
110  _aux = RooCubicSplineKnot( centres.begin(), centres.end() );
111  vector<double> errs;
112  init(name, values, errs, 0, constCoeffs);
113 }
114 
115 //_____________________________________________________________________________
116 RooCubicSplinePdf::RooCubicSplinePdf(const char* name, const char* title,
117  RooRealVar& x, const TH1* hist, double smooth,
118  bool constCoeffs) :
119  RooAbsPdf(name, title),
120  _x("x", "Dependent", this, x),
121  _coefList("coefficients","List of coefficients",this),
122  _aux(0,0)
123 {
124  // bin 0 is underflow, and bin nBins + 1 is overflow...
125  int nBins = hist->GetNbinsX();
126  vector<double> centres, values;
127  for (int i=0;i<nBins ;++i) {
128  centres.push_back(hist->GetBinCenter(1+i));
129  values.push_back(hist->GetBinContent(1+i));
130  }
131  _aux = RooCubicSplineKnot( centres.begin(), centres.end() );
132 
133  vector<double> errs;
134  if (smooth>0) for (int i=0;i<nBins ;++i) errs.push_back(hist->GetBinError(1+i));
135 
136  init(name, values, errs, smooth, constCoeffs);
137 }
138 //_____________________________________________________________________________
139 RooCubicSplinePdf::RooCubicSplinePdf(const char* name, const char* title,
140  RooRealVar& x, const TGraphErrors* graph, double smooth,
141  bool constCoeffs) :
142  RooAbsPdf(name, title),
143  _x("x", "Dependent", this, x),
144  _coefList("coefficients","List of coefficients",this),
145  _aux(0,0)
146 {
147  int nPoints = graph->GetN();
148  vector<double> centres, values;
149  for (int i=0;i<nPoints ;++i) {
150  double xx,yy;
151  graph->GetPoint(i,xx,yy);
152  centres.push_back(xx);
153  values.push_back(yy);
154  }
155  _aux = RooCubicSplineKnot( centres.begin(), centres.end() );
156 
157  vector<double> errs;
158  if (smooth>0) for (int i=0;i<nPoints ;++i) errs.push_back(graph->GetErrorY(i));
159 
160  init(name, values, errs, smooth, constCoeffs);
161 }
162 
163 //_____________________________________________________________________________
164 RooCubicSplinePdf::RooCubicSplinePdf(const char* name, const char* title,
165  RooRealVar& x, const char* knotBinningName,
166  const RooArgList& coefList):
167  RooAbsPdf(name, title),
168  _x("x", "Dependent", this, x),
169  _coefList("coefficients", "List of coefficients", this),
170  _aux(0,0)
171 {
172  // TODO: verify coefList is consistent with knots as specified by the knotBinningName binning
173  // should be N+2 coefficients for N bins...
174  const RooAbsBinning* binning = x.getBinningPtr(knotBinningName);
175  assert( binning!=0);
176  if ( coefList.getSize()!=2+binning->numBoundaries()) {
177  cout << TString::Format( "you have specified %d coefficients for %d knots. The differentce should 2!",coefList.getSize(),binning->numBoundaries()) << endl;
178  throw TString::Format( "you have specified %d coefficients for %d knots. The differentce should 2!",coefList.getSize(),binning->numBoundaries());
179  }
180  _coefList.add(coefList);
181 
182  Double_t* boundaries = binning->array();
183  _aux = RooCubicSplineKnot( boundaries, boundaries + binning->numBoundaries() );
184 }
185 
186 //_____________________________________________________________________________
187 RooCubicSplinePdf::RooCubicSplinePdf(const char* name, const char* title,
188  RooRealVar& x, const vector<double>& _knots,
189  const RooArgList& coefList):
190  RooAbsPdf(name, title),
191  _x("x", "Dependent", this, x),
192  _coefList("coefficients", "List of coefficients", this),
193  _aux(_knots.begin(), _knots.end())
194 {
195  assert(size_t(coefList.getSize()) == 2 + _knots.size());
196  _coefList.add(coefList);
197 }
198 
199 //_____________________________________________________________________________
201  RooAbsPdf(other, name),
202  _x("x", this, other._x),
203  _coefList("coefList", this, other._coefList),
204  _aux(other._aux)
205 {
206 }
207 
208 //_____________________________________________________________________________
210 {
211 }
212 
213 //_____________________________________________________________________________
215 {
216  return _aux.evaluate(_x,_coefList);
217 }
218 
219 //_____________________________________________________________________________
220 Int_t RooCubicSplinePdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
221 {
222  // No analytical calculation available (yet) of integrals over subranges
223  if (_x.min(rangeName)!=_aux.knots().front() || _x.max(rangeName)!=_aux.knots().back() ) return 0;
224  if (matchArgs(allVars, analVars, _x)) return 1;
225  return 0;
226 }
227 //_____________________________________________________________________________
228 Double_t RooCubicSplinePdf::analyticalIntegral(Int_t code, const char* /* rangeName */) const
229 {
230  if (code != 1) {
231  coutE(InputArguments) << "RooCubicSplinePdf::analyticalIntegral(" << GetName()
232  << "): argument \"code\" can only have value 1" << std::endl;
233  assert(code==1) ;
234  }
236 }
237 
238 //_____________________________________________________________________________
240  const RooGaussModelAcceptance::K_n& K, double _offset,
241  double* sc) const
242 {
243  RooCubicSplineKnot::S_edge s_jk( _aux.S_jk_edge( left, _coefList ), _offset );
244  return dM(0)*( s_jk(0,0) * K(0) * sc[0] + s_jk(0,1) * K(1) * sc[1] )
245  + dM(1)*( s_jk(1,0) * K(0) * sc[1] );
246 }
247 
248 //_____________________________________________________________________________
249 std::complex<double>
250 RooCubicSplinePdf::productAnalyticalIntegral(Double_t umin, Double_t umax,
251  Double_t scale, Double_t _offset,
252  const std::complex<double>& z) const
253 {
255  assert(knotSize()>1);
257  std::vector<M_n> M; M.reserve( knotSize() );
258  for (unsigned int i=0;i<knotSize();++i) {
259  double x = (u(i)-_offset)/scale ;
260  M.push_back( M_n( x, z ) );
261  }
262  double sc[4]; for (int i=0;i<4;++i) sc[i] = pow(scale,i);
263  double lo = scale*umin+_offset;
264  double hi = scale*umax+_offset;
265  std::complex<double> sum(0,0);
266  //TODO: verify we remain within [lo,hi]
267  assert(hi>=u(0)); // front only if hi>u(0)!!!
268  if (lo<u(0)) sum += gaussIntegralE(true, M.front()-M_n( umin,z), K, _offset, sc);
269  for (unsigned i=0;i<knotSize()-1 && u(i)<hi ;++i) {
270  if (u(i+1)<lo) continue;
271  // FIXME:TODO: we currently assume that u(0),u(knotSize()-1)] fully contained in [lo,hi]
272  assert(lo<=u(i));
273  assert(u(i+1)<=hi);
274  M_n dM = M[i+1]-M[i]; // take M[i] if lo<=u(i) else M_n(lo) ; take M[i+1] if u(i+1)<=hi else M_n(hi)
275  RooCubicSplineKnot::S_jk s_jk( _aux.S_jk_sum( i, _coefList ), _offset ); // pass sc into S_jk, remove from loop
276  for (int j=0;j<4;++j) for (int k=0;k<4-j;++k) sum += dM(j)*s_jk(j,k)*K(k)*sc[j+k];
277  }
278  assert(lo<=u(knotSize()-1));// back only if lo<u(knotsiwze()-1)!!!
279  if (hi>u(knotSize()-1)) sum += gaussIntegralE(false, M_n(umax,z)-M.back(), K, _offset, sc);
280  return sum;
281 }
282 
283 //_____________________________________________________________________________
284 Int_t RooCubicSplinePdf::getMaxVal(const RooArgSet& vars) const
285 {
286  // check that vars only contains _x...
287  return ( vars.getSize() == 1 && vars.contains( _x.arg() ) ) ? 1 : 0;
288 }
289 //_____________________________________________________________________________
290 Double_t RooCubicSplinePdf::maxVal(Int_t code) const
291 {
292  if (code != 1) {
293  coutE(InputArguments) << "RooCubicSplinePdf::maxVal(" << GetName()
294  << "): argument \"code\" can only have value 1" << std::endl;
295  assert(code==1) ;
296  }
297  RooFIter iter = _coefList.fwdIterator();
298  RooAbsReal *c(0);
299  double res = 0;
300  while((c=(RooAbsReal*)iter.next())) {
301  double x = fabs(c->getVal());
302  if (x>res) { res = x; }
303  }
304  return res;
305 }
RooCubicSplineKnot::S_jk S_jk_sum(int i, const RooArgList &b) const
unsigned knotSize() const
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName) const
void init(const char *name, const std::vector< double > &heights, const std::vector< double > &errors, double smooth, bool constCoeffs)
Int_t getMaxVal(const RooArgSet &vars) const
RooCubicSplineKnot::S_edge S_jk_edge(bool left, const RooArgList &b) const
std::complex< double > gaussIntegralE(bool left, const RooGaussModelAcceptance::M_n< 4U > &dM, const RooGaussModelAcceptance::K_n &K, double offset, double *sc) const
std::complex< double > productAnalyticalIntegral(Double_t umin, Double_t umax, Double_t scale, Double_t offset, const std::complex< double > &z) const
double analyticalIntegral(const RooArgList &b) const
RooCubicSplineKnot _aux
double u(int i) const
const std::vector< double > & knots() const
double evaluate(double _u, const RooArgList &b) const
Double_t maxVal(Int_t code) const
Double_t analyticalIntegral(Int_t code, const char *rangeName) const
RooListProxy _coefList
Double_t evaluate() const