Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
PdfProposal.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Authors: Kevin Belasco 17/06/2009
3// Authors: Kyle Cranmer 17/06/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////////////////////////////////////////////////////////////////////////////////
13
14/** \class RooStats::PdfProposal
15 \ingroup Roostats
16
17PdfProposal is a concrete implementation of the ProposalFunction interface.
18It proposes points across the parameter space in the distribution of the
19given PDF.
20
21To make Propose(xPrime, x) dependent on x, configure with
22PdfProposal::AddMapping(varToUpdate, valueToUse). For example, suppose we have:
23
24~~~{.cpp}
25// our parameter
26RooRealVar p("p", "p", 5, 0, 10);
27
28// create mean and sigma for gaussian proposal function
29RooRealVar meanP("meanP", "meanP", 0, 10);
30RooRealVar sigma("sigma", "sigma", 1, 0, 5);
31RooGaussian pGaussian("pGaussian", "pGaussian", p, meanP, sigma);
32
33// configure proposal function
34PdfProposal pdfProposal(pGaussian);
35pdfProposal.AddMapping(meanP, p); // each call of Propose(xPrime, x), meanP in
36 // the proposal function will be updated to
37 // the value of p in x. this will center the
38 // proposal function about x's p when
39 // proposing for xPrime
40
41// To improve performance, PdfProposal has the ability to cache a specified
42// number of proposals. If you don't call this function, the default cache size
43// is 1, which can be slow.
44pdfProposal.SetCacheSize(desiredCacheSize);
45~~~
46
47PdfProposal currently uses a fixed cache size. Adaptive caching methods are in the works
48for future versions.
49*/
50
51
52#include "Rtypes.h"
53
56#include "RooArgSet.h"
57#include "RooDataSet.h"
58#include "RooAbsPdf.h"
59#include "RooMsgService.h"
60#include "RooRealVar.h"
61
62#include <map>
63
64
65using namespace RooFit;
66using namespace RooStats;
67
68////////////////////////////////////////////////////////////////////////////////
69/// By default, PdfProposal does NOT own the PDF that serves as the
70/// proposal density function
71
72PdfProposal::PdfProposal() = default;
73
74////////////////////////////////////////////////////////////////////////////////
75/// By default, PdfProposal does NOT own the PDF that serves as the
76/// proposal density function
77
79
80////////////////////////////////////////////////////////////////////////////////
81/// determine whether these two RooArgSets represent the same point
82
84{
85 if (x1.equals(x2)) {
86 for (auto const *r : static_range_cast<RooRealVar *>(x1)) {
87 if (r->getVal() != x2.getRealValue(r->GetName())) {
88 return false;
89 }
90 }
91 return true;
92 }
93 return false;
94}
95
96////////////////////////////////////////////////////////////////////////////////
97/// Populate xPrime with a new proposed point
98
100{
101 if (fLastX.empty()) {
102 // fLastX not yet initialized
104 // generate initial cache
106 if (!fMap.empty()) {
107 for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
108 fIt->first->setVal(fIt->second->getVal(&x));
109 }
110 fCache = std::unique_ptr<RooDataSet>{fPdf->generate(xPrime, fCacheSize)};
111 }
112
113 bool moved = false;
114 if (!fMap.empty()) {
115 moved = !Equals(fLastX, x);
116
117 // if we've moved, set the values of the variables in the PDF to the
118 // corresponding values of the variables in x, according to the
119 // mappings (i.e. let the variables in x set the given values for the
120 // PDF that will generate xPrime)
121 if (moved) {
122 // update the pdf parameters
124
125 for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
126 fIt->first->setVal(fIt->second->getVal(&x));
127
128 // save the new x in fLastX
130 }
131 }
132
133 // generate new cache if necessary
134 if (moved || fCachePosition >= fCacheSize) {
135 fCache = std::unique_ptr<RooDataSet>{fPdf->generate(xPrime, fCacheSize)};
136 fCachePosition = 0;
137 }
138
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// Determine whether or not the proposal density is symmetric for
146/// points x1 and x2 - that is, whether the propabilty of reaching x2
147/// from x1 is equal to the probability of reaching x1 from x2
148
150{
151 // kbelasco: is there a better way to do this?
152 return false;
153}
154
155////////////////////////////////////////////////////////////////////////////////
156/// Return the probability of proposing the point x1 given the starting
157/// point x2
158
160{
162 for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
163 fIt->first->setVal(fIt->second->getVal(&x2));
164 std::unique_ptr<RooArgSet> temp{fPdf->getObservables(x1)};
165 RooStats::SetParameters(&x1, temp.get());
166 return fPdf->getVal(&x1); // could just as well use x2
167}
168
169////////////////////////////////////////////////////////////////////////////////
170/// specify a mapping between a parameter of the proposal function and
171/// a parameter of interest. this mapping is used to set the value of
172/// proposalParam equal to the value of update to determine the
173/// proposal function.
174/// proposalParam is a parameter of the proposal function that must
175/// be set to the value of update (from the current point) in order to
176/// propose a new point.
177
179{
180 fMaster.add(*update.getParameters(static_cast<RooAbsData const*>(nullptr)));
181 if (update.getParameters(static_cast<RooAbsData const*>(nullptr))->empty())
183 fMap.insert(std::pair<RooRealVar*, RooAbsReal*>(&proposalParam, &update));
184}
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooFit::OwningPtr< RooDataSet > generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={})
See RooAbsPdf::generate(const RooArgSet&,const RooCmdArg&,const RooCmdArg&,const RooCmdArg&,...
Definition RooAbsPdf.h:57
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Variable that can be changed from the outside.
Definition RooRealVar.h:37
std::map< RooRealVar *, RooAbsReal * >::iterator fIt
map of values in pdf to update
RooArgSet fMaster
the cached proposal data set
bool IsSymmetric(RooArgSet &x1, RooArgSet &x2) override
Determine whether or not the proposal density is symmetric for points x1 and x2 - that is,...
virtual void AddMapping(RooRealVar &proposalParam, RooAbsReal &update)
specify a mapping between a parameter of the proposal function and a parameter of interest.
virtual bool Equals(RooArgSet &x1, RooArgSet &x2)
whether we own the proposal density function
RooArgSet fLastX
pdf iterator
void Propose(RooArgSet &xPrime, RooArgSet &x) override
Populate xPrime with a new proposed point.
Int_t fCacheSize
the last point we were at
Int_t fCachePosition
how many points to generate each time
std::map< RooRealVar *, RooAbsReal * > fMap
the proposal density function
PdfProposal()
By default, PdfProposal does NOT own the PDF that serves as the proposal density function.
std::unique_ptr< RooDataSet > fCache
our position in the cached proposal data set
double GetProposalDensity(RooArgSet &x1, RooArgSet &x2) override
Return the probability of proposing the point x1 given the starting point x2.
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)