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 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. 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#ifndef ROO_KEYS
17#define ROO_KEYS
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21
22class RooRealVar;
23
24class RooKeysPdf : public RooAbsPdf {
25public:
30 RooKeysPdf() ;
31 RooKeysPdf(const char *name, const char *title,
33 double rho=1);
34 RooKeysPdf(const char *name, const char *title,
36 double rho=1);
37 RooKeysPdf(const RooKeysPdf& other, const char* name=nullptr);
38 TObject* clone(const char* newname) const override {return new RooKeysPdf(*this,newname); }
39 ~RooKeysPdf() override;
40
42 const char* rangeName = nullptr) const override;
43 double analyticalIntegral(Int_t code, const char* rangeName = nullptr) const override;
44 Int_t getMaxVal(const RooArgSet& vars) const override;
45 double maxVal(Int_t code) const override;
46
48
49protected:
50
52 double evaluate() const override;
53
54private:
55 // how far you have to go out in a Gaussian until it is smaller than the
56 // machine precision
57 static const double _nSigma; //!
58
60 double *_dataPts = nullptr; //[_nEvents]
61 double *_dataWgts = nullptr; //[_nEvents]
62 double *_weights = nullptr; //[_nEvents]
63 double _sumWgt = 0.0;
64
65 constexpr static int _nPoints{1000};
67
68 double g(double x,double sigma) const;
69
70 bool _mirrorLeft = false;
71 bool _mirrorRight = false;
72 bool _asymLeft = false;
73 bool _asymRight = false;
74
75 // cached info on variable
77 double _lo, _hi, _binWidth;
78 double _rho;
79
80 ClassDefOverride(RooKeysPdf,2) // One-dimensional non-parametric kernel estimation p.d.f.
81};
82
83#endif
#define g(i)
Definition RSha256.hxx:105
char Char_t
Definition RtypesCore.h:37
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
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
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Container class to hold unbinned data.
Definition RooDataSet.h:33
Class RooKeysPdf implements a one-dimensional kernel estimation p.d.f which model the distribution of...
Definition RooKeysPdf.h:24
double _binWidth
Definition RooKeysPdf.h:77
double _sumWgt
Definition RooKeysPdf.h:63
static constexpr int _nPoints
Definition RooKeysPdf.h:65
double * _dataWgts
Definition RooKeysPdf.h:61
RooKeysPdf()
coverity[UNINIT_CTOR]
double _lookupTable[_nPoints+1]
Definition RooKeysPdf.h:66
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:78
double _hi
Definition RooKeysPdf.h:77
Char_t _varName[128]
Definition RooKeysPdf.h:76
Int_t _nEvents
Definition RooKeysPdf.h:59
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:51
TObject * clone(const char *newname) const override
Definition RooKeysPdf.h:38
double * _weights
Definition RooKeysPdf.h:62
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:70
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
~RooKeysPdf() override
double _lo
Definition RooKeysPdf.h:77
@ MirrorAsymRight
Definition RooKeysPdf.h:28
@ MirrorLeftAsymRight
Definition RooKeysPdf.h:28
@ MirrorAsymLeftRight
Definition RooKeysPdf.h:27
double * _dataPts
Definition RooKeysPdf.h:60
bool _asymRight
Definition RooKeysPdf.h:73
static const double _nSigma
Definition RooKeysPdf.h:57
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.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41
const Double_t sigma
Double_t x[n]
Definition legend1.C:17