Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooParametricStepFunction.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooParametricStepFunction.h,v 1.5 2007/05/11 09:13:07 verkerke Exp $
5 * Authors: *
6 * Aaron Roodman, Stanford Linear Accelerator Center, Stanford University *
7 * *
8 * Copyright (c) 2000-2005, Stanford University. All rights reserved. *
9 *
10 * *
11 * Redistribution and use in source and binary forms, *
12 * with or without modification, are permitted according to the terms *
13 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14 *****************************************************************************/
15#ifndef ROO_PARAMETRIC_STEP_FUNCTION
16#define ROO_PARAMETRIC_STEP_FUNCTION
17
18#include "TArrayD.h"
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21#include "RooListProxy.h"
22
23class RooRealVar;
24class RooArgList ;
25
27public:
28
30
31 RooParametricStepFunction(const char *name, const char *title,
32 RooAbsReal& x, const RooArgList& coefList, TArrayD& limits, Int_t nBins=1) ;
33
34 RooParametricStepFunction(const RooParametricStepFunction& other, const char* name = nullptr);
35 TObject* clone(const char* newname) const override { return new RooParametricStepFunction(*this, newname); }
36
37 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
38 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
39 Int_t getnBins() const { return _nBins; }
40 double* getLimits() { return _limits.GetArray(); }
41
42protected:
43
44 double lastBinValue() const ;
45
50
51 double evaluate() const override;
52
53 ClassDefOverride(RooParametricStepFunction,1) // Parametric Step Function Pdf
54};
55
56#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
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
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
The Parametric Step Function PDF is a binned distribution whose parameters are the heights of each bi...
TObject * clone(const char *newname) const override
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
Array of doubles (64 bits per element).
Definition TArrayD.h:27
const Double_t * GetArray() const
Definition TArrayD.h:43
Mother of all ROOT objects.
Definition TObject.h:41
Double_t x[n]
Definition legend1.C:17