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: " <<
48 size << ". Cache size unchanged." << std::endl;
49 }
50
51 virtual void SetUpdateProposalParameters(bool updateParams)
52 { fUseUpdates = updateParams; }
53
54 virtual void SetVariables(RooArgList& vars)
55 { fVars = &vars; }
56
57 virtual void SetVariables(const RooArgList& vars)
58 { fVars = new RooArgList(vars); fOwnsVars = true; }
59
60 /// set what fraction of the proposal density function should come from
61 /// a uniform proposal distribution
62 virtual void SetUniformFraction(double uniFrac) { fUniFrac = uniFrac; }
63
64 /// set what fraction of the proposal density function should come from
65 /// the bank of clues
66 virtual void SetCluesFraction(double cluesFrac) { fCluesFrac = cluesFrac; }
67
68 /// set the covariance matrix to use for a multi-variate Gaussian proposal
69 virtual void SetCovMatrix(const TMatrixDSym& covMatrix)
70 { fCovMatrix = new TMatrixDSym(covMatrix); }
71
72 /// set what divisor we will use when dividing the range of a variable to
73 /// determine the width of the proposal function for each dimension
74 /// e.g. divisor = 6 for sigma = 1/6th
75 virtual void SetWidthRangeDivisor(double divisor)
76 { if (divisor > 0.) fSigmaRangeDivisor = divisor; }
77
78 /// set the option string to pass to the RooNDKeysPdf constructor
79 /// if the bank of clues pdf is being automatically generated by this
80 /// ProposalHelper
81 virtual void SetCluesOptions(const Option_t* options)
82 { if (options != nullptr) fCluesOptions = options; }
83
84 virtual void SetVariables(RooArgSet& vars)
85 {
86 RooArgList* argList = new RooArgList(vars);
87 SetVariables(*argList);
88 fOwnsVars = true;
89 }
90
91 ~ProposalHelper() override
92 {
93 if (fOwnsPdfProp) delete fPdfProp;
94 if (fOwnsPdf) delete fPdf;
95 if (fOwnsCluesPdf) delete fCluesPdf;
96 if (fOwnsVars) delete fVars;
97 delete fCovMatrix;
98 delete fUniformPdf;
99 }
100
101 protected:
102 RooAbsPdf* fPdf; ///< the main proposal density function
103 RooAbsPdf* fCluesPdf; ///< proposal dens. func. with clues for certain points
104 RooAbsPdf* fUniformPdf; ///< uniform proposal dens. func.
105 RooDataSet* fClues; ///< data set of clues
106 TMatrixDSym* fCovMatrix; ///< covariance matrix for multi var gaussian pdf
107 PdfProposal* fPdfProp; ///< the PdfProposal we are (probably) going to return
108 RooArgList* fVars; ///< the RooRealVars to generate proposals for
109 Int_t fCacheSize; ///< for generating proposals from PDFs
110 double fSigmaRangeDivisor; ///< range divisor to get sigma for each variable
111 double fUniFrac; ///< what fraction of the PDF integral is uniform
112 double fCluesFrac; ///< what fraction of the PDF integral comes from clues
113 bool fOwnsPdfProp; ///< whether we own the PdfProposal; equivalent to:
114 ///< !(whether we have returned it in GetProposalFunction)
115 bool fOwnsPdf; ///< whether we created (and own) the main pdf
116 bool fOwnsCluesPdf; ///< whether we created (and own) the clues pdf
117 bool fOwnsVars; ///< whether we own fVars
118 bool fUseUpdates; ///< whether to set updates for proposal params in PdfProposal
119 const Option_t* fCluesOptions; ///< option string for clues RooNDKeysPdf
120
121 void CreatePdf();
122 void CreateCluesPdf();
123 void CreateUniformPdf();
124 void CreateCovMatrix(RooArgList& xVec);
125
127 };
128}
129#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:341
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:55
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
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:
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 Asimov.h:19