Logo ROOT   6.12/07
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 
15 UniformProposal is a concrete implementation of the ProposalFunction interface
16 for use with a Markov Chain Monte Carlo algorithm. This proposal function is
17 a uniformly random distribution over the parameter space. The proposal
18 ignores the current point when it proposes a new point. The proposal
19 function is symmetric, though it may not cause a MetropolisHastings run to
20 converge as quickly as other proposal functions.
21 
22 */
23 
24 #include "Rtypes.h"
25 
26 #include "RooStats/RooStatsUtils.h"
28 #include "RooArgSet.h"
29 #include "RooMsgService.h"
30 #include "RooRealVar.h"
31 #include "TIterator.h"
32 
33 using namespace std;
34 
36 
37 using namespace RooFit;
38 using namespace RooStats;
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 /// Populate xPrime with a new proposed point
42 
43 void 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 
55 Bool_t UniformProposal::IsSymmetric(RooArgSet& /* x1 */ , RooArgSet& /* x2 */)
56 {
57  return true;
58 }
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Return the probability of proposing the point x1 given the starting
62 /// point x2
63 
64 Double_t UniformProposal::GetProposalDensity(RooArgSet& /* x1 */,
65  RooArgSet& x2)
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 }
virtual Double_t getMin(const char *name=0) const
TIterator * createIterator(Bool_t dir=kIterForward) const
virtual Double_t getMax(const char *name=0) const
UniformProposal is a concrete implementation of the ProposalFunction interface for use with a Markov ...
bool Bool_t
Definition: RtypesCore.h:59
STL namespace.
Iterator abstract base class.
Definition: TIterator.h:30
static const double x2[5]
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
void RandomizeCollection(RooAbsCollection &set, Bool_t randomizeConstants=kTRUE)
Definition: RooStatsUtils.h:99
Namespace for the RooStats classes.
Definition: Asimov.h:20
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
virtual TObject * Next()=0