]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaAnalysis/AliAnalysisTaskTRDmon.cxx
Moveing AliNeutralTrackParam to STEERBase to avoid dependence of libAOD on libESD
[u/mrichter/AliRoot.git] / TRD / qaAnalysis / AliAnalysisTaskTRDmon.cxx
1 //**************************************************************************
2 //* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 //*                                                                        *
4 //* Author: The ALICE Off-line Project.                                    *
5 //* Contributors are mentioned in the code where appropriate.              *
6 //*                                                                        *
7 //* Permission to use, copy, modify and distribute this software and its   *
8 //* documentation strictly for non-commercial purposes is hereby granted   *
9 //* without fee, provided that the above copyright notice appears in all   *
10 //* copies and that both the copyright notice and this permission notice   *
11 //* appear in the supporting documentation. The authors make no claims     *
12 //* about the suitability of this software for any purpose. It is          *
13 //* provided "as is" without express or implied warranty.                  *
14 //**************************************************************************/
15 //////////////////////////////////////////////////
16 // Analysis task used for TRD monitoring        //
17 // Input:   TChain of ESDs                      //
18 // Output:  TObjArray of histograms             //
19 //////////////////////////////////////////////////
20
21 #include <TObjArray.h>
22 #include <TObject.h>
23 #include <TParticle.h>
24 #include <TH1F.h>
25 #include <TH2F.h>
26 #include <TH3F.h>
27 #include <TProfile.h>
28 #include <TProfile2D.h>
29 #include <TTree.h>
30 #include <TChain.h>
31 #include <TMath.h>
32 #include <TLegend.h>
33 #include <TCanvas.h>
34 #include <TPad.h>
35 #include <TLatex.h>
36
37 #include "AliAnalysisTask.h"
38 #include "AliAnalysisManager.h"
39 #include "AliLog.h"
40 #include "AliESDEvent.h"
41 #include "AliESDfriend.h"
42 #include "AliESDInputHandler.h"
43 #include "AliESDtrack.h"
44 #include "AliESDfriendTrack.h"
45 #include "AliTRDtrackV1.h"
46 #include "AliTRDseedV1.h"
47 #include "AliTRDcluster.h"
48 #include "AliTrackPointArray.h"
49
50 #include <cstdio>
51
52 #ifndef ALIANALYSISTASKTRDMON_H
53 #include "AliAnalysisTaskTRDmon.h"
54 #endif
55
56 ClassImp(AliAnalysisTaskTRDmon)
57
58 AliAnalysisTaskTRDmon::AliAnalysisTaskTRDmon(const Char_t *name):
59   AliAnalysisTask(name, ""),
60   fESD(0x0),
61   fESDfriend(0x0),
62   fEventTriggerName("CINT1B-ABCE-NOPF-ALL"),
63   fIsCollisionEvent(kTRUE),
64   fOutStorage(0x0),
65   fHzvert1(0x0),
66   fHzvert2(0x0),
67   fHntracks(0x0),
68   fHntracks2(0x0),
69   fHntracks3(0x0),
70   fHdca(0x0),
71   fHdcaz(0x0),
72   fHpt(0x0),
73   fHpt2(0x0),
74   fHpt3(0x0),
75   fHpt3n(0x0),
76   fHpt4(0x0),
77   fHpt4n(0x0),
78   fHtheta(0x0),
79   fHphi(0x0),
80   fHtpccl(0x0),
81   fHtpccl2(0x0),
82   fHdedxp(0x0),
83   fHetaphi(0x0),
84   fHetancl(0x0),
85   fHphincl(0x0),
86   fHtrdtr(0x0),
87   fHtrdtr2(0x0),
88   fHph(0x0),
89   fHph2d(0x0),
90   fHncltrkl(0x0),
91   fHntrkl(0x0),
92   fHntrklVsP(0x0),
93   fHsm(0x0),
94   fHetantr(0x0),
95   fHphintr(0x0),
96   fHcltime(0x0),
97   fHcldiff(0x0),
98   fHxyA(0x0),
99   fHxyC(0x0),
100   fHPropagXY(0x0),
101   fHPropagRZ(0x0),
102   fHnFriendTracks(0x0),
103   fHnCalibObjects(0x0),
104   fHTpcEtaPhiLocaln(0x0),
105   fHTrdEtaPhiLocaln(0x0),
106   fHTpcEtaPhiLocalp(0x0),
107   fHTrdEtaPhiLocalp(0x0),
108   fHTpcEtaPhiSagitan(0x0),
109   fHTrdEtaPhiSagitan(0x0),
110   fHTpcEtaPhiSagitap(0x0),
111   fHTrdEtaPhiSagitap(0x0),
112   fHTpcEtaPhiDeltan(0x0),
113   fHTrdEtaPhiDeltan(0x0),
114   fHTpcEtaPhiDeltap(0x0),
115   fHTrdEtaPhiDeltap(0x0),
116   fHTpcEtaPhin(0x0),
117   fHTrdEtaPhin(0x0),
118   fHTpcEtaPhip(0x0),
119   fHTrdEtaPhip(0x0),
120   fHTpcRef3Dpos(0x0),
121   fHTpcRef3Dneg(0x0),
122   fHTrdRef3Dpos(0x0),
123   fHTrdRef3Dneg(0x0),
124   fHTrdEtaPhiLocalNTracklets(0x0),
125   fHTrdEtaPhiSagitaNTracklets(0x0),
126   fHTrdEtaPhiDeltaNTracklets(0x0),
127   fHTrdEtaPhiNTracklets(0x0)
128 {
129   //
130   // constructor
131   //
132   DefineInput(0, TChain::Class());
133   DefineOutput(0, TObjArray::Class());
134 }
135
136 AliAnalysisTaskTRDmon::AliAnalysisTaskTRDmon(const AliAnalysisTaskTRDmon &task) :
137   AliAnalysisTask(task),
138   fESD(0x0),
139   fESDfriend(0x0),
140   fEventTriggerName(task.GetTriggerName()),
141   fIsCollisionEvent(task.GetIsCollisionEvent()),
142   fOutStorage(0x0),
143   fHzvert1(0x0),
144   fHzvert2(0x0),
145   fHntracks(0x0),
146   fHntracks2(0x0),
147   fHntracks3(0x0),
148   fHdca(0x0),
149   fHdcaz(0x0),
150   fHpt(0x0),
151   fHpt2(0x0),
152   fHpt3(0x0),
153   fHpt3n(0x0),
154   fHpt4(0x0),
155   fHpt4n(0x0),
156   fHtheta(0x0),
157   fHphi(0x0),
158   fHtpccl(0x0),
159   fHtpccl2(0x0),
160   fHdedxp(0x0),
161   fHetaphi(0x0),
162   fHetancl(0x0),
163   fHphincl(0x0),
164   fHtrdtr(0x0),
165   fHtrdtr2(0x0),
166   fHph(0x0),
167   fHph2d(0x0),
168   fHncltrkl(0x0),
169   fHntrkl(0x0),
170   fHntrklVsP(0x0),
171   fHsm(0x0),
172   fHetantr(0x0),
173   fHphintr(0x0),
174   fHcltime(0x0),
175   fHcldiff(0x0),
176   fHxyA(0x0),
177   fHxyC(0x0),
178   fHPropagXY(0x0),
179   fHPropagRZ(0x0),
180   fHnFriendTracks(0x0),
181   fHnCalibObjects(0x0),
182   fHTpcEtaPhiLocaln(0x0),
183   fHTrdEtaPhiLocaln(0x0),
184   fHTpcEtaPhiLocalp(0x0),
185   fHTrdEtaPhiLocalp(0x0),
186   fHTpcEtaPhiSagitan(0x0),
187   fHTrdEtaPhiSagitan(0x0),
188   fHTpcEtaPhiSagitap(0x0),
189   fHTrdEtaPhiSagitap(0x0),
190   fHTpcEtaPhiDeltan(0x0),
191   fHTrdEtaPhiDeltan(0x0),
192   fHTpcEtaPhiDeltap(0x0),
193   fHTrdEtaPhiDeltap(0x0),
194   fHTpcEtaPhin(0x0),
195   fHTrdEtaPhin(0x0),
196   fHTpcEtaPhip(0x0),
197   fHTrdEtaPhip(0x0),
198   fHTpcRef3Dpos(0x0),
199   fHTpcRef3Dneg(0x0),
200   fHTrdRef3Dpos(0x0),
201   fHTrdRef3Dneg(0x0),
202   fHTrdEtaPhiLocalNTracklets(0x0),
203   fHTrdEtaPhiSagitaNTracklets(0x0),
204   fHTrdEtaPhiDeltaNTracklets(0x0),
205   fHTrdEtaPhiNTracklets(0x0)
206 {
207   //
208   // copy constructor
209   //
210 }
211
212 AliAnalysisTaskTRDmon& AliAnalysisTaskTRDmon::operator=(const AliAnalysisTaskTRDmon &task) 
213 {
214   //
215   // Assignment operator
216   //
217   if(&task == this) return *this;
218   AliAnalysisTask::operator=(task);
219   fEventTriggerName = task.GetTriggerName();
220   fIsCollisionEvent = task.GetIsCollisionEvent();
221
222   return *this;
223 }
224
225 void AliAnalysisTaskTRDmon::ConnectInputData(Option_t *){
226   //
227   // Connect input data. Overloads AliAnalysisTask::ConnectInputData()
228   //
229   TTree* tree = dynamic_cast<TTree*> (GetInputData(0)); 
230   if (!tree) Printf("ERROR: Could not read chain from input slot 0"); 
231   else tree->SetBranchStatus("ESDfriend*", kTRUE); 
232   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
233   if(!esdH){
234     AliError("ERROR: Could not get ESDInputHandler");
235   } else {
236     fESD = esdH->GetEvent();
237     fESDfriend = dynamic_cast<AliESDfriend*> (fESD->FindListObject("AliESDfriend"));
238     
239   }
240   
241 }
242
243 void AliAnalysisTaskTRDmon::CreateOutputObjects(){
244   //
245   // Creates the output object list, initializes histograms, adds histograms to the output list
246   //
247
248   fOutStorage = new TObjArray(40);
249   fHzvert1 = new TH1F("fHzvert1", "Zvertex", 100, -25., 25.);
250   fHzvert2 = new TH1F("fHzvert2", "Zvertex TPC", 100, -25., 25.);
251   fHdca = new TH1F("fHdca", "DCA", 100, -100., 100.);
252   fHdcaz = new TH1F("fHdcaz", "DCA", 100, -100., 100.);
253   fHntracks = new TH1F("fHntracks", "Ntracks", 100, 0., 100.);
254   fHntracks2 = new TH1F("fHntracks2", "Ntracks dca.lt.1", 100, 0., 100.);
255   fHntracks3 = new TH1F("fHntracks3", "Ntracks dca.lt.1, NclTPC", 100, 0., 100.);
256   Float_t binLimits[33] = {0.0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9,
257                            1.0, 1.1, 1.2, 1.3, 1.4, 
258                            1.6, 1.8, 2.0, 2.2, 2.4, 2.6, 2.8, 3.0,
259                            3.4, 3.8, 4.2, 4.6, 5.0, 6.0, 7.0, 8.0, 9.0, 10.0};
260   fHpt = new TH1F("fHpt", "Pt", 32, binLimits);
261   fHpt2 = new TH1F("fHpt2", "Pt dca.lt.1", 32, binLimits);
262   fHpt2->Sumw2();
263   fHpt3 = new TH1F("fHpt3", "Pt dca.lt.1, NclTPC +", 32, binLimits);
264   fHpt3->Sumw2();
265   fHpt3n = new TH1F("fHpt3n", "Pt dca.lt.1, NclTPC -", 32, binLimits);
266   fHpt3n->Sumw2();
267   fHpt4 = new TH1F("fHpt4", "Pt dca.lt.1, NclTPC & TRD +", 32, binLimits);
268   fHpt4->Sumw2();
269   fHpt4n = new TH1F("fHpt4n", "Pt dca.lt.1, NclTPC & TRD -", 32, binLimits);
270   fHpt4n->Sumw2();
271   fHtheta = new TH1F("fHtheta", "theta", 220,.5,2.7);
272   fHphi = new TH1F("fHphi", "phi", 157,0,6.28);
273   fHtpccl = new TH1F("fHtpccl", "nClTPC", 160, 0., 160.);
274   fHtpccl2 = new TH1F("fHtpccl2", "nClTPC p.gt.1", 160, 0., 160.);
275   
276   //Float_t binx[21]={.1,.13,.16,.22,.3,.4,.55,.7,1,1.3,1.6,2.2,3,4,5.5,7,10,13,16,22,30};
277   
278   fHdedxp = new TH2F("fHdedxp", "dE/dx vs. p TPC", 100, 0.1,10.1, 120, 0,600.);
279   // fHdedxp->SetBit(TH1::kCanRebin); 
280   fHetaphi = new TH2F("fHetaphi", "eta-phi TPC", 50, -1, 1, 157, 0, 6.28);
281   fHetancl = new TH2F("fHetancl", "eta-Ncl TPC",50, -1, 1, 160, 0, 160.);
282   fHphincl = new TH2F("fHphincl", "phi-Ncl TPC",157, 0, 6.28, 160, 0, 160.);
283   fHtrdtr = new TH1F("fHtrdtr", "nTRDlayers", 7, -.5, 6.5);
284   fHtrdtr2 = new TH1F("fHtrdtr2", "nTRDlayers p.gt.1", 7, -.5, 6.5);
285   fHph = new TProfile("fHph", "ph", 30,0,30);
286   fHph2d = new TH2F("fHph2d", "PH2d TRD", 30, 0, 30, 200, 0, 1000.);
287   fHncltrkl = new TH1F("fHncltrkl", "Nclusters/tracklet", 50,0,50.);
288   fHntrkl = new TH1F("fHntrkl", "Ntracklets", 7,-.5,6.5);
289   Double_t binsP[19] = {0.0, 0.2, 0.4, 0.6, 0.8, 1.0, 1.2, 1.4, 1.7, 2.0,
290                         2.5, 3.0, 3.5, 4.0, 5.0, 6.0, 7.0, 9.0, 12.0};
291   fHntrklVsP = new TH2F("fHntrklVsP", "Ntracklets vs momentum", 18, binsP, 7, -0.5, 6.5);
292   fHsm = new TH1F("fHsm", "SM dep", 18,-.5,17.5);
293   fHetantr = new TH2F("fHetantr", "PH2d TRD", 50, -1, 1, 7, -.5, 6.5);
294   fHphintr = new TH2F("fHphintr", "PH2d TRD", 157, 0,6.28, 7, -.5, 6.5);
295   fHcltime = new TH1F("fHcltime", "Cluster-time", 30,0,30);
296   fHcldiff = new TH1F("fHcldiff", "Cluster time diff", 30,0,30);
297   fHxyA = new TH2F("fHxyA","x-y side A",800,-400,400,800,-400,400);
298   fHxyC = new TH2F("fHxyC","x-y side C",800,-400,400,800,-400,400);
299   fHPropagXY = new TH2F("fHPropagXY","x-y propagated point to TRD",1000,-500,500,1000,-500,500);
300   fHPropagRZ = new TH2F("fHPropagRZ","r-z propagated point to TRD",700, -350., 350.,500,0.,500.);
301   fHnFriendTracks = new TH1F("fHnFriendTracks", "NFriendTracks", 100, 0., 100.);
302   fHnCalibObjects = new TH1F("fHnCalibObjects", "NCalibObjects", 100, 0., 100.);
303   fHTpcEtaPhiLocaln = new TH2F("fHTpcEtaPhiLocaln", "(#eta - #varphi_{local}) distribution of neg TPC tracks", 
304                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
305   fHTrdEtaPhiLocaln = new TH2F("fHTrdEtaPhiLocaln", "(#eta - #varphi_{local}) distribution of neg TRD tracks", 
306                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
307   fHTpcEtaPhiLocalp = new TH2F("fHTpcEtaPhiLocalp", "(#eta - #varphi_{local}) distribution of pos TPC tracks", 
308                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
309   fHTrdEtaPhiLocalp = new TH2F("fHTrdEtaPhiLocalp", "(#eta - #varphi_{local}) distribution of pos TRD tracks", 
310                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
311   fHTpcEtaPhiSagitan = new TH2F("fHTpcEtaPhiSagitan", "(#eta - #varphi_{det}) distribution of neg TPC tracks", 
312                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
313   fHTrdEtaPhiSagitan = new TH2F("fHTrdEtaPhiSagitan", "(#eta - #varphi_{det}) distribution of neg TRD tracks", 
314                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
315   fHTpcEtaPhiSagitap = new TH2F("fHTpcEtaPhiSagitap", "(#eta - #varphi_{det}) distribution of pos TPC tracks", 
316                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
317   fHTrdEtaPhiSagitap = new TH2F("fHTrdEtaPhiSagitap", "(#eta - #varphi_{det}) distribution of pos TRD tracks", 
318                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
319   fHTpcEtaPhiDeltan = new TH2F("fHTpcEtaPhiDeltan", "(#eta - #Delta#varphi) distribution of neg TPC tracks", 
320                           50, -1.0, 1.0, 50, -0.4*TMath::Pi(), 0.4*TMath::Pi());
321   fHTrdEtaPhiDeltan = new TH2F("fHTrdEtaPhiDeltan", "(#eta - #Delta#varphi) distribution of neg TRD tracks", 
322                           50, -1.0, 1.0, 50, -0.4*TMath::Pi(), 0.4*TMath::Pi());
323   fHTpcEtaPhiDeltap = new TH2F("fHTpcEtaPhiDeltap", "(#eta - #Delta#varphi) distribution of pos TPC tracks", 
324                           50, -1.0, 1.0, 50, -0.4*TMath::Pi(), 0.4*TMath::Pi());
325   fHTrdEtaPhiDeltap = new TH2F("fHTrdEtaPhiDeltap", "(#eta - #Delta#varphi) distribution of pos TRD tracks", 
326                           50, -1.0, 1.0, 50, -0.4*TMath::Pi(), 0.4*TMath::Pi());
327   fHTpcEtaPhin = new TH2F("fHTpcEtaPhin", "(#eta - #varphi) distribution of neg TPC tracks", 
328                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
329   fHTrdEtaPhin = new TH2F("fHTrdEtaPhin", "(#eta - #varphi) distribution of neg TRD tracks", 
330                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
331   fHTpcEtaPhip = new TH2F("fHTpcEtaPhip", "(#eta - #varphi) distribution of pos TPC tracks", 
332                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
333   fHTrdEtaPhip = new TH2F("fHTrdEtaPhip", "(#eta - #varphi) distribution of pos TRD tracks", 
334                           50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
335   Float_t etaBinLimits[51];
336   for(Int_t i=0; i<51; i++) etaBinLimits[i] = -1.0 + i*2.0/50.;
337   Float_t phiBinLimits[151];
338   for(Int_t i=0; i<151; i++) phiBinLimits[i] = -1.1*TMath::Pi() + i*2.2*TMath::Pi()/150.;
339   fHTpcRef3Dpos = new TH3F("fHTpcRef3Dpos", "TPC positive ref. tracks", 
340                            50, etaBinLimits, 150, phiBinLimits, 32, binLimits);
341   fHTpcRef3Dneg = new TH3F("fHTpcRef3Dneg", "TPC negative ref. tracks",
342                            50, etaBinLimits, 150, phiBinLimits, 32, binLimits);
343   fHTrdRef3Dpos = new TH3F("fHTrdRef3Dpos", "TRD positive tracks",
344                            50, etaBinLimits, 150, phiBinLimits, 32, binLimits);
345   fHTrdRef3Dneg = new TH3F("fHTrdRef3Dneg", "TRD negative tracks",
346                            50, etaBinLimits, 150, phiBinLimits, 32, binLimits);
347   fHTrdEtaPhiLocalNTracklets = new TProfile2D("fHTrdEtaPhiLocalNTracklets", "(#eta - #varphi_{local}) distribution of <ntracklets>", 
348                                               50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
349   fHTrdEtaPhiSagitaNTracklets = new TProfile2D("fHTrdEtaPhiSagitaNTracklets", "(#eta - #varphi_{det}) distribution of <ntracklets>", 
350                                               50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
351   fHTrdEtaPhiDeltaNTracklets = new TProfile2D("fHTrdEtaPhiDeltaNTracklets", "(#eta - #Delta#varphi) distribution of <ntracklets>", 
352                                               50, -1.0, 1.0, 50, -0.4*TMath::Pi(), 0.4*TMath::Pi());
353   fHTrdEtaPhiNTracklets = new TProfile2D("fHTrdEtaPhiNTracklets", "(#eta - #varphi) distribution of <ntracklets>", 
354                                          50, -1.0, 1.0, 150, -1.1*TMath::Pi(), 1.1*TMath::Pi());
355   
356   fOutStorage->Add(fHzvert1);
357   fOutStorage->Add(fHzvert2);
358   fOutStorage->Add(fHntracks);
359   fOutStorage->Add(fHntracks2);
360   fOutStorage->Add(fHntracks3);
361   fOutStorage->Add(fHdca);
362   fOutStorage->Add(fHdcaz);
363   fOutStorage->Add(fHpt);
364   fOutStorage->Add(fHpt2);
365   fOutStorage->Add(fHpt3);
366   fOutStorage->Add(fHpt3n);
367   fOutStorage->Add(fHpt4);
368   fOutStorage->Add(fHpt4n);
369   fOutStorage->Add(fHtheta);
370   fOutStorage->Add(fHphi);
371   fOutStorage->Add(fHtpccl);
372   fOutStorage->Add(fHtpccl2);
373   fOutStorage->Add(fHdedxp);
374   fOutStorage->Add(fHetaphi);
375   fOutStorage->Add(fHetancl);
376   fOutStorage->Add(fHphincl);
377   fOutStorage->Add(fHtrdtr);
378   fOutStorage->Add(fHtrdtr2);
379   fOutStorage->Add(fHph);
380   fOutStorage->Add(fHph2d);
381   fOutStorage->Add(fHncltrkl);
382   fOutStorage->Add(fHntrkl);
383   fOutStorage->Add(fHntrklVsP);
384   fOutStorage->Add(fHetantr);
385   fOutStorage->Add(fHphintr);
386   fOutStorage->Add(fHcltime);
387   fOutStorage->Add(fHcldiff);
388   fOutStorage->Add(fHxyA);
389   fOutStorage->Add(fHxyC);
390   fOutStorage->Add(fHPropagXY);
391   fOutStorage->Add(fHPropagRZ);
392   fOutStorage->Add(fHnFriendTracks);
393   fOutStorage->Add(fHnCalibObjects);
394   fOutStorage->Add(fHTpcEtaPhiLocaln);
395   fOutStorage->Add(fHTpcEtaPhiLocalp);
396   fOutStorage->Add(fHTrdEtaPhiLocaln);
397   fOutStorage->Add(fHTrdEtaPhiLocalp);
398   fOutStorage->Add(fHTpcEtaPhiSagitan);
399   fOutStorage->Add(fHTpcEtaPhiSagitap);
400   fOutStorage->Add(fHTrdEtaPhiSagitan);
401   fOutStorage->Add(fHTrdEtaPhiSagitap);
402   fOutStorage->Add(fHTpcEtaPhiDeltan);
403   fOutStorage->Add(fHTpcEtaPhiDeltap);
404   fOutStorage->Add(fHTrdEtaPhiDeltan);
405   fOutStorage->Add(fHTrdEtaPhiDeltap);
406   fOutStorage->Add(fHTpcEtaPhin);
407   fOutStorage->Add(fHTpcEtaPhip);
408   fOutStorage->Add(fHTrdEtaPhin);
409   fOutStorage->Add(fHTrdEtaPhip);
410   fOutStorage->Add(fHTpcRef3Dpos);
411   fOutStorage->Add(fHTpcRef3Dneg);
412   fOutStorage->Add(fHTrdRef3Dpos);
413   fOutStorage->Add(fHTrdRef3Dneg);
414   fOutStorage->Add(fHTrdEtaPhiLocalNTracklets);
415   fOutStorage->Add(fHTrdEtaPhiSagitaNTracklets);
416   fOutStorage->Add(fHTrdEtaPhiDeltaNTracklets);
417   fOutStorage->Add(fHTrdEtaPhiNTracklets);
418 }
419
420 void AliAnalysisTaskTRDmon::Exec(Option_t *){
421   //
422   // Implementation of the Exec() function
423   //
424   if(!fESD){
425           AliError("ESD Event missing");
426           return;
427   }
428   //  const Float_t kCutDCAxy = 40.;      // cm
429   const Float_t kCutDCAxy = 3; //40.;      // cm
430   const Float_t kCutDCAz = 15.;      // cm
431   const Int_t kNclTPCmin = 100;      // Nclusters TPC
432   const Float_t kPmin = 0.2;      // min. pt
433   
434   // event trigger cut
435   if(fEventTriggerName.Data()[0]!='\0' && !fESD->IsTriggerClassFired(fEventTriggerName.Data())) {
436     PostData(0, fOutStorage);
437     return; 
438   }
439
440   fESDfriend = dynamic_cast<AliESDfriend*> (fESD->FindListObject("AliESDfriend"));
441   AliESDtrack *fTrack = 0x0;
442   AliESDfriendTrack *fFriendTrack = 0x0;
443   AliTRDtrackV1 *fTRDtrack=0x0; //TRD "calibration" object
444
445   //const AliESDVertex *EvVertex = fESD->GetPrimaryVertexTPC();
446   const AliESDVertex *evVertex = fESD->GetPrimaryVertex();
447   Float_t zvert1 = evVertex->GetZv();
448   Float_t nvtxContr = evVertex->GetNContributors();
449
450   // if the required trigger is a collision trigger then apply event vertex cut
451   if(fIsCollisionEvent && (TMath::Abs(zvert1)>15. || TMath::Abs(zvert1)<1.0e-10 || nvtxContr<=1)) {
452     PostData(0, fOutStorage);
453     return;
454   }
455   
456   Float_t zvert2 = fESD->GetPrimaryVertex()->GetZv();
457   
458   fHzvert1->Fill(zvert1);
459   fHzvert2->Fill(zvert2);
460
461  
462   Int_t ntracks = fESD->GetNumberOfTracks();
463   fHntracks->Fill((Float_t)ntracks);
464   
465   Int_t ntracksDca=0;
466   Int_t ntracksNcl=0;
467   Int_t nFriendTracks = 0;
468   Int_t nCalibObjects = 0;
469   
470   for(Int_t itrk = ntracks; itrk--;){ 
471     fTrack = fESD->GetTrack(itrk);
472     if(!fTrack) continue;
473     //printf("%d, Label %d\n",itrk,fTrack->GetLabel()); // ?
474     if(!(fTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
475     if(fTrack->GetKinkIndex(0) > 0) continue; //rejects kinks
476     
477     Double_t pt = fTrack->Pt();
478     Double_t p = fTrack->GetP();
479     
480     if ( pt < kPmin ) continue;
481     
482     Float_t theta=fTrack->Theta();
483     Float_t phi=fTrack->Phi();
484     Float_t eta=-TMath::Log(TMath::Tan(theta/2.));
485     
486     if ( TMath::Abs(eta)>0.9 || fTrack->GetTPCNcls()<10) continue;
487     
488     Int_t nTPCclus = fTrack->GetTPCNcls();
489     fHetancl->Fill(eta,(Float_t)nTPCclus);          
490     fHphincl->Fill(phi,(Float_t)nTPCclus);          
491     fHpt->Fill(pt);
492     
493     Float_t imPar[2], cov[3];
494     fTrack->GetImpactParameters(imPar, cov);
495     //Float_t dca = TMath::Sqrt(imPar[0]*imPar[0] + imPar[1]*imPar[1]);
496     Float_t dca =imPar[0]; //transverse
497     Float_t dcaz =imPar[1]; //z
498     fHdca->Fill(dca);
499     fHdcaz->Fill(dcaz);
500     
501     if (fIsCollisionEvent && TMath::Abs(dca) > kCutDCAxy ) continue;    
502     if (fIsCollisionEvent && TMath::Abs(dcaz) > kCutDCAz ) continue;    
503     
504     fHpt2->Fill(p);
505     ntracksDca++;
506     
507     Double_t dedxTPC=fTrack->GetTPCsignal();
508     
509     fHtpccl->Fill((Float_t)nTPCclus); //?
510     if ( pt > 1 ) fHtpccl2->Fill((Float_t)nTPCclus);
511     
512     // track points for x-y histograms
513     fFriendTrack = fESDfriend->GetTrack(itrk); //this works
514     Bool_t hasFriends = (fFriendTrack ? kTRUE : kFALSE);
515     if(hasFriends) nFriendTracks++;
516     const AliTrackPointArray *array = (hasFriends ? fFriendTrack->GetTrackPointArray() : 0);
517     if (array) {
518       Int_t npoints = array->GetNPoints();
519       for (int ipoint=0; ipoint<npoints; ipoint++) {
520         AliTrackPoint tp;
521         array->GetPoint(tp, ipoint);
522         if (tp.GetZ()<0) fHxyC->Fill(tp.GetX(),tp.GetY());
523         if (tp.GetZ()>0) fHxyA->Fill(tp.GetX(),tp.GetY());
524       }
525     }
526     
527     if ( nTPCclus < kNclTPCmin ) continue; 
528     
529     //    Bool_t propagGood = fTrack->PropagateTo(298., fESD->GetMagneticField());
530     Double_t localCoord[3] = {0., 0., 0.};
531     Double_t localMom[3]   = {0., 0., 0.};
532     Bool_t localCoordGood = fTrack->GetXYZAt(298., fESD->GetMagneticField(), localCoord);
533     if(localCoordGood) {
534       fHPropagXY->Fill(localCoord[0], localCoord[1]);
535       fHPropagRZ->Fill(localCoord[2], TMath::Sqrt(localCoord[0]*localCoord[0]+localCoord[1]*localCoord[1]));
536     }
537     Bool_t localMomGood   = fTrack->GetPxPyPzAt(298., fESD->GetMagneticField(), localMom);
538     Double_t localPhi = (localMomGood ? TMath::ATan2(localMom[1], localMom[0]) : 0.0);
539     Double_t localSagitaPhi = (localCoordGood ? TMath::ATan2(localCoord[1], localCoord[0]) : 0.0);    
540     //    if(localPhi<0.) localPhi = TMath::Abs(localPhi)+TMath::Pi();   
541     //    if(localSagitaPhi<0.) localSagitaPhi = TMath::Abs(localSagitaPhi)+TMath::Pi();     
542
543     if (fTrack->Charge()>0) {
544       fHpt3->Fill(p);
545       if(localCoordGood && localMomGood) {
546         fHTpcEtaPhiLocalp->Fill(eta,localPhi);
547         fHTpcEtaPhiSagitap->Fill(eta,localSagitaPhi);
548         fHTpcEtaPhiDeltap->Fill(eta,localPhi-localSagitaPhi);
549         fHTpcRef3Dpos->Fill(eta, localSagitaPhi, p);
550       }
551       fHTpcEtaPhip->Fill(eta,phi-TMath::Pi());
552     }
553     if (fTrack->Charge()<0) {
554       fHpt3n->Fill(p);
555       if(localCoordGood && localMomGood) {
556         fHTpcEtaPhiLocaln->Fill(eta,localPhi);
557         fHTpcEtaPhiSagitan->Fill(eta,localSagitaPhi);
558         fHTpcEtaPhiDeltan->Fill(eta,localPhi-localSagitaPhi);
559         fHTpcRef3Dneg->Fill(eta, localSagitaPhi, p);
560       }
561       fHTpcEtaPhin->Fill(eta,phi-TMath::Pi());
562     }
563     ntracksNcl++;
564     fHetaphi->Fill(eta,phi);
565     fHdedxp->Fill(p,dedxTPC);
566     Int_t nTRDtrkl = fTrack->GetTRDntracklets();  
567     
568     if ( nTRDtrkl < 1 ) continue;
569     if (fTrack->Charge()>0) {
570       fHpt4->Fill(p);
571       if(localCoordGood && localMomGood) {
572         fHTrdEtaPhiLocalp->Fill(eta,localPhi);
573         fHTrdEtaPhiSagitap->Fill(eta,localSagitaPhi);
574         fHTrdEtaPhiDeltap->Fill(eta,localPhi-localSagitaPhi);
575         fHTrdRef3Dpos->Fill(eta, localSagitaPhi, p);
576       }
577       fHTrdEtaPhip->Fill(eta,phi-TMath::Pi());
578     }
579     if (fTrack->Charge()<0) {
580       fHpt4n->Fill(p);
581       if(localCoordGood && localMomGood) {
582         fHTrdEtaPhiLocaln->Fill(eta,localPhi);
583         fHTrdEtaPhiSagitan->Fill(eta,localSagitaPhi);
584         fHTrdEtaPhiDeltan->Fill(eta,localPhi-localSagitaPhi);
585         fHTrdRef3Dneg->Fill(eta, localSagitaPhi, p);
586       }
587       fHTrdEtaPhin->Fill(eta,phi-TMath::Pi());
588     }
589     fHetantr->Fill(eta,(Float_t)nTRDtrkl);          
590     fHphintr->Fill(phi,(Float_t)nTRDtrkl);          
591     if(localCoordGood && localMomGood) {
592       fHTrdEtaPhiLocalNTracklets->Fill(eta, localPhi, (Float_t)nTRDtrkl);
593       fHTrdEtaPhiSagitaNTracklets->Fill(eta, localSagitaPhi, (Float_t)nTRDtrkl);
594       fHTrdEtaPhiDeltaNTracklets->Fill(eta, localPhi-localSagitaPhi, (Float_t)nTRDtrkl);
595     }
596     fHTrdEtaPhiNTracklets->Fill(eta, phi-TMath::Pi(), (Float_t)nTRDtrkl);
597     fHtrdtr->Fill((Float_t)nTRDtrkl);
598     if ( pt > 1 ) fHtrdtr2->Fill((Float_t)nTRDtrkl);
599     fHtheta->Fill(theta);
600     fHphi->Fill(phi);
601     
602     // start looking into TRD ...
603     TObject *o = 0x0;
604     Int_t ical = 0; 
605     fTRDtrack = 0x0;
606     while(hasFriends && (o = fFriendTrack->GetCalibObject(ical++))){
607       //while((o = (const_cast<AliESDfriendTrack*>(fFriendTrack))->GetCalibObject(ical++))) { //...this is fishy...
608       //printf("Object type: %s\n", o->IsA()->GetName());
609       if(strcmp(o->IsA()->GetName(),"AliTRDtrackV1") != 0){
610         continue;
611       }
612       fTRDtrack = dynamic_cast<AliTRDtrackV1 *>(o);
613       break;
614     }
615     
616     Bool_t hasTRDfriends = (fTRDtrack ? kTRUE : kFALSE);
617     //    if(!fTRDtrack) continue;
618     if(hasTRDfriends) nCalibObjects++;
619     Int_t nTracklets = 0;
620     Int_t tbinPrev=0;
621     
622     //    Int_t ntls = fTRDtrack->GetNumberOfTracklets();
623     Int_t ntls = fTrack->GetTRDntracklets();
624     //Int_t ncls = track->GetNumberOfClusters(); //per track
625     if(hasTRDfriends) {
626       for(Int_t itl = 0; itl < 6; itl++){ // TRD layers
627         AliTRDseedV1 * tracklet = fTRDtrack->GetTracklet(itl);
628         if(!tracklet || !tracklet->IsOK()) continue;
629         if(!tracklet) continue;
630         nTracklets++;
631         Int_t nclsTracklet = 0;
632         for(Int_t itb = 0; itb < AliTRDseedV1::kNtb; itb++){ // timebins
633           AliTRDcluster *c = tracklet->GetClusters(itb);
634           AliTRDcluster *c1 = tracklet->GetClusters(itb + AliTRDseedV1::kNtb);  // second pad row
635           //AliTRDcluster *c = tracklet->GetdQdl(itb);
636           if(!(c || c1)) continue;
637           AliTRDcluster *cptr = c ? c : c1;
638           Int_t tcal = cptr->GetLocalTimeBin();
639           fHcldiff->Fill((Float_t)(tcal-tbinPrev)); 
640           tbinPrev=tcal;
641           Int_t idet = cptr->GetDetector();
642           //If (itrk<100 && itb==0) cout << "...Detector  "<<idet<< endl;
643           Int_t sm = idet/30;
644           fHsm->Fill((Float_t)sm); 
645           idet -= sm * 30;
646           float sig = 0;
647           if(c) sig += TMath::Abs(c->GetQ()); //
648           if(c1) sig += TMath::Abs(c1->GetQ());
649           
650           fHph->Fill((Float_t)tcal, sig);  
651           fHph2d->Fill((Float_t)tcal, sig);  
652           fHcltime->Fill((Float_t)tcal); 
653           nclsTracklet++;
654         }
655         fHncltrkl->Fill(nclsTracklet);
656       } //tracklets
657     } // end if(hasTRDfriends)
658     //fHntrkl->Fill((Float_t)nTracklets);
659     fHntrkl->Fill((Float_t)ntls);
660     fHntrklVsP->Fill(p,(Float_t)ntls);
661   } //tracks
662   
663   fHntracks2->Fill((Float_t)ntracksDca);
664   fHntracks3->Fill((Float_t)ntracksNcl);
665   fHnFriendTracks->Fill(nFriendTracks);
666   fHnCalibObjects->Fill(nCalibObjects);
667   
668   PostData(0, fOutStorage);
669 }
670 void AliAnalysisTaskTRDmon::Terminate(Option_t *)
671 {    
672   //
673   // Do nothing
674   //
675 }
676