ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/cbrown/AnalysisFramework/Plotting/Modules/SUSYScan.C
Revision: 1.24
Committed: Wed Aug 31 11:38:26 2011 UTC (13 years, 8 months ago) by buchmann
Content type: text/plain
Branch: MAIN
Changes since 1.23: +13 -4 lines
Log Message:
Updated the way an mSUGRA scan is identified (via title instead of filename); added an ipoint map so that any possible holes in a map can be easily traced back to the individual point

File Contents

# Content
1 #include <iostream>
2 #include <vector>
3 #include <sys/stat.h>
4 #include <fstream>
5 #include <pthread.h>
6
7 #include <TCut.h>
8 #include <TROOT.h>
9 #include <TCanvas.h>
10 #include <TMath.h>
11 #include <TColor.h>
12 #include <TPaveText.h>
13 #include <TRandom.h>
14 #include <TH1.h>
15 #include <TH2.h>
16 #include <TF1.h>
17 #include <TSQLResult.h>
18 #include <TProfile.h>
19
20 //#include "TTbar_stuff.C"
21 using namespace std;
22
23 using namespace PlottingSetup;
24
25 void susy_scan_axis_labeling(TH2F *histo) {
26 histo->GetXaxis()->SetTitle("#Chi_{2}^{0}-LSP");
27 histo->GetXaxis()->CenterTitle();
28 histo->GetYaxis()->SetTitle("m_{#tilde{q}}");
29 histo->GetYaxis()->CenterTitle();
30 }
31
32 bool isThreadActive=false;
33
34 struct limit_args {
35 float mceff;
36 float toterr;
37 int ibin;
38 string mcjzb;
39 vector<float> sigmas;
40 };
41
42 void* compute_one_upper_limit_wrapper(void* data) {
43 isThreadActive=true;
44 limit_args* args = (limit_args*) data;
45 std::cout << "Thread to compute limits has been started" << std::endl;
46 args->sigmas=compute_one_upper_limit(args->mceff,args->toterr,args->ibin,args->mcjzb);
47 std::cout << "Thread to compute limits has finished" << std::endl;
48 isThreadActive=false;
49 return NULL; // oder in C++: return 0;// Damit kann man Werte zurückgeben
50 }
51
52 void do_limit_wrapper(float mceff,float toterr,int ibin,string mcjzb,vector<float> &sigmas) {
53 pthread_t limitthread;
54 limit_args limargs={mceff,toterr,ibin,mcjzb,sigmas};
55 pthread_create( &limitthread, NULL, compute_one_upper_limit_wrapper, (void*) &limargs);
56 int counter=0;
57 sleep(1); //waiting a second for the process to become active
58 while(counter<limitpatience*60 && isThreadActive) {
59 std::cout << "Limits are being calculated; Checking round " << counter << " ( corresponds to " << seconds_to_time(counter) << " ) , patience will end in " << seconds_to_time(60*60*limitpatience-counter) << std::endl;
60 counter++;
61 sleep(1);
62 }
63
64 if(!isThreadActive) {
65 cout << "Thread finished sucessfully" << endl;
66 pthread_join( limitthread, NULL );
67 std::cout<< __FUNCTION__ << " : going to save sigmas " << std::endl;
68 sigmas=limargs.sigmas;
69 } else {
70 pthread_cancel(limitthread);
71 std::cout << "DID NOT TERMINATE IN TIME - ABORTED!" << std::endl;
72 sigmas.push_back(-1);sigmas.push_back(-1);sigmas.push_back(-1);
73 }
74 }
75
76 void prepare_scan_axis(TH2 *variablemap,bool ismSUGRA) {
77 variablemap->GetXaxis()->SetTitle("m_{glu}");
78 variablemap->GetYaxis()->SetTitle("m_{LSP}");
79
80 if(ismSUGRA) {
81 variablemap->GetXaxis()->SetTitle("m_{0}");
82 variablemap->GetYaxis()->SetTitle("m_{1/2}");
83 }
84
85 variablemap->GetXaxis()->CenterTitle();
86 variablemap->GetYaxis()->CenterTitle();
87 }
88
89 void scan_SUSY_parameter_space(string mcjzb,string datajzb,vector<float> jzb_cut,bool requireZ, float peakerror, int ibin,float njobs=-1, float jobnumber=-1, bool systematicsonly=false,bool efficiencyonly=false) {
90 bool runninglocally=true;
91 if(njobs>-1&&jobnumber>-1) {
92 runninglocally=false;
93 dout << "Running on the GRID (this is job " << jobnumber << "/" << njobs << ") for a jzb cut at " << jzb_cut[ibin] << endl;
94 } else {
95 dout << "Running locally " << endl;
96 dout << "This will take a really, really long time - if you want to see the results THIS week try running the LimitWorkerScript on the grid (DistributedModelCalculation/Limits/)" << endl;
97 }
98
99 cout << "FILE USED: " << (scansample.collection)[0].filename << endl;
100 cout << "FILE USED (name): " << (scansample.collection)[0].samplename << endl;
101 jzbSel=jzb_cut[ibin];
102 geqleq="geq";
103 automatized=true;
104 mcjzbexpression=mcjzb;
105
106 string massgluname="MassGlu";
107 string massLSPname="MassLSP";
108
109 string prefix="SMS_";
110 // up to here, everything is set up for SMS; now we need to switch stuff around if we're dealing with an mSUGRA scan!
111 bool ismSUGRA=false;
112 if(Contains((scansample.collection)[0].samplename,"mSUGRA")) ismSUGRA=true;
113 if(ismSUGRA) {
114 massgluname="M0"; // this is the "x axis" in the limit plot (like the gluino in the SMS case)
115 massLSPname="M12"; // this is the "y axis" in the limit plot (like the LSP in the SMS case)
116 mglustart=m0start;
117 mgluend=m0end;
118 mglustep=m0step;
119 mLSPstart=m12start;
120 mLSPend=m12end;
121 mLSPstep=m12step;
122 prefix="mSUGRA_";
123 cout << "mSUGRA scan has been set up." << endl;
124 } else {
125 cout << "SMS scan has been set up." << endl;
126 }
127
128 Int_t MyPalette[100];
129 Double_t r[] = {0., 0.0, 1.0, 1.0, 1.0};
130 Double_t g[] = {0., 0.0, 0.0, 1.0, 1.0};
131 Double_t b[] = {0., 1.0, 0.0, 0.0, 1.0};
132 Double_t stop[] = {0., .25, .50, .75, 1.0};
133 Int_t FI = TColor::CreateGradientColorTable(5, stop, r, g, b, 110);
134 for (int i=0;i<100;i++) MyPalette[i] = FI+i;
135
136 gStyle->SetPalette(100, MyPalette);
137
138 TH2F *limitmap = new TH2F((prefix+"limitmap"+any2string(jzbSel)).c_str(), "",(mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
139 TH2F *exclmap = new TH2F((prefix+"exclusionmap"+any2string(jzbSel)).c_str(), "",(mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
140 TH2F *sysjesmap = new TH2F((prefix+"sysjes"+any2string(jzbSel)).c_str(), "",(mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
141 TH2F *sysjsumap = new TH2F((prefix+"sysjsu"+any2string(jzbSel)).c_str(), "",(mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
142 TH2F *sysresmap = new TH2F((prefix+"sysresmap"+any2string(jzbSel)).c_str(),"",(mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
143 TH2F *efficiencymap = new TH2F((prefix+"efficiencymap").c_str(),"",(mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
144 TH2F *Neventsmap = new TH2F((prefix+"Neventsmap").c_str(),"", (mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
145 TH2F *ipointmap = new TH2F((prefix+"ipointmap").c_str(),"", (mgluend-mglustart)/mglustep+1,mglustart-0.5*mglustep,mgluend+0.5*mglustep,(mLSPend-mLSPstart)/mLSPstep+1,mLSPstart-0.5*mLSPstep,mLSPend+0.5*mLSPstep);
146 float rightmargin=gStyle->GetPadRightMargin();
147 gStyle->SetPadRightMargin(0.15);
148
149 TCanvas *limcanvas = new TCanvas("limcanvas","Limit canvas");
150
151 int Npoints=0;
152 for(int mglu=mglustart;mglu<=mgluend;mglu+=mglustep) {
153 for (int mlsp=mLSPstart;mlsp<=mLSPend&&mlsp<=mglu;mlsp+=mLSPstep) Npoints++;
154 }
155
156 int ipoint=-1;
157 for(int mglu=mglustart;mglu<=mgluend;mglu+=mglustep) {
158 for (int mlsp=mLSPstart;mlsp<=mLSPend&&mlsp<=mglu;mlsp+=mLSPstep)
159 {
160 ipoint++;
161 if(!runninglocally&&!do_this_point(ipoint,Npoints,jobnumber,njobs)) continue;
162 float result=-987,resulterr=-987;
163 stringstream addcut;
164 addcut << "(TMath::Abs("<<massgluname<<"-"<<mglu<<")<5)&&(TMath::Abs("<<massLSPname<<"-"<<mlsp<<")<5)";
165 (scansample.collection)[0].events->Draw("eventNum",addcut.str().c_str(),"goff");
166 float nevents = (scansample.collection)[0].events->GetSelectedRows();
167 vector<vector<float> > systematics;
168 if(nevents==0) {
169 dout << "This point ("<<ipoint<<") with configuration ("<<massgluname<<"="<<mglu<<" , "<<massLSPname<<"="<<mlsp << ") contains no events and will be skipped."<< endl;
170 continue;
171 } else {
172 dout << "OK! This point ("<<ipoint<<") with configuration ("<<massgluname<<"="<<mglu<<" , "<<massLSPname<<"="<<mlsp << ") contains " << nevents << " and will therefore not be skipped."<< endl;
173 }
174 if(nevents!=0&&efficiencyonly) {
175 MCefficiency((scansample.collection)[0].events,result,resulterr,mcjzb,requireZ,nevents,addcut.str());
176 efficiencymap->Fill(mglu,mlsp,result);
177 }
178 Neventsmap->Fill(mglu,mlsp,nevents);
179 ipointmap->Fill(mglu,mlsp,ipoint);
180 if(efficiencyonly) continue;
181
182 do_systematics_for_one_file((scansample.collection)[0].events,nevents,"SUSY SCAN", systematics,mcjzb,datajzb,peakerror,requireZ, addcut.str(),ismSUGRA);
183 float JZBcutat = systematics[0][0];
184 float mceff = systematics[0][1];
185 float mcefferr = systematics[0][2];//MC stat error
186 float toterr = systematics[0][4];
187 float sys_jes = systematics[0][5]; // Jet Energy Scale
188 float sys_jsu = systematics[0][6]; // JZB scale uncertainty
189 float sys_res = systematics[0][7]; // resolution
190 efficiencymap->Fill(mglu,mlsp,mceff);
191 efficiencymap->SetBinError(efficiencymap->FindBin(mglu,mlsp),mcefferr);//blublu
192
193 if(mceff!=mceff||toterr!=toterr||mceff<0) {
194 dout << "Limits can't be calculated for this configuration (mglu="<<mglu<<" , mlsp="<<mlsp << ") as either the efficiency or its error are not positive numbers! (mceff="<<mceff<<" and toterr="<<toterr<<")"<< endl;
195 continue;
196 } else {
197 if(!systematicsonly&&!efficiencyonly) {
198 dout << "Calculating limit now for "<<massgluname<<"="<<mglu<<" , "<<massLSPname<<"="<<mlsp <<endl;
199 vector<float> sigmas;
200 do_limit_wrapper(mceff,toterr,ibin,mcjzb,sigmas);
201 cout << "back in " << __FUNCTION__ << endl;
202 if(sigmas[0]>-0.5) { // negative sigmas are the error signature of do_limit_wrapper, so we want to exclude them.
203 limitmap->Fill(mglu,mlsp,sigmas[0]);
204 sysjesmap->Fill(mglu,mlsp,sys_jes);
205 sysjsumap->Fill(mglu,mlsp,sys_jsu);
206 sysresmap->Fill(mglu,mlsp,sys_res);
207 efficiencymap->Fill(mglu,mlsp,result);
208 dout << "A limit has been added at " << sigmas[0] << " for m_{glu}="<<mglu << " and m_{lsp}="<<mlsp<<endl;
209 } //end of if sigma is positive
210 //end of not systematics only condition
211 }
212 if(systematicsonly) {
213 sysjesmap->Fill(mglu,mlsp,sys_jes);
214 sysjsumap->Fill(mglu,mlsp,sys_jsu);
215 sysresmap->Fill(mglu,mlsp,sys_res);
216 efficiencymap->Fill(mglu,mlsp,result);
217 }
218 }//efficiency is valid
219 }
220 }
221
222 prepare_scan_axis(limitmap,ismSUGRA);
223 prepare_scan_axis(sysjesmap,ismSUGRA);
224 prepare_scan_axis(sysjsumap,ismSUGRA);
225 prepare_scan_axis(sysresmap,ismSUGRA);
226
227 if(!systematicsonly&&!efficiencyonly) {
228 limcanvas->cd();
229 limitmap->Draw("COLZ");
230 string titletobewritten="Limits in LSP-Glu plane";
231 if(ismSUGRA) titletobewritten="Limits in m_{1/2}-m_{0} plane";
232 TText *title = write_title(titletobewritten);
233 title->Draw();
234 if(runninglocally) {
235 CompleteSave(limcanvas,"SUSYScan/Limits_JZB_geq"+any2string(jzb_cut[ibin]));
236 } else {
237 TFile *outputfile;
238 if(systematicsonly) outputfile=new TFile(("output/DistributedSystematics_job"+string(any2string(jobnumber))+"_of_"+string(any2string(njobs))+".root").c_str(),"UPDATE");//needs to be "UPDATE" as we can get to this point for different JZB cuts and don't want to erase the previous data :-)
239 else outputfile=new TFile(("output/DistributedLimits_job"+string(any2string(jobnumber))+"_of_"+string(any2string(njobs))+".root").c_str(),"UPDATE");//needs to be "UPDATE" as we can get to this point for
240 limitmap->Write();
241 sysjesmap->Write();
242 sysjsumap->Write();
243 efficiencymap->Write();
244 Neventsmap->Write();
245 ipointmap->Write();
246 outputfile->Close();
247 }
248 }
249 if(systematicsonly) { // systematics only :
250 limcanvas->cd();
251 sysjesmap->Draw("COLZ");
252 string titletobewritten="Jet Energy scale in LSP-Glu plane";
253 if(ismSUGRA) titletobewritten="Limits in m_{1/2}-m_{0} plane";
254 TText *title = write_title(titletobewritten);
255 title->Draw();
256 if(runninglocally) {
257 CompleteSave(limcanvas,"SUSYScan/JES_geq"+any2string(jzb_cut[ibin]));
258 }
259 sysjsumap->Draw("COLZ");
260 titletobewritten="JZB Scale Uncertainty in LSP-Glu plane";
261 if(ismSUGRA) titletobewritten="JZB Scale Uncertainty in m_{1/2}-m_{0} plane";
262 TText *title2 = write_title(titletobewritten);
263 title2->Draw();
264 if(runninglocally) {
265 CompleteSave(limcanvas,"SUSYScan/JSU_geq"+any2string(jzb_cut[ibin]));
266 }
267 sysresmap->Draw("COLZ");
268 titletobewritten="Resolution in LSP-Glu plane";
269 if(ismSUGRA) titletobewritten="Resolution in m_{1/2}-m_{0} plane";
270 TText *title3 = write_title(titletobewritten);
271 title3->Draw();
272 if(runninglocally) {
273 CompleteSave(limcanvas,"SUSYScan/Resolution_geq"+any2string(jzb_cut[ibin]));
274 }
275
276 if(!runninglocally) {
277 TFile *outputfile=new TFile(("output/DistributedSystematics_job"+string(any2string(jobnumber))+"_of_"+string(any2string(njobs))+".root").c_str(),"UPDATE");
278 sysjesmap->Write();
279 sysjsumap->Write();
280 sysresmap->Write();
281 efficiencymap->Write();
282 Neventsmap->Write();
283 ipointmap->Write();
284 outputfile->Close();
285 }
286 }//end of systematics only
287 if(efficiencyonly) {
288 limcanvas->cd();
289 efficiencymap->Draw("COLZ");
290 string titletobewritten="Efficiencies in LSP-Glu plane";
291 if(ismSUGRA) titletobewritten="Efficiencies in m_{1/2}-m_{0} plane";
292 TText *title = write_title(titletobewritten);
293 title->Draw();
294 if(runninglocally) {
295 CompleteSave(limcanvas,"SUSYScan/Efficiency_geq"+any2string(jzb_cut[ibin]));
296 }
297 limcanvas->cd();
298 sysjesmap->Draw("COLZ");
299 titletobewritten="N(events) in LSP-Glu plane";
300 if(ismSUGRA) titletobewritten="N(events) in m_{1/2}-m_{0} plane";
301 TText *title2 = write_title(titletobewritten);
302 title2->Draw();
303 if(runninglocally) {
304 CompleteSave(limcanvas,"SUSYScan/Nevents_geq"+any2string(jzb_cut[ibin]));
305 }
306 if(!runninglocally) {
307 TFile *outputfile=new TFile(("output/DistributedSystematics_job"+string(any2string(jobnumber))+"_of_"+string(any2string(njobs))+".root").c_str(),"UPDATE");
308 ipointmap->Write();
309 Neventsmap->Write();
310 efficiencymap->Write();
311 outputfile->Close();
312 }
313 }//end of efficiencies only
314 delete limcanvas;
315 }
316
317 void scan_SUSY_parameter_space(string mcjzb,string datajzb,vector<float> jzb_cut,bool requireZ, float peakerror, float njobs=-1, float jobnumber=-1, bool systonly=false, bool effonly=false) {
318 dout << "Starting the SUSY scan now with all " << jzb_cut.size() << " bin(s)" << endl;
319 for(int ibin=0;ibin<jzb_cut.size();ibin++) {
320 scan_SUSY_parameter_space(mcjzb,datajzb,jzb_cut,requireZ, peakerror, ibin, njobs, jobnumber,systonly,effonly);
321 }
322 }