]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TestAnalisys.C
Extra header added to the list
[u/mrichter/AliRoot.git] / TPC / TestAnalisys.C
1 /// \file TestAnalisys.C
2 ///
3 /// ~~~{.cxx}
4 /// .L AliGenInfo.C+
5 /// .L TestAnalisys.C+
6 /// AddChains(868);    // AddChains(runNumber);
7 ///  Select();          // make default selection of data
8 ///  MakePictures("pic868");
9 /// ~~~
10
11 #include "TFile.h"
12 #include "TTree.h"
13 #include "TChain.h"
14 #include "TString.h"
15 #include "TLegend.h"
16 #include "TStyle.h"
17 #include "TCut.h"
18 #include "TF1.h"
19 #include "TProfile.h"
20 #include "TProfile2D.h"
21 #include "TH1F.h"
22 #include "TH2F.h"
23 #include "TPad.h"
24 #include "TCanvas.h"
25 #include "TSystem.h"
26 #include "TParticle.h"
27 #include "TRandom.h"
28 #include "TEventList.h"
29
30 #include "AliTrackReference.h"
31 #include "AliTPCParam.h"
32 #include "AliDetector.h"
33 #include "AliStack.h" 
34 #include "AliGenInfo.h"
35
36
37
38
39
40
41 void Select();                // make default selection 
42 void SelectLaser();           // make default selection  for laser tracks
43
44 void AddChains(Int_t run);    // add all the trees with selected run number to the chain
45 void MakePictures(char *dirname);           // make default pictures
46
47 void PRFYZ(TCut cut0, TCut cut1,  char * description);  
48 void PRFZZ(TCut cut0, TCut cut1,  char * description);  
49 void P5Z(TCut cut0, TCut cut1,  char * description);  
50 void P3Z(TCut cut0, TCut cut1,  char * description);  
51 void ResYZ(TCut cut0, TCut cut1,  char * description);
52 void SysYX(TCut cut0,  char * description);
53 void SysZX(TCut cut0,  char * description);
54
55 TProfile * ProfileMaxRow(TCut cut0, char *name, Int_t max);
56 TProfile * ProfileMaxPhi(TCut cut0, char *name, Int_t max);
57 TProfile * ProfileMaxZ(TCut cut0, char *name, Int_t max);
58 TProfile * ProfileQRow(TCut cut0, char *name, Int_t max);
59 TProfile * ProfileQPhi(TCut cut0, char *name, Int_t max);
60 TProfile * ProfileQZ(TCut cut0, char *name, Int_t max);
61 TCanvas *  NoiseSector(TCut cut0,  char * description, Int_t maxrow, Int_t maxpad);
62
63 //
64 // global variables
65 //
66 TChain chaincl("Tracks","Tracks");   // tpc tracks and clusters
67 TChain chaincl2("Tracks","Tracks");   // tpc tracks and clusters
68 TChain chainSignal("SignalB","SignalB");   // signals over threshold 50
69
70 TChain chainFit("Fit","Fit");        // fitted signals with fit parameters
71 TChain chainPed("Fit","Fit");        // fitted pedestal with noise
72 TString runDesc="Run ";              // run descriptor
73 //
74 AliComparisonDraw comp;
75 AliComparisonDraw compCl;
76 AliComparisonDraw compF;
77 AliComparisonDraw compP;
78 //
79 //
80 // selection of data for analysis
81 //
82 TEventList * listTracks    = new TEventList("listTracks","listTracks");
83 TEventList * listFitS      = new TEventList("listFitS","listFitS");
84 TEventList * listFitPed    = new TEventList("listFitPed","listFitPed");
85
86
87
88
89 void MakePictures(char *dirname){
90   /// Define Uli Style
91
92   gROOT->SetStyle("Plain");
93   gStyle->SetFillColor(10);
94   gStyle->SetPadColor(10);
95   gStyle->SetCanvasColor(10);
96   gStyle->SetStatColor(10);
97   
98   gStyle->SetPalette(1,0);
99   gStyle->SetNumberContours(50);
100   //
101   const Int_t kMinCl = 200;
102   char chshell[100];
103   sprintf(chshell,"mkdir %s", dirname);
104   gSystem->Exec(chshell); 
105   sprintf(chshell,"cd %s", dirname);
106   gSystem->Exec(chshell);
107   //
108   //
109   //
110   TCanvas c(dirname,dirname);
111   for (Int_t isector=0; isector<36; isector++){
112     char chcut1[100];
113     char chcut2[100];
114     char chdesc[100];
115     TProfile * prof;
116     sprintf(chshell,"Cl.fX>0&&Cl.fDetector==%d",isector);
117     Int_t ncl = comp.fTree->Draw("Cl.fY",chshell);
118     if (ncl<kMinCl) continue;
119     printf("MakePictures sector %d\n",isector);
120     //
121     // charge pictures
122     //
123     //
124     // charge row
125     //
126     c.cd();
127     sprintf(chdesc,"%s Sector %d IROC",runDesc.Data(), isector);
128     sprintf(chcut1,"Cl.fDetector==%d", isector);
129     prof = ProfileQRow(chcut1, chdesc, 70);
130     sprintf(chshell,"%s/qrow_sec%dIROC.eps", dirname,isector);
131     prof->Draw();
132     c.Update();
133     c.Print(chshell);
134     prof = ProfileMaxRow(chcut1, chdesc, 70);
135     sprintf(chshell,"%s/maxrow_sec%dIROC.eps", dirname,isector);
136     prof->Draw();
137     c.Update();
138     c.Print(chshell);
139     //
140     sprintf(chdesc,"%s Sector %d OROC",runDesc.Data(), isector);
141     sprintf(chcut1,"Cl.fDetector==%d", isector+36);
142     prof = ProfileQRow(chcut1, chdesc, 100);
143     sprintf(chshell,"%s/qrow_sec%dOROC.eps", dirname,isector);
144     prof->Draw();
145     c.Update();
146     c.Print(chshell);
147     prof = ProfileMaxRow(chcut1, chdesc, 100);
148     sprintf(chshell,"%s/maxrow_sec%dOROC.eps", dirname,isector);
149     prof->Draw();
150     c.Update();
151     c.Print(chshell);
152     //
153     // charge phi
154     //
155     sprintf(chdesc,"%s Sector %d IROC",runDesc.Data(), isector);
156     sprintf(chcut1,"Cl.fDetector==%d", isector);
157     prof = ProfileMaxPhi(chcut1, chdesc,20);
158     sprintf(chshell,"%s/qphi_sec%dIROC.eps", dirname,isector);
159     prof->Draw();
160     c.Update();
161     c.Print(chshell);
162     //
163     sprintf(chdesc,"%s Sector %d OROC",runDesc.Data(), isector);
164     sprintf(chcut1,"Cl.fDetector==%d", isector+36);
165     prof = ProfileMaxPhi(chcut1, chdesc,20);
166     sprintf(chshell,"%s/qphi_sec%dOROC.eps", dirname,isector);
167     prof->Draw();
168     c.Update();
169     c.Print(chshell);
170     //
171     //   charge z
172     //
173     c.cd();
174     sprintf(chdesc,"%s Sector %d IROC",runDesc.Data(), isector);
175     sprintf(chcut1,"Cl.fDetector==%d", isector);
176     prof = ProfileQZ(chcut1, chdesc,20);
177     sprintf(chshell,"%s/qz_sec%dIROC.eps", dirname,isector);
178     //    prof->Draw();
179     c.Update();
180     c.Print(chshell);
181     //
182     c.cd();
183     sprintf(chdesc,"%s Sector %d OROC",runDesc.Data(), isector);
184     sprintf(chcut1,"Cl.fDetector==%d", isector+36);
185     prof = ProfileQZ(chcut1, chdesc,20);
186     sprintf(chshell,"%s/qz_sec%dOROC.eps", dirname,isector);
187     //prof->Draw();
188     c.Update();
189     c.Print(chshell);
190     //
191     // Picture noise
192     //
193     sprintf(chdesc,"%s Sector %d IROC",runDesc.Data(), isector);
194     sprintf(chcut1,"Sector==%d", isector);
195     TCanvas *cnoise = NoiseSector(chcut1, chdesc,70,70);
196     sprintf(chshell,"%s/noise_sec%dIROC.eps", dirname,isector);
197     cnoise->Print(chshell);
198     sprintf(chdesc,"%s Sector %d OROC",runDesc.Data(), isector);
199     sprintf(chcut1,"Sector==%d", isector+36);
200     cnoise = NoiseSector(chcut1, chdesc,70,70);
201     sprintf(chshell,"%s/noise_sec%dOROC.eps", dirname,isector);
202     cnoise->Print(chshell);
203     //
204     // systematic
205     //
206     c.cd();
207     sprintf(chdesc,"%s Sector %d",runDesc.Data(), isector);
208     sprintf(chcut1,"Cl.fDetector==%d||Cl.fDetector==%d", isector, isector+36);
209     SysYX(chcut1,chdesc);
210     sprintf(chshell,"%s/deltayx_sec%d.eps", dirname,isector);
211     c.Print(chshell);
212     c.cd();
213     SysZX(chcut1,chdesc);
214     sprintf(chshell,"%s/deltazx_sec%d.eps", dirname,isector);
215     c.Print(chshell);
216     
217     //
218     // picture prf
219     //  
220     if (ncl<500) continue;  //not enough statistic
221     //
222     sprintf(chdesc,"%s Sector %d",runDesc.Data(), isector);
223     sprintf(chcut1,"Cl.fDetector==%d", isector);
224     sprintf(chcut2,"Cl.fDetector==%d", isector+36);
225     c.cd();
226     PRFYZ(chcut1, chcut2,chdesc);
227     sprintf(chshell,"%s/prfyz_sec%d.eps", dirname,isector);
228     c.Print(chshell);
229     sprintf(chcut1,"Sector==%d", isector);
230     sprintf(chcut2,"Sector==%d", isector+36);
231     PRFZZ(chcut1, chcut2,chdesc);
232     sprintf(chshell,"%s/prfzz_sec%d.eps", dirname,isector);
233     c.Print(chshell);
234     //
235     // y resolution
236     //
237     sprintf(chdesc,"%s Sector %d",runDesc.Data(), isector);
238     sprintf(chcut1,"Cl.fDetector==%d", isector);
239     sprintf(chcut2,"Cl.fDetector==%d", isector+36);
240     c.cd();
241     ResYZ(chcut1, chcut2,chdesc);
242     sprintf(chshell,"%s/resyz_sec%d.eps", dirname,isector);
243     c.Print(chshell);
244     //
245   }
246 }
247
248
249
250
251
252
253 void AddChains(Int_t run){
254   /// add files to the chains + check consistency
255
256   ifstream in0;
257   ifstream in1;
258   ifstream in2;
259   ifstream in3;
260   ifstream in4;
261   TString sfile;
262   char strcl[100];
263   runDesc+=run;
264   // TPC tracks
265   //
266   sprintf(strcl,"ls  *%d*/TPCtracks.root > files.txt", run);
267   gSystem->Exec(strcl);
268   in0.open("files.txt");
269   for (;in0>>sfile;){
270     if (sfile.Length()==0) break;
271     printf("%s\n",sfile.Data());
272     TFile f(sfile.Data());
273     TTree * tree = (TTree*)f.Get("Tracks");
274     if (tree){      
275       f.Close();
276       chaincl.Add(sfile.Data());
277     }
278   }
279   //
280   // Fitted signals
281   sprintf(strcl,"ls  *%d*/FitSignal.root > files.txt", run);
282   gSystem->Exec(strcl);
283   in1.open("files.txt");
284   for (;in1>>sfile;){
285     if (sfile.Length()==0) break;
286     printf("%s\n",sfile.Data()); 
287     TFile f(sfile.Data());
288     TTree * tree =(TTree*)f.Get("Fit");
289     if (tree){      
290       f.Close();
291       chainFit.Add(sfile.Data());
292     }
293   }
294   //
295   // Fitted pedestal
296   sprintf(strcl,"ls  *%d*/TPCsignal.root > files.txt", run);
297   gSystem->Exec(strcl);
298   in2.open("files.txt");
299   for (;in2>>sfile;){
300     if (sfile.Length()==0) break;
301     printf("%s\n",sfile.Data());
302     TFile f(sfile.Data());
303     TTree * tree =(TTree*)f.Get("Fit");
304     if (tree){      
305       f.Close();
306       chainPed.Add(sfile.Data());
307     }
308     //    chainPed.Add(sfile.Data());
309   }
310   //
311   // Random signals
312   sprintf(strcl,"ls  *%d*/TPCsignal.root > files.txt", run);
313   gSystem->Exec(strcl);
314   in4.open("files.txt");
315   for (;in4>>sfile;){
316     if (sfile.Length()==0) break;
317     printf("%s\n",sfile.Data());
318     TFile f(sfile.Data());
319     TTree * tree =(TTree*)f.Get("SignalB");
320     if (tree){      
321       f.Close();
322       chainSignal.Add(sfile.Data());
323     }
324     //    chainPed.Add(sfile.Data());
325   }
326   //
327   // Rec points trees
328   //
329   printf("\n IMPORT REC points");
330   sprintf(strcl,"ls  *%d*/*RecPoints* > files.txt", run);
331   gSystem->Exec(strcl);
332   in3.open("files.txt");
333   for (;in3>>sfile;){
334     if (sfile.Length()==0) break;
335     printf("%s\n",sfile.Data());    
336     TFile fcl(sfile.Data());
337     char tname[100];
338     sprintf(tname,"%s/%s/TreeR",sfile.Data(),fcl.GetListOfKeys()->At(0)->GetName());
339     chaincl2.Add(tname);
340     //    chainPed.Add(sfile.Data());
341   }
342
343   comp.fTree = &chaincl;
344   compF.fTree = &chainFit;
345   compP.fTree = &chainPed;
346 }
347
348 void Select(){
349   /// base cut on the tracks
350
351   comp.fTree->Draw(">>listTracks","Etrack.fTPCncls>30&&abs(Etrack.fIp.fP[4])<1");
352   comp.fTree->SetEventList(listTracks);
353   //
354   compF.fTree->Draw(">>listFitS","p2>0&&p2<5&&p1<900&&p0<10000&&p4<1&&p4>0&&p5<p3&&chi2<150");
355   compF.fTree->SetEventList(listFitS);
356 }
357
358 void SelectLaser(){
359   /// base cut on the tracks
360
361   comp.fTree->Draw(">>listTracks","Etrack.fTPCncls>20&&abs(Etrack.fIp.fP[4])<1&&abs(Etrack.fIp.fP[3])<0.01");
362   comp.fTree->SetEventList(listTracks);
363   //
364   compF.fTree->Draw(">>listFitS","p2>0&&p2<5&&p1<900&&p0<10000&&p4<1&&p4>0&&p5<p3&&chi2<150");
365   compF.fTree->SetEventList(listFitS);
366   //
367   // make default aliases
368   //
369   //  laser z beam
370   comp.fTree->SetAlias("lz0","abs(Etrack.fIp.fP[1]-20)<5");
371   comp.fTree->SetAlias("lz1","abs(Etrack.fIp.fP[1]-70)<20");
372   comp.fTree->SetAlias("lz2","abs(Etrack.fIp.fP[1]-150)<20");
373   comp.fTree->SetAlias("lz3","abs(Etrack.fIp.fP[1]-210)<20");
374
375 }
376
377
378
379 void PRFYZ(TCut cut0, TCut cut1,  char * description){
380   /// plot Pad response function as funtion of drift z
381
382   TF1 * f1 = new TF1("fdiff","sqrt([0]*[0]+(250-x)*[1]*[1])");
383   f1->SetParameter(1,0.2);
384   f1->SetParameter(0,0.2);
385   comp.DrawXY("abs(Cl.fZ)","sqrt(Cl.fSigmaY2)","abs(Track.fTrackPoints.GetAngleY())<0.05","Track.fTrackPoints.fTX>0"+cut0,5,10,240,-0,1);
386   TH1F * prfInnerY = (TH1F*)comp.fMean->Clone();
387
388   comp.DrawXY("abs(Cl.fZ)","sqrt(Cl.fSigmaY2)","abs(Track.fTrackPoints.GetAngleY())<0.05","Track.fTrackPoints.fTX>0"+cut1,5,10,240,-0,1);
389   TH1F * prfOuterY = (TH1F*)comp.fMean->Clone();
390   //
391   //
392   prfOuterY->SetMinimum(0);
393   prfOuterY->SetMarkerStyle(23);
394   prfInnerY->SetMarkerStyle(24);
395   prfOuterY->SetXTitle("Z position (cm)");
396   prfOuterY->SetYTitle("PRF width (cm)");
397   char chouter[100];
398   char chinner[100];
399   prfOuterY->Fit(f1);
400   sprintf(chouter,"Outer sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
401   prfInnerY->Fit(f1);
402   sprintf(chinner,"Inner sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
403   prfOuterY->Draw();
404   prfInnerY->Draw("same");
405   TString desc = description;
406   TLegend *legend = new TLegend(0.25,0.12,0.85,0.35, desc+"\nTPC cluster shape Fit: #sigma = #sqrt{p_{0}^{2}+(z_{d}-z)p_{1}^{2}}");
407   legend->SetBorderSize(1);
408   legend->AddEntry(prfOuterY,chouter);
409   legend->AddEntry(prfInnerY,chinner);
410   legend->Draw();
411 }
412
413
414
415 void PRFZZ(TCut cut0, TCut cut1,  char * description){
416   TF1 * f1 = new TF1("fdiff","sqrt([0]*[0]+x*[1]*[1])");
417   f1->SetParameter(1,0.2);
418   f1->SetParameter(0,0.2);
419   compF.DrawXY("p1*0.285","p2*0.285","p2>0",cut0,8,20,250,-0,2);
420   TH1F * prfInnerY = (TH1F*)compF.fMean->Clone();
421   compF.DrawXY("p1*0.285","p2*0.285","p2>0",cut1,8,20,250,-0,2);
422   TH1F * prfOuterY = (TH1F*)compF.fMean->Clone();
423   //
424   //
425   prfOuterY->SetMinimum(0);
426   prfOuterY->SetMarkerStyle(23);
427   prfInnerY->SetMarkerStyle(24);
428   prfOuterY->SetXTitle("Drift length(cm)");
429   prfOuterY->SetYTitle("Z Sigma (cm)");
430   char chouter[100];
431   char chinner[100];
432   prfOuterY->Fit(f1);
433   sprintf(chouter,"Outer sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
434   prfInnerY->Fit(f1);
435   sprintf(chinner,"Inner sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
436   prfOuterY->Draw();
437   prfInnerY->Draw("same");
438   TString desc = description;
439   TLegend *legend = new TLegend(0.25,0.12,0.85,0.35, desc+"TPC signal shape Fit: #sigma = #sqrt{p_{0}^{2}+(z)p_{1}^{2}}");
440   legend->SetBorderSize(1);
441   legend->AddEntry(prfOuterY,chouter);
442   legend->AddEntry(prfInnerY,chinner);
443   legend->Draw();
444 }
445
446
447 void ResYZ(TCut cut0, TCut cut1,  char * description){
448   /// resolution in y coordinate as function of z
449
450   TF1 * f1 = new TF1("fdiff","sqrt([0]*[0]+(250-x)*[1]*[1])");
451   f1->SetParameter(1,0.2);
452   f1->SetParameter(0,0.2);
453   comp.DrawXY("abs(Cl.fZ)","Track.fTrackPoints.GetY()-Cl.GetY()","abs(Track.fTrackPoints.GetAngleY())<0.05","Track.fTrackPoints.fTX>0"+cut0,5,10,240,-0.5,0.5);
454   TH1F * prfInnerY = (TH1F*)comp.fRes->Clone();
455
456   comp.DrawXY("abs(Cl.fZ)","Track.fTrackPoints.GetY()-Cl.GetY()","abs(Track.fTrackPoints.GetAngleY())<0.05","Track.fTrackPoints.fTX>0"+cut1,5,10,240,-0.5,0.5);
457   TH1F * prfOuterY = (TH1F*)comp.fRes->Clone();
458   //
459   //
460   prfOuterY->SetMinimum(0);
461   prfOuterY->SetMaximum(0.15);  
462   prfOuterY->SetMarkerStyle(23);
463   prfInnerY->SetMarkerStyle(24);
464   prfOuterY->SetXTitle("Z position (cm)");
465   prfOuterY->SetYTitle("Y resolution (cm)");
466   char chouter[100];
467   char chinner[100];
468   prfOuterY->Fit(f1);
469   sprintf(chouter,"Outer sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
470   prfInnerY->Fit(f1);
471   sprintf(chinner,"Inner sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
472   prfOuterY->Draw();
473   prfInnerY->Draw("same");
474   TString desc = description;
475   TLegend *legend = new TLegend(0.25,0.12,0.85,0.35, desc+"TPC cluster resolution: #sigma = #sqrt{p_{0}^{2}+(z_{d}-z)p_{1}^{2}}");
476   legend->SetBorderSize(1);
477   legend->AddEntry(prfOuterY,chouter);
478   legend->AddEntry(prfInnerY,chinner);
479   legend->Draw();
480 }
481
482 void SysYX(TCut cut0,  char * description){
483   ///
484
485   TProfile * profA = new TProfile("profY","profY",70,89,250);
486   comp.fTree->Draw("Cl.fY-Track.fTrackPoints.GetY():Track.fTrackPoints.GetX()>>profY","abs(Cl.fY-Track.fTrackPoints.GetY())<1&&Track.fTrackPoints.fTX>10"+cut0,"prof");
487   profA->SetXTitle("Local X (cm)");
488   profA->SetYTitle("Mean #Delta Y (cm)");  
489   TLegend *legend = new TLegend(0.55,0.25,0.85,0.30, description);
490   legend->Draw();
491 }
492
493 void SysZX(TCut cut0,  char * description){
494   ///
495
496   TProfile * profA = new TProfile("profZ","profZ",70,89,250);
497   comp.fTree->Draw("abs(Cl.fZ)-abs(Track.fTrackPoints.GetZ()):Track.fTrackPoints.GetX()>>profZ","abs(abs(Cl.fZ)-abs(Track.fTrackPoints.GetZ()))<1&&Track.fTrackPoints.fTX>10"+cut0,"prof");
498   profA->SetXTitle("Local X (cm)");
499   profA->SetYTitle("Mean #Delta Z (cm)");  
500   TLegend *legend = new TLegend(0.55,0.25,0.85,0.30, description);
501   legend->Draw(); 
502 }
503
504 TProfile * ProfileMaxRow(TCut cut0, char *name, Int_t max){ 
505   /// make profile histrogram of amplitudes
506
507   TProfile *profA = new TProfile(name,name,max,0,max-1);
508   char expr[100];
509   sprintf(expr,"Cl.fMax:Cl.fRow>>%s",name);
510   comp.fTree->Draw(expr,"abs(Cl.fZ)>0&&Cl.fMax<500"+cut0,"prof");
511   profA->SetXTitle("Pad Row");
512   profA->SetYTitle("Amplitude at maxima (ADC)");
513   return profA;
514 }
515
516 TProfile * ProfileMaxPhi(TCut cut0, char *name, Int_t max){ 
517   /// make profile histrogram of amplitudes
518
519   TProfile *profA = new TProfile(name,name,max,-0.14,0.14);
520   char expr[100];
521   sprintf(expr,"Cl.fMax:Cl.fY/Cl.fX>>%s",name);
522   comp.fTree->Draw(expr,"abs(Cl.fZ)>0&&Cl.fMax<500"+cut0,"prof");
523   profA->SetXTitle("Local #phi(rad)");
524   profA->SetYTitle("Amplitude at maxima (ADC)");
525   return profA; 
526 }
527
528 TProfile * ProfileQRow(TCut cut0, char *name, Int_t max){ 
529   /// make profile histrogram of amplitudes
530
531   TProfile *profA = new TProfile(name,name,max,0,max-1);
532   char expr[100];
533   sprintf(expr,"Cl.fQ:Cl.fRow>>%s",name);
534   comp.fTree->Draw(expr,"abs(Cl.fZ)>0&&Cl.fMax<500"+cut0,"prof");
535   profA->SetXTitle("Pad Row");
536   profA->SetYTitle("Total charge(ADC)");
537   return profA;
538 }
539
540 TProfile * ProfileQPhi(TCut cut0, char *name, Int_t max){ 
541   /// make profile histrogram of amplitudes
542
543   TProfile *profA = new TProfile(name,name,max,-0.14,0.14);
544   char expr[100];
545   sprintf(expr,"Cl.fQ:Cl.fY/Cl.fX>>%s",name);
546   comp.fTree->Draw(expr,"abs(Cl.fZ)>0&&Cl.fMax<500"+cut0,"prof");
547   profA->SetXTitle("Local #phi(rad)");
548   profA->SetYTitle("Total charge (ADC)");
549   return profA;
550 }
551
552 TProfile * ProfileQZ(TCut cut0, char *name, Int_t max){ 
553   /// make profile histrogram of amplitudes
554
555   TF1 * fline = new TF1("fline","[0]+[1]*[0]*(250-x)");
556   TF1 * f1 = new TF1("f1","[0]*exp(-[1]*(250-x))");
557   TProfile *profA = new TProfile(name,name,max,0,250);
558   char expr[100];
559   sprintf(expr,"Cl.fQ:abs(Cl.fZ)>>%s",name);
560   comp.fTree->Draw(expr,"abs(Cl.fZ)>0&&Cl.fMax<500"+cut0,"prof");
561   profA->SetXTitle("Z position (cm)"); 
562   profA->SetYTitle("Amplitude (ADC)");
563   char chc[100];
564   profA->Fit(fline);
565   f1->SetParameter(0,fline->GetParameter(0));
566   f1->SetParameter(1,fline->GetParameter(1));
567   profA->Fit(f1);
568   sprintf(chc,"Exponential fit params: p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
569   printf("%s",chc);
570   TLegend *legend = new TLegend(0.25,0.12,0.85,0.25, chc);
571   legend->Draw();
572   return profA;
573 }
574
575 TProfile * ProfileMaxZ(TCut cut0, char *name, Int_t max){ 
576   /// make profile histrogram of amplitudes
577
578   TF1 * f1 = new TF1("f1","[0]+[1]*[0]*(250-x)");
579   TProfile *profA = new TProfile(name,name,max,0,250);
580   char expr[100];
581   sprintf(expr,"Cl.fMax:abs(Cl.fZ)>>%s",name);
582   comp.fTree->Draw(expr,"abs(Cl.fZ)>0&&Cl.fMax<500"+cut0,"prof");
583   profA->SetXTitle("Z position (cm)"); 
584   profA->SetYTitle("Amplitude at maxima (ADC)");
585   char chc[100];
586   profA->Fit(f1);
587   sprintf(chc,"p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
588   TLegend *legend = new TLegend(0.25,0.12,0.85,0.35, chc);
589   legend->Draw();
590   return profA;
591 }
592
593
594 void P3Z(TCut cut0, TCut cut1,  char * description){
595   /// first exponenent as function of z drift
596
597   TF1 * f1 = new TF1("fdiff","[0]+[1]/[0]*x");
598   f1->SetParameter(1,0.2);   
599   f1->SetParameter(0,0.2);
600   compF.DrawXY("p1*0.285","p3","Max>250",cut0,5,20,250,-0,2);
601   TH1F * prfInnerY = (TH1F*)compF.fMean->Clone();
602   compF.DrawXY("p1*0.285","p3","Max>250",cut1,5,20,250,-0,2);
603   TH1F * prfOuterY = (TH1F*)compF.fMean->Clone();
604   //
605   //
606   prfOuterY->SetMinimum(0);
607   prfOuterY->SetMaximum(1);
608   prfOuterY->SetMarkerStyle(23);
609   prfInnerY->SetMarkerStyle(24);
610   prfOuterY->SetXTitle("Drift length (cm)");
611   prfOuterY->SetYTitle("Lambda 0 (Time Bin)");
612   char chouter[100];
613   char chinner[100];
614   prfOuterY->Fit(f1);
615   sprintf(chouter,"Outer sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
616   prfInnerY->Fit(f1);
617   sprintf(chinner,"Inner sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
618   prfOuterY->Draw();
619   prfInnerY->Draw("same");
620   TString desc = description;
621   TLegend *legend = new TLegend(0.25,0.12,0.85,0.35, desc+"TPC cluster shape fit Parameter Lambda0 - P3:");
622   legend->SetBorderSize(1);
623   legend->AddEntry(prfOuterY,chouter);
624   legend->AddEntry(prfInnerY,chinner);
625   legend->Draw();
626 }
627
628
629 void P5Z(TCut cut0, TCut cut1,  char * description){
630   /// second exponenent as function of z drift
631
632   TF1 * f1 = new TF1("fdiff","[0]+[1]/[0]*x");
633   f1->SetParameter(1,0.2);
634   f1->SetParameter(0,0.2);
635   compF.DrawXY("p1*0.285","p5","Max>250",cut0,5,20,250,-0,0.2);
636   TH1F * prfInnerY = (TH1F*)compF.fMean->Clone();
637   compF.DrawXY("p1*0.285","p5","Max>250",cut1,5,20,250,-0,0.2);
638   TH1F * prfOuterY = (TH1F*)compF.fMean->Clone();
639   //
640   //
641   prfOuterY->SetMinimum(0);
642   prfOuterY->SetMaximum(0.15);
643   prfOuterY->SetMarkerStyle(23);
644   prfInnerY->SetMarkerStyle(24);
645   prfOuterY->SetXTitle("Drift length (Time Bin)");
646   prfOuterY->SetYTitle("Lambda1 (Time bin)");
647   char chouter[100];
648   char chinner[100];
649   prfOuterY->Fit(f1);
650   sprintf(chouter,"Outer sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
651   prfInnerY->Fit(f1);
652   sprintf(chinner,"Inner sector : p_{0} = %f  p_{1} = %f",f1->GetParameter(0),f1->GetParameter(1));
653   prfOuterY->Draw();
654   prfInnerY->Draw("same");
655    TString desc = description;
656   TLegend *legend = new TLegend(0.25,0.12,0.85,0.35, desc+"TPC cluster shape fit Parameter Lambda1 - P5");
657   legend->SetBorderSize(1);
658   legend->AddEntry(prfOuterY,chouter);
659   legend->AddEntry(prfInnerY,chinner);
660   legend->Draw();
661 }
662
663 TCanvas *  NoiseSector(TCut cut0,  char * description, Int_t maxrow, Int_t maxpad){
664   /// draw plots of the noise
665
666   TCanvas * c = new TCanvas;
667   c->Divide(2,1);
668   c->Draw();
669   c->cd(1);
670   compP.fTree->Draw("GSigma","GSigma<5"+cut0);
671   c->cd(2);
672   Float_t rand  = gRandom->Rndm();
673   char name[100];
674   sprintf(name,"prof%f",rand);
675   TProfile2D * prof= new TProfile2D(name,name,maxrow, 0, maxrow-1, 2*maxpad,-maxpad,maxpad);
676   char expr[100];
677   sprintf(expr,"GSigma:RPad:Row>>%s",name);
678   prof->SetXTitle("Pad row");
679   prof->SetYTitle("Pad number");
680   compP.fTree->Draw(expr,cut0,"profcolz");
681   c->cd(1);
682   TString desc = description;
683   TLegend *legend = new TLegend(0.25,0.30,0.85,0.85, desc+"Noise map");
684   legend->Draw();
685   
686   return c;
687 }
688
689