Logo ROOT   6.08/07
Reference Guide
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 
15 
16 #ifndef ROOT_Rtypes
17 #include "Rtypes.h"
18 #endif
19 
20 #ifndef ROOSTATS_PdfProposal
21 #include "RooStats/PdfProposal.h"
22 #endif
23 #ifndef RooStats_RooStatsUtils
24 #include "RooStats/RooStatsUtils.h"
25 #endif
26 #ifndef ROO_ARG_SET
27 #include "RooArgSet.h"
28 #endif
29 #ifndef ROO_DATA_SET
30 #include "RooDataSet.h"
31 #endif
32 #ifndef ROO_ABS_PDF
33 #include "RooAbsPdf.h"
34 #endif
35 #ifndef ROO_MSG_SERVICE
36 #include "RooMsgService.h"
37 #endif
38 #ifndef ROO_REAL_VAR
39 #include "RooRealVar.h"
40 #endif
41 #ifndef ROOT_TIterator
42 #include "TIterator.h"
43 #endif
44 
45 #include <map>
46 
48 
49 using namespace RooFit;
50 using namespace RooStats;
51 using namespace std;
52 
53 // By default, PdfProposal does NOT own the PDF that serves as the
54 // proposal density function
55 PdfProposal::PdfProposal() : ProposalFunction()
56 {
57  fPdf = NULL;
58  fOwnsPdf = kFALSE;
59  fCacheSize = 1;
60  fCachePosition = 0;
61  fCache = NULL;
62 }
63 
64 // By default, PdfProposal does NOT own the PDF that serves as the
65 // proposal density function
67 {
68  fPdf = &pdf;
69  fOwnsPdf = kFALSE;
70  fCacheSize = 1;
71  fCachePosition = 0;
72  fCache = NULL;
73 }
74 
76 {
77  if (x1.equals(x2)) {
78  TIterator* it = x1.createIterator();
79  RooRealVar* r;
80  while ((r = (RooRealVar*)it->Next()) != NULL)
81  if (r->getVal() != x2.getRealValue(r->GetName())) {
82  delete it;
83  return kFALSE;
84  }
85  delete it;
86  return kTRUE;
87  }
88  return kFALSE;
89 }
90 
91 // Populate xPrime with a new proposed point
93 {
94  if (fLastX.getSize() == 0) {
95  // fLastX not yet initialized
96  fLastX.addClone(x);
97  // generate initial cache
99  if (fMap.size() > 0) {
100  for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
101  fIt->first->setVal(fIt->second->getVal(&x));
102  }
103  fCache = fPdf->generate(xPrime, fCacheSize);
104  }
105 
106  Bool_t moved = false;
107  if (fMap.size() > 0) {
108  moved = !Equals(fLastX, x);
109 
110  // if we've moved, set the values of the variables in the PDF to the
111  // corresponding values of the variables in x, according to the
112  // mappings (i.e. let the variables in x set the given values for the
113  // PDF that will generate xPrime)
114  if (moved) {
115  // update the pdf parameters
117 
118  for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
119  fIt->first->setVal(fIt->second->getVal(&x));
120 
121  // save the new x in fLastX
123  }
124  }
125 
126  // generate new cache if necessary
127  if (moved || fCachePosition >= fCacheSize) {
128  delete fCache;
129  fCache = fPdf->generate(xPrime, fCacheSize);
130  fCachePosition = 0;
131  }
132 
133  const RooArgSet* proposal = fCache->get(fCachePosition);
134  fCachePosition++;
135  RooStats::SetParameters(proposal, &xPrime);
136 }
137 
138 // Determine whether or not the proposal density is symmetric for
139 // points x1 and x2 - that is, whether the probabilty of reaching x2
140 // from x1 is equal to the probability of reaching x1 from x2
142 {
143  // kbelasco: is there a better way to do this?
144  return false;
145 }
146 
147 // Return the probability of proposing the point x1 given the starting
148 // point x2
150 {
152  for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
153  fIt->first->setVal(fIt->second->getVal(&x2));
154  RooArgSet* temp = fPdf->getObservables(x1);
155  RooStats::SetParameters(&x1, temp);
156  delete temp;
157  return fPdf->getVal(&x1); // could just as well use x2
158 }
159 
161 {
163  if (update.getParameters((RooAbsData*)NULL)->getSize() == 0)
164  fMaster.add(update);
165  fMap.insert(pair<RooRealVar*, RooAbsReal*>(&proposalParam, &update));
166 }
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
TIterator * createIterator(Bool_t dir=kIterForward) const
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:86
RooArgSet * getObservables(const RooArgSet &set, Bool_t valueOnly=kTRUE) const
Definition: RooAbsArg.h:194
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)
Definition: RooStatsUtils.h:69
virtual Bool_t Equals(RooArgSet &x1, RooArgSet &x2)
whether we own the proposal density function
Definition: PdfProposal.cxx:75
bool Bool_t
Definition: RtypesCore.h:59
Bool_t equals(const RooAbsCollection &otherColl) const
Check if this and other collection have identically named contents.
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
Iterator abstract base class.
Definition: TIterator.h:32
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
virtual Double_t GetProposalDensity(RooArgSet &x1, RooArgSet &x2)
Return the probability of proposing the point x1 given the starting point x2.
virtual void AddMapping(RooRealVar &proposalParam, RooAbsReal &update)
specify a mapping between a parameter of the proposal function and a parameter of interest...
Int_t fCacheSize
the last point we were at
Definition: PdfProposal.h:163
std::map< RooRealVar *, RooAbsReal * >::iterator fIt
map of values in pdf to update
Definition: PdfProposal.h:161
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual void addClone(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling addOwned() for each element in the source...
Definition: RooArgSet.h:94
Int_t getSize() const
TRandom2 r(17)
RooDataSet * fCache
our position in the cached proposal data set
Definition: PdfProposal.h:165
RooArgSet fLastX
pdf iterator
Definition: PdfProposal.h:162
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event &#39;index&#39;.
Bool_t fOwnsPdf
pointers to master variables needed for updates
Definition: PdfProposal.h:167
Namespace for the RooStats classes.
Definition: Asimov.h:20
virtual void Propose(RooArgSet &xPrime, RooArgSet &x)
Populate xPrime with a new proposed point.
Definition: PdfProposal.cxx:92
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooArgSet * getParameters(const RooAbsData *data, Bool_t stripDisconnected=kTRUE) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don&#39;t match any of...
Definition: RooAbsArg.cxx:560
RooArgSet fMaster
the cached proposal data set
Definition: PdfProposal.h:166
PdfProposal is a concrete implementation of the ProposalFunction interface.
Definition: PdfProposal.h:85
std::map< RooRealVar *, RooAbsReal * > fMap
the proposal density function
Definition: PdfProposal.h:160
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Definition: RooArgSet.cxx:527
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none())
Generate a new dataset containing the specified variables with events sampled from our distribution...
Definition: RooAbsPdf.cxx:1702
virtual TObject * Next()=0
#define NULL
Definition: Rtypes.h:82
virtual Bool_t IsSymmetric(RooArgSet &x1, RooArgSet &x2)
Determine whether or not the proposal density is symmetric for points x1 and x2 - that is...
const Bool_t kTRUE
Definition: Rtypes.h:91
Int_t fCachePosition
how many points to generate each time
Definition: PdfProposal.h:164