Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ParamHistFunc.h
Go to the documentation of this file.
1// @(#)root/roostats:$Id: cranmer $
2// Author: George Lewis, Kyle Cranmer
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11
12#ifndef ROO_PARAMHISTFUNC
13#define ROO_PARAMHISTFUNC
14
15#include "RooAbsReal.h"
16#include "RooListProxy.h"
17#include "RooObjCacheManager.h"
18#include "RooDataHist.h"
19
20// Forward Declarations
21class RooRealVar;
22class RooWorkspace;
23
24class ParamHistFunc : public RooAbsReal {
25public:
26
28 ParamHistFunc(const char *name, const char *title, const RooArgList& vars, const RooArgList& paramSet );
29 ParamHistFunc(const char *name, const char *title, const RooArgList& vars, const RooArgList& paramSet, const TH1* hist );
30
31 ParamHistFunc(const ParamHistFunc& other, const char *name = nullptr);
32 TObject* clone(const char* newname) const override { return new ParamHistFunc(*this, newname); }
33
34 const RooArgList& paramList() const { return _paramSet ; }
35
36 Int_t numBins() const { return _dataSet.numEntries(); } // Number of bins (called numEntries in RooDataHist)
37
38 void setParamConst( Int_t, bool=true );
39 void setConstant(bool constant);
40
41 void setShape(TH1* shape);
42
43 RooAbsReal& getParameter() const ;
44 RooAbsReal& getParameter( Int_t masterIdx ) const ;
45
46 const RooArgSet* get(Int_t masterIdx) const { return _dataSet.get( masterIdx ) ; }
47 const RooArgSet* get(const RooArgSet& coord) const { return _dataSet.get( coord ) ; }
48
49 RooDataHist const& dataHist() const { return _dataSet; }
50
51 double binVolume() const { return _dataSet.binVolume(); }
52
53 bool forceAnalyticalInt(const RooAbsArg&) const override { return true ; }
54
55 Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet,const char* rangeName=nullptr) const override;
56 double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=nullptr) const override;
57
58 static RooArgList createParamSet(RooWorkspace& w, const std::string&, const RooArgList& Vars);
59 static RooArgList createParamSet(RooWorkspace& w, const std::string&, const RooArgList& Vars, double, double);
60 static RooArgList createParamSet(const std::string&, Int_t, double, double);
61
62 std::list<double>* binBoundaries(RooAbsRealLValue& /*obs*/, double /*xlo*/, double /*xhi*/) const override;
63 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override;
64 bool isBinnedDistribution(const RooArgSet& obs) const override { return _dataVars.overlaps(obs); }
65
66 const RooArgList& dataVars() const { return _dataVars; }
67
68protected:
69
71 public:
74 ret.add(_lowIntList);
75 ret.add(_highIntList);
76 return ret ;
77 }
81 // will want std::vector<RooRealVar*> for low and high also
82 } ;
83 mutable RooObjCacheManager _normIntMgr ; ///<! The integration cache manager
84
85 RooListProxy _dataVars; ///< The RooRealVars
86 RooListProxy _paramSet ; ///< interpolation parameters
87
89 struct NumBins {
91 NumBins(int nx, int ny, int nz) : x{nx}, y{ny}, z{nz}, xy{x*y}, xz{x*z}, yz{y*z}, xyz{xy*z} {}
92 int x = 0;
93 int y = 0;
94 int z = 0;
95 int xy = 0;
96 int xz = 0;
97 int yz = 0;
98 int xyz = 0;
99 };
102
103 Int_t getCurrentBin() const;
104 Int_t addParamSet( const RooArgList& params );
105 static Int_t GetNumBins( const RooArgSet& vars );
106 double evaluate() const override;
107 void doEval(RooFit::EvalContext &) const override;
108
109 private:
110 static NumBins getNumBinsPerDim(RooArgSet const& vars);
111
113};
114
115#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
char name[80]
Definition TGX11.cxx:110
RooArgList containedArgs(Action) override
A class which maps the current values of a RooRealVar (or a set of RooRealVars) to one of a number of...
static NumBins getNumBinsPerDim(RooArgSet const &vars)
const RooArgSet * get(Int_t masterIdx) const
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
RooDataHist _dataSet
static Int_t GetNumBins(const RooArgSet &vars)
TObject * clone(const char *newname) const override
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
void setConstant(bool constant)
const RooArgList & paramList() const
Int_t getCurrentBin() const
Get the index of the gamma parameter associated with the current bin.
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise that all integrals can be handled internally.
RooObjCacheManager _normIntMgr
! The integration cache manager
Int_t numBins() const
double evaluate() const override
Find the bin corresponding to the current value of the observable, and evaluate the associated parame...
const RooArgSet * get(const RooArgSet &coord) const
void doEval(RooFit::EvalContext &) const override
Find all bins corresponding to the values of the observables in evalData, and evaluate the associated...
RooAbsReal & getParameter() const
bool forceAnalyticalInt(const RooAbsArg &) const override
Int_t addParamSet(const RooArgList &params)
NumBins _numBinsPerDim
static RooArgList createParamSet(RooWorkspace &w, const std::string &, const RooArgList &Vars)
Create the list of RooRealVar parameters which represent the height of the histogram bins.
void setParamConst(Int_t, bool=true)
const RooArgList & dataVars() const
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
RooListProxy _paramSet
interpolation parameters
RooDataHist const & dataHist() const
void setShape(TH1 *shape)
double binVolume() const
RooListProxy _dataVars
The RooRealVars.
bool isBinnedDistribution(const RooArgSet &obs) const override
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
Abstract base class for objects to be stored in RooAbsCache cache manager objects.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
bool overlaps(Iterator_t otherCollBegin, Iterator_t otherCollEnd) const
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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:24
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
double binVolume(std::size_t i) const
Return bin volume of i-th bin.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
Implementation of a RooCacheManager<RooAbsCacheElement> that specializes in the storage of cache elem...
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Persistable container for RooFit projects.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
Mother of all ROOT objects.
Definition TObject.h:41
NumBins(int nx, int ny, int nz)