54 const Double_t RooKeysPdf::_nSigma = std::sqrt(-2. *
55 std::log(std::numeric_limits<Double_t>::epsilon()));
60 RooKeysPdf::RooKeysPdf() : _nEvents(0), _dataPts(0), _dataWgts(0), _weights(0), _sumWgt(0),
61 _mirrorLeft(kFALSE), _mirrorRight(kFALSE),
62 _asymLeft(kFALSE), _asymRight(kFALSE)
70 RooKeysPdf::RooKeysPdf(
const char *name,
const char *title,
71 RooAbsReal& x, RooDataSet& data,
72 Mirror mirror, Double_t rho) :
73 RooAbsPdf(name,title),
74 _x(
"x",
"observable",this,x),
79 _mirrorLeft(mirror==MirrorLeft || mirror==MirrorBoth || mirror==MirrorLeftAsymRight),
80 _mirrorRight(mirror==MirrorRight || mirror==MirrorBoth || mirror==MirrorAsymLeftRight),
81 _asymLeft(mirror==MirrorAsymLeft || mirror==MirrorAsymLeftRight || mirror==MirrorAsymBoth),
82 _asymRight(mirror==MirrorAsymRight || mirror==MirrorLeftAsymRight || mirror==MirrorAsymBoth),
85 snprintf(_varName, 128,
"%s", x.GetName());
86 RooAbsRealLValue& real= (RooRealVar&)(_x.arg());
89 _binWidth = (_hi-_lo)/(_nPoints-1);
99 RooKeysPdf::RooKeysPdf(
const char *name,
const char *title,
100 RooAbsReal& xpdf, RooRealVar& xdata, RooDataSet& data,
101 Mirror mirror, Double_t rho) :
102 RooAbsPdf(name,title),
103 _x(
"x",
"Observable",this,xpdf),
108 _mirrorLeft(mirror==MirrorLeft || mirror==MirrorBoth || mirror==MirrorLeftAsymRight),
109 _mirrorRight(mirror==MirrorRight || mirror==MirrorBoth || mirror==MirrorAsymLeftRight),
110 _asymLeft(mirror==MirrorAsymLeft || mirror==MirrorAsymLeftRight || mirror==MirrorAsymBoth),
111 _asymRight(mirror==MirrorAsymRight || mirror==MirrorLeftAsymRight || mirror==MirrorAsymBoth),
114 snprintf(_varName, 128,
"%s", xdata.GetName());
115 RooAbsRealLValue& real= (RooRealVar&)(xdata);
118 _binWidth = (_hi-_lo)/(_nPoints-1);
127 RooKeysPdf::RooKeysPdf(
const RooKeysPdf& other,
const char* name):
128 RooAbsPdf(other,name), _x(
"x",this,other._x), _nEvents(other._nEvents),
129 _dataPts(0), _dataWgts(0), _weights(0), _sumWgt(0),
130 _mirrorLeft( other._mirrorLeft ), _mirrorRight( other._mirrorRight ),
131 _asymLeft(other._asymLeft), _asymRight(other._asymRight),
134 snprintf(_varName, 128,
"%s", other._varName );
137 _binWidth = other._binWidth;
148 for (Int_t i= 0; i<_nPoints+1; i++)
149 _lookupTable[i]= other._lookupTable[i];
156 RooKeysPdf::~RooKeysPdf() {
174 inline bool operator()(
const struct Data& a,
const struct Data& b)
const
175 {
return a.x < b.x; }
178 void RooKeysPdf::LoadDataSet( RooDataSet& data) {
183 std::vector<Data> tmp;
184 tmp.reserve((1 + _mirrorLeft + _mirrorRight) * data.numEntries());
185 Double_t x0 = 0., x1 = 0., x2 = 0.;
188 RooRealVar& real = (RooRealVar&)(data.get()->operator[](_varName));
189 for (Int_t i = 0; i < data.numEntries(); ++i) {
191 const Double_t x = real.getVal();
192 const Double_t w = data.weight();
196 _sumWgt += Double_t(1 + _mirrorLeft + _mirrorRight) * w;
211 std::sort(tmp.begin(), tmp.end(), cmp());
214 _nEvents = tmp.size();
215 _dataPts =
new Double_t[_nEvents];
216 _dataWgts =
new Double_t[_nEvents];
217 for (
unsigned i = 0; i < tmp.size(); ++i) {
218 _dataPts[i] = tmp[i].x;
219 _dataWgts[i] = tmp[i].w;
223 std::vector<Data> tmp2;
227 Double_t meanv=x1/x0;
228 Double_t sigmav=std::sqrt(x2/x0-meanv*meanv);
229 Double_t h=std::pow(Double_t(4)/Double_t(3),0.2)*std::pow(_nEvents,-0.2)*_rho;
230 Double_t hmin=h*sigmav*std::sqrt(2.)/10;
231 Double_t norm=h*std::sqrt(sigmav)/(2.0*std::sqrt(3.0));
233 _weights=
new Double_t[_nEvents];
234 for(Int_t j=0;j<_nEvents;++j) {
235 _weights[j]=norm/std::sqrt(g(_dataPts[j],h*sigmav));
236 if (_weights[j]<hmin) _weights[j]=hmin;
243 for (Int_t i=0;i<_nPoints+1;++i) _lookupTable[i] = 0.;
244 for(Int_t j=0;j<_nEvents;++j) {
245 const Double_t xlo = std::min(_hi,
246 std::max(_lo, _dataPts[j] - _nSigma * _weights[j]));
247 const Double_t xhi = std::max(_lo,
248 std::min(_hi, _dataPts[j] + _nSigma * _weights[j]));
249 if (xlo >= xhi)
continue;
250 const Double_t chi2incr = _binWidth / _weights[j] / std::sqrt(2.);
251 const Double_t weightratio = _dataWgts[j] / _weights[j];
252 const Int_t binlo =
static_cast<Int_t
>(std::floor((xlo - _lo) / _binWidth));
253 const Int_t binhi =
static_cast<Int_t
>(_nPoints - std::floor((_hi - xhi) / _binWidth));
254 const Double_t x = (Double_t(_nPoints - binlo) * _lo +
255 Double_t(binlo) * _hi) / Double_t(_nPoints);
256 Double_t chi = (x - _dataPts[j]) / _weights[j] / std::sqrt(2.);
257 for (Int_t k = binlo; k <= binhi; ++k, chi += chi2incr) {
258 _lookupTable[k] += weightratio * std::exp(- chi * chi);
262 for(Int_t j=0;j<_nEvents;++j) {
263 const Double_t xlo = std::min(_hi,
264 std::max(_lo, 2. * _lo - _dataPts[j] + _nSigma * _weights[j]));
265 const Double_t xhi = std::max(_lo,
266 std::min(_hi, 2. * _lo - _dataPts[j] - _nSigma * _weights[j]));
267 if (xlo >= xhi)
continue;
268 const Double_t chi2incr = _binWidth / _weights[j] / std::sqrt(2.);
269 const Double_t weightratio = _dataWgts[j] / _weights[j];
270 const Int_t binlo =
static_cast<Int_t
>(std::floor((xlo - _lo) / _binWidth));
271 const Int_t binhi =
static_cast<Int_t
>(_nPoints - std::floor((_hi - xhi) / _binWidth));
272 const Double_t x = (Double_t(_nPoints - binlo) * _lo +
273 Double_t(binlo) * _hi) / Double_t(_nPoints);
274 Double_t chi = (x - (2. * _lo - _dataPts[j])) / _weights[j] / std::sqrt(2.);
275 for (Int_t k = binlo; k <= binhi; ++k, chi += chi2incr) {
276 _lookupTable[k] -= weightratio * std::exp(- chi * chi);
281 for(Int_t j=0;j<_nEvents;++j) {
282 const Double_t xlo = std::min(_hi,
283 std::max(_lo, 2. * _hi - _dataPts[j] + _nSigma * _weights[j]));
284 const Double_t xhi = std::max(_lo,
285 std::min(_hi, 2. * _hi - _dataPts[j] - _nSigma * _weights[j]));
286 if (xlo >= xhi)
continue;
287 const Double_t chi2incr = _binWidth / _weights[j] / std::sqrt(2.);
288 const Double_t weightratio = _dataWgts[j] / _weights[j];
289 const Int_t binlo =
static_cast<Int_t
>(std::floor((xlo - _lo) / _binWidth));
290 const Int_t binhi =
static_cast<Int_t
>(_nPoints - std::floor((_hi - xhi) / _binWidth));
291 const Double_t x = (Double_t(_nPoints - binlo) * _lo +
292 Double_t(binlo) * _hi) / Double_t(_nPoints);
293 Double_t chi = (x - (2. * _hi - _dataPts[j])) / _weights[j] / std::sqrt(2.);
294 for (Int_t k = binlo; k <= binhi; ++k, chi += chi2incr) {
295 _lookupTable[k] -= weightratio * std::exp(- chi * chi);
299 static const Double_t sqrt2pi(std::sqrt(2*TMath::Pi()));
300 for (Int_t i=0;i<_nPoints+1;++i)
301 _lookupTable[i] /= sqrt2pi * _sumWgt;
306 Double_t RooKeysPdf::evaluate()
const {
307 Int_t i = (Int_t)floor((Double_t(_x)-_lo)/_binWidth);
320 Double_t dx = (Double_t(_x)-(_lo+i*_binWidth))/_binWidth;
324 Double_t ret = (_lookupTable[i]+dx*(_lookupTable[i+1]-_lookupTable[i]));
329 Int_t RooKeysPdf::getAnalyticalIntegral(
330 RooArgSet& allVars, RooArgSet& analVars,
const char* )
const
332 if (matchArgs(allVars, analVars, _x))
return 1;
336 Double_t RooKeysPdf::analyticalIntegral(Int_t code,
const char* rangeName)
const
338 R__ASSERT(1 == code);
341 const Double_t xmin = std::max(_lo, _x.min(rangeName));
342 const Double_t xmax = std::min(_hi, _x.max(rangeName));
343 const Int_t imin = (Int_t)floor((xmin - _lo) / _binWidth);
344 const Int_t imax = std::min((Int_t)floor((xmax - _lo) / _binWidth),
349 sum += _lookupTable[imin + 1] + _lookupTable[imax];
350 for (Int_t i = imin + 2; i < imax; ++i)
351 sum += 2. * _lookupTable[i];
352 sum *= _binWidth * 0.5;
354 const Double_t dxmin = (xmin - (_lo + imin * _binWidth)) / _binWidth;
355 const Double_t dxmax = (xmax - (_lo + imax * _binWidth)) / _binWidth;
358 sum += _binWidth * (1. - dxmin) * 0.5 * (_lookupTable[imin + 1] +
359 _lookupTable[imin] + dxmin *
360 (_lookupTable[imin + 1] - _lookupTable[imin]));
362 sum += _binWidth * dxmax * 0.5 * (_lookupTable[imax] +
363 _lookupTable[imax] + dxmax *
364 (_lookupTable[imax + 1] - _lookupTable[imax]));
365 }
else if (imin == imax) {
367 sum += _binWidth * (dxmax - dxmin) * 0.5 * (
368 _lookupTable[imin] + dxmin *
369 (_lookupTable[imin + 1] - _lookupTable[imin]) +
370 _lookupTable[imax] + dxmax *
371 (_lookupTable[imax + 1] - _lookupTable[imax]));
376 Int_t RooKeysPdf::getMaxVal(
const RooArgSet& vars)
const
378 if (vars.contains(*_x.absArg()))
return 1;
382 Double_t RooKeysPdf::maxVal(Int_t code)
const
384 R__ASSERT(1 == code);
385 Double_t max = -std::numeric_limits<Double_t>::max();
386 for (Int_t i = 0; i <= _nPoints; ++i)
387 if (max < _lookupTable[i]) max = _lookupTable[i];
393 Double_t RooKeysPdf::g(Double_t x,Double_t sigmav)
const {
397 Double_t* it = std::lower_bound(_dataPts, _dataPts + _nEvents,
398 x - _nSigma * sigmav);
399 if (it >= (_dataPts + _nEvents))
return 0.;
400 Double_t* iend = std::upper_bound(it, _dataPts + _nEvents,
401 x + _nSigma * sigmav);
402 for ( ; it < iend; ++it) {
403 const Double_t r = (x - *it) / sigmav;
404 y += std::exp(-0.5 * r * r);
407 static const Double_t sqrt2pi(std::sqrt(2*TMath::Pi()));
408 return y/(sigmav*sqrt2pi*_nEvents);