Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooKeysPdf.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooKeysPdf.h,v 1.10 2007/05/11 09:13:07 verkerke Exp $
5 * Authors: *
6 * GR, Gerhard Raven, UC San Diego, raven@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
9 * *
10 * Copyright (c) 2000-2005, Regents of the University of California *
11 * and Stanford University. All rights reserved. *
12 * *
13 * Redistribution and use in source and binary forms, *
14 * with or without modification, are permitted according to the terms *
15 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
16 *****************************************************************************/
17#ifndef ROO_KEYS
18#define ROO_KEYS
19
20#include "RooAbsPdf.h"
21#include "RooRealProxy.h"
22
23class RooRealVar;
24
25class RooKeysPdf : public RooAbsPdf {
26public:
31 RooKeysPdf() ;
32 RooKeysPdf(const char *name, const char *title,
34 double rho=1);
35 RooKeysPdf(const char *name, const char *title,
37 double rho=1);
38 RooKeysPdf(const RooKeysPdf& other, const char* name=nullptr);
39 TObject* clone(const char* newname) const override {return new RooKeysPdf(*this,newname); }
40 ~RooKeysPdf() override;
41
43 const char* rangeName = nullptr) const override;
44 double analyticalIntegral(Int_t code, const char* rangeName = nullptr) const override;
45 Int_t getMaxVal(const RooArgSet& vars) const override;
46 double maxVal(Int_t code) const override;
47
49
50protected:
51
53 double evaluate() const override;
54
55private:
56 // how far you have to go out in a Gaussian until it is smaller than the
57 // machine precision
58 static const double _nSigma; //!
59
61 double *_dataPts; //[_nEvents]
62 double *_dataWgts; //[_nEvents]
63 double *_weights; //[_nEvents]
64 double _sumWgt ;
65
66 constexpr static int _nPoints{1000};
68
69 double g(double x,double sigma) const;
70
73
74 // cached info on variable
76 double _lo, _hi, _binWidth;
77 double _rho;
78
79 ClassDefOverride(RooKeysPdf,2) // One-dimensional non-parametric kernel estimation p.d.f.
80};
81
82#endif
#define g(i)
Definition RSha256.hxx:105
char Char_t
Definition RtypesCore.h:37
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
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
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
Class RooKeysPdf implements a one-dimensional kernel estimation p.d.f which model the distribution of...
Definition RooKeysPdf.h:25
double _binWidth
Definition RooKeysPdf.h:76
double _sumWgt
Definition RooKeysPdf.h:64
static constexpr int _nPoints
Definition RooKeysPdf.h:66
double * _dataWgts
Definition RooKeysPdf.h:62
RooKeysPdf()
coverity[UNINIT_CTOR]
double _lookupTable[_nPoints+1]
Definition RooKeysPdf.h:67
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
bool _mirrorRight
Definition RooKeysPdf.h:71
double _rho
Definition RooKeysPdf.h:77
double _hi
Definition RooKeysPdf.h:76
Char_t _varName[128]
Definition RooKeysPdf.h:75
Int_t _nEvents
Definition RooKeysPdf.h:60
void LoadDataSet(RooDataSet &data)
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooRealProxy _x
Definition RooKeysPdf.h:52
TObject * clone(const char *newname) const override
Definition RooKeysPdf.h:39
double * _weights
Definition RooKeysPdf.h:63
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
bool _mirrorLeft
Definition RooKeysPdf.h:71
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
~RooKeysPdf() override
double _lo
Definition RooKeysPdf.h:76
@ MirrorAsymRight
Definition RooKeysPdf.h:29
@ MirrorLeftAsymRight
Definition RooKeysPdf.h:29
@ MirrorAsymLeftRight
Definition RooKeysPdf.h:28
double * _dataPts
Definition RooKeysPdf.h:61
bool _asymRight
Definition RooKeysPdf.h:72
static const double _nSigma
Definition RooKeysPdf.h:58
bool _asymLeft
Definition RooKeysPdf.h:72
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise capability to determine maximum value of function for given set of observables.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
Mother of all ROOT objects.
Definition TObject.h:41
const Double_t sigma
Double_t x[n]
Definition legend1.C:17