]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AnalyzeESDtracks.C
Analysis of central electrode signal added (Marian)
[u/mrichter/AliRoot.git] / TPC / AnalyzeESDtracks.C
1 //
2 // Process ESD tracks  - 
3 // Extract TPC tracks  - write them to tree
4 //
5 /*
6   .L AnalyzeESDtracks.C+
7   .L AliGenInfo.C+
8   AnalyzeESDtracks(567);   // process tracks
9   // Tracks are written to the file "TPCtracks.root"
10   // Now yo can analyze it
11   TFile fesd("AliESDs.root");
12   TTree * treeE = (TTree*)fesd.Get("esdTree");
13   TFile f("TPCtracks.root")
14   TTree * tree =(TTree*)f.Get("Tracks");
15   AliComparisonDraw comp;
16   comp->fTree = tree;
17
18 */
19
20 // 
21
22 /*
23 .L AnalyzeESDtracks.C+
24 .L AliGenInfo.C+
25
26 //AnalyzeESDtracks(567);
27
28
29 TFile fesd("AliESDs.root");
30 TTree * treeE = (TTree*)fesd.Get("esdTree");
31 TFile f("TPCtracks.root")
32 TTree * tree =(TTree*)f.Get("Tracks");
33 AliComparisonDraw comp;
34 comp->fTree = tree;
35
36 TFile fs("TPCsignal.root");
37 TTree *treeB =(TTree*)fs.Get("SignalB");
38 TTree *treeN =(TTree*)fs.Get("SignalN");
39 TTree *treeS =(TTree*)fs.Get("Signal");
40 TTree *treef =(TTree*)fs.Get("Fit");
41
42
43 FitSignals(treeB,"Max-Median>100&&RMS06<2.5")
44 TFile ffit("FitSignal.root");
45 TTree * treeF = (TTree*)ffit->Get("Fit");
46
47
48 TChain chaincl("TreeR","TreeR")
49 chaincl.Add("TPC.RecPoints.root/Event0/TreeR")
50 chaincl.Add("TPC.RecPoints1.root/Event1/TreeR")
51 chaincl.Add("TPC.RecPoints2.root/Event2/TreeR")
52 chaincl.Add("TPC.RecPoints3.root/Event3/TreeR")
53
54 */
55
56
57
58 #include "TObject.h"
59 #include "TFile.h"
60 #include "TTree.h"
61 #include "TChain.h"
62 #include "TBranch.h"
63 #include "TTreeStream.h"
64 #include "TEventList.h"
65 #include "TCut.h"
66 #include "TFitter.h"
67 #include  "TMatrixD.h"
68 #include  "TRobustEstimator.h"
69
70 #include  "AliLog.h"
71 #include  "AliMagF.h"
72
73 #include  "AliESD.h"
74 #include  "AliESDfriend.h"
75 #include  "AliESDtrack.h"
76 #include  "AliTracker.h"
77 #include  "AliTPCseed.h"
78 #include  "AliTPCclusterMI.h"
79 #include  "AliTPCParamSR.h"
80 #include  "AliTPCROC.h"
81
82
83 #include  "TTreeStream.h"
84 #include  "TF1.h"
85 #include  "TGraph.h"
86 #include  "AliSignalProcesor.h"
87 #include  "TCanvas.h"
88
89
90 void FitSignals(TTree * treeB, TCut cut="Max-Median>150&&RMS06<2&&abs(Median-Mean09)<0.5", Int_t maxS=100);
91 void LaserCalib(TTreeSRedirector & cstream, TTree * chain, Float_t tmin, Float_t tmax, Float_t fraction=0.8);
92
93 void AnalyzeESDtracks(Int_t run){
94   //
95   // output redirect 
96   TTreeSRedirector *  pcstream = new TTreeSRedirector("TPCtracks.root");
97   TTreeSRedirector &cstream = *pcstream;
98   //
99   // dummy magnetic field
100   AliMagF mag("aaa","aaa",1,1,10);
101   AliTracker::SetFieldMap(&mag,kTRUE);
102   TFile f("AliESDs.root");
103   TTree * tree =(TTree*)f.Get("esdTree"); 
104   AliESD * esd =0;
105   tree->SetBranchAddress("ESD",&esd);
106   AliESDfriend *evf=0;
107   tree->AddFriend("esdFriendTree","AliESDfriends.root");
108   tree->SetBranchAddress("ESDfriend",&evf);
109
110   Int_t nevents = tree->GetEntries();
111   TClonesArray *clusters = new TClonesArray("AliTPCclusterMI",160);
112   for (Int_t irow=0; irow<160; irow++){
113     new ((*clusters)[irow])  AliTPCclusterMI;   // iitial dummy clusters
114   }
115   
116   for (Int_t ievent=0; ievent<nevents; ievent++){
117     tree->GetEntry(ievent);
118     if (!esd) continue;
119     if (!evf) continue;
120     esd->SetESDfriend(evf); //Attach the friend to the ESD      
121     for (Int_t itrack =0; itrack<esd->GetNumberOfTracks(); itrack++){
122       // Int_t itrack = 0;      
123       if (esd->GetTrack(itrack)->GetFriendTrack()==0) continue;
124       AliESDtrack * etrack = esd->GetTrack(itrack);
125       AliESDfriendTrack * ftrack = (AliESDfriendTrack *)esd->GetTrack(itrack)->GetFriendTrack();
126       AliTPCseed * seed =  (AliTPCseed*)(ftrack->GetCalibObject(0));
127       if (!seed) continue;  
128       for (Int_t irow=0; irow<160; irow++){
129         if (seed->GetClusterFast(irow)){
130           AliTPCclusterMI * cl = new ((*clusters)[irow])  AliTPCclusterMI(*(seed->GetClusterFast(irow)));
131           cl->SetLabel(itrack,0);
132         }
133         else{
134           AliTPCclusterMI * cl = (AliTPCclusterMI*)clusters->At(irow);
135           cl->SetX(0); cl->SetY(0); cl->SetZ(0); cl->SetQ(0); cl->SetLabel(-1,0);
136         }
137       }
138       Float_t dEdx = seed->GetdEdx();
139       Float_t dEdxI = seed->CookdEdx(0.05,0.6,0,77);
140       Float_t dEdxO = seed->CookdEdx(0.05,0.6,78,155);
141       Int_t ncl = seed->GetNumberOfClusters();
142       cstream<<"Tracks"<<
143         "Run="<<run<<
144         "Ncl="<<ncl<<
145         "Event="<<ievent<<
146         "dEdx="<<dEdx<<
147         "dEdxI="<<dEdxI<<
148         "dEdxO="<<dEdxO<<
149         "Track.="<<seed<<
150         "Etrack.="<<etrack<<
151         "Cl.="<<clusters<<
152         "\n";
153     }  
154   }
155   delete pcstream;
156   //
157   // Fit signal part
158   //
159   TFile fs("TPCsignal.root");
160   TTree *treeB =(TTree*)fs.Get("SignalB");
161   //
162   // Fit central electrode part
163   //
164   TTreeSRedirector * pcestream = new  TTreeSRedirector("TimeRoot.root");
165   TTree * treece = (TTree*)fs.Get("Signalce");
166   if (tree) {
167     LaserCalib(*pcestream, treece, 900,1000, 0.8);
168     delete pcestream;
169   }
170   FitSignals(treeB,"Max-Median>150&&RMS06<1.0&&RMS09<1.5&&abs(Median-Mean09)<0.2&&abs(Mean06-Mean09)<0.2",200);
171   //
172 }
173
174
175 void FitSignals(TTree * treeB, TCut cut, Int_t max){
176   AliSignalProcesor proc;
177   TF1 * f1 = proc.GetAsymGauss();
178   TTreeSRedirector cstream("FitSignal.root");
179   TFile *f = cstream.GetFile();
180
181   char lname[100];
182   sprintf(lname,"Fit%s", cut.GetTitle());
183   TEventList *list = new TEventList(lname,lname);
184   sprintf(lname,">>Fit%s", cut.GetTitle());
185   treeB->Draw(lname,cut);
186   treeB->SetEventList(list);
187   Int_t nFits=0;
188   for (Int_t ievent=0; ievent<list->GetN(); ievent++){
189     if (nFits>max) break;
190     if (nFits%50==0) printf("%d\n",nFits);
191     char ename[100];
192     sprintf(ename,"Fit%d", ievent);
193     Double_t nsample = treeB->Draw("Graph.fY-Mean09:Graph.fX","","",1,ievent);
194     Double_t * signal  = treeB->GetV1();
195     Double_t * time  = treeB->GetV2();
196     Double_t maxpos =0;
197     Double_t max = 0;
198     for (Int_t ipos = 0; ipos<nsample; ipos++){
199       if (signal[ipos]>max){
200         max    = signal[ipos];
201         maxpos = ipos;
202       }
203     }
204
205     Int_t first = TMath::Max(maxpos-10,0.);
206     Int_t last  = TMath::Min(maxpos+60, nsample);
207     //
208     f->cd();
209     TH1F his(ename,ename,last-first,first,last);
210     for (Int_t ipos=0; ipos<last-first; ipos++){
211       his.SetBinContent(ipos+1,signal[ipos+first]);
212     }
213     treeB->Draw("Sector:Row:Pad","","",1,ievent);
214     Double_t sector = treeB->GetV1()[0];
215     Double_t row    = treeB->GetV2()[0];
216     Double_t pad    = treeB->GetV3()[0];
217     //    TGraph  graph(last-first,&time[first],&signal[first]);
218     f1->SetParameters(0.75*max,maxpos,1.1,0.8,0.25,0.2);
219     //    TH1F * his = (TH1F*)graph.GetHistogram();
220     his.Fit(f1,"q");
221     his.Write(ename);
222     gPad->Clear();
223     his.Draw();
224     gPad->Update();
225     Double_t params[6];
226     for (Int_t ipar=0; ipar<6; ipar++) params[ipar] = f1->GetParameters()[ipar];
227     Double_t chi2 = TFitter::GetFitter()->Chisquare(6,params);
228     TMatrixD cov(6,6);
229     cov.SetMatrixArray(TFitter::GetFitter()->GetCovarianceMatrix());
230     //
231     // tail cancellation
232     //
233     Double_t x0[1000];
234     Double_t x1[1000];
235     Double_t x2[1000];
236     for (Int_t ipos=0; ipos<last-first; ipos++){
237       x0[ipos] = signal[ipos+first];
238     }
239     proc.TailCancelationALTRO1(x0,x1,0.85*0.339,0.09,last-first);
240     proc.TailCancelationALTRO1(x1,x2,0.85,0.789,last-first);
241     //
242     sprintf(ename,"Cancel1_%d", ievent);
243     TH1F his1(ename,ename,last-first,first,last);
244     for (Int_t ipos=0; ipos<last-first; ipos++){
245       his1.SetBinContent(ipos+1,x1[ipos]);
246     }
247     his1.Write(ename);
248     sprintf(ename,"Cancel2_%d", ievent);
249     TH1F his2(ename,ename,last-first,first,last);
250     for (Int_t ipos=0; ipos<last-first; ipos++){
251       his2.SetBinContent(ipos+1,x1[ipos]);
252     }
253     f1->SetParameters(0.75*max,maxpos,1.1,0.8,0.25,0.2);
254     his2.Fit(f1,"q");
255     his2.Write(ename);
256     Double_t params2[6];
257     for (Int_t ipar=0; ipar<6; ipar++) params2[ipar] = f1->GetParameters()[ipar];
258     Double_t chi22 = TFitter::GetFitter()->Chisquare(6,params2);    
259     TMatrixD cov2(6,6);
260     cov2.SetMatrixArray(TFitter::GetFitter()->GetCovarianceMatrix());
261
262     TGraph gr0(last-first, &time[first],x0);
263     TGraph gr1(last-first, &time[first],x1);
264     TGraph gr2(last-first, &time[first],x2);
265     //
266     cstream<<"Fit"<<
267       "Sector="<<sector<<
268       "Row="<<row<<
269       "Pad="<<pad<<
270       "First="<<first<<
271       "Max="<<max<<
272       "MaxPos="<<maxpos<<
273       "chi2="<<chi2<<
274       "chi22="<<chi22<<
275       "Cov="<<&cov<<
276       "Cov2="<<&cov2<<
277       "gr0.="<<&gr0<<
278       "gr1.="<<&gr1<<
279       "gr2.="<<&gr2<<
280       "p0="<<params[0]<<
281       "p1="<<params[1]<<
282       "p2="<<params[2]<<
283       "p3="<<params[3]<<
284       "p4="<<params[4]<<
285       "p5="<<params[5]<<
286       "p02="<<params2[0]<<
287       "p12="<<params2[1]<<
288       "p22="<<params2[2]<<
289       "p32="<<params2[3]<<
290       "p42="<<params2[4]<<
291       "p52="<<params2[5]<<
292       "\n";
293     //    delete his;
294     nFits++;
295   }
296
297 }
298
299
300 void LaserCalib(TTreeSRedirector & cstream, TTree * chain, Float_t tmin, Float_t tmax, Float_t fraction){
301   //
302   //
303   //
304   AliTPCParamSR param;
305   param.Update();
306   TFile fce("TPCsignal.root");
307   TTree   * treece =(TTree*)fce.Get("Signalce");
308   if (chain) treece=chain;
309   //
310   TBranch * brsector  = treece->GetBranch("Sector");
311   TBranch * brpad     = treece->GetBranch("Pad");
312   TBranch * brrow     = treece->GetBranch("Row");
313   //
314   TBranch * brtime    = treece->GetBranch("Time");
315   TBranch * brrms     = treece->GetBranch("RMS06");
316   TBranch * brmax     = treece->GetBranch("Max");
317   TBranch * brsum     = treece->GetBranch("Qsum");
318
319   Int_t sector, pad, row=0;
320   Double_t time=0, rms=0, qMax=0, qSum=0;
321   brsector->SetAddress(&sector);
322   brrow->SetAddress(&row);
323   brpad->SetAddress(&pad);
324   
325   brtime->SetAddress(&time);
326   brrms->SetAddress(&rms);
327   brmax->SetAddress(&qMax);
328   brsum->SetAddress(&qSum);
329
330
331   brsector->GetEntry(0);
332   //
333   Int_t firstSector  = sector;
334   Int_t lastSector   = sector;
335   Int_t fentry = 0;
336   Int_t lentry = 0;
337   //
338   // find time offset for differnt events
339   //
340   Int_t count = 0;
341   Double_t padTimes[500000];
342   TRobustEstimator restim;
343   Double_t meanS[72], sigmaS[72];
344   Int_t   firstS[72], lastS[72];
345   Double_t   sectorsID[72];
346   for (Int_t isector=0;isector<72; isector++){
347     firstS[isector]=-1; 
348     lastS[isector] =-1;
349   };
350
351   for (Int_t ientry=0; ientry<treece->GetEntriesFast(); ientry++){
352     treece->GetEvent(ientry);
353     //
354     if (sector!=lastSector && sector==firstSector){
355       //if (sector!=lastSector){
356       lentry = ientry;
357       printf("\nEvent\t%d\tFirst\t%d\tLast\t%d\t%d\n",count, fentry, lentry, lentry-fentry);
358       //
359       //
360       Int_t ngood=0;      
361       for (Int_t ientry2=fentry; ientry2<lentry; ientry2++){
362         //      brtime->GetEvent(ientry2);
363         //  brsector->GetEvent(ientry2);
364         treece->GetEvent(ientry2);
365         if (time>tmin&&time<tmax){
366           padTimes[ngood]=time;
367           ngood++;      
368           if (firstS[sector]<0)  firstS[sector]= ngood;
369           if (firstS[sector]>=0) lastS[sector] = ngood;
370         }       
371       }
372       //
373       //
374       Double_t mean,sigma;
375       restim.EvaluateUni(ngood,padTimes,mean, sigma,int(ngood*fraction));
376       printf("Event\t%d\t%f\t%f\n",count, mean, sigma);
377       for (Int_t isector=0; isector<72; isector++){
378         sectorsID[isector]=sector;
379         if (firstS[isector]>=0 &&lastS[isector]>=0 && lastS[isector]>firstS[isector] ){
380           Int_t ngoodS = lastS[isector]-firstS[isector];
381           restim.EvaluateUni(ngoodS, &padTimes[firstS[isector]],meanS[isector], 
382                              sigmaS[isector],int(ngoodS*fraction));
383         }
384       }
385       TGraph  graphM(72,sectorsID,meanS);
386       TGraph  graphS(72,sectorsID,sigmaS);
387       cstream<<"TimeS"<<
388         "Event="<<count<<
389         "GM.="<<&graphM<<
390         "GS.="<<&graphS<<
391         "\n";
392       
393
394       for (Int_t ientry2=fentry; ientry2<lentry-1; ientry2++){
395         treece->GetEvent(ientry2);
396         Double_t x      = param.GetPadRowRadii(sector,row);
397         Int_t    maxpad = AliTPCROC::Instance()->GetNPads(sector,row);
398         Double_t y = (pad - 2.5 - 0.5*maxpad)*param.GetPadPitchWidth(sector);
399         Double_t alpha = TMath::DegToRad()*(10.+20.*(sector%18));       
400         Double_t gx = x*TMath::Cos(alpha)-y*TMath::Sin(alpha);
401         Double_t gy = y*TMath::Cos(alpha)+x*TMath::Sin(alpha);
402         
403         Int_t npadS = lastS[sector]-firstS[sector];
404         cstream<<"Time"<<
405           "Event="<<count<<
406           "x="<<x<<
407           "y="<<y<<
408           "gx="<<gx<<
409           "gy="<<gy<<
410           "Sector="<<sector<<
411           "Row="<<row<<
412           "Pad="<<pad<<
413           "Time="<<time<<
414           "RMS="<<rms<<
415           "Time0="<<mean<<
416           "Sigma0="<<sigma<< 
417           "TimeS0="<<meanS[sector]<<
418           "SigmaS0="<<sigmaS[sector]<<
419           "npad0="<<ngood<<
420           "npadS="<<npadS<<
421           "Max="<<qMax<<
422           "Sum="<<qSum<<
423           "\n";
424       }
425       treece->GetEvent(ientry);
426       fentry = ientry;
427       count++;      
428       for (Int_t isector=0;isector<72; isector++){
429         firstS[isector]=-1; 
430         lastS[isector] =-1;
431       }
432     }
433     lastSector=sector;
434   }
435 }
436
437
438
439
440 TChain *MakeChainCL(Int_t first, Int_t last){
441   TChain *chaincl = new TChain("TreeR","TreeR");
442   //
443   char fname[100];
444   for (Int_t i=first;i<last; i++){
445     if (i>0) sprintf(fname,"TPC.RecPoints%d.root/Event%d/TreeR",i,i);
446     if (i==0) sprintf(fname,"TPC.RecPoints.root/Event%d/TreeR",i);
447     chaincl->Add(fname);
448   }
449   return chaincl;
450 }
451
452 TTree* GetTree(Int_t ievent){
453   char fname[100];
454   char tname[100];
455   if (ievent>0) sprintf(fname,"TPC.RecPoints%d.root",ievent);
456   if (ievent==0) sprintf(fname,"TPC.RecPoints.root");
457   sprintf(tname,"Event%d/TreeR",ievent);
458   TFile * f  = new TFile(fname);
459   TTree * tree = (TTree*)f->Get(tname);
460   return tree;
461
462 }