PdfProposal.cxx

Go to the documentation of this file.
00001 // @(#)root/roostats:$Id: PdfProposal.cxx 34109 2010-06-24 15:00:16Z moneta $
00002 // Authors: Kevin Belasco        17/06/2009
00003 // Authors: Kyle Cranmer         17/06/2009
00004 /*************************************************************************
00005  * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //_________________________________________________
00013 /*
00014 BEGIN_HTML
00015 <p>
00016 PdfProposal is a concrete implementation of the ProposalFunction interface.
00017 It proposes points across the parameter space in the distribution of the
00018 given PDF.
00019 </p>
00020 <p>
00021 To make Propose(xPrime, x) dependent on x, configure with
00022 PdfProposal::AddMapping(varToUpdate, valueToUse).  For example, suppose we have:
00023 </p>
00024 
00025 <p>
00026 // our parameter
00027 RooRealVar p("p", "p", 5, 0, 10);
00028 
00029 // create mean and sigma for gaussian proposal function
00030 RooRealVar meanP("meanP", "meanP", 0, 10);
00031 RooRealVar sigma("sigma", "sigma", 1, 0, 5);
00032 RooGaussian pGaussian("pGaussian", "pGaussian", p, meanP, sigma);
00033 
00034 // configure proposal function
00035 PdfProposal pdfProposal(pGaussian);
00036 pdfProposal.AddMapping(meanP, p); // each call of Propose(xPrime, x), meanP in
00037                                   // the proposal function will be updated to
00038                                   // the value of p in x.  this will center the
00039                                   // proposal function about x's p when
00040                                   // proposing for xPrime
00041 
00042 // To improve performance, PdfProposal has the ability to cache a specified
00043 // number of proposals. If you don't call this function, the default cache size
00044 // is 1, which can be slow.
00045 pdfProposal.SetCacheSize(desiredCacheSize);
00046 </p>
00047 
00048 <p>
00049 PdfProposal currently uses a fixed cache size. Adaptive caching methods are in the works
00050 for future versions.
00051 </p>
00052 
00053 END_HTML
00054 */
00055 //_________________________________________________
00056 
00057 #ifndef ROOT_Rtypes
00058 #include "Rtypes.h"
00059 #endif
00060 
00061 #ifndef ROOSTATS_PdfProposal
00062 #include "RooStats/PdfProposal.h"
00063 #endif
00064 #ifndef RooStats_RooStatsUtils
00065 #include "RooStats/RooStatsUtils.h"
00066 #endif
00067 #ifndef ROO_ARG_SET
00068 #include "RooArgSet.h"
00069 #endif
00070 #ifndef ROO_DATA_SET
00071 #include "RooDataSet.h"
00072 #endif
00073 #ifndef ROO_ABS_PDF
00074 #include "RooAbsPdf.h"
00075 #endif
00076 #ifndef ROO_MSG_SERVICE
00077 #include "RooMsgService.h"
00078 #endif
00079 #ifndef ROO_REAL_VAR
00080 #include "RooRealVar.h"
00081 #endif
00082 #ifndef ROOT_TIterator
00083 #include "TIterator.h"
00084 #endif
00085 
00086 #include <map>
00087 
00088 ClassImp(RooStats::PdfProposal);
00089 
00090 using namespace RooFit;
00091 using namespace RooStats;
00092 
00093 // By default, PdfProposal does NOT own the PDF that serves as the
00094 // proposal density function
00095 PdfProposal::PdfProposal() : ProposalFunction()
00096 {
00097    fPdf = NULL;
00098    fOwnsPdf = kFALSE;
00099    fCacheSize = 1;
00100    fCachePosition = 0;
00101    fCache = NULL;
00102 }
00103 
00104 // By default, PdfProposal does NOT own the PDF that serves as the
00105 // proposal density function
00106 PdfProposal::PdfProposal(RooAbsPdf& pdf) : ProposalFunction()
00107 {
00108    fPdf = &pdf;
00109    fOwnsPdf = kFALSE;
00110    fCacheSize = 1;
00111    fCachePosition = 0;
00112    fCache = NULL;
00113 }
00114 
00115 Bool_t PdfProposal::Equals(RooArgSet& x1, RooArgSet& x2)
00116 {
00117    if (x1.equals(x2)) {
00118       TIterator* it = x1.createIterator();
00119       RooRealVar* r;
00120       while ((r = (RooRealVar*)it->Next()) != NULL)
00121          if (r->getVal() != x2.getRealValue(r->GetName())) {
00122             delete it;
00123             return kFALSE;
00124          }
00125       delete it;
00126       return kTRUE;
00127    }
00128    return kFALSE;
00129 }
00130 
00131 // Populate xPrime with a new proposed point
00132 void PdfProposal::Propose(RooArgSet& xPrime, RooArgSet& x)
00133 {
00134    if (fLastX.getSize() == 0) {
00135       // fLastX not yet initialized
00136       fLastX.addClone(x);
00137       // generate initial cache
00138       RooStats::SetParameters(&x, &fMaster);
00139       if (fMap.size() > 0) {
00140          for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
00141             fIt->first->setVal(fIt->second->getVal(&x));
00142       }
00143       fCache = fPdf->generate(xPrime, fCacheSize);
00144    }
00145 
00146    Bool_t moved = false;
00147    if (fMap.size() > 0) {
00148       moved = !Equals(fLastX, x);
00149 
00150       // if we've moved, set the values of the variables in the PDF to the
00151       // corresponding values of the variables in x, according to the
00152       // mappings (i.e. let the variables in x set the given values for the
00153       // PDF that will generate xPrime)
00154       if (moved) {
00155          // update the pdf parameters
00156          RooStats::SetParameters(&x, &fMaster);
00157 
00158          for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
00159             fIt->first->setVal(fIt->second->getVal(&x));
00160 
00161          // save the new x in fLastX
00162          RooStats::SetParameters(&x, &fLastX);
00163       }
00164    }
00165 
00166    // generate new cache if necessary
00167    if (moved || fCachePosition >= fCacheSize) {
00168       delete fCache;
00169       fCache = fPdf->generate(xPrime, fCacheSize);
00170       fCachePosition = 0;
00171    }
00172 
00173    const RooArgSet* proposal = fCache->get(fCachePosition);
00174    fCachePosition++;
00175    RooStats::SetParameters(proposal, &xPrime);
00176 }
00177 
00178 // Determine whether or not the proposal density is symmetric for
00179 // points x1 and x2 - that is, whether the probabilty of reaching x2
00180 // from x1 is equal to the probability of reaching x1 from x2
00181 Bool_t PdfProposal::IsSymmetric(RooArgSet& /* x1 */, RooArgSet& /* x2 */)
00182 {
00183    // kbelasco: is there a better way to do this?
00184    return false;
00185 }
00186 
00187 // Return the probability of proposing the point x1 given the starting
00188 // point x2
00189 Double_t PdfProposal::GetProposalDensity(RooArgSet& x1, RooArgSet& x2)
00190 {
00191    RooStats::SetParameters(&x2, &fMaster);
00192    for (fIt = fMap.begin(); fIt != fMap.end(); fIt++)
00193       fIt->first->setVal(fIt->second->getVal(&x2));
00194    RooArgSet* temp = fPdf->getObservables(x1);
00195    RooStats::SetParameters(&x1, temp);
00196    delete temp;
00197    return fPdf->getVal(&x1); // could just as well use x2
00198 }
00199 
00200 void PdfProposal::AddMapping(RooRealVar& proposalParam, RooAbsReal& update)
00201 {
00202    fMaster.add(*update.getParameters((RooAbsData*)NULL));
00203    if (update.getParameters((RooAbsData*)NULL)->getSize() == 0)
00204       fMaster.add(update);
00205    fMap.insert(pair<RooRealVar*, RooAbsReal*>(&proposalParam, &update));
00206 }

Generated on Tue Jul 5 15:14:36 2011 for ROOT_528-00b_version by  doxygen 1.5.1