00001 // @(#)root/roostats:$Id: UniformProposal.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 UniformProposal is a concrete implementation of the ProposalFunction interface 00017 for use with a Markov Chain Monte Carlo algorithm. This proposal function is 00018 a uniformly random distribution over the parameter space. The proposal 00019 ignores the current point when it proposes a new point. The proposal 00020 function is symmetric, though it may not cause a MetropolisHastings run to 00021 converge as quickly as other proposal functions. 00022 </p> 00023 END_HTML 00024 */ 00025 // 00026 00027 #ifndef ROOT_Rtypes 00028 #include "Rtypes.h" 00029 #endif 00030 00031 #ifndef RooStats_RooStatsUtils 00032 #include "RooStats/RooStatsUtils.h" 00033 #endif 00034 #ifndef ROOSTATS_UniformProposal 00035 #include "RooStats/UniformProposal.h" 00036 #endif 00037 #ifndef ROO_ARG_SET 00038 #include "RooArgSet.h" 00039 #endif 00040 #ifndef ROO_MSG_SERVICE 00041 #include "RooMsgService.h" 00042 #endif 00043 #ifndef ROO_REAL_VAR 00044 #include "RooRealVar.h" 00045 #endif 00046 #ifndef ROOT_TIterator 00047 #include "TIterator.h" 00048 #endif 00049 00050 ClassImp(RooStats::UniformProposal); 00051 00052 using namespace RooFit; 00053 using namespace RooStats; 00054 00055 // Populate xPrime with a new proposed point 00056 void UniformProposal::Propose(RooArgSet& xPrime, RooArgSet& /* x */) 00057 { 00058 // kbelasco: remember xPrime and x have not been checked for containing 00059 // only RooRealVars 00060 RooStats::RandomizeCollection(xPrime); 00061 } 00062 00063 // Determine whether or not the proposal density is symmetric for 00064 // points x1 and x2 - that is, whether the probabilty of reaching x2 00065 // from x1 is equal to the probability of reaching x1 from x2 00066 Bool_t UniformProposal::IsSymmetric(RooArgSet& /* x1 */ , RooArgSet& /* x2 */) 00067 { 00068 return true; 00069 } 00070 00071 // Return the probability of proposing the point x1 given the starting 00072 // point x2 00073 Double_t UniformProposal::GetProposalDensity(RooArgSet& /* x1 */, 00074 RooArgSet& x2) 00075 { 00076 // For a uniform proposal, all points have equal probability and the 00077 // value of the proposal density function is: 00078 // 1 / (N-dimensional volume of interval) 00079 Double_t volume = 1.0; 00080 TIterator* it = x2.createIterator(); 00081 RooRealVar* var; 00082 while ((var = (RooRealVar*)it->Next()) != NULL) 00083 volume *= (var->getMax() - var->getMin()); 00084 delete it; 00085 return 1.0 / volume; 00086 }