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 "AliGeomManager.h"
19 #include <TGeoGlobalMagField.h>
20 #include "AliESDInputHandlerRP.h"
21 /**************************************************************************
22 * Copyright(c) 1998-2010, ALICE Experiment at CERN, All rights reserved. *
24 * Author: The ALICE Off-line Project. *
25 * Contributors are mentioned in the code where appropriate. *
27 * Permission to use, copy, modify and distribute this software and its *
28 * documentation strictly for non-commercial purposes is hereby granted *
29 * without fee, provided that the above copyright notice appears in all *
30 * copies and that both the copyright notice and this permission notice *
31 * appear in the supporting documentation. The authors make no claims *
32 * about the suitability of this software for any purpose. It is *
33 * provided "as is" without express or implied warranty. *
34 **************************************************************************/
36 //*************************************************************************
37 // Implementation of class AliAnalysiTaskSDDRP
38 // AliAnalysisTaskSE to extract from ESD + ESDfreinds + ITS rec points
39 // performance plots for SDD detector
41 // Author: F. Prino, prino@to.infn.it
42 //*************************************************************************
45 #include "AliAnalysisTaskSDDRP.h"
47 ClassImp(AliAnalysisTaskSDDRP)
48 //______________________________________________________________________________
49 AliAnalysisTaskSDDRP::AliAnalysisTaskSDDRP() : AliAnalysisTaskSE("SDD RecPoints"),
73 fDriftTimeTPNoExtra(0),
78 fUseITSsaTracks(kFALSE),
86 DefineOutput(1, TList::Class());
90 //___________________________________________________________________________
91 AliAnalysisTaskSDDRP::~AliAnalysisTaskSDDRP(){
100 //___________________________________________________________________________
102 void AliAnalysisTaskSDDRP::UserCreateOutputObjects() {
104 fOutput = new TList();
106 fOutput->SetName("OutputHistos");
108 fHistNEvents = new TH1F("hNEvents", "Number of processed events",3,-1.5,1.5);
109 fHistNEvents->Sumw2();
110 fHistNEvents->SetMinimum(0);
111 fOutput->Add(fHistNEvents);
115 fHistAllPMod = new TH1F("hAllPmod","Crossing Tracks vs. Module",260,239.5,499.5);
116 fHistAllPMod->Sumw2();
117 fHistAllPMod->SetMinimum(0);
118 fOutput->Add(fHistAllPMod);
120 fHistGoodPMod = new TH1F("hGoodPmod","PointsAssocToTrack per Module",260,239.5,499.5);
121 fHistGoodPMod->Sumw2();
122 fHistGoodPMod->SetMinimum(0);
123 fOutput->Add(fHistGoodPMod);
125 fHistBadRegMod = new TH1F("hBadRegmod","Tracks in BadRegion per Module",260,239.5,499.5);
126 fHistBadRegMod->Sumw2();
127 fHistBadRegMod->SetMinimum(0);
128 fOutput->Add(fHistBadRegMod);
130 fHistMissPMod = new TH1F("hMissPmod","Missing Points per Module",260,239.5,499.5);
131 fHistMissPMod->Sumw2();
132 fHistMissPMod->SetMinimum(0);
133 fOutput->Add(fHistMissPMod);
135 fHistSkippedMod = new TH1F("hSkippedmod","Tracks in Skipped Module",260,239.5,499.5);
136 fHistSkippedMod->Sumw2();
137 fHistSkippedMod->SetMinimum(0);
138 fOutput->Add(fHistSkippedMod);
140 fHistOutAccMod = new TH1F("hOutAccmod","Tracks outside zAcc per Module",260,239.5,499.5);
141 fHistOutAccMod->Sumw2();
142 fHistOutAccMod->SetMinimum(0);
143 fOutput->Add(fHistOutAccMod);
145 fHistNoRefitMod = new TH1F("hNoRefitmod","Points rejected in refit per Module",260,239.5,499.5);
146 fHistNoRefitMod->Sumw2();
147 fHistNoRefitMod->SetMinimum(0);
148 fOutput->Add(fHistNoRefitMod);
150 fHistdEdxL3VsP=new TH2F("hdEdxL3VsP","dE/dx vs. p lay3",40,0.,2.,100,0.,500.);
151 fHistdEdxL3VsP->Sumw2();
152 fHistdEdxL3VsP->SetMinimum(0);
153 fOutput->Add(fHistdEdxL3VsP);
155 fHistdEdxL4VsP=new TH2F("hdEdxL4VsP","dE/dx vs. p lay4",40,0.,2.,100,0.,500);
156 fHistdEdxL4VsP->Sumw2();
157 fHistdEdxL4VsP->SetMinimum(0);
158 fOutput->Add(fHistdEdxL4VsP);
160 fHistdEdxVsMod=new TH2F("hdEdxVsMod","dE/dx vs. mod",260,239.5,499.5,100,0.,500.);
161 fHistdEdxVsMod->Sumw2();
162 fHistdEdxVsMod->SetMinimum(0);
163 fOutput->Add(fHistdEdxVsMod);
166 fRecPMod = new TH1F("hRPMod","Rec Points per Module",260,239.5,499.5);
168 fRecPMod->SetMinimum(0);
169 fOutput->Add(fRecPMod);
171 fTrackPMod = new TH1F("hTPMod","Track Points per Module",260,239.5,499.5);
173 fTrackPMod->SetMinimum(0);
174 fOutput->Add(fTrackPMod);
176 fGoodAnMod = new TH1F("hGAMod","Good Anodes per Module",260,239.5,499.5);
177 fOutput->Add(fGoodAnMod);
181 fRecPLadLay3 = new TH1F("hRPLad3","Rec Points per Ladder Layer 3",14,-0.5,13.5);
182 fRecPLadLay3->Sumw2();
183 fRecPLadLay3->SetMinimum(0);
184 fOutput->Add(fRecPLadLay3);
186 fRecPLadLay4 = new TH1F("hRPLad4","Rec Points per Ladder Layer 4",22,-0.5,21.5);
187 fRecPLadLay4->Sumw2();
188 fRecPLadLay4->SetMinimum(0);
189 fOutput->Add(fRecPLadLay4);
191 fTrackPLadLay3 = new TH1F("hTPLad3","Track Points per Ladder Layer 3",14,-0.5,13.5);
192 fTrackPLadLay3->Sumw2();
193 fTrackPLadLay3->SetMinimum(0);
194 fOutput->Add(fTrackPLadLay3);
196 fTrackPLadLay4 = new TH1F("hTPLad4","Track Points per Ladder Layer 4",22,-0.5,21.5);
197 fTrackPLadLay4->Sumw2();
198 fTrackPLadLay4->SetMinimum(0);
199 fOutput->Add(fTrackPLadLay4);
201 fGoodAnLadLay3 = new TH1F("hGALad3","Good Anodes per Ladder Layer 3",14,-0.5,13.5);
202 fOutput->Add(fGoodAnLadLay3);
204 fGoodAnLadLay4 = new TH1F("hGALad4","Good Anodes per Ladder Layer 4",22,-0.5,21.5);
205 fOutput->Add(fGoodAnLadLay4);
207 fDriftTimeRP=new TH1F("hDrTimRP","Drift Time from Rec Points (ns)",640,0.,6400.);
208 fDriftTimeRP->Sumw2();
209 fDriftTimeRP->SetMinimum(0.);
210 fOutput->Add(fDriftTimeRP);
212 fDriftTimeTPAll=new TH1F("hDrTimTPAll","Drift Time from Track Points (ns)",640,0.,6400.);
213 fDriftTimeTPAll->Sumw2();
214 fDriftTimeTPAll->SetMinimum(0.);
215 fOutput->Add(fDriftTimeTPAll);
217 fDriftTimeTPNoExtra=new TH1F("hDrTimTPNoExtra","Drift Time from Track Points (ns)",640,0.,6400.);
218 fDriftTimeTPNoExtra->Sumw2();
219 fDriftTimeTPNoExtra->SetMinimum(0.);
220 fOutput->Add(fDriftTimeTPNoExtra);
222 fDriftTimeTPExtra=new TH1F("hDrTimTPExtra","Drift Time from Track Points (ns)",640,0.,6400.);
223 fDriftTimeTPExtra->Sumw2();
224 fDriftTimeTPExtra->SetMinimum(0.);
225 fOutput->Add(fDriftTimeTPExtra);
227 for(Int_t it=0; it<8; it++){
228 fSignalTime[it]=new TH1F(Form("hSigTimeInt%d",it),Form("hSigTimeInt%d",it),100,0.,300.);
229 fSignalTime[it]->Sumw2();
230 fSignalTime[it]->SetMinimum(0);
231 fOutput->Add(fSignalTime[it]);
234 //______________________________________________________________________________
235 void AliAnalysisTaskSDDRP::UserExec(Option_t *)
238 fESD = (AliESDEvent*) (InputEvent());
239 fESDfriend = static_cast<AliESDfriend*>(fESD->FindListObject("AliESDfriend"));
242 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
245 printf("ERROR: Could not get ESDInputHandler\n");
248 fESD = esdH->GetEvent();
251 printf("AliAnalysisTaskSDDRP::Exec(): bad ESD\n");
254 // fESDfriend = static_cast<AliESDfriend*>(fESD->FindListObject("AliESDfriend"));
258 printf("AliAnalysisTaskSDDRP::Exec(): bad ESDfriend\n");
263 PostData(1, fOutput);
264 fESD->SetESDfriend(fESDfriend);
269 AliCDBManager* man = AliCDBManager::Instance();
270 man->SetDefaultStorage("raw://");
271 man->SetRun(fESD->GetRunNumber());
274 AliCDBEntry* eR=(AliCDBEntry*)man->Get("ITS/Calib/RespSDD");
275 fResp=(AliITSresponseSDD*)eR->GetObject();
277 AliCDBEntry* eC=(AliCDBEntry*)man->Get("ITS/Calib/CalibSDD");
278 TObjArray* calsdd=(TObjArray*)eC->GetObject();
279 Int_t countGood3[14];
280 Int_t countGood4[22];
281 Int_t countGoodMod[260];
282 for(Int_t ilad=0;ilad<14;ilad++) countGood3[ilad]=0;
283 for(Int_t ilad=0;ilad<22;ilad++) countGood4[ilad]=0;
284 for(Int_t imod=0;imod<260;imod++) countGoodMod[imod]=0;
285 for(Int_t imod=0;imod<260;imod++){
286 AliITSCalibrationSDD* cal=(AliITSCalibrationSDD*)calsdd->At(imod);
287 if(cal->IsBad()) continue;
288 Int_t modId=imod+AliITSgeomTGeo::GetModuleIndex(3,1,1);
290 AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
291 if(!CheckModule(lay,lad,det)) continue;
292 for(Int_t ian=0; ian<512; ian++){
293 if(cal->IsBadChannel(ian)) continue;
294 countGoodMod[imod]++;
295 if(lay==3) countGood3[lad-1]++;
296 else if(lay==4) countGood4[lad-1]++;
299 for(Int_t imod=0;imod<260;imod++) fGoodAnMod->SetBinContent(imod+1,countGoodMod[imod]);
300 fGoodAnMod->SetMinimum(0);
301 for(Int_t ilad=0;ilad<14;ilad++) fGoodAnLadLay3->SetBinContent(ilad+1,countGood3[ilad]);
302 fGoodAnLadLay3->SetMinimum(0);
303 for(Int_t ilad=0;ilad<22;ilad++) fGoodAnLadLay4->SetBinContent(ilad+1,countGood4[ilad]);
304 fGoodAnLadLay4->SetMinimum(0);
308 fHistNEvents->Fill(0);
310 if(!fESD->IsTriggerClassFired("CINT1B-ABCE-NOPF-ALL")) return;
311 fHistNEvents->Fill(1);
313 const AliTrackPointArray *array = 0;
314 Int_t ntracks = fESD->GetNumberOfTracks();
315 for (Int_t itrack=0; itrack < ntracks; itrack++) {
316 AliESDtrack * track = fESD->GetTrack(itrack);
317 if (!track) continue;
321 if(track->GetNcls(1)>0) accept=kFALSE;
323 if(track->GetNcls(1)<fMinTPCpts) accept=kFALSE;
325 if(track->GetNcls(0) < fMinITSpts) accept=kFALSE;
326 Int_t trstatus=track->GetStatus();
327 if(!(trstatus&AliESDtrack::kITSrefit)) accept=kFALSE;
328 if(!accept) continue;
331 track->GetITSdEdxSamples(dedx);
332 Float_t mom=track->P();
335 for(Int_t iLay=2; iLay<=3; iLay++){
336 Bool_t ok=track->GetITSModuleIndexInfo(iLay,iMod,status,xloc,zloc);
339 fHistAllPMod->Fill(iMod);
341 fHistGoodPMod->Fill(iMod);
342 if(mom>fMinPfordEdx) fHistdEdxVsMod->Fill(iMod,dedx[iLay-2]);
343 if(iLay==2) fHistdEdxL3VsP->Fill(mom,dedx[0]);
344 else fHistdEdxL4VsP->Fill(mom,dedx[1]);
346 else if(status==2) fHistBadRegMod->Fill(iMod);
347 else if(status==3) fHistSkippedMod->Fill(iMod);
348 else if(status==4) fHistOutAccMod->Fill(iMod);
349 else if(status==5) fHistMissPMod->Fill(iMod);
350 else if(status==6) fHistNoRefitMod->Fill(iMod);
355 array = track->GetTrackPointArray();
357 for(Int_t ipt=0; ipt<array->GetNPoints(); ipt++) {
360 array->GetPoint(point,ipt);
361 Int_t volId = point.GetVolumeID();
362 Int_t layerId = AliGeomManager::VolUIDToLayer(volId,modId);
363 if(layerId<3 || layerId>4) continue;
364 modId+=AliITSgeomTGeo::GetModuleIndex(layerId,1,1);
366 AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
367 if(!CheckModule(lay,lad,det)) continue;
368 fTrackPMod->Fill(modId);
369 fDriftTimeTPAll->Fill(point.GetDriftTime());
370 if(point.IsExtra()) fDriftTimeTPExtra->Fill(point.GetDriftTime());
371 else fDriftTimeTPNoExtra->Fill(point.GetDriftTime());
372 Float_t dtime=point.GetDriftTime()-fResp->GetTimeZero(modId);
373 Int_t theBin=int(dtime/6500.*8.);
375 fTrackPLadLay3->Fill(lad-1);
376 if(dedx[0]>0. && track->P()>fMinPfordEdx) fSignalTime[theBin]->Fill(dedx[0]);
379 fTrackPLadLay4->Fill(lad-1);
380 if(dedx[1]>0.&& track->P()>fMinPfordEdx) fSignalTime[theBin]->Fill(dedx[1]);
385 AliESDInputHandlerRP *hand = dynamic_cast<AliESDInputHandlerRP*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
386 TTree* tR = hand->GetTreeR("ITS");
388 TClonesArray *ITSrec= new TClonesArray("AliITSRecPoint");
389 TBranch *branch =tR->GetBranch("ITSRecPoints");
390 branch->SetAddress(&ITSrec);
391 for (Int_t modId=240; modId<500; modId++){
393 AliITSgeomTGeo::GetModuleId(modId,lay,lad,det);
394 if(!CheckModule(lay,lad,det)) continue;
395 branch->GetEvent(modId);
396 Int_t nrecp = ITSrec->GetEntries();
397 fRecPMod->Fill(modId,nrecp);
398 if(lay==3) fRecPLadLay3->Fill(lad-1,nrecp);
399 if(lay==4) fRecPLadLay4->Fill(lad-1,nrecp);
400 for(Int_t irec=0;irec<nrecp;irec++) {
401 AliITSRecPoint *recp = (AliITSRecPoint*)ITSrec->At(irec);
402 fDriftTimeRP->Fill(recp->GetDriftTime());
411 //______________________________________________________________________________
412 Bool_t AliAnalysisTaskSDDRP::CheckModule(Int_t lay, Int_t lad, Int_t det) const{
415 if(lad==3 && det==5) return kFALSE; // 1500 V
416 if(lad==3 && det==6) return kFALSE; // 0 V
417 if(lad==3 && det==7) return kFALSE; // 1500 V
418 if(lad==4 && det==1) return kFALSE; // 0 V
419 if(lad==4 && det==2) return kFALSE; // 1500 V
420 if(lad==7 && det==5) return kFALSE; // 0 MV
421 if(lad==9 && det==3) return kFALSE; // 1500 V
422 if(lad==9 && det==4) return kFALSE; // 0 V
423 if(lad==9 && det==5) return kFALSE; // 1500 V
424 if(lad==11 && det==6) return kFALSE; // 1500 V
425 if(lad==11 && det==7) return kFALSE; // 0 V
426 if(lad==11 && det==8) return kFALSE; // 1500 V
427 if(lad==18 && det==5) return kFALSE; // 1500 V
428 if(lad==18 && det==6) return kFALSE; // 0 V
429 if(lad==18 && det==7) return kFALSE; // 1500 V
430 if(lad==22 && det==1) return kFALSE; // 0 V
431 if(lad==22 && det==2) return kFALSE; // 1500 V
434 if(lad==4 && det==4) return kFALSE; // 1500 V
435 if(lad==3) return kFALSE; // swapped in geometry
440 //______________________________________________________________________________
441 void AliAnalysisTaskSDDRP::Terminate(Option_t */*option*/)
443 // Terminate analysis
444 fOutput = dynamic_cast<TList*> (GetOutputData(1));
446 printf("ERROR: fOutput not available\n");
449 fHistNEvents= dynamic_cast<TH1F*>(fOutput->FindObject("hNEvents"));
451 fHistAllPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hAllPMod"));
452 fHistGoodPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hGoodPMod"));
453 fHistBadRegMod= dynamic_cast<TH1F*>(fOutput->FindObject("hBadRegMod"));
454 fHistMissPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hMissPMod"));
455 fHistSkippedMod= dynamic_cast<TH1F*>(fOutput->FindObject("hSkippedMod"));
456 fHistOutAccMod= dynamic_cast<TH1F*>(fOutput->FindObject("hOutAccMod"));
457 fHistNoRefitMod= dynamic_cast<TH1F*>(fOutput->FindObject("hNoRefitMod"));
459 fHistdEdxL3VsP= dynamic_cast<TH2F*>(fOutput->FindObject("hdEdxL3VsP"));
460 fHistdEdxL4VsP= dynamic_cast<TH2F*>(fOutput->FindObject("hdEdxL4VsP"));
461 fHistdEdxVsMod= dynamic_cast<TH2F*>(fOutput->FindObject("hdEdxVsMod"));
463 fRecPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hRPMod"));
464 fTrackPMod= dynamic_cast<TH1F*>(fOutput->FindObject("hTPMod"));
465 fGoodAnMod= dynamic_cast<TH1F*>(fOutput->FindObject("hGAMod"));
467 fRecPLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hRPLad3"));
468 fRecPLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hRPLad4"));
469 fTrackPLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hTPLad3"));
470 fTrackPLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hTPLad4"));
471 fGoodAnLadLay3= dynamic_cast<TH1F*>(fOutput->FindObject("hGALad3"));
472 fGoodAnLadLay4= dynamic_cast<TH1F*>(fOutput->FindObject("hGALad4"));
474 fDriftTimeRP= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimRP"));
475 fDriftTimeTPAll= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimTPAll"));
476 fDriftTimeTPNoExtra= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimTPNoExtra"));
477 fDriftTimeTPExtra= dynamic_cast<TH1F*>(fOutput->FindObject("hDrTimTPExtra"));
479 for(Int_t it=0; it<8; it++){
480 fSignalTime[it]= dynamic_cast<TH1F*>(fOutput->FindObject(Form("hSigTimeInt%d",it)));