MINT2
Public Member Functions | Private Member Functions | Private Attributes | List of all members
RooCubicSplinePdf Class Reference

#include <RooCubicSplinePdf.h>

Inheritance diagram for RooCubicSplinePdf:

Public Member Functions

 RooCubicSplinePdf ()
 
 RooCubicSplinePdf (const char *name, const char *title, RooRealVar &x, const std::vector< double > &knots, const std::vector< double > &values, const std::vector< double > &errors=std::vector< double >(), double smooth=0, bool constCoeffs=true)
 
 RooCubicSplinePdf (const char *name, const char *title, RooRealVar &x, const TGraph *graph, bool constCoeffs=true)
 
 RooCubicSplinePdf (const char *name, const char *title, RooRealVar &x, const TH1 *hist, double smooth=0, bool constCoeffs=true)
 
 RooCubicSplinePdf (const char *name, const char *title, RooRealVar &x, const TGraphErrors *graph, double smooth=0, bool constCoeffs=true)
 
 RooCubicSplinePdf (const char *name, const char *title, RooRealVar &x, const char *knotBinningName, const RooArgList &coefList)
 
 RooCubicSplinePdf (const char *name, const char *title, RooRealVar &x, const std::vector< double > &knots, const RooArgList &coefList)
 
 ~RooCubicSplinePdf ()
 
 RooCubicSplinePdf (const RooCubicSplinePdf &other, const char *name=0)
 
TObject * clone (const char *newname) const
 
Int_t getAnalyticalIntegral (RooArgSet &allVars, RooArgSet &analVars, const char *rangeName) const
 
Double_t analyticalIntegral (Int_t code, const char *rangeName) const
 
Int_t getMaxVal (const RooArgSet &vars) const
 
Double_t maxVal (Int_t code) const
 
std::complex< double > productAnalyticalIntegral (Double_t umin, Double_t umax, Double_t scale, Double_t offset, const std::complex< double > &z) const
 
unsigned knotSize () const
 
double u (int i) const
 
const std::vector< double > & knots () const
 
const RooArgList & coefficients () const
 

Private Member Functions

void init (const char *name, const std::vector< double > &heights, const std::vector< double > &errors, double smooth, bool constCoeffs)
 
Double_t evaluate () const
 
std::complex< double > gaussIntegralE (bool left, const RooGaussModelAcceptance::M_n< 4U > &dM, const RooGaussModelAcceptance::K_n &K, double offset, double *sc) const
 

Private Attributes

RooRealProxy _x
 
RooListProxy _coefList
 
RooCubicSplineKnot _aux
 

Detailed Description

Definition at line 28 of file RooCubicSplinePdf.h.

Constructor & Destructor Documentation

◆ RooCubicSplinePdf() [1/8]

RooCubicSplinePdf::RooCubicSplinePdf ( )

Definition at line 75 of file RooCubicSplinePdf.cpp.

76  : _aux(0,0)
77 {
78 }
RooCubicSplineKnot _aux

◆ RooCubicSplinePdf() [2/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const char *  name,
const char *  title,
RooRealVar &  x,
const std::vector< double > &  knots,
const std::vector< double > &  values,
const std::vector< double > &  errors = std::vector< double >(),
double  smooth = 0,
bool  constCoeffs = true 
)

◆ RooCubicSplinePdf() [3/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const char *  name,
const char *  title,
RooRealVar &  x,
const TGraph *  graph,
bool  constCoeffs = true 
)

Definition at line 95 of file RooCubicSplinePdf.cpp.

96  :
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 }
void init(const char *name, const std::vector< double > &heights, const std::vector< double > &errors, double smooth, bool constCoeffs)
RooCubicSplineKnot _aux
RooListProxy _coefList

◆ RooCubicSplinePdf() [4/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const char *  name,
const char *  title,
RooRealVar &  x,
const TH1 *  hist,
double  smooth = 0,
bool  constCoeffs = true 
)

Definition at line 116 of file RooCubicSplinePdf.cpp.

118  :
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 }
void init(const char *name, const std::vector< double > &heights, const std::vector< double > &errors, double smooth, bool constCoeffs)
RooCubicSplineKnot _aux
RooListProxy _coefList

◆ RooCubicSplinePdf() [5/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const char *  name,
const char *  title,
RooRealVar &  x,
const TGraphErrors *  graph,
double  smooth = 0,
bool  constCoeffs = true 
)

Definition at line 139 of file RooCubicSplinePdf.cpp.

141  :
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 }
void init(const char *name, const std::vector< double > &heights, const std::vector< double > &errors, double smooth, bool constCoeffs)
RooCubicSplineKnot _aux
RooListProxy _coefList

◆ RooCubicSplinePdf() [6/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const char *  name,
const char *  title,
RooRealVar &  x,
const char *  knotBinningName,
const RooArgList &  coefList 
)

Definition at line 164 of file RooCubicSplinePdf.cpp.

166  :
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 }
RooCubicSplineKnot _aux
RooListProxy _coefList

◆ RooCubicSplinePdf() [7/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const char *  name,
const char *  title,
RooRealVar &  x,
const std::vector< double > &  knots,
const RooArgList &  coefList 
)

◆ ~RooCubicSplinePdf()

RooCubicSplinePdf::~RooCubicSplinePdf ( )

Definition at line 209 of file RooCubicSplinePdf.cpp.

210 {
211 }

◆ RooCubicSplinePdf() [8/8]

