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
38
39using namespace RooFit;
40using namespace RooStats;
41using std::endl;
42
43//static const double DEFAULT_UNI_FRAC = 0.10;
44static const double DEFAULT_CLUES_FRAC = 0.20;
45//static const double SIGMA_RANGE_DIVISOR = 6;
46static const double SIGMA_RANGE_DIVISOR = 5;
47//static const Int_t DEFAULT_CACHE_SIZE = 100;
48//static const Option_t* CLUES_OPTIONS = "a";
49
50////////////////////////////////////////////////////////////////////////////////
51
53
54////////////////////////////////////////////////////////////////////////////////
55
57{
58 if (fPdf == nullptr)
59 CreatePdf();
60 RooArgList components;
61 RooArgList coeffs;
62 if (fCluesPdf == nullptr)
64 if (fCluesPdf != nullptr) {
65 if (fCluesFrac < 0) {
67 }
68 oocoutI(nullptr,InputArguments) << "added clues from dataset " << fClues->GetName() << " with fraction " << fCluesFrac << std::endl;
69 components.add(*fCluesPdf);
70 coeffs.add(RooConst(fCluesFrac));
71 }
72 if (fUniFrac > 0.) {
74 components.add(*fUniformPdf);
75 coeffs.add(RooConst(fUniFrac));
76 }
77 components.add(*fPdf);
78 RooAddPdf* addPdf = new RooAddPdf("proposalFunction", "Proposal Density",
79 components, coeffs);
80 fPdfProp->SetPdf(*addPdf);
81 fPdfProp->SetOwnsPdf(true);
82 if (fCacheSize > 0)
83 fPdfProp->SetCacheSize(fCacheSize);
84 fOwnsPdfProp = false;
85 return fPdfProp;
86}
87
88////////////////////////////////////////////////////////////////////////////////
89
91{
92 if (fVars == nullptr) {
93 coutE(InputArguments) << "ProposalHelper::CreatePdf(): " <<
94 "Variables to create proposal function for are not set." << std::endl;
95 return;
96 }
97 RooArgList xVec{};
98 RooArgList muVec{};
99 RooRealVar* clone;
100 for (auto *r : static_range_cast<RooRealVar *> (*fVars)){
101 xVec.add(*r);
102 TString cloneName = TString::Format("%s%s", "mu__", r->GetName());
103 clone = static_cast<RooRealVar*>(r->clone(cloneName.Data()));
104 muVec.add(*clone);
105 if (fUseUpdates)
106 fPdfProp->AddMapping(*clone, *r);
107 }
108 if (fCovMatrix == nullptr)
109 CreateCovMatrix(xVec);
110 fPdf = new RooMultiVarGaussian("mvg", "MVG Proposal", xVec, muVec,
111 *fCovMatrix);
112}
113
114////////////////////////////////////////////////////////////////////////////////
115
117{
118 Int_t size = xVec.size();
120 RooRealVar* r;
121 for (Int_t i = 0; i < size; i++) {
122 r = static_cast<RooRealVar*>(xVec.at(i));
123 double range = r->getMax() - r->getMin();
124 (*fCovMatrix)(i,i) = range / fSigmaRangeDivisor;
125 }
126}
127
128////////////////////////////////////////////////////////////////////////////////
129
131{
132 if (fClues != nullptr) {
133 if (fCluesOptions == nullptr) {
134 fCluesPdf = new RooNDKeysPdf("cluesPdf", "Clues PDF", *fVars, *fClues);
135 } else {
136 fCluesPdf = new RooNDKeysPdf("cluesPdf", "Clues PDF", *fVars, *fClues, fCluesOptions);
137 }
138 }
139}
140
141////////////////////////////////////////////////////////////////////////////////
142
144{
145 fUniformPdf = new RooUniform("uniform", "Uniform Proposal PDF",
146 RooArgSet(*fVars));
147}
ROOT::R::TRInterface & r
Definition Object.C:4
static const double SIGMA_RANGE_DIVISOR
static const double DEFAULT_CLUES_FRAC
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define oocoutI(o, a)
#define coutE(a)
int Int_t
Signed integer 4 bytes (int).
Definition RtypesCore.h:59
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:24
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
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
Basic string class.
Definition TString.h:138
const char * Data() const
Definition TString.h:384
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:2385
RooConstVar & RooConst(double val)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:72
@ InputArguments
Namespace for the RooStats classes.
Definition CodegenImpl.h:66