Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
UniformProposal.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/** \class RooStats::UniformProposal
13 \ingroup Roostats
14
15UniformProposal is a concrete implementation of the ProposalFunction interface
16for use with a Markov Chain Monte Carlo algorithm. This proposal function is
17a uniformly random distribution over the parameter space. The proposal
18ignores the current point when it proposes a new point. The proposal
19function is symmetric, though it may not cause a MetropolisHastings run to
20converge as quickly as other proposal functions.
21
22*/
23
24#include "Rtypes.h"
25
28#include "RooArgSet.h"
29#include "RooMsgService.h"
30#include "RooRealVar.h"
31
32using namespace std;
33
35
36using namespace RooFit;
37using namespace RooStats;
38
39////////////////////////////////////////////////////////////////////////////////
40/// Populate xPrime with a new proposed point
41
43{
44 // kbelasco: remember xPrime and x have not been checked for containing
45 // only RooRealVars
47}
48
49////////////////////////////////////////////////////////////////////////////////
50/// Determine whether or not the proposal density is symmetric for
51/// points x1 and x2 - that is, whether the probability of reaching x2
52/// from x1 is equal to the probability of reaching x1 from x2
53
55{
56 return true;
57}
58
59////////////////////////////////////////////////////////////////////////////////
60/// Return the probability of proposing the point x1 given the starting
61/// point x2
62
65{
66 // For a uniform proposal, all points have equal probability and the
67 // value of the proposal density function is:
68 // 1 / (N-dimensional volume of interval)
69 double volume = 1.0;
70 for (auto const *var : static_range_cast<RooRealVar *> (x2))
71 volume *= (var->getMax() - var->getMin());
72 return 1.0 / volume;
73}
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char x2
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
UniformProposal is a concrete implementation of the ProposalFunction interface for use with a Markov ...
void Propose(RooArgSet &xPrime, RooArgSet &x) override
Populate xPrime with a new proposed point.
double GetProposalDensity(RooArgSet &x1, RooArgSet &x2) override
Return the probability of proposing the point x1 given the starting point x2.
bool IsSymmetric(RooArgSet &x1, RooArgSet &x2) override
Determine whether or not the proposal density is symmetric for points x1 and x2 - that is,...
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition Common.h:18
Namespace for the RooStats classes.
Definition Asimov.h:19
void RandomizeCollection(RooAbsCollection &set, bool randomizeConstants=true)
assuming all values in set are RooRealVars, randomize their values