RooCubicSplinePdf::RooCubicSplinePdf ( const RooCubicSplinePdf other,
const char *  name = 0 
)

Definition at line 200 of file RooCubicSplinePdf.cpp.

200  :
201  RooAbsPdf(other, name),
202  _x("x", this, other._x),
203  _coefList("coefList", this, other._coefList),
204  _aux(other._aux)
205 {
206 }
RooCubicSplineKnot _aux
RooListProxy _coefList

Member Function Documentation

◆ analyticalIntegral()

Double_t RooCubicSplinePdf::analyticalIntegral ( Int_t  code,
const char *  rangeName 
) const

Definition at line 228 of file RooCubicSplinePdf.cpp.

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 }
double analyticalIntegral(const RooArgList &b) const
RooCubicSplineKnot _aux
RooListProxy _coefList

◆ clone()

TObject* RooCubicSplinePdf::clone ( const char *  newname) const
inline

Definition at line 53 of file RooCubicSplinePdf.h.

53 { return new RooCubicSplinePdf(*this, newname); }

◆ coefficients()

const RooArgList& RooCubicSplinePdf::coefficients ( ) const
inline

Definition at line 69 of file RooCubicSplinePdf.h.

69 { return _coefList; }
RooListProxy _coefList

◆ evaluate()

Double_t RooCubicSplinePdf::evaluate ( ) const
private

Definition at line 214 of file RooCubicSplinePdf.cpp.

215 {
216  return _aux.evaluate(_x,_coefList);
217 }
RooCubicSplineKnot _aux
double evaluate(double _u, const RooArgList &b) const
RooListProxy _coefList

◆ gaussIntegralE()

complex< double > RooCubicSplinePdf::gaussIntegralE ( bool  left,
const RooGaussModelAcceptance::M_n< 4U > &  dM,
const RooGaussModelAcceptance::K_n K,
double  offset,
double *  sc 
) const
private

Definition at line 239 of file RooCubicSplinePdf.cpp.

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 }
RooCubicSplineKnot::S_edge S_jk_edge(bool left, const RooArgList &b) const
RooCubicSplineKnot _aux
RooListProxy _coefList

◆ getAnalyticalIntegral()

Int_t RooCubicSplinePdf::getAnalyticalIntegral ( RooArgSet &  allVars,
RooArgSet &  analVars,
const char *  rangeName 
) const

Definition at line 220 of file RooCubicSplinePdf.cpp.

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 }
RooCubicSplineKnot _aux
const std::vector< double > & knots() const

◆ getMaxVal()

Int_t RooCubicSplinePdf::getMaxVal ( const RooArgSet &  vars) const

Definition at line 284 of file RooCubicSplinePdf.cpp.

285 {
286  // check that vars only contains _x...
287  return ( vars.getSize() == 1 && vars.contains( _x.arg() ) ) ? 1 : 0;
288 }

◆ init()

void RooCubicSplinePdf::init ( const char *  name,
const std::vector< double > &  heights,
const std::vector< double > &  errors,
double  smooth,
bool  constCoeffs 
)
private

Definition at line 49 of file RooCubicSplinePdf.cpp.

52  {
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 }
void smooth(std::vector< double > &y, const std::vector< double > &dy, double lambda) const
RooCubicSplineKnot _aux
void computeCoefficients(std::vector< double > &y, BoundaryConditions bc=BoundaryConditions()) const
RooListProxy _coefList

◆ knots()

const std::vector<double>& RooCubicSplinePdf::knots ( ) const
inline

Definition at line 68 of file RooCubicSplinePdf.h.

68 { return _aux.knots(); }
RooCubicSplineKnot _aux
const std::vector< double > & knots() const

◆ knotSize()

unsigned RooCubicSplinePdf::knotSize ( ) const
inline

Definition at line 66 of file RooCubicSplinePdf.h.

66 { return _aux.size(); }
RooCubicSplineKnot _aux

◆ maxVal()

Double_t RooCubicSplinePdf::maxVal ( Int_t  code) const

Definition at line 290 of file RooCubicSplinePdf.cpp.

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 }
RooListProxy _coefList

◆ productAnalyticalIntegral()

std::complex< double > RooCubicSplinePdf::productAnalyticalIntegral ( Double_t  umin,
Double_t  umax,
Double_t  scale,
Double_t  offset,
const std::complex< double > &  z 
) const

Definition at line 250 of file RooCubicSplinePdf.cpp.

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 }
RooCubicSplineKnot::S_jk S_jk_sum(int i, const RooArgList &b) const
unsigned knotSize() const
std::complex< double > gaussIntegralE(bool left, const RooGaussModelAcceptance::M_n< 4U > &dM, const RooGaussModelAcceptance::K_n &K, double offset, double *sc) const
RooCubicSplineKnot _aux
double u(int i) const
RooListProxy _coefList

◆ u()

double RooCubicSplinePdf::u ( int  i) const
inline

Definition at line 67 of file RooCubicSplinePdf.h.

67 { return _aux.u(i); }
RooCubicSplineKnot _aux
double u(int i) const

Member Data Documentation

◆ _aux

RooCubicSplineKnot RooCubicSplinePdf::_aux
private

Definition at line 75 of file RooCubicSplinePdf.h.

◆ _coefList

RooListProxy RooCubicSplinePdf::_coefList
private

Definition at line 74 of file RooCubicSplinePdf.h.

◆ _x

RooRealProxy RooCubicSplinePdf::_x
private

Definition at line 73 of file RooCubicSplinePdf.h.


The documentation for this class was generated from the following files: