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
33
34using namespace RooFit;
35using namespace RooStats;
36
37static const char* WEIGHT_NAME = "weight_MarkovChain_local_";
38static const char* NLL_NAME = "nll_MarkovChain_local_";
39static const char* DATASET_NAME = "dataset_MarkovChain_local_";
40static const char* DEFAULT_NAME = "_markov_chain";
41static const char* DEFAULT_TITLE = "Markov Chain";
42
44
46{
47 SetParameters(parameters);
48}
49
50MarkovChain::MarkovChain(const char *name, const char *title, RooArgSet &parameters) : TNamed(name, title)
51{
52 SetParameters(parameters);
53}
54
56{
57 delete fChain;
58 delete fParameters;
59 delete fDataEntry;
60
61 fParameters = new RooArgSet();
62 fParameters->addClone(parameters);
63
64 // kbelasco: consider setting fDataEntry = fChain->get()
65 // to see if that makes it possible to get values of variables without
66 // doing string comparison
67 RooRealVar nll(NLL_NAME, "-log Likelihood", 0);
68
69 fDataEntry = new RooArgSet();
70 fDataEntry->addClone(parameters);
71 fDataEntry->addClone(nll);
72 fNLL = static_cast<RooRealVar*>(fDataEntry->find(NLL_NAME));
73
75}
76
77void MarkovChain::Add(RooArgSet& entry, double nllValue, double weight)
78{
79 if (fParameters == nullptr)
80 SetParameters(entry);
82 fNLL->setVal(nllValue);
83 //kbelasco: this is stupid, but some things might require it, so be doubly sure
84 fChain->add(*fDataEntry, weight);
85 //fChain->add(*fDataEntry);
86}
87
89{
90 // Discards the first n accepted points.
91
92 if(fParameters == nullptr) SetParameters(*const_cast<RooArgSet*>(otherChain.Get()));
93 int counter = 0;
94 for( int i=0; i < otherChain.Size(); i++ ) {
95 RooArgSet* entry = const_cast<RooArgSet*>(otherChain.Get(i));
96 counter += 1;
97 if( counter > burnIn ) {
98 AddFast( *entry, otherChain.NLL(), otherChain.Weight() );
99 }
100 }
101}
102void MarkovChain::Add(MarkovChain& otherChain, double discardEntries)
103{
104 // Discards the first entries. This is different to the definition of
105 // burn-in used in the Bayesian calculator where the first n accepted
106 // terms from the proposal function are discarded.
107
108 if(fParameters == nullptr) SetParameters(*const_cast<RooArgSet*>(otherChain.Get()));
109 double counter = 0.0;
110 for( int i=0; i < otherChain.Size(); i++ ) {
111 RooArgSet* entry = const_cast<RooArgSet*>(otherChain.Get(i));
112 counter += otherChain.Weight();
113 if( counter > discardEntries ) {
114 AddFast( *entry, otherChain.NLL(), otherChain.Weight() );
115 }
116 }
117}
118
119void MarkovChain::AddFast(RooArgSet& entry, double nllValue, double weight)
120{
122 fNLL->setVal(nllValue);
123 //kbelasco: this is stupid, but some things might require it, so be doubly sure
124 fChain->addFast(*fDataEntry, weight);
125 //fChain->addFast(*fDataEntry);
126}
127
129{
130 RooArgSet args;
131 if (whichVars == nullptr) {
132 //args.add(*fParameters);
133 //args.add(*fNLL);
134 args.add(*fDataEntry);
135 } else {
136 args.add(*whichVars);
137 }
138
139 return RooFit::makeOwningPtr<RooDataSet>(std::unique_ptr<RooAbsData>{fChain->reduce(args)});
140}
141
143{
144 RooArgSet args;
145 if (whichVars == nullptr) {
146 args.add(*fParameters);
147 //args.add(*fNLL);
148 //args.add(*fDataEntry);
149 } else {
150 args.add(*whichVars);
151 }
152
153 std::unique_ptr<RooAbsData> data{fChain->reduce(args)};
154 return RooFit::makeOwningPtr(std::unique_ptr<RooDataHist>{static_cast<RooDataSet&>(*data).binnedClone()});
155}
156
158{
159 RooArgList axes;
160 if (whichVars == nullptr) {
161 axes.add(*fParameters);
162 } else {
163 axes.add(*whichVars);
164 }
165
166 Int_t dim = axes.size();
167 std::vector<double> min(dim);
168 std::vector<double> max(dim);
169 std::vector<Int_t> bins(dim);
170 std::vector<const char *> names(dim);
171 Int_t i = 0;
172 for (auto const *var : static_range_cast<RooRealVar *>(axes)) {
173 names[i] = var->GetName();
174 min[i] = var->getMin();
175 max[i] = var->getMax();
176 bins[i] = var->numBins();
177 ++i;
178 }
179
180 THnSparseF* sparseHist = new THnSparseF("posterior", "MCMC Posterior Histogram",
181 dim, &bins[0], &min[0], &max[0]);
182
183 // kbelasco: it appears we need to call Sumw2() just to get the
184 // histogram to keep a running total of the weight so that Getsumw doesn't
185 // just return 0
186 sparseHist->Sumw2();
187
188 // Fill histogram
190 const RooArgSet* entry;
191 std::vector<double> x(dim);
192 for ( i = 0; i < size; i++) {
193 entry = fChain->get(i);
194
195 for (Int_t ii = 0; ii < dim; ii++) {
196 //LM: doing this is probably quite slow
197 x[ii] = entry->getRealValue( names[ii]);
198 sparseHist->Fill(x.data(), fChain->weight());
199 }
200 }
201
202 return sparseHist;
203}
204
205double MarkovChain::NLL(Int_t i) const
206{
207 // kbelasco: how to do this?
208 //fChain->get(i);
209 //return fNLL->getVal();
210 return fChain->get(i)->getRealValue(NLL_NAME);
211}
212
213double MarkovChain::NLL() const
214{
215 // kbelasco: how to do this?
216 //fChain->get();
217 //return fNLL->getVal();
218 return fChain->get()->getRealValue(NLL_NAME);
219}
220
222{
223 return fChain->weight();
224}
225
227{
228 fChain->get(i);
229 return fChain->weight();
230}
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:382
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:221
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.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
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.
RooFit::OwningPtr< RooAbsData > reduce(const RooCmdArg &arg1, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Create a reduced copy of this dataset.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
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:24
Container class to hold unbinned data.
Definition RooDataSet.h:33
const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
RooFit::OwningPtr< RooDataHist > binnedClone(const char *newName=nullptr, const char *newTitle=nullptr) const
Return binned clone of this dataset.
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.
void add(const RooArgSet &row, double weight, double weightError)
Add one ore more rows of data.
double weight() const override
Return event weight of current event.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
Stores the steps in a Markov Chain of points.
Definition MarkovChain.h:26
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 RooFit::OwningPtr< RooDataHist > GetAsDataHist(RooArgSet *whichVars=nullptr) const
get this MarkovChain as a RooDataHist whose entries contain the values of whichVars.
RooArgSet * fDataEntry
Definition MarkovChain.h:95
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
virtual RooFit::OwningPtr< RooDataSet > GetAsDataSet(RooArgSet *whichVars=nullptr) const
get this MarkovChain as a RooDataSet whose entries contain the values of whichVars.
RooArgSet * fParameters
Definition MarkovChain.h:94
virtual THnSparse * GetAsSparseHist(RooAbsCollection *whichVars=nullptr) const
Get a clone of the markov chain on which this interval is based as a sparse histogram.
RooDataSet * fChain
Definition MarkovChain.h:96
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:47
virtual double Weight() const
get the weight of the current (last indexed) entry
virtual Int_t Size() const
get the number of steps in the chain
Definition MarkovChain.h:45
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:207
Efficient multidimensional histogram.
Definition THnSparse.h:37
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 CodegenImpl.h:64
OwningPtr< T > makeOwningPtr(std::unique_ptr< T > &&ptr)
Internal helper to turn a std::unique_ptr<T> into an OwningPtr.
Definition Config.h:40
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:35
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
void SetParameters(const RooArgSet *desiredVals, RooArgSet *paramsToChange)