Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ProposalHelper.h
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Authors: Kevin Belasco 7/22/2009
3// Authors: Kyle Cranmer 7/22/2009
4/*************************************************************************
5 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#ifndef RooStats_ProposalHelper
13#define RooStats_ProposalHelper
14
15#include "Rtypes.h"
19#include "RooArgSet.h"
20#include "RooMsgService.h"
21#include "RooRealVar.h"
22#include "TMatrixDSym.h"
23#include "TObject.h"
24
25
26
27namespace RooStats {
28
29 class ProposalHelper : public TObject {
30
31 public:
33
34 /// Set the PDF to be the proposal density function
35 virtual void SetPdf(RooAbsPdf& pdf) { fPdf = &pdf; }
36 /// Set the bank of clues to add to the current proposal density function
37 virtual void SetClues(RooDataSet& clues) { fClues = &clues; }
38
39 /// Get the ProposalFunction that we've been designing
41
42 virtual void SetCacheSize(Int_t size)
43 {
44 if (size > 0) {
46 } else {
47 coutE(Eval) << "Warning: Requested non-positive cache size: " << size << ". Cache size unchanged."
48 << std::endl;
49 }
50 }
51
52 virtual void SetUpdateProposalParameters(bool updateParams)
53 { fUseUpdates = updateParams; }
54
55 virtual void SetVariables(RooArgList& vars)
56 { fVars = &vars; }
57
58 virtual void SetVariables(const RooArgList& vars)
59 { fVars = new RooArgList(vars); fOwnsVars = true; }
60
61 /// set what fraction of the proposal density function should come from
62 /// a uniform proposal distribution
63 virtual void SetUniformFraction(double uniFrac) { fUniFrac = uniFrac; }
64
65 /// set what fraction of the proposal density function should come from
66 /// the bank of clues
67 virtual void SetCluesFraction(double cluesFrac) { fCluesFrac = cluesFrac; }
68
69 /// set the covariance matrix to use for a multi-variate Gaussian proposal
70 virtual void SetCovMatrix(const TMatrixDSym& covMatrix)
71 { fCovMatrix = new TMatrixDSym(covMatrix); }
72
73 /// set what divisor we will use when dividing the range of a variable to
74 /// determine the width of the proposal function for each dimension
75 /// e.g. divisor = 6 for sigma = 1/6th
76 virtual void SetWidthRangeDivisor(double divisor)
77 { if (divisor > 0.) fSigmaRangeDivisor = divisor; }
78
79 /// set the option string to pass to the RooNDKeysPdf constructor
80 /// if the bank of clues pdf is being automatically generated by this
81 /// ProposalHelper
82 virtual void SetCluesOptions(const Option_t* options)
83 { if (options != nullptr) fCluesOptions = options; }
84
85 virtual void SetVariables(RooArgSet& vars)
86 {
87 RooArgList* argList = new RooArgList(vars);
88 SetVariables(*argList);
89 fOwnsVars = true;
90 }
91
92 ~ProposalHelper() override
93 {
94 if (fOwnsPdfProp) delete fPdfProp;
95 if (fOwnsPdf) delete fPdf;
96 if (fOwnsCluesPdf) delete fCluesPdf;
97 if (fOwnsVars) delete fVars;
98 delete fCovMatrix;
99 delete fUniformPdf;
100 }
101
102 protected:
103 RooAbsPdf *fPdf = nullptr; ///< the main proposal density function
104 RooAbsPdf *fCluesPdf = nullptr; ///< proposal dens. func. with clues for certain points
105 RooAbsPdf *fUniformPdf = nullptr; ///< uniform proposal dens. func.
106 RooDataSet *fClues = nullptr; ///< data set of clues
107 TMatrixDSym *fCovMatrix = nullptr; ///< covariance matrix for multi var gaussian pdf
108 PdfProposal* fPdfProp = nullptr; ///< the PdfProposal we are (probably) going to return
109 RooArgList *fVars = nullptr; ///< the RooRealVars to generate proposals for
110 Int_t fCacheSize = -1; ///< for generating proposals from PDFs
111 double fSigmaRangeDivisor; ///< range divisor to get sigma for each variable
112 double fUniFrac = -1; ///< what fraction of the PDF integral is uniform
113 double fCluesFrac = -1; ///< what fraction of the PDF integral comes from clues
114 bool fOwnsPdfProp = true; ///< whether we own the PdfProposal; equivalent to:
115 ///< !(whether we have returned it in GetProposalFunction)
116 bool fOwnsPdf = false; ///< whether we created (and own) the main pdf
117 bool fOwnsCluesPdf = false; ///< whether we created (and own) the clues pdf
118 bool fOwnsVars = false; ///< whether we own fVars
119 bool fUseUpdates = false; ///< whether to set updates for proposal params in PdfProposal
120 const Option_t *fCluesOptions = nullptr; ///< option string for clues RooNDKeysPdf
121
122 void CreatePdf();
123 void CreateCluesPdf();
124 void CreateUniformPdf();
125 void CreateCovMatrix(RooArgList& xVec);
126
128 };
129}
130#endif
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutE(a)
const char Option_t
Definition RtypesCore.h:66
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
TMatrixTSym< Double_t > TMatrixDSym
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
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 unbinned data.
Definition RooDataSet.h:34
PdfProposal is a concrete implementation of the ProposalFunction interface.
Definition PdfProposal.h:30
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
virtual void SetUniformFraction(double uniFrac)
set what fraction of the proposal density function should come from a uniform proposal distribution
double fUniFrac
what fraction of the PDF integral is uniform
virtual void SetCovMatrix(const TMatrixDSym &covMatrix)
set the covariance matrix to use for a multi-variate Gaussian proposal
RooDataSet * fClues
data set of clues
bool fOwnsCluesPdf
whether we created (and own) the clues pdf
bool fUseUpdates
whether to set updates for proposal params in PdfProposal
void CreateCovMatrix(RooArgList &xVec)
bool fOwnsVars
whether we own fVars
virtual void SetVariables(RooArgSet &vars)
double fCluesFrac
what fraction of the PDF integral comes from clues
virtual void SetPdf(RooAbsPdf &pdf)
Set the PDF to be the proposal density function.
virtual ProposalFunction * GetProposalFunction()
Get the ProposalFunction that we've been designing.
const Option_t * fCluesOptions
option string for clues RooNDKeysPdf
virtual void SetVariables(RooArgList &vars)
virtual void SetVariables(const RooArgList &vars)
RooAbsPdf * fPdf
the main proposal density function
bool fOwnsPdfProp
whether we own the PdfProposal; equivalent to: !(whether we have returned it in GetProposalFunction)
Int_t fCacheSize
for generating proposals from PDFs
RooAbsPdf * fUniformPdf
uniform proposal dens. func.
virtual void SetClues(RooDataSet &clues)
Set the bank of clues to add to the current proposal density function.
RooArgList * fVars
the RooRealVars to generate proposals for
TMatrixDSym * fCovMatrix
covariance matrix for multi var gaussian pdf
PdfProposal * fPdfProp
the PdfProposal we are (probably) going to return
virtual void SetCacheSize(Int_t size)
virtual void SetWidthRangeDivisor(double divisor)
set what divisor we will use when dividing the range of a variable to determine the width of the prop...
virtual void SetUpdateProposalParameters(bool updateParams)
double fSigmaRangeDivisor
range divisor to get sigma for each variable
virtual void SetCluesOptions(const Option_t *options)
set the option string to pass to the RooNDKeysPdf constructor if the bank of clues pdf is being autom...
virtual void SetCluesFraction(double cluesFrac)
set what fraction of the proposal density function should come from the bank of clues
RooAbsPdf * fCluesPdf
proposal dens. func. with clues for certain points
bool fOwnsPdf
whether we created (and own) the main pdf
Mother of all ROOT objects.
Definition TObject.h:41
Namespace for the RooStats classes.
Definition CodegenImpl.h:58