Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
MarkovChain.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::MarkovChain
13 \ingroup Roostats
14
15Stores the steps in a Markov Chain of points. Allows user to access the
16weight and NLL value (if applicable) with which a point was added to the
17MarkovChain.
18
19*/
20
21#include "Rtypes.h"
22
23#include "TNamed.h"
25#include "RooDataSet.h"
26#include "RooArgSet.h"
27#include "RooRealVar.h"
29#include "RooDataHist.h"
30#include "THnSparse.h"
31
32using namespace std;
33
35
36using namespace RooFit;
37using namespace RooStats;
38
39static const char* WEIGHT_NAME = "weight_MarkovChain_local_";
40static const char* NLL_NAME = "nll_MarkovChain_local_";
41static const char* DATASET_NAME = "dataset_MarkovChain_local_";
42static const char* DEFAULT_NAME = "_markov_chain";
43static const char* DEFAULT_TITLE = "Markov Chain";
44
47{
48 fParameters = nullptr;
49 fDataEntry = nullptr;
50 fChain = nullptr;
51 fNLL = nullptr;
52}
53
56{
57 fParameters = nullptr;
58 fDataEntry = nullptr;
59 fChain = nullptr;
60 fNLL = nullptr;
61 SetParameters(parameters);
62}
63
64MarkovChain::MarkovChain(const char* name, const char* title,
65 RooArgSet& parameters) : TNamed(name, title)
66{
67 fParameters = nullptr;
68 fDataEntry = nullptr;
69 fChain = nullptr;
70 fNLL = nullptr;
71 SetParameters(parameters);
72}
73
75{
76 delete fChain;
77 delete fParameters;
78 delete fDataEntry;
79
80 fParameters = new RooArgSet();
81 fParameters->addClone(parameters);
82
83 // kbelasco: consider setting fDataEntry = fChain->get()
84 // to see if that makes it possible to get values of variables without
85 // doing string comparison
86 RooRealVar nll(NLL_NAME, "-log Likelihood", 0);
87
88 fDataEntry = new RooArgSet();
89 fDataEntry->addClone(parameters);
90 fDataEntry->addClone(nll);
92
94}
95
96void MarkovChain::Add(RooArgSet& entry, double nllValue, double weight)
97{
98 if (fParameters == nullptr)
99 SetParameters(entry);
101 fNLL->setVal(nllValue);
102 //kbelasco: this is stupid, but some things might require it, so be doubly sure
103 fChain->add(*fDataEntry, weight);
104 //fChain->add(*fDataEntry);
105}
106
108{
109 // Discards the first n accepted points.
110
111 if(fParameters == nullptr) SetParameters(*(RooArgSet*)otherChain.Get());
112 int counter = 0;
113 for( int i=0; i < otherChain.Size(); i++ ) {
114 RooArgSet* entry = (RooArgSet*)otherChain.Get(i);
115 counter += 1;
116 if( counter > burnIn ) {
117 AddFast( *entry, otherChain.NLL(), otherChain.Weight() );
118 }
119 }
120}
121void MarkovChain::Add(MarkovChain& otherChain, double discardEntries)
122{
123 // Discards the first entries. This is different to the definition of
124 // burn-in used in the Bayesian calculator where the first n accepted
125 // terms from the proposal function are discarded.
126
127 if(fParameters == nullptr) SetParameters(*(RooArgSet*)otherChain.Get());
128 double counter = 0.0;
129 for( int i=0; i < otherChain.Size(); i++ ) {
130 RooArgSet* entry = (RooArgSet*)otherChain.Get(i);
131 counter += otherChain.Weight();
132 if( counter > discardEntries ) {
133 AddFast( *entry, otherChain.NLL(), otherChain.Weight() );
134 }
135 }
136}
137
138void MarkovChain::AddFast(RooArgSet& entry, double nllValue, double weight)
139{
141 fNLL->setVal(nllValue);
142 //kbelasco: this is stupid, but some things might require it, so be doubly sure
143 fChain->addFast(*fDataEntry, weight);
144 //fChain->addFast(*fDataEntry);
145}
146
148{
149 RooArgSet args;
150 if (whichVars == nullptr) {
151 //args.add(*fParameters);
152 //args.add(*fNLL);
153 args.add(*fDataEntry);
154 } else {
155 args.add(*whichVars);
156 }
157
159 //data = dynamic_cast<RooDataSet*>(fChain->reduce(args));
160 data = (RooDataSet*)fChain->reduce(args);
161
162 return data;
163}
164
166 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
167 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
168{
170 data = (RooDataSet*)fChain->reduce(arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8);
171 return data;
172}
173
175{
176 RooArgSet args;
177 if (whichVars == nullptr) {
178 args.add(*fParameters);
179 //args.add(*fNLL);
180 //args.add(*fDataEntry);
181 } else {
182 args.add(*whichVars);
183 }
184
186 RooDataHist* hist = data->binnedClone();
187 delete data;
188
189 return hist;
190}
191
193 const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5,
194 const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
195{
197 RooDataHist* hist;
198 data = (RooDataSet*)fChain->reduce(arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8);
199 hist = data->binnedClone();
200 delete data;
201
202 return hist;
203}
204
206{
207 RooArgList axes;
208 if (whichVars == nullptr)
209 axes.add(*fParameters);
210 else
211 axes.add(*whichVars);
212
213 Int_t dim = axes.getSize();
214 std::vector<double> min(dim);
215 std::vector<double> max(dim);
216 std::vector<Int_t> bins(dim);
217 std::vector<const char *> names(dim);
218 Int_t i = 0;
219 for (auto const *var : static_range_cast<RooRealVar *>(axes)) {
220 names[i] = var->GetName();
221 min[i] = var->getMin();
222 max[i] = var->getMax();
223 bins[i] = var->numBins();
224 ++i;
225 }
226
227 THnSparseF* sparseHist = new THnSparseF("posterior", "MCMC Posterior Histogram",
228 dim, &bins[0], &min[0], &max[0]);
229
230 // kbelasco: it appears we need to call Sumw2() just to get the
231 // histogram to keep a running total of the weight so that Getsumw doesn't
232 // just return 0
233 sparseHist->Sumw2();
234
235 // Fill histogram
237 const RooArgSet* entry;
238 std::vector<double> x(dim);
239 for ( i = 0; i < size; i++) {
240 entry = fChain->get(i);
241
242 for (Int_t ii = 0; ii < dim; ii++) {
243 //LM: doing this is probably quite slow
244 x[ii] = entry->getRealValue( names[ii]);
245 sparseHist->Fill(x.data(), fChain->weight());
246 }
247 }
248
249 return sparseHist;
250}
251
252double MarkovChain::NLL(Int_t i) const
253{
254 // kbelasco: how to do this?
255 //fChain->get(i);
256 //return fNLL->getVal();
257 return fChain->get(i)->getRealValue(NLL_NAME);
258}
259
260double MarkovChain::NLL() const
261{
262 // kbelasco: how to do this?
263 //fChain->get();
264 //return fNLL->getVal();
265 return fChain->get()->getRealValue(NLL_NAME);
266}
267
269{
270 return fChain->weight();
271}
272
274{
275 fChain->get(i);
276 return fChain->weight();
277}
static const char * DATASET_NAME
static const char * NLL_NAME
static const char * DEFAULT_TITLE
static const char * WEIGHT_NAME
static const char * DEFAULT_NAME
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
THnSparseT< TArrayF > THnSparseF
Definition THnSparse.h:220
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
double getRealValue(const char *name, double defVal=0.0, bool verbose=false) const
Get value of a RooAbsReal stored in set with given name.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
RooAbsData * reduce(const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg())
Create a reduced copy of this dataset.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooCmdArg is a named container for two doubles, two integers two object points and three string point...
Definition RooCmdArg.h:26
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:57
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
void add(const RooArgSet &row, double weight=1.0, double weightError=0.0) override
Add one ore more rows of data.
virtual void addFast(const RooArgSet &row, double weight=1.0, double weightError=0.0)
Add a data point, with its coordinates specified in the 'data' argset, to the data set.
double weight() const override
Return event weight of current event.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setVal(double value) override
Set value of variable to 'value'.
Stores the steps in a Markov Chain of points.
Definition MarkovChain.h:30
virtual RooDataHist * GetAsDataHist(RooArgSet *whichVars=nullptr) const
get this MarkovChain as a RooDataHist whose entries contain the values of whichVars.
virtual void AddFast(RooArgSet &entry, double nllValue, double weight=1.0)
add an entry to the chain ONLY IF you have constructed with parameters or called SetParameters
virtual double NLL() const
get the NLL value of the current (last indexed) entry
virtual void AddWithBurnIn(MarkovChain &otherChain, Int_t burnIn=0)
add another markov chain
virtual void Add(RooArgSet &entry, double nllValue, double weight=1.0)
safely add an entry to the chain
virtual double NLL(Int_t i) const
get the NLL value of entry at position i
RooArgSet * fParameters
virtual THnSparse * GetAsSparseHist(RooAbsCollection *whichVars=nullptr) const
Get a clone of the markov chain on which this interval is based as a sparse histogram.
virtual void SetParameters(RooArgSet &parameters)
set which of your parameters this chain should store
virtual const RooArgSet * Get(Int_t i) const
get the entry at position i
Definition MarkovChain.h:51
virtual double Weight() const
get the weight of the current (last indexed) entry
virtual RooDataSet * GetAsDataSet(RooArgSet *whichVars=nullptr) const
get this MarkovChain as a RooDataSet whose entries contain the values of whichVars.
virtual Int_t Size() const
get the number of steps in the chain
Definition MarkovChain.h:49
Long64_t Fill(const Double_t *x, Double_t w=1.)
Definition THnBase.h:149
Templated implementation of the abstract base THnSparse.
Definition THnSparse.h:206
Efficient multidimensional histogram.
Definition THnSparse.h:36
void Sumw2() override
Enable calculation of errors.
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
Double_t x[n]
Definition legend1.C:17
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 SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)