Logo ROOT  
Reference Guide
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#include "TIterator.h"
32
33using namespace std;
34
36
37using namespace RooFit;
38using namespace RooStats;
39
40////////////////////////////////////////////////////////////////////////////////
41/// Populate xPrime with a new proposed point
42
43void UniformProposal::Propose(RooArgSet& xPrime, RooArgSet& /* x */)
44{
45 // kbelasco: remember xPrime and x have not been checked for containing
46 // only RooRealVars
48}
49
50////////////////////////////////////////////////////////////////////////////////
51/// Determine whether or not the proposal density is symmetric for
52/// points x1 and x2 - that is, whether the probability of reaching x2
53/// from x1 is equal to the probability of reaching x1 from x2
54
56{
57 return true;
58}
59
60////////////////////////////////////////////////////////////////////////////////
61/// Return the probability of proposing the point x1 given the starting
62/// point x2
63
66{
67 // For a uniform proposal, all points have equal probability and the
68 // value of the proposal density function is:
69 // 1 / (N-dimensional volume of interval)
70 Double_t volume = 1.0;
71 TIterator* it = x2.createIterator();
72 RooRealVar* var;
73 while ((var = (RooRealVar*)it->Next()) != NULL)
74 volume *= (var->getMax() - var->getMin());
75 delete it;
76 return 1.0 / volume;
77}
static const double x2[5]
#define ClassImp(name)
Definition: Rtypes.h:361
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
UniformProposal is a concrete implementation of the ProposalFunction interface for use with a Markov ...
virtual Bool_t IsSymmetric(RooArgSet &x1, RooArgSet &x2)
Determine whether or not the proposal density is symmetric for points x1 and x2 - that is,...
virtual Double_t GetProposalDensity(RooArgSet &x1, RooArgSet &x2)
Return the probability of proposing the point x1 given the starting point x2.
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Namespace for the RooStats classes.
Definition: Asimov.h:19
void RandomizeCollection(RooAbsCollection &set, Bool_t randomizeConstants=kTRUE)