Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ProposalFunction.h
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
15#ifndef ROOSTATS_ProposalFunction
16#define ROOSTATS_ProposalFunction
17
18#include "Rtypes.h"
19
20#include "RooArgSet.h"
21#include "RooMsgService.h"
22#include "RooRealVar.h"
23
24
25namespace RooStats {
26
27/** \class ProposalFunction
28 \ingroup Roostats
29ProposalFunction is an interface for all proposal functions that would be used
30with a Markov Chain Monte Carlo algorithm.
31Given a current point in the parameter space it proposes a new point.
32Proposal functions may or may not be symmetric, in the sense that the
33probability to propose X1 given we are at X2
34need not be the same as the probability to propose X2 given that we are at X1.
35In this case, the IsSymmetric method
36should return false, and the Metropolis algorithm will need to take into account
37the proposal density to maintain detailed balance.
38*/
39
40
41 class ProposalFunction : public TObject {
42
43 public:
44 ///Default constructor
46
47 /// Populate xPrime with the new proposed point,
48 /// possibly based on the current point x
49 virtual void Propose(RooArgSet& xPrime, RooArgSet& x) = 0;
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 virtual bool IsSymmetric(RooArgSet& x1, RooArgSet& x2) = 0;
55
56 /// Return the probability of proposing the point x1 given the starting
57 /// point x2
58 virtual double GetProposalDensity(RooArgSet& x1, RooArgSet& x2) = 0;
59
60 /// Check the parameters for which the ProposalFunction will
61 /// propose values to make sure they are all RooRealVars
62 /// Return true if all objects are RooRealVars, false otherwise
63 virtual bool CheckParameters(RooArgSet& params)
64 {
65 for (auto *obj : params){
66 if (!dynamic_cast<RooRealVar*>(obj)) {
67 coutE(Eval) << "Error when checking parameters in"
68 << "ProposalFunction: "
69 << "Object \"" << obj->GetName() << "\" not of type "
70 << "RooRealVar" << std::endl;
71 return false;
72 }
73 }
74 // Made it here, so all parameters are RooRealVars
75 return true;
76 }
77
78 protected:
79 ClassDefOverride(ProposalFunction,1) /// Interface for the proposal function used with Markov Chain Monte Carlo
80 };
81}
82
83#endif
#define coutE(a)
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
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
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
virtual void Propose(RooArgSet &xPrime, RooArgSet &x)=0
Populate xPrime with the new proposed point, possibly based on the current point x.
virtual double GetProposalDensity(RooArgSet &x1, RooArgSet &x2)=0
Return the probability of proposing the point x1 given the starting point x2.
ProposalFunction()
Default constructor.
virtual bool IsSymmetric(RooArgSet &x1, RooArgSet &x2)=0
Determine whether or not the proposal density is symmetric for points x1 and x2 - that is,...
virtual bool CheckParameters(RooArgSet &params)
Check the parameters for which the ProposalFunction will propose values to make sure they are all Roo...
Mother of all ROOT objects.
Definition TObject.h:41
Double_t x[n]
Definition legend1.C:17
Namespace for the RooStats classes.
Definition CodegenImpl.h:58