]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AnalyzeESDtracks.C
Raw2SDigits implemented (A.M., M. Masera)
[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   AnalyzESDtracks(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 //AnalyzESDtracks(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
69 #include  "AliLog.h"
70 #include  "AliMagF.h"
71
72 #include  "AliESD.h"
73 #include  "AliESDfriend.h"
74 #include  "AliESDtrack.h"
75 #include  "AliTracker.h"
76 #include  "AliTPCseed.h"
77 #include  "AliTPCclusterMI.h"
78 #include  "TTreeStream.h"
79 #include  "TF1.h"
80 #include  "TGraph.h"
81 #include  "AliSignalProcesor.h"
82 #include  "TCanvas.h"
83
84
85 void FitSignals(TTree * treeB, TCut cut="Max-Median>150&&RMS06<2&&abs(Median-Mean09)<0.5");
86
87
88 void AnalyzeESDtracks(Int_t run){
89   //
90   // output redirect 
91   TTreeSRedirector  cstream("TPCtracks.root");
92   //
93   // dummy magnetic field
94   AliMagF mag("aaa","aaa",1,1,10);
95   AliTracker::SetFieldMap(&mag,kTRUE);
96   TFile f("AliESDs.root");
97   TTree * tree =(TTree*)f.Get("esdTree"); 
98   AliESD * esd =0;
99   tree->SetBranchAddress("ESD",&esd);
100   AliESDfriend *evf=0;
101   tree->AddFriend("esdFriendTree","AliESDfriends.root");
102   tree->SetBranchAddress("ESDfriend",&evf);
103
104   Int_t nevents = tree->GetEntries();
105   TClonesArray *clusters = new TClonesArray("AliTPCclusterMI",160);
106   for (Int_t irow=0; irow<160; irow++){
107     new ((*clusters)[irow])  AliTPCclusterMI;   // iitial dummy clusters
108   }
109   
110   for (Int_t ievent=0; ievent<nevents; ievent++){
111     tree->GetEntry(ievent);
112     if (!esd) continue;
113     if (!evf) continue;
114     esd->SetESDfriend(evf); //Attach the friend to the ESD      
115     for (Int_t itrack =0; itrack<esd->GetNumberOfTracks(); itrack++){
116       // Int_t itrack = 0;      
117       if (esd->GetTrack(itrack)->GetFriendTrack()==0) continue;
118       AliESDtrack * etrack = esd->GetTrack(itrack);
119       AliESDfriendTrack * ftrack = (AliESDfriendTrack *)esd->GetTrack(itrack)->GetFriendTrack();
120       AliTPCseed * seed =  (AliTPCseed*)(ftrack->GetCalibObject(0));
121       if (!seed) continue;  
122       for (Int_t irow=0; irow<160; irow++){
123         if (seed->GetClusterFast(irow)){
124           AliTPCclusterMI * cl = new ((*clusters)[irow])  AliTPCclusterMI(*(seed->GetClusterFast(irow)));
125           cl->SetLabel(itrack,0);
126         }
127         else{
128           AliTPCclusterMI * cl = (AliTPCclusterMI*)clusters->At(irow);
129           cl->SetX(0); cl->SetY(0); cl->SetZ(0); cl->SetQ(0); cl->SetLabel(-1,0);
130         }
131       }
132       Float_t dEdx = seed->GetdEdx();
133       Float_t dEdxI = seed->CookdEdx(0.05,0.6,0,77);
134       Float_t dEdxO = seed->CookdEdx(0.05,0.6,78,155);
135       Int_t ncl = seed->GetNumberOfClusters();
136       cstream<<"Tracks"<<
137         "Run="<<run<<
138         "Ncl="<<ncl<<
139         "Event="<<ievent<<
140         "dEdx="<<dEdx<<
141         "dEdxI="<<dEdxI<<
142         "dEdxO="<<dEdxO<<
143         "Track.="<<seed<<
144         "Etrack.="<<etrack<<
145         "Cl.="<<clusters<<
146         "\n";
147     }  
148   }
149   //
150   // Fit signal part
151   //
152   TFile fs("TPCsignal.root");
153   TTree *treeB =(TTree*)fs.Get("SignalB");
154   FitSignals(treeB,"Max-Median>150&&RMS06<2.5");
155 }
156
157
158 void FitSignals(TTree * treeB, TCut cut){
159   AliSignalProcesor proc;
160   TF1 * f1 = proc.GetAsymGauss();
161   TTreeSRedirector cstream("FitSignal.root");
162   TFile *f = cstream.GetFile();
163
164   char lname[100];
165   sprintf(lname,"Fit%s", cut.GetTitle());
166   TEventList *list = new TEventList(lname,lname);
167   sprintf(lname,">>Fit%s", cut.GetTitle());
168   treeB->Draw(lname,cut);
169   treeB->SetEventList(list);
170   for (Int_t ievent=0; ievent<list->GetN(); ievent++){
171     char ename[100];
172     sprintf(ename,"Fit%d", ievent);
173     Double_t nsample = treeB->Draw("Graph.fY-Mean09:Graph.fX","","",1,ievent);
174     Double_t * signal  = treeB->GetV1();
175     Double_t * time  = treeB->GetV2();
176     Double_t maxpos =0;
177     Double_t max = 0;
178     for (Int_t ipos = 0; ipos<nsample; ipos++){
179       if (signal[ipos]>max){
180         max    = signal[ipos];
181         maxpos = ipos;
182       }
183     }
184
185     Int_t first = TMath::Max(maxpos-10,0.);
186     Int_t last  = TMath::Min(maxpos+60, nsample);
187     //
188     f->cd();
189     TH1F his(ename,ename,last-first,first,last);
190     for (Int_t ipos=0; ipos<last-first; ipos++){
191       his.SetBinContent(ipos+1,signal[ipos+first]);
192     }
193     treeB->Draw("Sector:Row:Pad","","",1,ievent);
194     Double_t sector = treeB->GetV1()[0];
195     Double_t row    = treeB->GetV2()[0];
196     Double_t pad    = treeB->GetV3()[0];
197     //    TGraph  graph(last-first,&time[first],&signal[first]);
198     f1->SetParameters(0.75*max,maxpos,1.1,0.8,0.25,0.2);
199     //    TH1F * his = (TH1F*)graph.GetHistogram();
200     his.Fit(f1,"q");
201     his.Write(ename);
202     gPad->Clear();
203     his.Draw();
204     gPad->Update();
205     Double_t params[6];
206     for (Int_t ipar=0; ipar<6; ipar++) params[ipar] = f1->GetParameters()[ipar];
207     Double_t chi2 = TFitter::GetFitter()->Chisquare(6,params);
208     TMatrixD cov(6,6);
209     cov.SetMatrixArray(TFitter::GetFitter()->GetCovarianceMatrix());
210     //
211     // tail cancellation
212     //
213     Double_t x0[1000];
214     Double_t x1[1000];
215     Double_t x2[1000];
216     for (Int_t ipos=0; ipos<last-first; ipos++){
217       x0[ipos] = signal[ipos+first];
218     }
219     proc.TailCancelationALTRO1(x0,x1,0.85*0.339,0.09,last-first);
220     proc.TailCancelationALTRO1(x1,x2,0.85,0.789,last-first);
221     //
222     sprintf(ename,"Cancel1_%d", ievent);
223     TH1F his1(ename,ename,last-first,first,last);
224     for (Int_t ipos=0; ipos<last-first; ipos++){
225       his1.SetBinContent(ipos+1,x1[ipos]);
226     }
227     his1.Write(ename);
228     sprintf(ename,"Cancel2_%d", ievent);
229     TH1F his2(ename,ename,last-first,first,last);
230     for (Int_t ipos=0; ipos<last-first; ipos++){
231       his2.SetBinContent(ipos+1,x1[ipos]);
232     }
233     f1->SetParameters(0.75*max,maxpos,1.1,0.8,0.25,0.2);
234     his2.Fit(f1,"q");
235     his2.Write(ename);
236     Double_t params2[6];
237     for (Int_t ipar=0; ipar<6; ipar++) params2[ipar] = f1->GetParameters()[ipar];
238     Double_t chi22 = TFitter::GetFitter()->Chisquare(6,params2);    
239     TMatrixD cov2(6,6);
240     cov2.SetMatrixArray(TFitter::GetFitter()->GetCovarianceMatrix());
241
242     TGraph gr0(last-first, &time[first],x0);
243     TGraph gr1(last-first, &time[first],x1);
244     TGraph gr2(last-first, &time[first],x2);
245     //
246     cstream<<"Fit"<<
247       "Sector="<<sector<<
248       "Row="<<row<<
249       "Pad="<<pad<<
250       "First="<<first<<
251       "Max="<<max<<
252       "MaxPos="<<maxpos<<
253       "chi2="<<chi2<<
254       "chi22="<<chi22<<
255       "Cov="<<&cov<<
256       "Cov2="<<&cov2<<
257       "gr0.="<<&gr0<<
258       "gr1.="<<&gr1<<
259       "gr2.="<<&gr2<<
260       "p0="<<params[0]<<
261       "p1="<<params[1]<<
262       "p2="<<params[2]<<
263       "p3="<<params[3]<<
264       "p4="<<params[4]<<
265       "p5="<<params[5]<<
266       "p02="<<params2[0]<<
267       "p12="<<params2[1]<<
268       "p22="<<params2[2]<<
269       "p32="<<params2[3]<<
270       "p42="<<params2[4]<<
271       "p52="<<params2[5]<<
272       "\n";
273     //    delete his;
274   }
275
276 }
277
278
279
280 TChain *MakeChainCL(Int_t first, Int_t last){
281   TChain *chaincl = new TChain("TreeR","TreeR");
282   //
283   char fname[100];
284   for (Int_t i=first;i<last; i++){
285     if (i>0) sprintf(fname,"TPC.RecPoints%d.root/Event%d/TreeR",i,i);
286     if (i==0) sprintf(fname,"TPC.RecPoints.root/Event%d/TreeR",i);
287     chaincl->Add(fname);
288   }
289   return chaincl;
290 }
291
292 TTree* GetTree(Int_t ievent){
293   char fname[100];
294   char tname[100];
295   if (ievent>0) sprintf(fname,"TPC.RecPoints%d.root",ievent);
296   if (ievent==0) sprintf(fname,"TPC.RecPoints.root");
297   sprintf(tname,"Event%d/TreeR",ievent);
298   TFile * f  = new TFile(fname);
299   TTree * tree = (TTree*)f->Get(tname);
300   return tree;
301
302 }