]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/ITS/AliAnalysisTaskSDDRP.cxx
MFT track shit tool added
[u/mrichter/AliRoot.git] / PWGPP / ITS / AliAnalysisTaskSDDRP.cxx
1 #include "AliAnalysisTask.h"
2 #include "AliAnalysisManager.h"
3 #include "AliAnalysisDataContainer.h"
4 #include "AliITSRecPoint.h"
5 #include "AliESDEvent.h"
6 #include "AliTrackPointArray.h"
7 #include "AliITSgeomTGeo.h"
8 #include "AliESDfriend.h"
9 #include "AliCDBManager.h"
10 #include "AliCDBEntry.h"
11 #include "AliITSCalibrationSDD.h"
12 #include "AliITSresponseSDD.h"
13 #include "AliTriggerConfiguration.h"
14 #include "AliGeomManager.h"
15 #include <TSystem.h>
16 #include <TTree.h>
17 #include <TH1F.h>
18 #include <TH2F.h>
19 #include <TChain.h>
20 #include <TGeoGlobalMagField.h>
21 #include "AliESDInputHandlerRP.h"
22 /**************************************************************************
23  * Copyright(c) 1998-2010, ALICE Experiment at CERN, All rights reserved. *
24  *                                                                        *
25  * Author: The ALICE Off-line Project.                                    *
26  * Contributors are mentioned in the code where appropriate.              *
27  *                                                                        *
28  * Permission to use, copy, modify and distribute this software and its   *
29  * documentation strictly for non-commercial purposes is hereby granted   *
30  * without fee, provided that the above copyright notice appears in all   *
31  * copies and that both the copyright notice and this permission notice   *
32  * appear in the supporting documentation. The authors make no claims     *
33  * about the suitability of this software for any purpose. It is          *
34  * provided "as is" without express or implied warranty.                  *
35  **************************************************************************/
36
37 //*************************************************************************
38 // Implementation of class AliAnalysiTaskSDDRP
39 // AliAnalysisTaskSE to extract from ESD + ESDfreinds + ITS rec points
40 // performance plots for SDD detector
41 //
42 // Author: F. Prino, prino@to.infn.it
43 //*************************************************************************
44
45
46 #include "AliAnalysisTaskSDDRP.h"
47
48 ClassImp(AliAnalysisTaskSDDRP)
49 //______________________________________________________________________________
50 AliAnalysisTaskSDDRP::AliAnalysisTaskSDDRP() : AliAnalysisTaskSE("SDD RecPoints"), 
51   fOutput(0),
52   fHistNEvents(0),
53   fHistCluInLay(0),
54   fHistAllPMod(0),
55   fHistGoodPMod(0),
56   fHistBadRegMod(0),
57   fHistMissPMod(0),
58   fHistSkippedMod(0),
59   fHistOutAccMod(0),
60   fHistNoRefitMod(0),
61   fHistAllPXloc(0),
62   fHistGoodPXloc(0),
63   fHistBadRegXloc(0),
64   fHistMissPXloc(0),
65   fHistAllPZloc(0),
66   fHistGoodPZloc(0),
67   fHistBadRegZloc(0),
68   fHistMissPZloc(0),
69   fHistdEdxL3VsP(0),
70   fHistdEdxL4VsP(0),
71   fHistdEdxVsMod(0),
72   fRecPMod(0),
73   fTrackPMod(0),
74   fGoodAnMod(0),
75   fRecPLadLay3(0),
76   fRecPLadLay4(0),
77   fTrackPLadLay3(0),
78   fTrackPLadLay4(0),
79   fGoodAnLadLay3(0),
80   fGoodAnLadLay4(0),
81   fDriftTimeRP(0),
82   fDriftTimeTPAll(0),
83   fDriftTimeTPNoExtra(0),
84   fDriftTimeTPExtra(0),
85   fCluSizAnVsTime(0),
86   fCluSizTbVsTime(0),
87   fResp(0),
88   fTrigConfig(0),
89   fUseITSsaTracks(kFALSE),
90   fMinITSpts(3),
91   fMinTPCpts(70),
92   fMinPfordEdx(0.5),
93   fTriggerClass(""),
94   fOnlyEventsWithSDD(kTRUE),
95   fExcludeBadMod(kFALSE)
96 {
97   //
98   DefineOutput(1, TList::Class());
99 }
100
101
102 //___________________________________________________________________________
103 AliAnalysisTaskSDDRP::~AliAnalysisTaskSDDRP(){
104   //
105   if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
106     delete fOutput;
107     fOutput = 0;
108   }
109   
110 }
111
112
113 //___________________________________________________________________________
114
115 void AliAnalysisTaskSDDRP::UserCreateOutputObjects() {
116
117   fOutput = new TList();
118   fOutput->SetOwner();
119   fOutput->SetName("OutputHistos");
120
121   fHistNEvents = new TH1F("hNEvents", "Number of processed events",5,-1.5,3.5);
122   fHistNEvents->GetXaxis()->SetBinLabel(1,"All events");
123   fHistNEvents->GetXaxis()->SetBinLabel(2,"Selected triggers"); 
124   fHistNEvents->GetXaxis()->SetBinLabel(3,"Without SDD"); 
125   fHistNEvents->GetXaxis()->SetBinLabel(4,"With SDD"); 
126   fHistNEvents->GetXaxis()->SetBinLabel(5,"Analyzed events"); 
127   fHistNEvents->Sumw2();
128   fHistNEvents->SetMinimum(0);
129   fOutput->Add(fHistNEvents);
130
131   fHistCluInLay = new TH1F("hCluInLay","hCluInLay",7,-1.5,5.5);
132   fHistCluInLay->Sumw2();
133   fHistCluInLay->SetMinimum(0);
134   fOutput->Add(fHistCluInLay);
135
136   // -- Module histos
137
138   fHistAllPMod  = new TH1F("hAllPmod","Crossing Tracks vs. Module",260,239.5,499.5);
139   fHistAllPMod->Sumw2();
140   fHistAllPMod->SetMinimum(0);
141   fOutput->Add(fHistAllPMod);
142
143   fHistGoodPMod  = new TH1F("hGoodPmod","PointsAssocToTrack per Module",260,239.5,499.5);
144   fHistGoodPMod->Sumw2();
145   fHistGoodPMod->SetMinimum(0);
146   fOutput->Add(fHistGoodPMod);
147
148   fHistBadRegMod  = new TH1F("hBadRegmod","Tracks in BadRegion per Module",260,239.5,499.5);
149   fHistBadRegMod->Sumw2();
150   fHistBadRegMod->SetMinimum(0);
151   fOutput->Add(fHistBadRegMod);
152
153   fHistMissPMod  = new TH1F("hMissPmod","Missing Points per Module",260,239.5,499.5);
154   fHistMissPMod->Sumw2();
155   fHistMissPMod->SetMinimum(0);
156   fOutput->Add(fHistMissPMod);
157
158   fHistSkippedMod  = new TH1F("hSkippedmod","Tracks in Skipped Module",260,239.5,499.5);
159   fHistSkippedMod->Sumw2();
160   fHistSkippedMod->SetMinimum(0);
161   fOutput->Add(fHistSkippedMod);
162
163   fHistOutAccMod  = new TH1F("hOutAccmod","Tracks outside zAcc per Module",260,239.5,499.5);
164   fHistOutAccMod->Sumw2();
165   fHistOutAccMod->SetMinimum(0);
166   fOutput->Add(fHistOutAccMod);
167
168   fHistNoRefitMod  = new TH1F("hNoRefitmod","Points rejected in refit per Module",260,239.5,499.5);
169   fHistNoRefitMod->Sumw2();
170   fHistNoRefitMod->SetMinimum(0);
171   fOutput->Add(fHistNoRefitMod);
172
173
174
175   fRecPMod = new TH1F("hRPMod","Rec Points per Module",260,239.5,499.5);
176   fRecPMod->Sumw2();
177   fRecPMod->SetMinimum(0);
178   fOutput->Add(fRecPMod);
179
180   fTrackPMod = new TH1F("hTPMod","Track Points per Module",260,239.5,499.5);
181   fTrackPMod->Sumw2();
182   fTrackPMod->SetMinimum(0);
183   fOutput->Add(fTrackPMod);
184
185   fGoodAnMod = new TH1F("hGAMod","Good Anodes per Module",260,239.5,499.5);
186   fOutput->Add(fGoodAnMod);
187
188   // -- Local coordinates
189
190   fHistAllPXloc  = new TH1F("hAllPxloc","Crossing Tracks vs. Xloc",75, -3.75, 3.75);
191   fHistAllPXloc->Sumw2();
192   fHistAllPXloc->SetMinimum(0);
193   fOutput->Add(fHistAllPXloc);
194
195   fHistGoodPXloc  = new TH1F("hGoodPxloc","PointsAssocToTrack vs. Xloc",75, -3.75, 3.75);
196   fHistGoodPXloc->Sumw2();
197   fHistGoodPXloc->SetMinimum(0);
198   fOutput->Add(fHistGoodPXloc);
199
200   fHistBadRegXloc  = new TH1F("hBadRegxloc","Tracks in BadRegion vs. Xloc",75, -3.75, 3.75);
201   fHistBadRegXloc->Sumw2();
202   fHistBadRegXloc->SetMinimum(0);
203   fOutput->Add(fHistBadRegXloc);
204
205   fHistMissPXloc  = new TH1F("hMissPxloc","Missing Points vs. Xloc",75, -3.75, 3.75);
206   fHistMissPXloc->Sumw2();
207   fHistMissPXloc->SetMinimum(0);
208   fOutput->Add(fHistMissPXloc);
209
210   fHistAllPZloc  = new TH1F("hAllPzloc","Crossing Tracks vs. Zloc",77, -3.85, 3.85);
211   fHistAllPZloc->Sumw2();
212   fHistAllPZloc->SetMinimum(0);
213   fOutput->Add(fHistAllPZloc);
214
215   fHistGoodPZloc  = new TH1F("hGoodPzloc","PointsAssocToTrack vs. Zloc",77, -3.85, 3.85);
216   fHistGoodPZloc->Sumw2();
217   fHistGoodPZloc->SetMinimum(0);
218   fOutput->Add(fHistGoodPZloc);
219
220   fHistBadRegZloc  = new TH1F("hBadRegzloc","Tracks in BadRegion vs. Zloc",77, -3.85, 3.85);
221   fHistBadRegZloc->Sumw2();
222   fHistBadRegZloc->SetMinimum(0);
223   fOutput->Add(fHistBadRegZloc);
224
225   fHistMissPZloc  = new TH1F("hMissPzloc","Missing Points vs. Zloc",77, -3.85, 3.85);
226   fHistMissPZloc->Sumw2();
227   fHistMissPZloc->SetMinimum(0);
228   fOutput->Add(fHistMissPZloc);
229
230   // -- Ladder histos
231
232   fRecPLadLay3 = new TH1F("hRPLad3","Rec Points per Ladder Layer 3",14,-0.5,13.5);
233   fRecPLadLay3->Sumw2();
234   fRecPLadLay3->SetMinimum(0);
235   fOutput->Add(fRecPLadLay3);
236
237   fRecPLadLay4 = new TH1F("hRPLad4","Rec Points per Ladder Layer 4",22,-0.5,21.5);
238   fRecPLadLay4->Sumw2();
239   fRecPLadLay4->SetMinimum(0);
240   fOutput->Add(fRecPLadLay4);
241
242   fTrackPLadLay3 = new TH1F("hTPLad3","Track Points per Ladder Layer 3",14,-0.5,13.5);
243   fTrackPLadLay3->Sumw2();
244   fTrackPLadLay3->SetMinimum(0);
245   fOutput->Add(fTrackPLadLay3);
246
247   fTrackPLadLay4 = new TH1F("hTPLad4","Track Points per Ladder Layer 4",22,-0.5,21.5);
248   fTrackPLadLay4->Sumw2();
249   fTrackPLadLay4->SetMinimum(0);
250   fOutput->Add(fTrackPLadLay4);
251
252   fGoodAnLadLay3 = new TH1F("hGALad3","Good Anodes per Ladder Layer 3",14,-0.5,13.5);
253   fOutput->Add(fGoodAnLadLay3);
254
255   fGoodAnLadLay4 = new TH1F("hGALad4","Good Anodes per Ladder Layer 4",22,-0.5,21.5);
256   fOutput->Add(fGoodAnLadLay4);
257
258   fDriftTimeRP=new TH1F("hDrTimRP","Drift Time from Rec Points (ns)",640,0.,6400.);
259   fDriftTimeRP->Sumw2();
260   fDriftTimeRP->SetMinimum(0.);
261   fOutput->Add(fDriftTimeRP);
262
263   fDriftTimeTPAll=new TH1F("hDrTimTPAll","Drift Time from Track Points (ns)",640,0.,6400.);
264   fDriftTimeTPAll->Sumw2();
265   fDriftTimeTPAll->SetMinimum(0.);
266   fOutput->Add(fDriftTimeTPAll);
267
268   fDriftTimeTPNoExtra=new TH1F("hDrTimTPNoExtra","Drift Time from Track Points (ns)",640,0.,6400.);
269   fDriftTimeTPNoExtra->Sumw2();
270   fDriftTimeTPNoExtra->SetMinimum(0.);
271   fOutput->Add(fDriftTimeTPNoExtra);
272
273   fDriftTimeTPExtra=new TH1F("hDrTimTPExtra","Drift Time from Track Points (ns)",640,0.,6400.);
274   fDriftTimeTPExtra->Sumw2();
275   fDriftTimeTPExtra->SetMinimum(0.);
276   fOutput->Add(fDriftTimeTPExtra);
277
278   // dE/dx histos
279
280   fHistdEdxL3VsP=new TH2F("hdEdxL3VsP","dE/dx vs. p lay3",40,0.,2.,100,0.,500.);
281   fHistdEdxL3VsP->Sumw2();
282   fHistdEdxL3VsP->SetMinimum(0);
283   fOutput->Add(fHistdEdxL3VsP);
284
285   fHistdEdxL4VsP=new TH2F("hdEdxL4VsP","dE/dx vs. p lay4",40,0.,2.,100,0.,500);
286   fHistdEdxL4VsP->Sumw2();
287   fHistdEdxL4VsP->SetMinimum(0);
288   fOutput->Add(fHistdEdxL4VsP);
289
290   fHistdEdxVsMod=new TH2F("hdEdxVsMod","dE/dx vs. mod",260,239.5,499.5,100,0.,500.);
291   fHistdEdxVsMod->Sumw2();
292   fHistdEdxVsMod->SetMinimum(0);
293   fOutput->Add(fHistdEdxVsMod);
294
295   for(Int_t it=0; it<8; it++){
296     fSignalTime[it]=new TH1F(Form("hSigTimeInt%d",it),Form("hSigTimeInt%d",it),100,0.,300.);
297     fSignalTime[it]->Sumw2();
298     fSignalTime[it]->SetMinimum(0);
299     fOutput->Add(fSignalTime[it]);
300   }
301
302   // cluster size histos
303   
304   fCluSizAnVsTime = new TH2F("hCluSizAn","hCluSizAn",40,0.,6400.,15,-0.5,14.5);
305   fCluSizAnVsTime->Sumw2();
306   fCluSizAnVsTime->SetMinimum(0);
307   fOutput->Add(fCluSizAnVsTime);
308
309   fCluSizTbVsTime = new TH2F("hCluSizTb","hCluSizTb",40,0.,6400.,15,-0.5,14.5);
310   fCluSizTbVsTime->Sumw2();
311   fCluSizTbVsTime->SetMinimum(0);
312   fOutput->Add(fCluSizTbVsTime);
313
314
315     // Read dead channels from OCDB
316
317 // == Not allowed in QA train. This is set by AliTaskCDBconnect. (A.G. 14/10/2011)
318 // If needed for independent use of the task, protect via a flag that is OFF by default.
319 /*
320   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
321   if (!mgr) AliFatal("No analysis manager");
322
323   AliCDBManager* man = AliCDBManager::Instance();
324   man->SetDefaultStorage("raw://");
325   Int_t nrun=mgr->GetRunFromPath();
326   man->SetRun(nrun);
327 */  
328     
329   PostData(1, fOutput);
330
331 }
332 //______________________________________________________________________________
333 void AliAnalysisTaskSDDRP::UserExec(Option_t *)
334 {
335   //
336   static Bool_t initCalib = kFALSE;
337   AliESDEvent *esd = (AliESDEvent*) (InputEvent());
338   if(!esd) {
339     printf("AliAnalysisTaskSDDRP::Exec(): bad ESD\n");
340     return;
341   } 
342
343
344   if(!ESDfriend()) {
345     printf("AliAnalysisTaskSDDRP::Exec(): bad ESDfriend\n");
346     return;
347   }
348   
349 // Code below moved from UserCreateOutputObjects where the run for OCDB may not 
350 // be yet properly set. Make sure this is called only once. (A.G. 14/10/2011)
351
352 /************/
353   if (!initCalib) {
354     AliCDBManager* man = AliCDBManager::Instance();
355     if (!man) {
356        AliFatal("CDB not set but needed by AliAnalysisTaskSDDRP");
357        return;
358     }   
359     AliCDBEntry* eT=(AliCDBEntry*)man->Get("GRP/CTP/Config");
360     if(eT){
361       eT->PrintId();
362       eT->PrintMetaData();
363       fTrigConfig=(AliTriggerConfiguration*)eT->GetObject();
364     }
365     if(!eT || !fTrigConfig){
366       AliError("Cannot retrieve CDB entry for GRP/CTP/Config");
367       return;      
368     }
369     AliCDBEntry* eR=(AliCDBEntry*)man->Get("ITS/Calib/RespSDD");
370     if (eR) {
371       eR->PrintId();
372       eR->PrintMetaData();
373       fResp=(AliITSresponseSDD*)eR->GetObject();
374     }else{
375       AliError("Cannot retrieve CDB entry for ITS/Calib/RespSDD");
376       return;
377     }
378     AliCDBEntry* eC=(AliCDBEntry*)man->Get("ITS/Calib/CalibSDD");
379     if (!eC) {
380        AliError("Cannot retrieve CDB entry for ITS/Calib/CalibSDD");
381        return;
382     }   
383     Int_t countGood3[14];
384     Int_t countGood4[22];
385     Int_t countGoodMod[260];
386     eC->PrintId();
387     eC->PrintMetaData();
388     TObjArray* calsdd=(TObjArray*)eC->GetObject();
389     for(Int_t ilad=0;ilad<14;ilad++) countGood3[ilad]=0;
390     for(Int_t ilad=0;ilad<22;ilad++) countGood4[ilad]=0;
391     for(Int_t imod=0;imod<260;imod++) countGoodMod[imod]=0;
392     for(Int_t imod=0;imod<260;imod++){
393       AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)calsdd->At(imod);
394       if(cal->IsBad()) continue;
395            Int_t modId=imod+AliITSgeomTGeo::GetModuleIndex(3,1,1);
396       Int_t lay,lad,det;
397       AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
398       if(fExcludeBadMod && !CheckModule(lay,lad,det)) continue;
399       for(Int_t ian=0; ian<512; ian++){
400         if(cal->IsBadChannel(ian)) continue;
401         countGoodMod[imod]++;
402         if(lay==3) countGood3[lad-1]++;
403         else if(lay==4) countGood4[lad-1]++;
404       }
405     }
406   
407     for(Int_t imod=0;imod<260;imod++) fGoodAnMod->SetBinContent(imod+1,countGoodMod[imod]);
408     fGoodAnMod->SetMinimum(0);
409     for(Int_t ilad=0;ilad<14;ilad++) fGoodAnLadLay3->SetBinContent(ilad+1,countGood3[ilad]);
410     fGoodAnLadLay3->SetMinimum(0);    
411     for(Int_t ilad=0;ilad<22;ilad++) fGoodAnLadLay4->SetBinContent(ilad+1,countGood4[ilad]);
412     fGoodAnLadLay4->SetMinimum(0);
413     initCalib = kTRUE;
414   }  
415 /************/
416   
417   PostData(1, fOutput);
418
419   fHistNEvents->Fill(-1);
420
421   TString firedTriggerClasses=esd->GetFiredTriggerClasses();
422   if(!firedTriggerClasses.Contains(fTriggerClass.Data())) return;
423   fHistNEvents->Fill(0);
424
425   Bool_t sddOK=esd->IsDetectorInTriggerCluster("ITSSDD",fTrigConfig);
426   if(!sddOK) fHistNEvents->Fill(1.);
427   else fHistNEvents->Fill(2.);
428   if(fOnlyEventsWithSDD && !sddOK) return; 
429   else fHistNEvents->Fill(3.);
430   
431   const AliTrackPointArray *array = 0;
432   Int_t ntracks = esd->GetNumberOfTracks();
433   for (Int_t itrack=0; itrack < ntracks; itrack++) {
434     AliESDtrack * track = esd->GetTrack(itrack);
435     if (!track) continue;
436
437     Bool_t accept=kTRUE;
438     if(fUseITSsaTracks){ 
439       if(track->GetNcls(1)>0) accept=kFALSE;
440     }else{
441       if(track->GetNcls(1)<fMinTPCpts) accept=kFALSE;
442     }
443     if(track->GetNcls(0) < fMinITSpts) accept=kFALSE;    
444     Int_t trstatus=track->GetStatus();
445     if(!(trstatus&AliESDtrack::kITSrefit)) accept=kFALSE;
446     if(!accept) continue;
447
448     fHistCluInLay->Fill(-1.); // bin -1 counts accepted tracks
449     UChar_t clumap=track->GetITSClusterMap();
450     for(Int_t iBit=0; iBit<6; iBit++){
451       if(clumap&(1<<iBit)) fHistCluInLay->Fill(iBit);
452     }
453
454
455     Double_t dedx[4];
456     track->GetITSdEdxSamples(dedx);
457     Float_t mom=track->P();
458     Int_t iMod,status;
459     Float_t xloc,zloc;
460     for(Int_t iLay=2; iLay<=3; iLay++){
461       Bool_t ok=track->GetITSModuleIndexInfo(iLay,iMod,status,xloc,zloc);
462       if(ok){
463         iMod+=240;
464         fHistAllPMod->Fill(iMod);
465         fHistAllPXloc->Fill(xloc);
466         fHistAllPZloc->Fill(zloc);
467         if(status==1){
468           fHistGoodPMod->Fill(iMod);
469           fHistGoodPXloc->Fill(xloc);
470           fHistGoodPZloc->Fill(zloc);
471           if(mom>fMinPfordEdx) fHistdEdxVsMod->Fill(iMod,dedx[iLay-2]);
472           if(iLay==2) fHistdEdxL3VsP->Fill(mom,dedx[0]);
473           else fHistdEdxL4VsP->Fill(mom,dedx[1]);
474         }
475         else if(status==2){ 
476           fHistBadRegMod->Fill(iMod);
477           fHistBadRegXloc->Fill(xloc);
478           fHistBadRegZloc->Fill(zloc);
479         }
480         else if(status==3) fHistSkippedMod->Fill(iMod);
481         else if(status==4) fHistOutAccMod->Fill(iMod);
482         else if(status==5){
483           fHistMissPMod->Fill(iMod);
484           fHistMissPXloc->Fill(xloc);
485           fHistMissPZloc->Fill(zloc);
486         }
487         else if(status==6) fHistNoRefitMod->Fill(iMod);
488       }
489     }
490
491
492     array = track->GetTrackPointArray();
493     if(!array) continue;
494     for(Int_t ipt=0; ipt<array->GetNPoints(); ipt++) {
495       AliTrackPoint point;
496       Int_t modId;
497       array->GetPoint(point,ipt);
498       Int_t volId = point.GetVolumeID();
499       Int_t layerId = AliGeomManager::VolUIDToLayer(volId,modId);
500       if(layerId<3 || layerId>4) continue;
501       modId+=AliITSgeomTGeo::GetModuleIndex(layerId,1,1);
502       Int_t lay,lad,det;
503       AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
504       if(fExcludeBadMod && !CheckModule(lay,lad,det)) continue;
505       fTrackPMod->Fill(modId);
506       fDriftTimeTPAll->Fill(point.GetDriftTime());
507       if(point.IsExtra()) fDriftTimeTPExtra->Fill(point.GetDriftTime());
508       else fDriftTimeTPNoExtra->Fill(point.GetDriftTime());
509       Float_t dtime=point.GetDriftTime()-fResp->GetTimeZero(modId);
510       Int_t cluTyp=point.GetClusterType();
511       Int_t clSizAn=(cluTyp>>8)&0xFF;
512       Int_t clSizTb=cluTyp&0xFF;
513       fCluSizAnVsTime->Fill(dtime,clSizAn);
514       fCluSizTbVsTime->Fill(dtime,clSizTb);
515       Int_t theBin=int(dtime/6500.*8.);
516       if(layerId==3){ 
517         fTrackPLadLay3->Fill(lad-1);      
518         if(dedx[0]>0. && track->P()>fMinPfordEdx) fSignalTime[theBin]->Fill(dedx[0]);
519       }
520       if(layerId==4){
521         fTrackPLadLay4->Fill(lad-1);
522         if(dedx[1]>0.&& track->P()>fMinPfordEdx) fSignalTime[theBin]->Fill(dedx[1]);
523       }
524     }
525   }
526
527   AliESDInputHandlerRP *hand = dynamic_cast<AliESDInputHandlerRP*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
528   TTree* tR = 0;
529   if (hand) tR = hand->GetTreeR("ITS");
530   if (tR){
531     TClonesArray *ITSrec= new TClonesArray("AliITSRecPoint");
532     TBranch *branch =tR->GetBranch("ITSRecPoints");
533     branch->SetAddress(&ITSrec);
534     for (Int_t modId=240; modId<500; modId++){
535       Int_t lay,lad,det;
536       AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
537       if(fExcludeBadMod && !CheckModule(lay,lad,det)) continue;
538       branch->GetEvent(modId);
539       Int_t nrecp = ITSrec->GetEntries();       
540       fRecPMod->Fill(modId,nrecp);        
541       if(lay==3) fRecPLadLay3->Fill(lad-1,nrecp);
542       if(lay==4) fRecPLadLay4->Fill(lad-1,nrecp);
543       for(Int_t irec=0;irec<nrecp;irec++) {
544         AliITSRecPoint *recp = (AliITSRecPoint*)ITSrec->At(irec);
545         fDriftTimeRP->Fill(recp->GetDriftTime());
546       }
547     }
548     ITSrec->Delete();
549     delete ITSrec;
550   }
551   PostData(1,fOutput);
552   
553 }
554 //______________________________________________________________________________
555 Bool_t AliAnalysisTaskSDDRP::CheckModule(Int_t lay, Int_t lad, Int_t det) const{
556   //
557   if(lay==4){
558     if(lad==3 && det==5) return kFALSE; // 1500 V
559     if(lad==3 && det==6) return kFALSE; // 0 V
560     if(lad==3 && det==7) return kFALSE; // 1500 V
561     if(lad==4 && det==1) return kFALSE; // 0 V
562     if(lad==4 && det==2) return kFALSE; // 1500 V
563     if(lad==7 && det==5) return kFALSE; // 0 MV
564     if(lad==9 && det==3) return kFALSE; // 1500 V
565     if(lad==9 && det==4) return kFALSE; // 0 V
566     if(lad==9 && det==5) return kFALSE; // 1500 V
567     if(lad==11 && det==6) return kFALSE; // 1500 V
568     if(lad==11 && det==7) return kFALSE; // 0 V
569     if(lad==11 && det==8) return kFALSE; // 1500 V
570     if(lad==18 && det==5) return kFALSE; // 1500 V
571     if(lad==18 && det==6) return kFALSE; // 0 V
572     if(lad==18 && det==7) return kFALSE; // 1500 V
573     if(lad==22 && det==1) return kFALSE; // 0 V
574     if(lad==22 && det==2) return kFALSE; // 1500 V
575   }
576   if(lay==3){
577     if(lad==4 && det==4) return kFALSE; // 1500 V 
578     if(lad==3) return kFALSE;  // swapped in geometry
579   }
580   return kTRUE;
581 }
582
583 //______________________________________________________________________________
584 void AliAnalysisTaskSDDRP::Terminate(Option_t */*option*/)
585 {
586   // Terminate analysis
587   fOutput = dynamic_cast<TList*> (GetOutputData(1));
588   if (!fOutput) {     
589     printf("ERROR: fOutput not available\n");
590     return;
591   }
592   fHistNEvents= dynamic_cast<TH1F*>(fOutput->FindObject("hNEvents"));
593   fHistCluInLay= dynamic_cast<TH1F*>(fOutput->FindObject("hCluInLay"));
594
595   fHistAllPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hAllPMod"));
596   fHistGoodPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hGoodPMod"));
597   fHistBadRegMod= dynamic_cast<TH1F*>(fOutput->FindObject("hBadRegMod"));
598   fHistMissPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hMissPMod"));
599   fHistSkippedMod= dynamic_cast<TH1F*>(fOutput->FindObject("hSkippedMod"));
600   fHistOutAccMod= dynamic_cast<TH1F*>(fOutput->FindObject("hOutAccMod"));
601   fHistNoRefitMod= dynamic_cast<TH1F*>(fOutput->FindObject("hNoRefitMod"));
602
603   fHistAllPXloc= dynamic_cast<TH1F*>(fOutput->FindObject("hAllPxloc"));
604   fHistGoodPXloc= dynamic_cast<TH1F*>(fOutput->FindObject("hGoodPxloc"));
605   fHistBadRegXloc= dynamic_cast<TH1F*>(fOutput->FindObject("hBadRegxloc"));
606   fHistMissPXloc= dynamic_cast<TH1F*>(fOutput->FindObject("hMissPxloc"));
607   fHistAllPZloc= dynamic_cast<TH1F*>(fOutput->FindObject("hAllPzloc"));
608   fHistGoodPZloc= dynamic_cast<TH1F*>(fOutput->FindObject("hGoodPzloc"));
609   fHistBadRegZloc= dynamic_cast<TH1F*>(fOutput->FindObject("fHistBadRegZloc"));
610   fHistMissPZloc= dynamic_cast<TH1F*>(fOutput->FindObject("hMissPzloc"));
611
612   fHistdEdxL3VsP= dynamic_cast<TH2F*>(fOutput->FindObject("hdEdxL3VsP"));
613   fHistdEdxL4VsP= dynamic_cast<TH2F*>(fOutput->FindObject("hdEdxL4VsP"));
614   fHistdEdxVsMod= dynamic_cast<TH2F*>(fOutput->FindObject("hdEdxVsMod"));
615
616   fRecPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hRPMod"));
617   fTrackPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hTPMod"));
618   fGoodAnMod= dynamic_cast<TH1F*>(fOutput->FindObject("hGAMod"));
619
620   fRecPLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hRPLad3"));
621   fRecPLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hRPLad4"));
622   fTrackPLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hTPLad3"));
623   fTrackPLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hTPLad4"));
624   fGoodAnLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hGALad3"));
625   fGoodAnLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hGALad4"));
626
627   fDriftTimeRP= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimRP"));
628   fDriftTimeTPAll= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimTPAll"));
629   fDriftTimeTPNoExtra= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimTPNoExtra"));
630   fDriftTimeTPExtra= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimTPExtra"));
631
632   for(Int_t it=0; it<8; it++){
633     fSignalTime[it]= dynamic_cast<TH1F*>(fOutput->FindObject(Form("hSigTimeInt%d",it)));
634   }
635   fCluSizAnVsTime= dynamic_cast<TH2F*>(fOutput->FindObject("hCluSizAn"));
636   fCluSizTbVsTime= dynamic_cast<TH2F*>(fOutput->FindObject("hCluSizTb"));
637
638   return;
639 }
640
641
642
643
644