00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #ifndef __CINT__
00016 #include "RooGlobalFunc.h"
00017 #endif
00018
00019 #include "RooProfileLL.h"
00020 #include "RooAbsPdf.h"
00021 #include "RooStats/HypoTestResult.h"
00022 #include "RooRealVar.h"
00023 #include "RooPlot.h"
00024 #include "RooDataSet.h"
00025 #include "RooTreeDataStore.h"
00026 #include "TTree.h"
00027 #include "TCanvas.h"
00028 #include "TLine.h"
00029 #include "TStopwatch.h"
00030
00031 #include "RooStats/ProfileLikelihoodCalculator.h"
00032 #include "RooStats/MCMCCalculator.h"
00033 #include "RooStats/UniformProposal.h"
00034 #include "RooStats/FeldmanCousins.h"
00035 #include "RooStats/NumberCountingPdfFactory.h"
00036 #include "RooStats/ConfInterval.h"
00037 #include "RooStats/PointSetInterval.h"
00038 #include "RooStats/LikelihoodInterval.h"
00039 #include "RooStats/LikelihoodIntervalPlot.h"
00040 #include "RooStats/RooStatsUtils.h"
00041 #include "RooStats/ModelConfig.h"
00042 #include "RooStats/MCMCInterval.h"
00043 #include "RooStats/MCMCIntervalPlot.h"
00044 #include "RooStats/ProposalFunction.h"
00045 #include "RooStats/ProposalHelper.h"
00046 #include "RooFitResult.h"
00047
00048
00049
00050 using namespace RooFit ;
00051 using namespace RooStats ;
00052
00053
00054 void rs101_limitexample()
00055 {
00056
00057
00058
00059
00060
00061 TStopwatch t;
00062 t.Start();
00063
00064
00065
00066
00067 RooWorkspace* wspace = new RooWorkspace();
00068 wspace->factory("Poisson::countingModel(obs[150,0,300], sum(s[50,0,120]*ratioSigEff[1.,0,2.],b[100,0,300]*ratioBkgEff[1.,0.,2.]))");
00069
00070
00071 wspace->factory("Gaussian::sigConstraint(1,ratioSigEff,0.05)");
00072 wspace->factory("Gaussian::bkgConstraint(1,ratioBkgEff,0.1)");
00073 wspace->factory("PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)");
00074 wspace->Print();
00075
00076 RooAbsPdf* modelWithConstraints = wspace->pdf("modelWithConstraints");
00077 RooRealVar* obs = wspace->var("obs");
00078 RooRealVar* s = wspace->var("s");
00079 RooRealVar* b = wspace->var("b");
00080 b->setConstant();
00081 RooRealVar* ratioSigEff = wspace->var("ratioSigEff");
00082 RooRealVar* ratioBkgEff = wspace->var("ratioBkgEff");
00083 RooArgSet constrainedParams(*ratioSigEff, *ratioBkgEff);
00084
00085
00086 obs->setVal(160.);
00087 RooDataSet* data = new RooDataSet("exampleData", "exampleData", RooArgSet(*obs));
00088 data->add(*obs);
00089
00090 RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);
00091
00092
00093 modelWithConstraints->fitTo(*data, RooFit::Constrain(RooArgSet(*ratioSigEff, *ratioBkgEff)));
00094
00095
00096 RooArgSet paramOfInterest(*s);
00097
00098 ModelConfig modelConfig(new RooWorkspace());
00099 modelConfig.SetPdf(*modelWithConstraints);
00100 modelConfig.SetParametersOfInterest(paramOfInterest);
00101
00102
00103
00104
00105 ProfileLikelihoodCalculator plc(*data, modelConfig);
00106 plc.SetTestSize(.05);
00107 ConfInterval* lrint = plc.GetInterval();
00108
00109
00110 TCanvas* dataCanvas = new TCanvas("dataCanvas");
00111 dataCanvas->Divide(2,1);
00112
00113 dataCanvas->cd(1);
00114 LikelihoodIntervalPlot plotInt((LikelihoodInterval*)lrint);
00115 plotInt.SetTitle("Profile Likelihood Ratio and Posterior for S");
00116 plotInt.Draw();
00117
00118
00119 FeldmanCousins fc(*data, modelConfig);
00120 fc.UseAdaptiveSampling(true);
00121 fc.FluctuateNumDataEntries(false);
00122 fc.SetNBins(100);
00123 fc.SetTestSize(.05);
00124
00125 ConfInterval* fcint = NULL;
00126 fcint = fc.GetInterval();
00127
00128 RooFitResult* fit = modelWithConstraints->fitTo(*data, Save(true));
00129
00130
00131
00132
00133 ProposalHelper ph;
00134 ph.SetVariables((RooArgSet&)fit->floatParsFinal());
00135 ph.SetCovMatrix(fit->covarianceMatrix());
00136 ph.SetUpdateProposalParameters(true);
00137 ph.SetCacheSize(100);
00138 ProposalFunction* pdfProp = ph.GetProposalFunction();
00139
00140 MCMCCalculator mc(*data, modelConfig);
00141 mc.SetNumIters(20000);
00142 mc.SetTestSize(.05);
00143 mc.SetNumBurnInSteps(40);
00144 mc.SetProposalFunction(*pdfProp);
00145 mc.SetLeftSideTailFraction(0.5);
00146 MCMCInterval* mcInt = (MCMCInterval*)mc.GetInterval();
00147
00148
00149
00150 cout << "Profile lower limit on s = " << ((LikelihoodInterval*) lrint)->LowerLimit(*s) << endl;
00151 cout << "Profile upper limit on s = " << ((LikelihoodInterval*) lrint)->UpperLimit(*s) << endl;
00152
00153
00154 if (fcint != NULL) {
00155 double fcul = ((PointSetInterval*) fcint)->UpperLimit(*s);
00156 double fcll = ((PointSetInterval*) fcint)->LowerLimit(*s);
00157 cout << "FC lower limit on s = " << fcll << endl;
00158 cout << "FC upper limit on s = " << fcul << endl;
00159 TLine* fcllLine = new TLine(fcll, 0, fcll, 1);
00160 TLine* fculLine = new TLine(fcul, 0, fcul, 1);
00161 fcllLine->SetLineColor(kRed);
00162 fculLine->SetLineColor(kRed);
00163 fcllLine->Draw("same");
00164 fculLine->Draw("same");
00165 dataCanvas->Update();
00166 }
00167
00168
00169 MCMCIntervalPlot mcPlot(*mcInt);
00170 mcPlot.SetLineColor(kMagenta);
00171 mcPlot.SetLineWidth(2);
00172 mcPlot.Draw("same");
00173
00174 double mcul = mcInt->UpperLimit(*s);
00175 double mcll = mcInt->LowerLimit(*s);
00176 cout << "MCMC lower limit on s = " << mcll << endl;
00177 cout << "MCMC upper limit on s = " << mcul << endl;
00178 cout << "MCMC Actual confidence level: "
00179 << mcInt->GetActualConfidenceLevel() << endl;
00180
00181
00182 dataCanvas->cd(2);
00183
00184 TTree& chain = ((RooTreeDataStore*) mcInt->GetChainAsDataSet()->store())->tree();
00185 chain.SetMarkerStyle(6);
00186 chain.SetMarkerColor(kRed);
00187 chain.Draw("s:ratioSigEff:ratioBkgEff","weight_MarkovChain_local_","box");
00188
00189
00190 TTree& parameterScan = ((RooTreeDataStore*) fc.GetPointsToScan()->store())->tree();
00191 parameterScan.SetMarkerStyle(24);
00192 parameterScan.Draw("s:ratioSigEff:ratioBkgEff","","same");
00193
00194 delete wspace;
00195 delete lrint;
00196 delete mcInt;
00197 delete fcint;
00198 delete data;
00199
00200
00201 t.Stop();
00202 t.Print();
00203 }