Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
ProposalHelper.h
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 #ifndef RooStats_ProposalHelper
13 #define RooStats_ProposalHelper
14 
15 #include "Rtypes.h"
18 #include "RooStats/PdfProposal.h"
19 #include "RooArgSet.h"
20 #include "RooMsgService.h"
21 #include "RooRealVar.h"
22 #include "TObject.h"
23 
24 
25 
26 namespace RooStats {
27 
28  class ProposalHelper : public TObject {
29 
30  public:
31  ProposalHelper();
32 
33  // Set the PDF to be the proposal density function
34  virtual void SetPdf(RooAbsPdf& pdf) { fPdf = &pdf; }
35  // Set the bank of clues to add to the current proposal density function
36  virtual void SetClues(RooDataSet& clues) { fClues = &clues; }
37 
38  // Get the ProposalFunction that we've been designing
39  virtual ProposalFunction* GetProposalFunction();
40 
41  virtual void SetCacheSize(Int_t size)
42  {
43  if (size > 0)
44  fCacheSize = size;
45  else
46  coutE(Eval) << "Warning: Requested non-positive cache size: " <<
47  size << ". Cache size unchanged." << std::endl;
48  }
49 
50  virtual void SetUpdateProposalParameters(Bool_t updateParams)
51  { fUseUpdates = updateParams; }
52 
53  virtual void SetVariables(RooArgList& vars)
54  { fVars = &vars; }
55 
56  virtual void SetVariables(const RooArgList& vars)
57  { fVars = new RooArgList(vars); fOwnsVars = kTRUE; }
58 
59  // set what fraction of the proposal density function should come from
60  // a uniform proposal distribution
61  virtual void SetUniformFraction(Double_t uniFrac) { fUniFrac = uniFrac; }
62 
63  // set what fraction of the proposal density function should come from
64  // the bank of clues
65  virtual void SetCluesFraction(Double_t cluesFrac) { fCluesFrac = cluesFrac; }
66 
67  // set the covariance matrix to use for a multi-variate Gaussian proposal
68  virtual void SetCovMatrix(const TMatrixDSym& covMatrix)
69  { fCovMatrix = new TMatrixDSym(covMatrix); }
70 
71  // set what divisor we will use when dividing the range of a variable to
72  // determine the width of the proposal function for each dimension
73  // e.g. divisor = 6 for sigma = 1/6th
74  virtual void SetWidthRangeDivisor(Double_t divisor)
75  { if (divisor > 0.) fSigmaRangeDivisor = divisor; }
76 
77  // set the option string to pass to the RooNDKeysPdf constructor
78  // if the bank of clues pdf is being automatically generated by this
79  // ProposalHelper
80  virtual void SetCluesOptions(const Option_t* options)
81  { if (options != NULL) fCluesOptions = options; }
82 
83  virtual void SetVariables(RooArgSet& vars)
84  {
85  RooArgList* argList = new RooArgList(vars);
86  SetVariables(*argList);
87  fOwnsVars = kTRUE;
88  }
89 
90  virtual ~ProposalHelper()
91  {
92  if (fOwnsPdfProp) delete fPdfProp;
93  if (fOwnsPdf) delete fPdf;
94  if (fOwnsCluesPdf) delete fCluesPdf;
95  if (fOwnsVars) delete fVars;
96  delete fCovMatrix;
97  delete fUniformPdf;
98  }
99 
100  protected:
101  RooAbsPdf* fPdf; // the main proposal density function
102  RooAbsPdf* fCluesPdf; // proposal dens. func. with clues for certain points
103  RooAbsPdf* fUniformPdf; // uniform proposal dens. func.
104  RooDataSet* fClues; // data set of clues
105  TMatrixDSym* fCovMatrix; // covariance matrix for multi var gaussian pdf
106  PdfProposal* fPdfProp; // the PdfProposal we are (probably) going to return
107  RooArgList* fVars; // the RooRealVars to generate proposals for
108  Int_t fCacheSize; // for generating proposals from PDFs
109  Double_t fSigmaRangeDivisor; // range divisor to get sigma for each variable
110  Double_t fUniFrac; // what fraction of the PDF integral is uniform
111  Double_t fCluesFrac; // what fraction of the PDF integral comes from clues
112  Bool_t fOwnsPdfProp; // whether we own the PdfProposal; equivalent to:
113  // !(whether we have returned it in GetProposalFunction)
114  Bool_t fOwnsPdf; // whether we created (and own) the main pdf
115  Bool_t fOwnsCluesPdf; // whether we created (and own) the clues pdf
116  Bool_t fOwnsVars; // whether we own fVars
117  Bool_t fUseUpdates; // whether to set updates for proposal params in PdfProposal
118  const Option_t* fCluesOptions; // option string for clues RooNDKeysPdf
119 
120  void CreatePdf();
121  void CreateCluesPdf();
122  void CreateUniformPdf();
123  void CreateCovMatrix(RooArgList& xVec);
124 
125  ClassDef(ProposalHelper,1)
126  };
127 }
128 #endif