1 |
#include <iostream>
|
2 |
#include <iomanip>
|
3 |
#include <sstream>
|
4 |
#include <cstdlib>
|
5 |
#include <cmath>
|
6 |
#include <algorithm>
|
7 |
|
8 |
|
9 |
#ifndef __CINT__
|
10 |
#include "RooGlobalFunc.h"
|
11 |
#endif
|
12 |
|
13 |
#include "RooProfileLL.h"
|
14 |
#include "RooAbsPdf.h"
|
15 |
#include "RooStats/HypoTestResult.h"
|
16 |
#include "RooRealVar.h"
|
17 |
#include "RooPlot.h"
|
18 |
#include "RooDataSet.h"
|
19 |
#include "RooTreeDataStore.h"
|
20 |
#include "TTree.h"
|
21 |
#include "TCanvas.h"
|
22 |
#include "TLine.h"
|
23 |
#include "TStopwatch.h"
|
24 |
|
25 |
#include "RooStats/ProfileLikelihoodCalculator.h"
|
26 |
#include "RooStats/MCMCCalculator.h"
|
27 |
#include "RooStats/UniformProposal.h"
|
28 |
#include "RooStats/FeldmanCousins.h"
|
29 |
#include "RooStats/NumberCountingPdfFactory.h"
|
30 |
#include "RooStats/ConfInterval.h"
|
31 |
#include "RooStats/PointSetInterval.h"
|
32 |
#include "RooStats/LikelihoodInterval.h"
|
33 |
#include "RooStats/LikelihoodIntervalPlot.h"
|
34 |
#include "RooStats/RooStatsUtils.h"
|
35 |
#include "RooStats/ModelConfig.h"
|
36 |
#include "RooStats/MCMCInterval.h"
|
37 |
#include "RooStats/MCMCIntervalPlot.h"
|
38 |
#include "RooStats/ProposalFunction.h"
|
39 |
#include "RooStats/ProposalHelper.h"
|
40 |
#include "RooFitResult.h"
|
41 |
|
42 |
#include <iostream>
|
43 |
#include <iomanip>
|
44 |
#include <cmath>
|
45 |
#include <fstream>
|
46 |
#include <time.h>
|
47 |
#include "ConfigFile.h"
|
48 |
|
49 |
// use this order for safety on library loading
|
50 |
using namespace RooFit ;
|
51 |
using namespace RooStats ;
|
52 |
|
53 |
|
54 |
|
55 |
double simpleProfile2(ConfigFile * config, string Type,
|
56 |
int dat, double bkg, double bkg_uncertainty, double sig, double
|
57 |
sig_uncertainty, double xsec, double ExpNsigLimit ) //absolute uncertainties!
|
58 |
{
|
59 |
TStopwatch t;
|
60 |
t.Start();
|
61 |
|
62 |
/////////////////////////////////////////
|
63 |
// The Model building stage
|
64 |
/////////////////////////////////////////
|
65 |
char * name = new char[1024];
|
66 |
RooWorkspace* wspace = new RooWorkspace();
|
67 |
sprintf(name,"Poisson::countingModel(obs[150,0,300],sum(s[%d,0,%d]*ratioSigEff[1.,0,2.],b[100,0,300]*ratioBkgEff[1.,0.,2.]))", (int)ExpNsigLimit,(int)(ExpNsigLimit*2.5));
|
68 |
wspace->factory(name); // counting model
|
69 |
//wspace->factory("Gaussian::sigConstraint(1,ratioSigEff,0.257)"); // 5% signal efficiency uncertainty
|
70 |
//wspace->factory("Gaussian::bkgConstraint(1,ratioBkgEff,0.136)"); // 10% background efficiency uncertainty
|
71 |
sprintf(name,"Gaussian::sigConstraint(1,ratioSigEff,%f)",sig_uncertainty/sig);
|
72 |
wspace->factory(name); // 5% signal efficiency uncertainty
|
73 |
sprintf(name,"Gaussian::bkgConstraint(1,ratioBkgEff,%f)",sig_uncertainty/sig);
|
74 |
wspace->factory(name); // 10% background efficiency uncertainty
|
75 |
wspace->factory("PROD::modelWithConstraints(countingModel,sigConstraint,bkgConstraint)"); // product of terms
|
76 |
wspace->Print();
|
77 |
|
78 |
RooAbsPdf* modelWithConstraints = wspace->pdf("modelWithConstraints"); // get the model
|
79 |
RooRealVar* obs = wspace->var("obs"); // get the observable
|
80 |
RooRealVar* s = wspace->var("s"); // get the signal we care about
|
81 |
RooRealVar* b = wspace->var("b"); // get the background and set it to a constant. Uncertainty included in ratioBkgEff
|
82 |
b->setConstant();
|
83 |
RooRealVar* ratioSigEff = wspace->var("ratioSigEff"); // get uncertaint parameter to constrain
|
84 |
RooRealVar* ratioBkgEff = wspace->var("ratioBkgEff"); // get uncertaint parameter to constrain
|
85 |
RooArgSet constrainedParams(*ratioSigEff, *ratioBkgEff); // need to constrain these in the fit (should change default behavior)
|
86 |
|
87 |
// Create an example dataset with 160 observed events
|
88 |
obs->setVal(dat);
|
89 |
b->setVal(bkg);
|
90 |
s->setVal(sig);
|
91 |
RooDataSet* data = new RooDataSet("exampleData", "exampleData", RooArgSet(*obs));
|
92 |
data->add(*obs);
|
93 |
|
94 |
RooArgSet all(*s, *ratioBkgEff, *ratioSigEff);
|
95 |
|
96 |
// not necessary
|
97 |
modelWithConstraints->fitTo(*data, RooFit::Constrain(RooArgSet(*ratioSigEff, *ratioBkgEff)));
|
98 |
|
99 |
// Now let's make some confidence intervals for s, our parameter of interest
|
100 |
RooArgSet paramOfInterest(*s);
|
101 |
|
102 |
ModelConfig modelConfig(new RooWorkspace());
|
103 |
modelConfig.SetPdf(*modelWithConstraints);
|
104 |
modelConfig.SetParametersOfInterest(paramOfInterest);
|
105 |
|
106 |
|
107 |
// First, let's use a Calculator based on the Profile Likelihood Ratio
|
108 |
//ProfileLikelihoodCalculator plc(*data, *modelWithConstraints, paramOfInterest);
|
109 |
ProfileLikelihoodCalculator plc(*data, modelConfig);
|
110 |
plc.SetTestSize(.05);
|
111 |
ConfInterval* lrint = plc.GetInterval(); // that was easy.
|
112 |
|
113 |
// Let's make a plot
|
114 |
TCanvas* dataCanvas = new TCanvas("dataCanvas");
|
115 |
dataCanvas->Divide(2,2);
|
116 |
|
117 |
dataCanvas->cd(3);
|
118 |
// --- Generate a toyMC sample from composite PDF ---
|
119 |
//RooDataSet *pseudodata = modelWithConstraints->generate(*b,2000) ;
|
120 |
// --- Perform extended ML fit of composite PDF to toy data ---
|
121 |
//modelWithConstraints.fitTo(*pseudodata,Extended());
|
122 |
|
123 |
// don't skip drawing, doesnt work without!!!
|
124 |
RooPlot* mesframe = b->frame();
|
125 |
modelWithConstraints->plotOn(mesframe);
|
126 |
mesframe->Draw();
|
127 |
|
128 |
dataCanvas->cd(1);
|
129 |
LikelihoodIntervalPlot plotInt((LikelihoodInterval*)lrint);
|
130 |
plotInt.SetTitle("Profile Likelihood Ratio and Posterior for S");
|
131 |
plotInt.Draw();
|
132 |
|
133 |
// Second, use a Calculator based on the Feldman Cousins technique
|
134 |
FeldmanCousins fc(*data, modelConfig);
|
135 |
fc.UseAdaptiveSampling(true);
|
136 |
fc.FluctuateNumDataEntries(false); // number counting analysis: dataset always has 1 entry with N events observed
|
137 |
fc.SetNBins(100); // number of points to test per parameter
|
138 |
fc.SetTestSize(.05);
|
139 |
// fc.SaveBeltToFile(true); // optional
|
140 |
ConfInterval* fcint = NULL;
|
141 |
fcint = fc.GetInterval(); // that was easy.
|
142 |
|
143 |
RooFitResult* fit = modelWithConstraints->fitTo(*data, Save(true));
|
144 |
|
145 |
// Third, use a Calculator based on Markov Chain monte carlo
|
146 |
// Before configuring the calculator, let's make a ProposalFunction
|
147 |
// that will achieve a high acceptance rate
|
148 |
ProposalHelper ph;
|
149 |
ph.SetVariables((RooArgSet&)fit->floatParsFinal());
|
150 |
ph.SetCovMatrix(fit->covarianceMatrix());
|
151 |
ph.SetUpdateProposalParameters(true);
|
152 |
ph.SetCacheSize(100);
|
153 |
ProposalFunction* pdfProp = ph.GetProposalFunction(); // that was easy
|
154 |
|
155 |
MCMCCalculator mc(*data, modelConfig);
|
156 |
mc.SetNumIters(20000); // steps to propose in the chain
|
157 |
mc.SetTestSize(.05); // 95% CL
|
158 |
mc.SetNumBurnInSteps(40); // ignore first N steps in chain as "burn in"
|
159 |
mc.SetProposalFunction(*pdfProp);
|
160 |
mc.SetLeftSideTailFraction(0.5); // find a "central" interval
|
161 |
MCMCInterval* mcInt = (MCMCInterval*)mc.GetInterval(); // that was easy
|
162 |
|
163 |
// Get Lower and Upper limits from Profile Calculator
|
164 |
double up = ((LikelihoodInterval*) lrint)->UpperLimit(*s);
|
165 |
cout << Type << ": d:"<<dat<<", b:"<<bkg<<"+-"<<bkg_uncertainty
|
166 |
<<"; s:"<<sig<<"+-"<<sig_uncertainty <<std::endl;
|
167 |
cout << "Profile lower limit on s = " << ((LikelihoodInterval*) lrint)->LowerLimit(*s) << endl;
|
168 |
cout << "Profile upper limit on s = " << up << endl;
|
169 |
config->add("RooSimpleProfile.signal."+Type+"UpperLimit", up);
|
170 |
config->add("RooSimpleProfile.xsec."+Type+"UpperLimit", up/sig * xsec);
|
171 |
|
172 |
// Get Lower and Upper limits from FeldmanCousins with profile construction
|
173 |
if (fcint != NULL) {
|
174 |
double fcul = ((PointSetInterval*) fcint)->UpperLimit(*s);
|
175 |
double fcll = ((PointSetInterval*) fcint)->LowerLimit(*s);
|
176 |
cout << "FC lower limit on s = " << fcll << endl;
|
177 |
cout << "FC upper limit on s = " << fcul << endl;
|
178 |
config->add("RooFC.signal."+Type+"UpperLimit", fcul);
|
179 |
config->add("RooFC.xsec."+Type+"UpperLimit", fcul/sig * xsec);
|
180 |
//TLine* fcllLine = new TLine(fcll, 0, fcll, 1);
|
181 |
//TLine* fculLine = new TLine(fcul, 0, fcul, 1);
|
182 |
//fcllLine->SetLineColor(kRed);
|
183 |
//fculLine->SetLineColor(kRed);
|
184 |
//fcllLine->Draw("same");
|
185 |
//fculLine->Draw("same");
|
186 |
dataCanvas->Update();
|
187 |
}
|
188 |
/*
|
189 |
|
190 |
// Plot MCMC interval and print some statistics
|
191 |
MCMCIntervalPlot mcPlot(*mcInt);
|
192 |
mcPlot.SetLineColor(kMagenta);
|
193 |
mcPlot.SetLineWidth(2);
|
194 |
mcPlot.Draw("same");
|
195 |
*/
|
196 |
double mcul = mcInt->UpperLimit(*s);
|
197 |
double mcll = mcInt->LowerLimit(*s);
|
198 |
cout << "MCMC lower limit on s = " << mcll << endl;
|
199 |
cout << "MCMC upper limit on s = " << mcul << endl;
|
200 |
cout << "MCMC Actual confidence level: "
|
201 |
<< mcInt->GetActualConfidenceLevel() << endl;
|
202 |
|
203 |
config->add("RooMCMC.signal."+Type+"UpperLimit", mcul);
|
204 |
config->add("RooMCMC.xsec."+Type+"UpperLimit", mcul/sig * xsec);
|
205 |
/*
|
206 |
// 3-d plot of the parameter points
|
207 |
dataCanvas->cd(2);
|
208 |
// also plot the points in the markov chain
|
209 |
TTree& chain = ((RooTreeDataStore*) mcInt->GetChainAsDataSet()->store())->tree();
|
210 |
chain.SetMarkerStyle(6);
|
211 |
chain.SetMarkerColor(kRed);
|
212 |
chain.Draw("s:ratioSigEff:ratioBkgEff","weight_MarkovChain_local_","box"); // 3-d box proporional to posterior
|
213 |
|
214 |
// the points used in the profile construction
|
215 |
TTree& parameterScan = ((RooTreeDataStore*) fc.GetPointsToScan()->store())->tree();
|
216 |
parameterScan.SetMarkerStyle(24);
|
217 |
parameterScan.Draw("s:ratioSigEff:ratioBkgEff","","same");
|
218 |
|
219 |
delete wspace;
|
220 |
delete lrint;
|
221 |
delete mcInt;
|
222 |
delete fcint;
|
223 |
delete data;
|
224 |
*/
|
225 |
/// print timing info
|
226 |
t.Stop();
|
227 |
t.Print();
|
228 |
|
229 |
return up;
|
230 |
}
|
231 |
|
232 |
void rooStat(int argc, char *argv[])
|
233 |
{
|
234 |
//myPDF();
|
235 |
//rooStatOrig();
|
236 |
if (argc<1){
|
237 |
std::cerr<<"Error: Expected '"<<argv[0]<<" <limit config>'!"<<std::endl;
|
238 |
exit(1);
|
239 |
}
|
240 |
|
241 |
for (int file=1; file<argc; ++file){
|
242 |
std::cout << "opening: "<<argv[file]<<std::endl;
|
243 |
ConfigFile config(argv[file]);
|
244 |
int data = config.read<int>("data");
|
245 |
double bkg = config.read<double>("background");
|
246 |
double bkgUncert = config.read<double>("background.uncertainty");
|
247 |
double sig = config.read<double>("signal");
|
248 |
double sigUncert = config.read<double>("signal.uncertainty");
|
249 |
double xsec = config.read<double>("Xsection",1.0);
|
250 |
double ExpNsigLimit = config.read<double>("ExpNsigLimit",100.0);
|
251 |
|
252 |
simpleProfile2(&config, "Obs", data, bkg, bkgUncert, sig, sigUncert, xsec, ExpNsigLimit);
|
253 |
simpleProfile2(&config, "Exp", (int)bkg, bkg, bkgUncert, sig, sigUncert, xsec, ExpNsigLimit);
|
254 |
|
255 |
//write stuff:
|
256 |
time_t rawtime;
|
257 |
struct tm * timeinfo;
|
258 |
time ( &rawtime );
|
259 |
timeinfo = localtime ( &rawtime );
|
260 |
ofstream ofile;
|
261 |
ofile.open (argv[file]);
|
262 |
ofile << config.getComment() << asctime (timeinfo)
|
263 |
<< config.getComment()<< "\n"
|
264 |
<< config;
|
265 |
}
|
266 |
}
|
267 |
|
268 |
|
269 |
int main(int argc, char *argv[])
|
270 |
{
|
271 |
rooStat(argc,argv);
|
272 |
}
|