ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/cbrown/AnalysisFramework/Plotting/Modules/LimitCalculation.C
Revision: 1.4
Committed: Wed Jul 20 12:26:11 2011 UTC (13 years, 9 months ago) by buchmann
Content type: text/plain
Branch: MAIN
Changes since 1.3: +7 -4 lines
Log Message:
Making limit calculation more customizable

File Contents

# Content
1 #include <iostream>
2 #include <vector>
3 #include <sys/stat.h>
4
5 #include <TCut.h>
6 #include <TROOT.h>
7 #include <TCanvas.h>
8 #include <TMath.h>
9 #include <TColor.h>
10 #include <TPaveText.h>
11 #include <TRandom.h>
12 #include <TH1.h>
13 #include <TH2.h>
14 #include <TF1.h>
15 #include <TSQLResult.h>
16 #include <TProfile.h>
17
18 //#include "TTbar_stuff.C"
19 using namespace std;
20
21 using namespace PlottingSetup;
22
23
24 void rediscover_the_top(string mcjzb, string datajzb) {
25 dout << "Hi! today we are going to (try to) rediscover the top!" << endl;
26 TCanvas *c3 = new TCanvas("c3","c3");
27 c3->SetLogy(1);
28 vector<float> binning;
29 //binning=allsamples.get_optimal_binsize(mcjzb,cutmass&&cutOSSF&&cutnJets,20,50,800);
30 /*
31 binning.push_back(50);
32 binning.push_back(100);
33 binning.push_back(150);
34 binning.push_back(200);
35 binning.push_back(500);
36
37
38 TH1F *dataprediction = allsamples.Draw("dataprediction", "-"+datajzb, binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data, luminosity);
39 TH1F *puresignal = allsamples.Draw("puresignal", datajzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,data, luminosity);
40 // TH1F *puresignal = allsamples.Draw("puresignal", mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,mc, luminosity,allsamples.FindSample("TTJets"));
41 TH1F *observed = allsamples.Draw("observed", datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data,luminosity);
42 /*
43 ofstream myfile;
44 TH1F *ratio = (TH1F*)observed->Clone();
45 ratio->Divide(dataprediction);
46 ratio->GetYaxis()->SetTitle("Ratio obs/pred");
47 ratio->Draw();
48 c3->SaveAs("testratio.png");
49 myfile.open ("ShapeFit_log.txt");
50 establish_upper_limits(observed,dataprediction,puresignal,"LM4",myfile);
51 myfile.close();
52 */
53
54
55 int nbins=100;
56 float low=0;
57 float hi=500;
58 TCanvas *c4 = new TCanvas("c4","c4",900,900);
59 c4->Divide(2,2);
60 c4->cd(1);
61 c4->cd(1)->SetLogy(1);
62 TH1F *datapredictiont = allsamples.Draw("datapredictiont", "-"+datajzb, nbins,low,hi, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data, luminosity);
63 TH1F *datapredictiono = allsamples.Draw("datapredictiono", "-"+datajzb, nbins,low,hi, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,data, luminosity);
64 datapredictiont->Add(datapredictiono,-1);
65 dout << "Second way of doing this !!!! Analytical shape to the left :-D" << endl;
66 vector<TF1*> functions = do_cb_fit_to_plot(datapredictiont,10);
67 datapredictiont->SetMarkerColor(kRed);
68 datapredictiont->SetLineColor(kRed);
69 datapredictiont->Draw();
70 functions[1]->Draw("same");
71 TText *title1 = write_title("Top Background Prediction (JZB<0, with osof subtr)");
72 title1->Draw();
73
74 c4->cd(2);
75 c4->cd(2)->SetLogy(1);
76 TH1F *observedt = allsamples.Draw("observedt", datajzb, nbins,low,hi, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data, luminosity);
77 observedt->Draw();
78 datapredictiont->Draw("histo,same");
79 functions[1]->Draw("same");
80 TText *title2 = write_title("Observed and predicted background");
81 title2->Draw();
82
83 c4->cd(3);
84 c4->cd(3)->SetLogy(1);
85 // TH1F *ratio = (TH1F*)observedt->Clone();
86
87 TH1F *analytical_background_prediction= new TH1F("analytical_background_prediction","",nbins,low,hi);
88 for(int i=0;i<=nbins;i++) {
89 analytical_background_prediction->SetBinContent(i+1,functions[1]->Eval(((hi-low)/((float)nbins))*(i+0.5)));
90 analytical_background_prediction->SetBinError(i+1,TMath::Sqrt(functions[1]->Eval(((hi-low)/((float)nbins))*(i+0.5))));
91 }
92 analytical_background_prediction->GetYaxis()->SetTitle("JZB [GeV]");
93 analytical_background_prediction->GetYaxis()->CenterTitle();
94 TH1F *analyticaldrawonly = (TH1F*)analytical_background_prediction->Clone();
95 analytical_background_prediction->SetFillColor(TColor::GetColor("#3399FF"));
96 analytical_background_prediction->SetMarkerSize(0);
97 analytical_background_prediction->Draw("e5");
98 analyticaldrawonly->Draw("histo,same");
99 functions[1]->Draw("same");
100 TText *title3 = write_title("Analytical bg pred histo");
101 title3->Draw();
102
103 c4->cd(4);
104 // c4->cd(4)->SetLogy(1);
105 vector<float> ratio_binning;
106 ratio_binning.push_back(0);
107 ratio_binning.push_back(5);
108 ratio_binning.push_back(10);
109 ratio_binning.push_back(20);
110 ratio_binning.push_back(50);
111 // ratio_binning.push_back(60);
112 /*
113 ratio_binning.push_back(51);
114 ratio_binning.push_back(52);
115 ratio_binning.push_back(53);
116 ratio_binning.push_back(54);
117 ratio_binning.push_back(55);
118 ratio_binning.push_back(56);
119 ratio_binning.push_back(57);
120 ratio_binning.push_back(58);
121 ratio_binning.push_back(59);
122 ratio_binning.push_back(60);
123 // ratio_binning.push_back(70);*/
124 // ratio_binning.push_back(80);
125 // ratio_binning.push_back(90);
126 ratio_binning.push_back(80);
127 // ratio_binning.push_back(110);
128 ratio_binning.push_back(500);
129
130 TH1F *observedtb = allsamples.Draw("observedtb", datajzb, ratio_binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data, luminosity);
131 TH1F *datapredictiontb = allsamples.Draw("datapredictiontb", "-"+datajzb, ratio_binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data, luminosity);
132 TH1F *datapredictiontbo = allsamples.Draw("datapredictiontbo", "-"+datajzb, ratio_binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,data, luminosity);
133 datapredictiontb->Add(datapredictiontbo,-1);
134 TH1F *analytical_background_predictionb = allsamples.Draw("analytical_background_predictionb",datajzb, ratio_binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets&&"mll<2",data, luminosity);
135 for(int i=0;i<=ratio_binning.size();i++) {
136 analytical_background_predictionb->SetBinContent(i+1,functions[1]->Eval(analytical_background_predictionb->GetBinCenter(i)));
137 analytical_background_predictionb->SetBinError(i+1,TMath::Sqrt(functions[1]->Eval(analytical_background_predictionb->GetBinCenter(i))));
138 }
139
140 TH1F *ratio = (TH1F*) observedtb->Clone();
141 ratio->Divide(datapredictiontb);
142
143 for (int i=0;i<=ratio_binning.size();i++) {
144 dout << observedtb->GetBinLowEdge(i+1) << ";"<<observedtb->GetBinContent(i+1) << ";" << datapredictiontb->GetBinContent(i+1) << " --> " << ratio->GetBinContent(i+1) << "+/-" << ratio->GetBinError(i+1) << endl;
145 }
146
147 // ratio->Divide(datapredictiontb);
148 // ratio->Divide(analytical_background_predictionb);
149 // TGraphAsymmErrors *JZBratio= histRatio(observedtb,analytical_background_predictionb,data,ratio_binning);
150 // ratio->Divide(analytical_background_prediction);
151 // ratio->Divide(datapredictiont);
152 // ratio->GetYaxis()->SetTitle("obs/pred");
153 // JZBratio->Draw("AP");
154 ratio->GetYaxis()->SetRangeUser(0,10);
155 ratio->Draw();
156 //analytical_background_predictionb->Draw();
157 // JZBratio->SetTitle("");
158 TText *title4 = write_title("Ratio of observed to predicted");
159 title4->Draw();
160
161 // CompleteSave(c4,"test/ttbar_discovery_dataprediction___analytical_function");
162 CompleteSave(c4,"test/ttbar_discovery_dataprediction__analytical__new_binning_one_huge_bin_from_80");
163
164
165
166
167
168 }
169
170 void calculate_upper_limits(string mcjzb, string datajzb) {
171 write_warning("calculate_upper_limits","Upper limit calculation temporarily deactivated");
172 // write_warning("calculate_upper_limits","Calculation of SUSY upper limits has been temporarily suspended in favor of top discovery");
173 // rediscover_the_top(mcjzb,datajzb);
174 /*
175 TCanvas *c3 = new TCanvas("c3","c3");
176 c3->SetLogy(1);
177 vector<float> binning;
178 //binning=allsamples.get_optimal_binsize(mcjzb,cutmass&&cutOSSF&&cutnJets,20,50,800);
179 binning.push_back(50);
180 binning.push_back(100);
181 binning.push_back(150);
182 binning.push_back(200);
183 binning.push_back(500);
184 TH1F *datapredictiona = allsamples.Draw("datapredictiona", "-"+datajzb, binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,mc, luminosity);
185 TH1F *datapredictionb = allsamples.Draw("datapredictionb", "-"+datajzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,mc, luminosity);
186 TH1F *datapredictionc = allsamples.Draw("datapredictionc", datajzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,mc, luminosity);
187 TH1F *dataprediction = (TH1F*)datapredictiona->Clone();
188 dataprediction->Add(datapredictionb,-1);
189 dataprediction->Add(datapredictionc);
190 TH1F *puresignal = allsamples.Draw("puresignal", mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
191 TH1F *signalpred = allsamples.Draw("signalpred", "-"+mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
192 TH1F *signalpredlo = allsamples.Draw("signalpredlo", "-"+mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
193 TH1F *signalpredro = allsamples.Draw("signalpredro", mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
194 TH1F *puredata = allsamples.Draw("puredata", datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data,luminosity);
195 signalpred->Add(signalpredlo,-1);
196 signalpred->Add(signalpredro);
197 puresignal->Add(signalpred,-1);//subtracting signal contamination
198 ofstream myfile;
199 myfile.open ("ShapeFit_log.txt");
200 establish_upper_limits(puredata,dataprediction,puresignal,"LM4",myfile);
201 myfile.close();
202 */
203 }
204
205 vector<float> compute_one_upper_limit(float mceff,float mcefferr, int ibin, string mcjzb, bool doobserved=false) {
206 float sigma95=0.0,sigma95A=0.0;
207 int nuisancemodel=1;
208 dout << "Now calling : CL95(" << luminosity << "," << lumiuncert*luminosity << "," << mceff << "," << mcefferr << "," << Npred[ibin] << "," << Nprederr[ibin] << "," << Nobs[ibin] << "," << false << "," << nuisancemodel<< ") " << endl;
209 sigma95 = CL95(luminosity, lumiuncert*luminosity, mceff, mcefferr, Npred[ibin], Nprederr[ibin], Nobs[ibin], false, nuisancemodel);
210 if(doobserved) {
211 dout << "Now calling : CLA(" << luminosity << "," << lumiuncert*luminosity << "," << mceff << "," << mcefferr << "," << Npred[ibin] << "," << Nprederr[ibin] << "," << nuisancemodel<< ") " << endl;
212 sigma95A = CLA(luminosity, lumiuncert*luminosity, mceff, mcefferr, Npred[ibin], Nprederr[ibin], nuisancemodel);
213 }
214 vector<float> sigmas;
215 sigmas.push_back(sigma95);
216 sigmas.push_back(sigma95A);
217 return sigmas;
218 }
219
220 void compute_upper_limits_from_counting_experiment(vector<vector<float> > uncertainties,vector<float> jzbcuts, string mcjzb, bool doobserved) {
221 dout << "Doing counting experiment ... " << endl;
222 vector<vector<string> > limits;
223 vector<vector<float> > vlimits;
224
225
226 for(int isample=0;isample<signalsamples.collection.size();isample++) {
227 vector<string> rows;
228 vector<float> vrows;
229 dout << "Considering sample " << signalsamples.collection[isample].samplename << endl;
230 rows.push_back(signalsamples.collection[isample].samplename);
231 for(int ibin=0;ibin<jzbcuts.size();ibin++) {
232 dout << "_________________________________________________________________________________" << endl;
233 float JZBcutat=uncertainties[isample*jzbcuts.size()+ibin][0];
234 float mceff=uncertainties[isample*jzbcuts.size()+ibin][1];
235 float staterr=uncertainties[isample*jzbcuts.size()+ibin][2];
236 float systerr=uncertainties[isample*jzbcuts.size()+ibin][3];
237 float toterr =uncertainties[isample*jzbcuts.size()+ibin][4];
238 float observed,null,result;
239 fill_result_histos(observed, null,null,null,null,null,null,null,mcjzb,JZBcutat,(int)5,result,(signalsamples.FindSample(signalsamples.collection[isample].filename)),signalsamples);
240 observed-=result;//this is the actual excess we see!
241 float expected=observed/luminosity;
242
243 dout << "Sample: " << signalsamples.collection[isample].samplename << ", JZB>"<<JZBcutat<< " : " << mceff << " +/- " << staterr << " (stat) +/- " << systerr << " (syst) --> toterr = " << toterr << endl;
244 vector<float> sigmas = compute_one_upper_limit(mceff,toterr,ibin,mcjzb,doobserved);
245
246 if(doobserved) {
247 rows.push_back(any2string(sigmas[0])+";"+any2string(sigmas[1])+";"+"("+any2string(expected)+")");
248 vrows.push_back(sigmas[0]);
249 vrows.push_back(sigmas[1]);
250 vrows.push_back(expected);
251 }
252 else {
253 rows.push_back(any2string(sigmas[0])+"("+any2string(expected)+")");
254 vrows.push_back(sigmas[0]);
255 vrows.push_back(expected);
256 }
257 }//end of bin loop
258 limits.push_back(rows);
259 vlimits.push_back(vrows);
260 }//end of sample loop
261 dout << endl << endl << "PAS table 3: " << endl << endl;
262 dout << "\t";
263 for (int irow=0;irow<jzbcuts.size();irow++) {
264 dout << jzbcuts[irow] << "\t";
265 }
266 dout << endl;
267 for(int irow=0;irow<limits.size();irow++) {
268 for(int ientry=0;ientry<limits[irow].size();ientry++) {
269 dout << limits[irow][ientry] << "\t";
270 }
271 dout << endl;
272 }
273
274 if(!doobserved) {
275 dout << endl << endl << "LIMITS: " << endl;
276 dout << "\t";
277 for (int irow=0;irow<jzbcuts.size();irow++) {
278 dout << jzbcuts[irow] << "\t";
279 }
280 dout << endl;
281 for(int irow=0;irow<limits.size();irow++) {
282 dout << limits[irow][0] << "\t";
283 for(int ientry=0;ientry<jzbcuts.size();ientry++) {
284 dout << Round(vlimits[irow][2*ientry] / vlimits[irow][2*ientry+1],3)<< "\t";
285 }
286 dout << endl;
287 }
288 }//do observed
289
290 dout << endl << endl << "Final selection efficiencies with total statistical and systematic errors, and corresponding observed and expected upper limits (UL) on ($\\sigma\\times$ BR $\\times$ acceptance) for the LM4 and LM8 scenarios, in the different regions. The last column contains the predicted ($\\sigma \\times $BR$\\times$ acceptance) at NLO obtained from Monte Carlo simulation." << endl;
291 dout << "Scenario \t Efficiency [%] \t Upper limits [pb] \t Prediction [pb]" << endl;
292 for(int icut=0;icut<jzbcuts.size();icut++) {
293 dout << "Region with JZB>" << jzbcuts[icut] << endl;
294 for(int isample=0;isample<signalsamples.collection.size();isample++) {
295 dout << limits[icut][0] << "\t" << Round(100*uncertainties[isample*jzbcuts.size()+icut][1],1) << "+/-" << Round(100*uncertainties[isample*jzbcuts.size()+icut][2],1) << " (stat) +/- " << Round(100*uncertainties[isample*jzbcuts.size()+icut][3],1) << " (syst) \t" << Round((vlimits[isample][2*icut]),3) << "\t" << Round(vlimits[isample][2*icut+1],3) << endl;
296 }
297 dout << endl;
298 }
299
300 write_warning("compute_upper_limits_from_counting_experiment","Still need to update the script");
301 }
302
303 void susy_scan_axis_labeling(TH2F *histo) {
304 histo->GetXaxis()->SetTitle("#Chi_{2}^{0}-LSP");
305 histo->GetXaxis()->CenterTitle();
306 histo->GetYaxis()->SetTitle("m_{#tilde{q}}");
307 histo->GetYaxis()->CenterTitle();
308 }
309
310 void scan_susy_space(string mcjzb, string datajzb) {
311 TCanvas *c3 = new TCanvas("c3","c3");
312 vector<float> binning;
313 binning=allsamples.get_optimal_binsize(mcjzb,cutmass&&cutOSSF&&cutnJets,20,50,800);
314 float arrbinning[binning.size()];
315 for(int i=0;i<binning.size();i++) arrbinning[i]=binning[i];
316 TH1F *puredata = allsamples.Draw("puredata", datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data,luminosity);
317 puredata->SetMarkerSize(DataMarkerSize);
318 TH1F *allbgs = allsamples.Draw("allbgs", "-"+datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data,luminosity);
319 TH1F *allbgsb = allsamples.Draw("allbgsb", "-"+datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,data,luminosity);
320 TH1F *allbgsc = allsamples.Draw("allbgsc", datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,data,luminosity);
321 allbgs->Add(allbgsb,-1);
322 allbgs->Add(allbgsc);
323 int ndata=puredata->Integral();
324 ofstream myfile;
325 myfile.open ("susyscan_log.txt");
326 TFile *susyscanfile = new TFile("/scratch/fronga/SMS/T5z_SqSqToQZQZ_38xFall10.root");
327 TTree *suevents = (TTree*)susyscanfile->Get("events");
328 TH2F *exclusionmap = new TH2F("exclusionmap","",20,0,500,20,0,1000);
329 TH2F *exclusionmap1s = new TH2F("exclusionmap1s","",20,0,500,20,0,1000);
330 TH2F *exclusionmap2s = new TH2F("exclusionmap2s","",20,0,500,20,0,1000);
331 TH2F *exclusionmap3s = new TH2F("exclusionmap3s","",20,0,500,20,0,1000);
332
333 susy_scan_axis_labeling(exclusionmap);
334 susy_scan_axis_labeling(exclusionmap1s);
335 susy_scan_axis_labeling(exclusionmap2s);
336 susy_scan_axis_labeling(exclusionmap3s);
337
338 Int_t MyPalette[100];
339 Double_t r[] = {0., 0.0, 1.0, 1.0, 1.0};
340 Double_t g[] = {0., 0.0, 0.0, 1.0, 1.0};
341 Double_t b[] = {0., 1.0, 0.0, 0.0, 1.0};
342 Double_t stop[] = {0., .25, .50, .75, 1.0};
343 Int_t FI = TColor::CreateGradientColorTable(5, stop, r, g, b, 100);
344 for (int i=0;i<100;i++) MyPalette[i] = FI+i;
345
346 gStyle->SetPalette(100, MyPalette);
347
348 for(int m23=50;m23<500;m23+=25) {
349 for (int m0=(2*(m23-50)+150);m0<=1000;m0+=50)
350 {
351 c3->cd();
352 stringstream drawcondition;
353 drawcondition << "pfJetGoodNum>=3&&(TMath::Abs(masses[0]-"<<m0<<")<10&&TMath::Abs(masses[2]-masses[3]-"<<m23<<")<10)&&mll>5&&id1==id2";
354 TH1F *puresignal = new TH1F("puresignal","puresignal",binning.size()-1,arrbinning);
355 TH1F *puresignall= new TH1F("puresignall","puresignal",binning.size()-1,arrbinning);
356 stringstream drawvar,drawvar2;
357 drawvar<<mcjzb<<">>puresignal";
358 drawvar2<<"-"<<mcjzb<<">>puresignall";
359 suevents->Draw(drawvar.str().c_str(),drawcondition.str().c_str());
360 suevents->Draw(drawvar2.str().c_str(),drawcondition.str().c_str());
361 if(puresignal->Integral()<60) {
362 delete puresignal;
363 continue;
364 }
365 puresignal->Add(puresignall,-1);//we need to correct for the signal contamination - we effectively only see (JZB>0)-(JZB<0) !!
366 puresignal->Scale(ndata/(20*puresignal->Integral()));//normalizing it to 5% of the data
367 stringstream saveas;
368 saveas<<"Model_Scan/m0_"<<m0<<"__m23_"<<m23;
369 dout << "PLEASE KEEP IN MIND THAT SIGNAL CONTAMINATION IS NOT REALLY TAKEN CARE OF YET DUE TO LOW STATISTICS! SHOULD BE SOMETHING LIKE THIS : "<< endl;
370 // TH1F *signalpredlo = allsamples.Draw("signalpredlo", "-"+mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
371 // TH1F *signalpredro = allsamples.Draw("signalpredro", mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSOF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
372 // TH1F *puredata = allsamples.Draw("puredata", datajzb,binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,data,luminosity);
373 // signalpred->Add(signalpredlo,-1);
374 // signalpred->Add(signalpredro);
375 // puresignal->Add(signalpred,-1);//subtracting signal contamination
376 //---------------------
377 // dout << "(m0,m23)=("<<m0<<","<<m23<<") contains " << puresignal->Integral() << endl;
378 // TH1F *puresignal = allsamples.Draw("puresignal",mcjzb, binning, "JZB [GeV]", "events", cutmass&&cutOSSF&&cutnJets,mc, luminosity,allsamples.FindSample("LM4"));
379 vector<float> results=establish_upper_limits(puredata,allbgs,puresignal,saveas.str(),myfile);
380 if(results.size()==0) {
381 delete puresignal;
382 continue;
383 }
384 exclusionmap->Fill(m23,m0,results[0]);
385 exclusionmap1s->Fill(m23,m0,results[1]);
386 exclusionmap2s->Fill(m23,m0,results[2]);
387 exclusionmap3s->Fill(m23,m0,results[3]);
388 delete puresignal;
389 dout << "(m0,m23)=("<<m0<<","<<m23<<") : 3 sigma at " << results[3] << endl;
390 }
391 }//end of model scan for loop
392
393 dout << "Exclusion Map contains" << exclusionmap->Integral() << " (integral) and entries: " << exclusionmap->GetEntries() << endl;
394 c3->cd();
395 exclusionmap->Draw("CONTZ");
396 CompleteSave(c3,"Model_Scan/CONT/Model_Scan_Mean_values");
397 exclusionmap->Draw("COLZ");
398 CompleteSave(c3,"Model_Scan/COL/Model_Scan_Mean_values");
399
400 exclusionmap1s->Draw("CONTZ");
401 CompleteSave(c3,"Model_Scan/CONT/Model_Scan_1sigma_values");
402 exclusionmap1s->Draw("COLZ");
403 CompleteSave(c3,"Model_Scan/COL/Model_Scan_1sigma_values");
404
405 exclusionmap2s->Draw("CONTZ");
406 CompleteSave(c3,"Model_Scan/CONT/Model_Scan_2sigma_values");
407 exclusionmap2s->Draw("COLZ");
408 CompleteSave(c3,"Model_Scan/COL/Model_Scan_2sigma_values");
409
410 exclusionmap3s->Draw("CONTZ");
411 CompleteSave(c3,"Model_Scan/CONT/Model_Scan_3sigma_values");
412 exclusionmap3s->Draw("COLZ");
413 CompleteSave(c3,"Model_Scan/COL/Model_Scan_3sigma_values");
414
415 TFile *exclusion_limits = new TFile("exclusion_limits.root","RECREATE");
416 exclusionmap->Write();
417 exclusionmap1s->Write();
418 exclusionmap2s->Write();
419 exclusionmap3s->Write();
420 exclusion_limits->Close();
421 susyscanfile->Close();
422
423 myfile.close();
424 }
425
426
427
428