Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ProposalHelper.cxx
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/** \class RooStats::ProposalHelper
13 \ingroup Roostats
14*/
15
16#include "Rtypes.h"
20#include "RooArgSet.h"
21#include "RooDataSet.h"
22#include "RooAbsPdf.h"
23#include "RooAddPdf.h"
24#include "RooNDKeysPdf.h"
25#include "RooUniform.h"
26#include "RooMsgService.h"
27#include "RooRealVar.h"
28#include "RooMultiVarGaussian.h"
29#include "RooConstVar.h"
30#include "TString.h"
31
32#include <map>
33
34namespace RooStats {
35 class ProposalFunction;
36}
37
39
40using namespace RooFit;
41using namespace RooStats;
42using std::endl;
43
44//static const double DEFAULT_UNI_FRAC = 0.10;
45static const double DEFAULT_CLUES_FRAC = 0.20;
46//static const double SIGMA_RANGE_DIVISOR = 6;
47static const double SIGMA_RANGE_DIVISOR = 5;
48//static const Int_t DEFAULT_CACHE_SIZE = 100;
49//static const Option_t* CLUES_OPTIONS = "a";
50
51////////////////////////////////////////////////////////////////////////////////
52
53ProposalHelper::ProposalHelper() : fPdfProp(new PdfProposal()), fSigmaRangeDivisor(SIGMA_RANGE_DIVISOR) {}
54
55////////////////////////////////////////////////////////////////////////////////
56
58{
59 if (fPdf == nullptr)
60 CreatePdf();
61 RooArgList components;
62 RooArgList coeffs;
63 if (fCluesPdf == nullptr)
65 if (fCluesPdf != nullptr) {
66 if (fCluesFrac < 0)
68 printf("added clues from dataset %s with fraction %g\n",
70 components.add(*fCluesPdf);
71 coeffs.add(RooConst(fCluesFrac));
72 }
73 if (fUniFrac > 0.) {
75 components.add(*fUniformPdf);
76 coeffs.add(RooConst(fUniFrac));
77 }
78 components.add(*fPdf);
79 RooAddPdf* addPdf = new RooAddPdf("proposalFunction", "Proposal Density",
80 components, coeffs);
81 fPdfProp->SetPdf(*addPdf);
82 fPdfProp->SetOwnsPdf(true);
83 if (fCacheSize > 0)
85 fOwnsPdfProp = false;
86 return fPdfProp;
87}
88
89////////////////////////////////////////////////////////////////////////////////
90
92{
93 if (fVars == nullptr) {
94 coutE(InputArguments) << "ProposalHelper::CreatePdf(): " <<
95 "Variables to create proposal function for are not set." << endl;
96 return;
97 }
98 RooArgList xVec{};
99 RooArgList muVec{};
101 for (auto *r : static_range_cast<RooRealVar *> (*fVars)){
102 xVec.add(*r);
103 TString cloneName = TString::Format("%s%s", "mu__", r->GetName());
104 clone = static_cast<RooRealVar*>(r->clone(cloneName.Data()));
105 muVec.add(*clone);
106 if (fUseUpdates)
108 }
109 if (fCovMatrix == nullptr)
110 CreateCovMatrix(xVec);
111 fPdf = new RooMultiVarGaussian("mvg", "MVG Proposal", xVec, muVec,
112 *fCovMatrix);
113}
114
115////////////////////////////////////////////////////////////////////////////////
116
118{
119 Int_t size = xVec.size();
121 RooRealVar* r;
122 for (Int_t i = 0; i < size; i++) {
123 r = static_cast<RooRealVar*>(xVec.at(i));
124 double range = r->getMax() - r->getMin();
125 (*fCovMatrix)(i,i) = range / fSigmaRangeDivisor;
126 }
127}
128
129////////////////////////////////////////////////////////////////////////////////
130
132{
133 if (fClues != nullptr) {
134 if (fCluesOptions == nullptr) {
135 fCluesPdf = new RooNDKeysPdf("cluesPdf", "Clues PDF", *fVars, *fClues);
136 } else {
137 fCluesPdf = new RooNDKeysPdf("cluesPdf", "Clues PDF", *fVars, *fClues, fCluesOptions);
138 }
139 }
140}
141
142////////////////////////////////////////////////////////////////////////////////
143
145{
146 fUniformPdf = new RooUniform("uniform", "Uniform Proposal PDF",
147 RooArgSet(*fVars));
148}
static const double SIGMA_RANGE_DIVISOR
static const double DEFAULT_CLUES_FRAC
TObject * clone(const char *newname) const override
Definition RooChi2Var.h:9
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
TMatrixTSym< Double_t > TMatrixDSym
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Multivariate Gaussian p.d.f.
Generic N-dimensional implementation of a kernel estimation p.d.f.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
PdfProposal is a concrete implementation of the ProposalFunction interface.
Definition PdfProposal.h:30
virtual void AddMapping(RooRealVar &proposalParam, RooAbsReal &update)
specify a mapping between a parameter of the proposal function and a parameter of interest.
virtual void SetCacheSize(Int_t size)
Set how many points to generate each time we propose from a new point Default (and minimum) is 1.
Definition PdfProposal.h:79
virtual void SetOwnsPdf(bool ownsPdf)
set whether we own the PDF that serves as the proposal density function By default,...
Definition PdfProposal.h:91
virtual void SetPdf(RooAbsPdf &pdf)
Set the PDF to be the proposal density function.
Definition PdfProposal.h:49
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
double fUniFrac
what fraction of the PDF integral is uniform
RooDataSet * fClues
data set of clues
bool fUseUpdates
whether to set updates for proposal params in PdfProposal
void CreateCovMatrix(RooArgList &xVec)
double fCluesFrac
what fraction of the PDF integral comes from clues
virtual ProposalFunction * GetProposalFunction()
Get the ProposalFunction that we've been designing.
const Option_t * fCluesOptions
option string for clues RooNDKeysPdf
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.
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
double fSigmaRangeDivisor
range divisor to get sigma for each variable
RooAbsPdf * fCluesPdf
proposal dens. func. with clues for certain points
Flat p.d.f.
Definition RooUniform.h:24
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
RooConstVar & RooConst(double val)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition JSONIO.h:26
@ InputArguments
Namespace for the RooStats classes.
Definition Asimov.h:19