df8f8070b4cb5c16771442b38e7500f79d0116c6
[u/mrichter/AliRoot.git] / PWG1 / ITS / AliAnalysisTaskSPD.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
17 //-----------------------------------------------------------------------
18 // Author : A. Mastroserio
19 //-----------------------------------------------------------------------
20
21
22 #ifndef ALIANALYSISTASKSPD_CXX
23 #define ALIANALYSISTASKSPD_CXX
24
25 #include "AliAnalysisTaskSPD.h"
26 #include "TClonesArray.h"
27 #include "TBranch.h"
28 #include "TH1I.h"
29 #include "TH2D.h"
30 #include "TFile.h"
31 #include "AliSPDUtils.h"
32 #include "AliESDEvent.h"
33 #include "TChain.h"
34 #include "AliLog.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() :
48     fSegSPD(0x0),
49     fOutput(0x0),
50     fRunNb(999),
51     fOCDBLocation("local://$ALICE_ROOT/OCDB")
52 {
53   //
54   //Default ctor
55   //
56 }
57 //___________________________________________________________________________
58 AliAnalysisTaskSPD::AliAnalysisTaskSPD(const Char_t* name) :
59   AliAnalysisTaskSE(name),
60   fSegSPD(0x0),
61   fOutput(0x0),
62   fRunNb(999),
63   fOCDBLocation("local://$ALICE_ROOT/OCDB")
64 {
65   //
66   // Constructor. Initialization of Inputs and Outputs
67   //
68   Info("AliAnalysisTaskSPD","Calling Constructor");
69
70   DefineOutput(1,TList::Class());
71   // 
72 }
73
74 //___________________________________________________________________________
75 AliAnalysisTaskSPD& AliAnalysisTaskSPD::operator=(const AliAnalysisTaskSPD& c) 
76 {
77   //
78   // Assignment operator
79   //
80   if (this!=&c) {
81     AliAnalysisTaskSE::operator=(c) ;
82     fSegSPD = c.fSegSPD ;
83     fOutput = c.fOutput ;
84     fRunNb = c.fRunNb;
85     fOCDBLocation = c.fOCDBLocation;
86   }
87   return *this;
88 }
89
90 //___________________________________________________________________________
91 AliAnalysisTaskSPD::AliAnalysisTaskSPD(const AliAnalysisTaskSPD& c) :
92   AliAnalysisTaskSE(c),
93   fSegSPD(c.fSegSPD),
94   fOutput(c.fOutput),
95   fRunNb(c.fRunNb),
96   fOCDBLocation(c.fOCDBLocation)
97 {
98   //
99   // Copy Constructor
100   //
101
102 }
103
104 //___________________________________________________________________________
105 AliAnalysisTaskSPD::~AliAnalysisTaskSPD() {
106   //
107   //destructor
108   //
109   Info("~AliAnalysisTaskSPD","Calling Destructor");
110   if (fSegSPD) delete fSegSPD ;
111    if (fOutput) {
112     delete fOutput;
113     fOutput = 0;
114   }
115 }
116 //___________________________________________________________________________
117 void AliAnalysisTaskSPD::UserCreateOutputObjects() {
118
119   Info("CreateOutputObjects","CreateOutputObjects of task %s", GetName());
120   if(fRunNb!=999){
121     LoadGeometryFromOCDB();
122   }
123   fSegSPD = new AliITSsegmentationSPD();
124
125   //slot #1
126   OpenFile(1);
127   fOutput = new TList();
128   //
129   // Booking rec points related histograms
130   //0
131   TH2D *hLocalMapL1 = new TH2D("hLocalMapL1"," Local coordinates  - Layer 1",660,-16.5,16.5,410,-0.5,40.5); // safe limits for local coordinates in a module : z = -4,4, x = -1,1;
132   hLocalMapL1->SetXTitle("direction Z [cm]");
133   hLocalMapL1->SetYTitle("direction X [cm]");
134   fOutput->AddLast(hLocalMapL1);
135   //1
136   TH2D *hLocalMapL2 = new TH2D("hLocalMapL2"," Local coordinates  - Layer 2",660,-16.5,16.5,810,-0.5,80.5);
137   hLocalMapL2->SetXTitle("direction Z [cm]");
138   hLocalMapL2->SetYTitle("direction X [cm]");
139   fOutput->AddLast(hLocalMapL2);
140   //2
141   TH1F *hClusterModYield = new TH1F("hClusterModYield","Cluster yield in modules",241,-0.5,240.5);
142   hClusterModYield->SetXTitle("module number");
143   fOutput->AddLast(hClusterModYield);
144   // 3
145   TH1F *hClusterYield = new TH1F("hClusterYield","Cluster yield per chip",1200,-0.5,1199.5) ;
146   hClusterYield->SetXTitle("chip key");
147   fOutput->AddLast(hClusterYield); 
148   //4
149   TH1F *hClusterYieldOnline = new TH1F("hClusterYieldOnline","Cluster yield per chip (online coord eq*60+hs*10+chip)",1200,-0.5,1199.5) ;
150   hClusterYieldOnline->SetXTitle("chip");
151   fOutput->AddLast(hClusterYieldOnline);
152   //5
153   TH1F *hFiredChip = new TH1F("hFiredChip","Fired chip (at least one cluster)",1200,-0.5,1199.5);
154   hFiredChip->SetXTitle("chip key");
155   fOutput->AddLast(hFiredChip);
156   //6
157   TH1F *hFOFiredChip = new TH1F("hFOFiredChip","FO Fired chip ",1200,-0.5,1199.5);
158   hFOFiredChip->SetXTitle("chip key");
159   fOutput->AddLast(hFOFiredChip);
160   //7
161   TH1F *hFOgood = new TH1F("hFOgood"," FO-good (at least one cluster) ",1200,-0.5,1199.5);
162   hFOgood->SetXTitle("chip key");
163   fOutput->AddLast(hFOgood);
164   //8
165   TH2F *hFOgoodPerBC = new TH2F("hFOgoodPerBCmod4"," FO-good signals in BCmod4 ",1200,-0.5,1199.5,4,-0.5,3.5);
166   hFOgoodPerBC->SetXTitle("chip key");
167   fOutput->AddLast(hFOgoodPerBC);
168   //9
169   TH2F *hFiredChipsPerBC = new TH2F("hFiredChipsPerBCmod4"," fired chips in BCmod4 ",1200,-0.5,1199.5,4,-0.5,3.5);
170   hFiredChipsPerBC->SetXTitle("chip key");
171   fOutput->AddLast(hFiredChipsPerBC);
172   //10
173   TH1F *hFOnoisy = new TH1F("hFOnoisy","FO-noisy (no cluster)",1200,-0.5,1199.5);
174   hFOnoisy->SetXTitle("chip key");
175   fOutput->AddLast(hFOnoisy);
176
177   //
178   // Booking ESD related histograms
179   //
180
181   // 11
182   TH1I *hTracklets = new TH1I("hNtracklets","Tracklet distribution",80,0,80);
183   hTracklets->SetXTitle("# Tracklets");
184   fOutput->AddLast(hTracklets);
185   //12
186   TH2F *hSPDphivsSPDeta= new TH2F("hSPDphivsSPDeta", "Tracklets - #varphi vs #eta",120,-3.,3,360,0.,2*TMath::Pi());
187   hSPDphivsSPDeta->GetXaxis()->SetTitle("Pseudorapidity #eta");
188   hSPDphivsSPDeta->GetYaxis()->SetTitle("#varphi [rad]");
189   fOutput->AddLast(hSPDphivsSPDeta);
190   //13
191   TH1F *hSPDphiZpos= new TH1F("hSPDphiZpos", "Tracklets - #varphi (Z>0)",360,0.,2*TMath::Pi());
192   hSPDphiZpos->SetXTitle("#varphi [rad]");
193   fOutput->AddLast(hSPDphiZpos);
194   //14
195   TH1F *hSPDphiZneg= new TH1F("hSPDphiZneg", "Tracklets - #varphi (Z<0)",360,0.,2*TMath::Pi());
196   hSPDphiZneg->SetXTitle("#varphi [rad]");
197   fOutput->AddLast(hSPDphiZneg);
198   //15
199   TH1F *hVertexZ = new TH1F("hVertexZ","Vertex Z distribution",300,-15,15);
200   hVertexZ->SetXTitle("Z Vertex [cm]");
201   fOutput->AddLast(hVertexZ); 
202   //16
203   TH1I *hEventsProcessed = new TH1I("hEventsProcessed","Number of processed events",1,0,1) ;
204   fOutput->AddLast(hEventsProcessed);
205 }
206 //_________________________________________________
207 void AliAnalysisTaskSPD::UserExec(Option_t *)
208 {
209   //
210   // Main loop function
211   //
212   Info("UserExec","") ;
213
214   AliESDInputHandlerRP *hand = dynamic_cast<AliESDInputHandlerRP*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
215   if(!hand) {
216     printf("No AliESDInputHandlerRP \n");  
217     return;
218   }
219
220   AliESDEvent *ESD = hand->GetEvent();
221   if(!ESD) {
222     printf("No AliESDEvent \n");
223     return;
224   }
225
226   Bool_t recP = kTRUE;
227   TTree * treeRP = hand->GetTreeR("ITS");
228   if(!treeRP) {
229     //AliWarning("No ITS RecPoints tree ");
230     recP=kFALSE;
231   }
232
233   
234
235   // ESD related histograms
236  
237    const AliESDVertex *vertex = ESD->GetPrimaryVertexSPD();
238    const AliMultiplicity *mult = ESD->GetMultiplicity();
239  
240    // Event selection
241    if(!vertex) return;
242    if(!vertex->GetStatus()) return;
243    if(vertex->GetNContributors() < 1) return;
244
245   ((TH1I*)fOutput->At(16))->Fill(0);
246    
247   ((TH1I*)fOutput->At(11))->Fill(mult->GetNumberOfTracklets());
248   UInt_t bc = (UInt_t)ESD->GetBunchCrossNumber();
249   for(Int_t iChipKey=0; iChipKey < 1200; iChipKey++){
250     if(mult->TestFiredChipMap(iChipKey)) {
251      ((TH1F*)fOutput->At(5))->Fill(iChipKey);
252      if(bc>0)((TH2F*)fOutput->At(9))->Fill(iChipKey,bc%4);   
253      }
254     if(mult->TestFastOrFiredChips(iChipKey)) ((TH1F*)fOutput->At(6))->Fill(iChipKey);
255     if(mult->TestFastOrFiredChips(iChipKey) && mult->TestFiredChipMap(iChipKey)) {
256       ((TH1F*)fOutput->At(7))->Fill(iChipKey);
257       if(bc>0) ((TH2F*)fOutput->At(8))->Fill(iChipKey,bc%4);
258      
259     }
260     if(mult->TestFastOrFiredChips(iChipKey) && !mult->TestFiredChipMap(iChipKey)) ((TH1F*)fOutput->At(10))->Fill(iChipKey);
261   }
262   
263   
264   Double_t vtxPos[3] = {999, 999, 999};
265   if(vertex){
266     vertex->GetXYZ(vtxPos);
267  
268     ((TH1F*)fOutput->At(15))->Fill(vtxPos[2]);
269  
270     for(Int_t iTracklet =0; iTracklet < mult->GetNumberOfTracklets(); iTracklet++){
271
272       Float_t phiTr= mult->GetPhi(iTracklet);
273       Float_t etaTr =mult->GetEta(iTracklet);
274
275       ((TH2F*)fOutput->At(12))->Fill(etaTr,phiTr);
276
277       // Z pos or Z neg
278       Float_t z = vtxPos[2] + 3.9 / TMath::Tan(2 * TMath::ATan(TMath::Exp(- etaTr)));
279       if(z>0) ((TH1F*)fOutput->At(13))->Fill(phiTr);
280       else ((TH1F*)fOutput->At(14))->Fill(phiTr);
281     }
282   } 
283   
284   
285   if(recP){
286   // RecPoints info
287
288   TClonesArray statITSrec("AliITSRecPoint");
289   TClonesArray *ITSCluster = &statITSrec;
290   TBranch* branch=treeRP->GetBranch("ITSRecPoints");
291   if(!branch) {
292     printf("NO treeRP branch available. Exiting...\n");
293     return;
294   }
295
296   branch->SetAddress(&ITSCluster);
297
298   for(Int_t iMod=0;iMod<240;iMod++){
299     branch->GetEvent(iMod);
300     Int_t nrecp = statITSrec.GetEntries();
301     if(nrecp>0) ((TH1F*)fOutput->At(2))->Fill(iMod,nrecp);
302   
303     for(Int_t irec=0;irec<nrecp;irec++) {
304       AliITSRecPoint *recp = (AliITSRecPoint*)statITSrec.At(irec);
305       Int_t lay=recp->GetLayer();
306       if(lay>1) continue;
307
308       // ----  Filling maps (local coordinates rearranged) -----
309       Float_t local[3]={-1,-1};
310       local[1]=recp->GetDetLocalX();
311       local[0]=recp->GetDetLocalZ();
312       Int_t eq = AliSPDUtils::GetOnlineEqIdFromOffline(iMod);
313       Int_t hs = AliSPDUtils::GetOnlineHSFromOffline(iMod);
314       Int_t row, col;
315       fSegSPD->LocalToDet(0.5,local[0],row,col);
316       Int_t chip = AliSPDUtils::GetOnlineChipFromOffline(iMod,col);
317      
318       Double_t locx, locz, equip;
319       Double_t corrlocz =0;
320       if(lay==0) corrlocz=local[0];
321       else corrlocz = -local[0];
322       // rearranging local geometry
323       if(eq<10) equip=eq;
324       else equip=eq-10;
325       if(eq<10){
326         if(chip<5) locz =corrlocz +8 +4;
327         else locz = corrlocz+4;
328       } else {
329         if(chip<5) locz = corrlocz-8-4;
330         else locz = corrlocz-4;
331       }
332       // filling maps
333       if(lay==0){
334         locx = (local[1]+1) + hs*2 +equip*4;
335         ((TH2D*)fOutput->At(0))->Fill(locz,locx);
336         
337       } else {
338         locx = (local[1]+1) + (hs-2)*2 +equip*8;
339         ((TH2D*)fOutput->At(1))->Fill(locz,locx);
340       }
341       // ---- End Filling maps (local coordinates rearranged) -----
342     
343       ((TH1F*)fOutput->At(3))->Fill(AliSPDUtils::GetOfflineChipKeyFromOnline(eq,hs,chip));
344       ((TH1F*)fOutput->At(4))->Fill(eq*60+hs*10+chip);
345    
346     }
347   }
348  }// end if rec points are available
349  
350   
351   
352   
353   /* PostData(0) is taken care of by AliAnalysisTaskSE */
354   PostData(1,fOutput) ;
355
356 }
357
358
359 //___________________________________________________________________________
360 void AliAnalysisTaskSPD::Terminate(Option_t*)
361 {
362   // The Terminate() function is the last function to be called during
363   // a query. It always runs on the client, it can be used to present
364   // the results graphically or save the results to file.
365
366   Info("Terminate","");
367   AliAnalysisTaskSE::Terminate();
368 }
369
370
371
372 //___________________________________________________________________________
373 void AliAnalysisTaskSPD::LoadGeometryFromOCDB(){
374   //method to get the gGeomanager
375   // it is called at the CreatedOutputObject stage
376   // to comply with the CAF environment
377   AliCDBManager *man = AliCDBManager::Instance();
378   man->SetDefaultStorage(fOCDBLocation.Data());
379   man->SetRun(fRunNb);  
380   AliCDBEntry* obj = man->Get(AliCDBPath("GRP", "Geometry", "Data"));
381   AliGeomManager::SetGeometry((TGeoManager*)obj->GetObject());
382   AliGeomManager::GetNalignable("ITS");
383   AliGeomManager::ApplyAlignObjsFromCDB("ITS"); 
384 }
385
386
387 #endif