Logo ROOT  
Reference Guide
ProposalHelper.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Authors: Kevin Belasco 7/22/2009
3// Authors: Kyle Cranmer 7/22/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::ProposalHelper
13 \ingroup Roostats
14*/
15
16#include "Rtypes.h"
20#include "RooArgSet.h"
21#include "RooDataSet.h"
22#include "RooAbsPdf.h"
23#include "RooAddPdf.h"
24#include "RooNDKeysPdf.h"
25#include "RooUniform.h"
26#include "RooMsgService.h"
27#include "RooRealVar.h"
28#include "TIterator.h"
29#include "RooMultiVarGaussian.h"
30#include "RooConstVar.h"
31#include "TString.h"
32
33#include <map>
34
35namespace RooStats {
36 class ProposalFunction;
37}
38
40
41using namespace RooFit;
42using namespace RooStats;
43using namespace std;
44
45//static const Double_t DEFAULT_UNI_FRAC = 0.10;
46static const Double_t DEFAULT_CLUES_FRAC = 0.20;
47//static const Double_t SIGMA_RANGE_DIVISOR = 6;
49//static const Int_t DEFAULT_CACHE_SIZE = 100;
50//static const Option_t* CLUES_OPTIONS = "a";
51
52////////////////////////////////////////////////////////////////////////////////
53
54ProposalHelper::ProposalHelper()
55{
56 fPdfProp = new PdfProposal();
57 fVars = NULL;
63 fPdf = NULL;
65 fCluesPdf = NULL;
66 fUniformPdf = NULL;
67 fClues = NULL;
68 fCovMatrix = NULL;
69 fCluesFrac = -1;
70 fUniFrac = -1;
71 fCacheSize = -1;
72 fCluesOptions = NULL;
73}
74
75////////////////////////////////////////////////////////////////////////////////
76
78{
79 if (fPdf == NULL)
80 CreatePdf();
81 // kbelasco: check here for memory leaks: does RooAddPdf make copies or
82 // take ownership of components, coeffs
83 RooArgList* components = new RooArgList();
84 RooArgList* coeffs = new RooArgList();
85 if (fCluesPdf == NULL)
87 if (fCluesPdf != NULL) {
88 if (fCluesFrac < 0)
90 printf("added clues from dataset %s with fraction %g\n",
92 components->add(*fCluesPdf);
93 coeffs->add(RooConst(fCluesFrac));
94 }
95 if (fUniFrac > 0.) {
97 components->add(*fUniformPdf);
98 coeffs->add(RooConst(fUniFrac));
99 }
100 components->add(*fPdf);
101 RooAddPdf* addPdf = new RooAddPdf("proposalFunction", "Proposal Density",
102 *components, *coeffs);
103 fPdfProp->SetPdf(*addPdf);
105 if (fCacheSize > 0)
108 return fPdfProp;
109}
110
111////////////////////////////////////////////////////////////////////////////////
112
114{
115 // kbelasco: check here for memory leaks:
116 // does RooMultiVarGaussian make copies of xVec and muVec?
117 // or should we delete them?
118 if (fVars == NULL) {
119 coutE(InputArguments) << "ProposalHelper::CreatePdf(): " <<
120 "Variables to create proposal function for are not set." << endl;
121 return;
122 }
123 RooArgList* xVec = new RooArgList();
124 RooArgList* muVec = new RooArgList();
126 RooRealVar* r;
127 RooRealVar* clone;
128 while ((r = (RooRealVar*)it->Next()) != NULL) {
129 xVec->add(*r);
130 TString cloneName = TString::Format("%s%s", "mu__", r->GetName());
131 clone = (RooRealVar*)r->clone(cloneName.Data());
132 muVec->add(*clone);
133 if (fUseUpdates)
134 fPdfProp->AddMapping(*clone, *r);
135 }
136 if (fCovMatrix == NULL)
137 CreateCovMatrix(*xVec);
138 fPdf = new RooMultiVarGaussian("mvg", "MVG Proposal", *xVec, *muVec,
139 *fCovMatrix);
140 delete xVec;
141 delete muVec;
142 delete it;
143}
144
145////////////////////////////////////////////////////////////////////////////////
146
148{
149 Int_t size = xVec.getSize();
150 fCovMatrix = new TMatrixDSym(size);
151 RooRealVar* r;
152 for (Int_t i = 0; i < size; i++) {
153 r = (RooRealVar*)xVec.at(i);
154 Double_t range = r->getMax() - r->getMin();
155 (*fCovMatrix)(i,i) = range / fSigmaRangeDivisor;
156 }
157}
158
159////////////////////////////////////////////////////////////////////////////////
160
162{
163 if (fClues != NULL) {
164 if (fCluesOptions == NULL)
165 fCluesPdf = new RooNDKeysPdf("cluesPdf", "Clues PDF", *fVars, *fClues);
166 else
167 fCluesPdf = new RooNDKeysPdf("cluesPdf", "Clues PDF", *fVars, *fClues,
169 }
170}
171
172////////////////////////////////////////////////////////////////////////////////
173
175{
176 fUniformPdf = new RooUniform("uniform", "Uniform Proposal PDF",
177 RooArgSet(*fVars));
178}
ROOT::R::TRInterface & r
Definition: Object.C:4
static const Double_t SIGMA_RANGE_DIVISOR
static const Double_t DEFAULT_CLUES_FRAC
#define coutE(a)
Definition: RooMsgService.h:34
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
TMatrixTSym< Double_t > TMatrixDSym
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
TIterator * createIterator(Bool_t dir=kIterForward) const R__SUGGEST_ALTERNATIVE("begin()
TIterator-style iteration over contained elements.
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Multivariate Gaussian p.d.f.
Generic N-dimensional implementation of a kernel estimation p.d.f.
Definition: RooNDKeysPdf.h:48
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
PdfProposal is a concrete implementation of the ProposalFunction interface.
Definition: PdfProposal.h:30
virtual void AddMapping(RooRealVar &proposalParam, RooAbsReal &update)
specify a mapping between a parameter of the proposal function and a parameter of interest.
virtual void SetCacheSize(Int_t size)
Set how many points to generate each time we propose from a new point Default (and minimum) is 1.
Definition: PdfProposal.h:80
virtual void SetPdf(RooAbsPdf &pdf)
Set the PDF to be the proposal density function.
Definition: PdfProposal.h:49
virtual void SetOwnsPdf(Bool_t ownsPdf)
set whether we own the PDF that serves as the proposal density function By default,...
Definition: PdfProposal.h:91
ProposalFunction is an interface for all proposal functions that would be used with a Markov Chain Mo...
void CreateCovMatrix(RooArgList &xVec)
virtual ProposalFunction * GetProposalFunction()
const Option_t * fCluesOptions
Flat p.d.f.
Definition: RooUniform.h:24
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * Next()=0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition: TString.cxx:2311
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
@ InputArguments
Definition: RooGlobalFunc.h:68
RooConstVar & RooConst(Double_t val)
Namespace for the RooStats classes.
Definition: Asimov.h:20