00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035 #ifndef ROOT_Rtypes
00036 #include "Rtypes.h"
00037 #endif
00038 #ifndef ROO_GLOBAL_FUNC
00039 #include "RooGlobalFunc.h"
00040 #endif
00041 #ifndef ROO_ABS_REAL
00042 #include "RooAbsReal.h"
00043 #endif
00044 #ifndef ROO_ARG_SET
00045 #include "RooArgSet.h"
00046 #endif
00047 #ifndef ROO_ARG_LIST
00048 #include "RooArgList.h"
00049 #endif
00050 #ifndef ROOSTATS_ModelConfig
00051 #include "RooStats/ModelConfig.h"
00052 #endif
00053 #ifndef RooStats_RooStatsUtils
00054 #include "RooStats/RooStatsUtils.h"
00055 #endif
00056 #ifndef ROOSTATS_MCMCCalculator
00057 #include "RooStats/MCMCCalculator.h"
00058 #endif
00059 #ifndef ROOSTATS_MetropolisHastings
00060 #include "RooStats/MetropolisHastings.h"
00061 #endif
00062 #ifndef ROOSTATS_MarkovChain
00063 #include "RooStats/MarkovChain.h"
00064 #endif
00065 #ifndef RooStats_MCMCInterval
00066 #include "RooStats/MCMCInterval.h"
00067 #endif
00068 #ifndef ROOT_TIterator
00069 #include "TIterator.h"
00070 #endif
00071 #ifndef ROOSTATS_UniformProposal
00072 #include "RooStats/UniformProposal.h"
00073 #endif
00074 #ifndef ROOSTATS_PdfProposal
00075 #include "RooStats/PdfProposal.h"
00076 #endif
00077 #ifndef ROO_PROD_PDF
00078 #include "RooProdPdf.h"
00079 #endif
00080
00081 ClassImp(RooStats::MCMCCalculator);
00082
00083 using namespace RooFit;
00084 using namespace RooStats;
00085
00086
00087 MCMCCalculator::MCMCCalculator() :
00088 fPropFunc(0),
00089 fPdf(0),
00090 fPriorPdf(0),
00091 fData(0),
00092 fAxes(0)
00093 {
00094 fNumIters = 0;
00095 fNumBurnInSteps = 0;
00096 fNumBins = 0;
00097 fUseKeys = kFALSE;
00098 fUseSparseHist = kFALSE;
00099 fSize = -1;
00100 fIntervalType = MCMCInterval::kShortest;
00101 fLeftSideTF = -1;
00102 fEpsilon = -1;
00103 fDelta = -1;
00104 }
00105
00106
00107
00108 MCMCCalculator::MCMCCalculator(RooAbsData& data, const ModelConfig & model) :
00109 fPropFunc(0),
00110 fData(&data),
00111 fAxes(0)
00112 {
00113 SetModel(model);
00114 SetupBasicUsage();
00115 }
00116
00117 void MCMCCalculator::SetModel(const ModelConfig & model) {
00118
00119 fPdf = model.GetPdf();
00120 fPriorPdf = model.GetPriorPdf();
00121 fPOI.removeAll();
00122 fNuisParams.removeAll();
00123 if (model.GetParametersOfInterest())
00124 fPOI.add(*model.GetParametersOfInterest());
00125 if (model.GetNuisanceParameters())
00126 fNuisParams.add(*model.GetNuisanceParameters());
00127 }
00128
00129
00130
00131
00132
00133 void MCMCCalculator::SetupBasicUsage()
00134 {
00135 fPropFunc = 0;
00136 fNumIters = 10000;
00137 fNumBurnInSteps = 40;
00138 fNumBins = 50;
00139 fUseKeys = kFALSE;
00140 fUseSparseHist = kFALSE;
00141 SetTestSize(0.05);
00142 fIntervalType = MCMCInterval::kShortest;
00143 fLeftSideTF = -1;
00144 fEpsilon = -1;
00145 fDelta = -1;
00146 }
00147
00148 void MCMCCalculator::SetLeftSideTailFraction(Double_t a)
00149 {
00150 if (a < 0 || a > 1) {
00151 coutE(InputArguments) << "MCMCCalculator::SetLeftSideTailFraction: "
00152 << "Fraction must be in the range [0, 1]. "
00153 << a << "is not allowed." << endl;
00154 return;
00155 }
00156
00157 fLeftSideTF = a;
00158 fIntervalType = MCMCInterval::kTailFraction;
00159 }
00160
00161 MCMCInterval* MCMCCalculator::GetInterval() const
00162 {
00163
00164
00165 if (!fData || !fPdf ) return 0;
00166 if (fPOI.getSize() == 0) return 0;
00167
00168 if (fSize < 0) {
00169 coutE(InputArguments) << "MCMCCalculator::GetInterval: "
00170 << "Test size/Confidence level not set. Returning NULL." << endl;
00171 return NULL;
00172 }
00173
00174
00175 bool useDefaultPropFunc = (fPropFunc == 0);
00176 bool usePriorPdf = (fPriorPdf != 0);
00177 if (useDefaultPropFunc) fPropFunc = new UniformProposal();
00178
00179
00180 RooAbsPdf * prodPdf = fPdf;
00181 if (usePriorPdf) {
00182 TString prodName = TString("product_") + TString(fPdf->GetName()) + TString("_") + TString(fPriorPdf->GetName() );
00183 prodPdf = new RooProdPdf(prodName,prodName,RooArgList(*fPdf,*fPriorPdf) );
00184 }
00185
00186 RooArgSet* constrainedParams = prodPdf->getParameters(*fData);
00187 RooAbsReal* nll = prodPdf->createNLL(*fData, Constrain(*constrainedParams));
00188 delete constrainedParams;
00189
00190 RooArgSet* params = nll->getParameters(*fData);
00191 RemoveConstantParameters(params);
00192 if (fNumBins > 0) {
00193 SetBins(*params, fNumBins);
00194 SetBins(fPOI, fNumBins);
00195 if (dynamic_cast<PdfProposal*>(fPropFunc)) {
00196 RooArgSet* proposalVars = ((PdfProposal*)fPropFunc)->GetPdf()->
00197 getParameters((RooAbsData*)NULL);
00198 SetBins(*proposalVars, fNumBins);
00199 }
00200 }
00201
00202 MetropolisHastings mh;
00203 mh.SetFunction(*nll);
00204 mh.SetType(MetropolisHastings::kLog);
00205 mh.SetSign(MetropolisHastings::kNegative);
00206 mh.SetParameters(*params);
00207 mh.SetProposalFunction(*fPropFunc);
00208 mh.SetNumIters(fNumIters);
00209
00210 MarkovChain* chain = mh.ConstructChain();
00211
00212 TString name = TString("MCMCInterval_") + TString(GetName() );
00213 MCMCInterval* interval = new MCMCInterval(name, fPOI, *chain);
00214 if (fAxes != NULL)
00215 interval->SetAxes(*fAxes);
00216 if (fNumBurnInSteps > 0)
00217 interval->SetNumBurnInSteps(fNumBurnInSteps);
00218 interval->SetUseKeys(fUseKeys);
00219 interval->SetUseSparseHist(fUseSparseHist);
00220 interval->SetIntervalType(fIntervalType);
00221 if (fIntervalType == MCMCInterval::kTailFraction)
00222 interval->SetLeftSideTailFraction(fLeftSideTF);
00223 if (fEpsilon >= 0)
00224 interval->SetEpsilon(fEpsilon);
00225 if (fDelta >= 0)
00226 interval->SetDelta(fDelta);
00227 interval->SetConfidenceLevel(1.0 - fSize);
00228
00229 if (useDefaultPropFunc) delete fPropFunc;
00230 if (usePriorPdf) delete prodPdf;
00231 delete nll;
00232 delete params;
00233
00234 return interval;
00235 }