1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 //-----------------------------------------------------------------------
18 // Author : A. Mastroserio
19 //-----------------------------------------------------------------------
22 #ifndef ALIANALYSISTASKSPD_CXX
23 #define ALIANALYSISTASKSPD_CXX
25 #include "AliAnalysisTaskSPD.h"
26 #include "TClonesArray.h"
31 #include "AliSPDUtils.h"
32 #include "AliESDEvent.h"
35 #include "AliESDInputHandlerRP.h"
36 #include "AliAnalysisManager.h"
37 #include "AliMultiplicity.h"
38 #include "AliCDBPath.h"
39 #include "AliCDBManager.h"
40 #include "AliCDBEntry.h"
41 #include "AliCDBStorage.h"
42 #include "AliGeomManager.h"
43 #include "AliITSRecPoint.h"
44 #include "AliITSsegmentationSPD.h"
45 ClassImp(AliAnalysisTaskSPD)
46 //__________________________________________________________________________
47 AliAnalysisTaskSPD::AliAnalysisTaskSPD() :
51 fOCDBLocation("local://$ALICE_ROOT/OCDB")
57 //___________________________________________________________________________
58 AliAnalysisTaskSPD::AliAnalysisTaskSPD(const Char_t* name) :
59 AliAnalysisTaskSE(name),
63 fOCDBLocation("local://$ALICE_ROOT/OCDB")
66 // Constructor. Initialization of Inputs and Outputs
68 Info("AliAnalysisTaskSPD","Calling Constructor");
70 DefineOutput(1,TList::Class());
74 //___________________________________________________________________________
75 AliAnalysisTaskSPD& AliAnalysisTaskSPD::operator=(const AliAnalysisTaskSPD& c)
78 // Assignment operator
81 AliAnalysisTaskSE::operator=(c) ;
85 fOCDBLocation = c.fOCDBLocation;
90 //___________________________________________________________________________
91 AliAnalysisTaskSPD::AliAnalysisTaskSPD(const AliAnalysisTaskSPD& c) :
96 fOCDBLocation(c.fOCDBLocation)
104 //___________________________________________________________________________
105 AliAnalysisTaskSPD::~AliAnalysisTaskSPD() {
110 Info("~AliAnalysisTaskSPD","Calling Destructor");
111 if (AliAnalysisManager::GetAnalysisManager()->IsProofMode()) return;
112 if (fSegSPD) delete fSegSPD ;
118 //___________________________________________________________________________
119 void AliAnalysisTaskSPD::UserCreateOutputObjects() {
121 Info("CreateOutputObjects","CreateOutputObjects of task %s", GetName());
123 LoadGeometryFromOCDB();
125 fSegSPD = new AliITSsegmentationSPD();
129 fOutput = new TList();
131 // Booking rec points related histograms
133 TH2D *hLocalMapL1 = new TH2D("hLocalMapL1"," Local coordinates - Layer 1",330,-16.5,16.5,205,-0.5,40.5); // safe limits for local coordinates in a module : z = -4,4, x = -1,1;
134 hLocalMapL1->SetXTitle("direction Z [cm]");
135 hLocalMapL1->SetYTitle("direction X [cm]");
136 fOutput->AddLast(hLocalMapL1);
138 TH2D *hLocalMapL2 = new TH2D("hLocalMapL2"," Local coordinates - Layer 2",330,-16.5,16.5,405,-0.5,80.5);
139 hLocalMapL2->SetXTitle("direction Z [cm]");
140 hLocalMapL2->SetYTitle("direction X [cm]");
141 fOutput->AddLast(hLocalMapL2);
143 TH1F *hClusterModYield = new TH1F("hClusterModYield","Cluster yield in modules",241,-0.5,240.5);
144 hClusterModYield->SetXTitle("module number");
145 fOutput->AddLast(hClusterModYield);
147 TH1F *hClusterYield = new TH1F("hClusterYield","Cluster yield per chip",1200,-0.5,1199.5) ;
148 hClusterYield->SetXTitle("chip key");
149 fOutput->AddLast(hClusterYield);
151 TH1F *hClusterYieldOnline = new TH1F("hClusterYieldOnline","Cluster yield per chip (online coord eq*60+hs*10+chip)",1200,-0.5,1199.5) ;
152 hClusterYieldOnline->SetXTitle("chip");
153 fOutput->AddLast(hClusterYieldOnline);
155 TH1F *hFiredChip = new TH1F("hFiredChip","Fired chip (at least one cluster)",1200,-0.5,1199.5);
156 hFiredChip->SetXTitle("chip key");
157 fOutput->AddLast(hFiredChip);
159 TH1F *hFOFiredChip = new TH1F("hFOFiredChip","FO Fired chip ",1200,-0.5,1199.5);
160 hFOFiredChip->SetXTitle("chip key");
161 fOutput->AddLast(hFOFiredChip);
163 TH1F *hFOgood = new TH1F("hFOgood"," FO-good (at least one cluster) ",1200,-0.5,1199.5);
164 hFOgood->SetXTitle("chip key");
165 fOutput->AddLast(hFOgood);
167 TH2F *hFOgoodPerBC = new TH2F("hFOgoodPerBCmod4"," FO-good signals in BCmod4 ",1200,-0.5,1199.5,4,-0.5,3.5);
168 hFOgoodPerBC->SetXTitle("chip key");
169 fOutput->AddLast(hFOgoodPerBC);
171 TH2F *hFiredChipsPerBC = new TH2F("hFiredChipsPerBCmod4"," fired chips in BCmod4 ",1200,-0.5,1199.5,4,-0.5,3.5);
172 hFiredChipsPerBC->SetXTitle("chip key");
173 fOutput->AddLast(hFiredChipsPerBC);
175 TH1F *hFOnoisy = new TH1F("hFOnoisy","FO-noisy (no cluster)",1200,-0.5,1199.5);
176 hFOnoisy->SetXTitle("chip key");
177 fOutput->AddLast(hFOnoisy);
180 // Booking ESD related histograms
184 TH1I *hTracklets = new TH1I("hNtracklets","Tracklet distribution",300,0,12000);
185 hTracklets->SetXTitle("# Tracklets");
186 fOutput->AddLast(hTracklets);
188 TH2F *hSPDphivsSPDeta= new TH2F("hSPDphivsSPDeta", "Tracklets - #varphi vs #eta",120,-3.,3,360,0.,2*TMath::Pi());
189 hSPDphivsSPDeta->GetXaxis()->SetTitle("Pseudorapidity #eta");
190 hSPDphivsSPDeta->GetYaxis()->SetTitle("#varphi [rad]");
191 fOutput->AddLast(hSPDphivsSPDeta);
193 TH1F *hSPDphiZpos= new TH1F("hSPDphiZpos", "Tracklets - #varphi (Z>0)",360,0.,2*TMath::Pi());
194 hSPDphiZpos->SetXTitle("#varphi [rad]");
195 fOutput->AddLast(hSPDphiZpos);
197 TH1F *hSPDphiZneg= new TH1F("hSPDphiZneg", "Tracklets - #varphi (Z<0)",360,0.,2*TMath::Pi());
198 hSPDphiZneg->SetXTitle("#varphi [rad]");
199 fOutput->AddLast(hSPDphiZneg);
201 TH1F *hVertexZ = new TH1F("hVertexZ","Vertex Z distribution",300,-15,15);
202 hVertexZ->SetXTitle("Z Vertex [cm]");
203 fOutput->AddLast(hVertexZ);
205 TH2F *hTracklVsClu1 = new TH2F("hTrackVsClu1","Tracklet Vs Clusters Layer 1",100,0,8000,300,0,12000);
206 hTracklVsClu1->SetXTitle("clusters SPD Layer 1");
207 hTracklVsClu1->SetYTitle("tracklets");
208 fOutput->AddLast(hTracklVsClu1);
210 TH2F *hTracklVsClu2 = new TH2F("hTrackVsClu2","Tracklet Vs Clusters Layer 2",100,0,8000,300,0,12000);
211 hTracklVsClu2->SetXTitle("clusters SPD Layer 2");
212 hTracklVsClu2->SetYTitle("tracklets");
213 fOutput->AddLast(hTracklVsClu2);
215 TH1I *hEventsProcessed = new TH1I("hEventsProcessed","Number of processed events",1,0,1) ;
216 fOutput->AddLast(hEventsProcessed);
220 //_________________________________________________
221 void AliAnalysisTaskSPD::UserExec(Option_t *)
224 // Main loop function
226 AliESDInputHandlerRP *hand = dynamic_cast<AliESDInputHandlerRP*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
228 printf("No AliESDInputHandlerRP \n");
232 AliESDEvent *ESD = hand->GetEvent();
234 printf("No AliESDEvent \n");
239 TTree * treeRP = hand->GetTreeR("ITS");
241 //AliWarning("No ITS RecPoints tree ");
247 // ESD related histograms
249 const AliESDVertex *vertex = ESD->GetPrimaryVertexSPD();
250 const AliMultiplicity *mult = ESD->GetMultiplicity();
254 if(!vertex->GetStatus()) return;
255 if(vertex->GetNContributors() < 1) return;
257 ((TH1I*)fOutput->At(18))->Fill(0);
259 ((TH1I*)fOutput->At(11))->Fill(mult->GetNumberOfTracklets());
260 UInt_t bc = (UInt_t)ESD->GetBunchCrossNumber();
261 for(Int_t iChipKey=0; iChipKey < 1200; iChipKey++){
262 if(mult->TestFiredChipMap(iChipKey)) {
263 ((TH1F*)fOutput->At(5))->Fill(iChipKey);
264 if(bc>0)((TH2F*)fOutput->At(9))->Fill(iChipKey,bc%4);
266 if(mult->TestFastOrFiredChips(iChipKey)) ((TH1F*)fOutput->At(6))->Fill(iChipKey);
267 if(mult->TestFastOrFiredChips(iChipKey) && mult->TestFiredChipMap(iChipKey)) {
268 ((TH1F*)fOutput->At(7))->Fill(iChipKey);
269 if(bc>0) ((TH2F*)fOutput->At(8))->Fill(iChipKey,bc%4);
272 if(mult->TestFastOrFiredChips(iChipKey) && !mult->TestFiredChipMap(iChipKey)) ((TH1F*)fOutput->At(10))->Fill(iChipKey);
276 Double_t vtxPos[3] = {999, 999, 999};
278 vertex->GetXYZ(vtxPos);
280 ((TH1F*)fOutput->At(15))->Fill(vtxPos[2]);
282 ((TH2F*)fOutput->At(16))->Fill(mult->GetNumberOfITSClusters(0),mult->GetNumberOfTracklets());
283 ((TH2F*)fOutput->At(17))->Fill(mult->GetNumberOfITSClusters(1),mult->GetNumberOfTracklets());
285 for(Int_t iTracklet =0; iTracklet < mult->GetNumberOfTracklets(); iTracklet++){
287 Float_t phiTr= mult->GetPhi(iTracklet);
288 Float_t etaTr =mult->GetEta(iTracklet);
290 ((TH2F*)fOutput->At(12))->Fill(etaTr,phiTr);
293 Float_t z = vtxPos[2] + 3.9 / TMath::Tan(2 * TMath::ATan(TMath::Exp(- etaTr)));
294 if(z>0) ((TH1F*)fOutput->At(13))->Fill(phiTr);
295 else ((TH1F*)fOutput->At(14))->Fill(phiTr);
303 TClonesArray statITSrec("AliITSRecPoint");
304 TClonesArray *ITSCluster = &statITSrec;
305 TBranch* branch=treeRP->GetBranch("ITSRecPoints");
307 printf("NO treeRP branch available. Exiting...\n");
311 branch->SetAddress(&ITSCluster);
313 for(Int_t iMod=0;iMod<240;iMod++){
314 branch->GetEvent(iMod);
315 Int_t nrecp = statITSrec.GetEntries();
316 if(nrecp>0) ((TH1F*)fOutput->At(2))->Fill(iMod,nrecp);
318 for(Int_t irec=0;irec<nrecp;irec++) {
319 AliITSRecPoint *recp = (AliITSRecPoint*)statITSrec.At(irec);
320 Int_t lay=recp->GetLayer();
323 // ---- Filling maps (local coordinates rearranged) -----
324 Float_t local[3]={-1,-1};
325 local[1]=recp->GetDetLocalX();
326 local[0]=recp->GetDetLocalZ();
327 Int_t eq = AliSPDUtils::GetOnlineEqIdFromOffline(iMod);
328 Int_t hs = AliSPDUtils::GetOnlineHSFromOffline(iMod);
330 fSegSPD->LocalToDet(0.5,local[0],row,col);
331 Int_t chip = AliSPDUtils::GetOnlineChipFromOffline(iMod,col);
333 Double_t locx, locz, equip;
334 Double_t corrlocz =0;
335 if(lay==0) corrlocz=local[0];
336 else corrlocz = -local[0];
337 // rearranging local geometry
341 if(chip<5) locz =corrlocz +8 +4;
342 else locz = corrlocz+4;
344 if(chip<5) locz = corrlocz-8-4;
345 else locz = corrlocz-4;
349 locx = (local[1]+1) + hs*2 +equip*4;
350 ((TH2D*)fOutput->At(0))->Fill(locz,locx);
353 locx = (local[1]+1) + (hs-2)*2 +equip*8;
354 ((TH2D*)fOutput->At(1))->Fill(locz,locx);
356 // ---- End Filling maps (local coordinates rearranged) -----
358 ((TH1F*)fOutput->At(3))->Fill(AliSPDUtils::GetOfflineChipKeyFromOnline(eq,hs,chip));
359 ((TH1F*)fOutput->At(4))->Fill(eq*60+hs*10+chip);
363 }// end if rec points are available
368 /* PostData(0) is taken care of by AliAnalysisTaskSE */
369 PostData(1,fOutput) ;
374 //___________________________________________________________________________
375 void AliAnalysisTaskSPD::Terminate(Option_t*)
377 // The Terminate() function is the last function to be called during
378 // a query. It always runs on the client, it can be used to present
379 // the results graphically or save the results to file.
381 AliAnalysisTaskSE::Terminate();
386 //___________________________________________________________________________
387 void AliAnalysisTaskSPD::LoadGeometryFromOCDB(){
388 //method to get the gGeomanager
389 // it is called at the CreatedOutputObject stage
390 // to comply with the CAF environment
391 AliCDBManager *man = AliCDBManager::Instance();
392 man->SetDefaultStorage(fOCDBLocation.Data());
394 AliCDBEntry* obj = man->Get(AliCDBPath("GRP", "Geometry", "Data"));
396 AliGeomManager::SetGeometry((TGeoManager*)obj->GetObject());
397 AliGeomManager::GetNalignable("ITS");
398 AliGeomManager::ApplyAlignObjsFromCDB("ITS");