]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/PHOSTasks/PHOS_TriggerQA/AliAnalysisTaskPHOSTriggerQA.cxx
PHOS trigger QA is added (B.Polishchuk)
[u/mrichter/AliRoot.git] / PWGGA / PHOSTasks / PHOS_TriggerQA / AliAnalysisTaskPHOSTriggerQA.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 #include "AliAnalysisManager.h"
17 #include "AliAnalysisTaskSE.h"
18 #include "AliAnalysisTaskPHOSTriggerQA.h"
19 #include "AliESDCaloCluster.h"
20 #include "AliPHOSGeometry.h"
21 #include "AliESDEvent.h"
22 #include "AliESDCaloCells.h"
23 #include "AliLog.h"
24 #include "TObjArray.h"
25 #include "TList.h"
26 #include "TH1.h"
27 #include "TH2.h"
28
29 // QA of PHOS Trigger data.
30 //...
31 // Author: Boris Polishchuk 
32 // Date  : 06.02.2012
33
34 ClassImp(AliAnalysisTaskPHOSTriggerQA)
35
36 //________________________________________________________________________
37 AliAnalysisTaskPHOSTriggerQA::AliAnalysisTaskPHOSTriggerQA() : AliAnalysisTaskSE(),
38   fOutputContainer(0),fPHOSGeo(0),fEventCounter(0)
39 {
40   //Default constructor.  
41   // Initialize the PHOS geometry 
42   fPHOSGeo = AliPHOSGeometry::GetInstance("IHEP") ;
43   
44 }
45
46 //________________________________________________________________________
47 AliAnalysisTaskPHOSTriggerQA::AliAnalysisTaskPHOSTriggerQA(const char *name) 
48 : AliAnalysisTaskSE(name),
49   fOutputContainer(0),fPHOSGeo(0),fEventCounter(0)
50 {
51   
52   // Output slots #0 write into a TH1 container
53   DefineOutput(1,TList::Class());
54
55   // Initialize the PHOS geometry
56   fPHOSGeo = AliPHOSGeometry::GetInstance("IHEP") ;
57
58 }
59
60 //________________________________________________________________________
61 void AliAnalysisTaskPHOSTriggerQA::UserCreateOutputObjects()
62 {
63   // Create histograms
64   // Called once
65   
66   // ESD histograms
67   if(fOutputContainer != NULL){
68     delete fOutputContainer;
69   }
70
71   fOutputContainer = new TList();
72   fOutputContainer->SetOwner(kTRUE);
73
74   //Bin 1: total number of processed events.
75   //Bin 2: number of events contained PHOS trigger digits.
76   fOutputContainer->Add(new TH1F("hNev","Number of events",10,0.,10.));
77   
78   char key[55],titl[55];
79   Int_t nCols = 56, nRows = 64;
80   
81   Int_t    nPtPhot  = 1000 ;
82   Double_t ptPhotMax =  100. ;
83
84   Int_t nTrMax = 1200;
85   Float_t trMax = 600.;
86
87   fOutputContainer->Add(new TH1F("hNtr","Number of fired 4x4 regions",nTrMax,0.,trMax));
88
89   for(Int_t sm=1; sm<4; sm++) {
90
91     snprintf(key,55,"hNtrSM%d",sm);
92     snprintf(titl,55,"Number of fired 4x4 regions in SM%d",sm);
93     fOutputContainer->Add(new TH1F(key,titl,nTrMax/3,0.,trMax/3));
94     
95     snprintf(key,55,"h4x4SM%d",sm);
96     snprintf(titl,55,"SM%d 4x4 occupancy",sm);
97     fOutputContainer->Add(new TH2F(key,titl,nRows,0.,nRows,nCols,0.,nCols));
98
99     snprintf(key,55,"hCluSM%d",sm);
100     snprintf(titl,55,"SM%d cluster occupancy",sm);
101     fOutputContainer->Add(new TH2F(key,titl,nRows,0.,nRows,nCols,0.,nCols));
102  
103     snprintf(key,55,"hCluTSM%d",sm);
104     snprintf(titl,55,"SM%d triggered cluster occupancy",sm);
105     fOutputContainer->Add(new TH2F(key,titl,nRows,0.,nRows,nCols,0.,nCols));
106    
107     snprintf(key,55,"hPhotAllSM%d",sm);
108     snprintf(titl,55,"SM%d cluster energy",sm);
109     fOutputContainer->Add(new TH1F(key,titl,nPtPhot,0.,ptPhotMax));
110
111     snprintf(key,55,"hPhotTrigSM%d",sm);
112     snprintf(titl,55,"SM%d triggered cluster energy",sm);
113     fOutputContainer->Add(new TH1F(key,titl,nPtPhot,0.,ptPhotMax));
114   }
115   
116   PostData(1, fOutputContainer);
117   
118 }
119
120 //________________________________________________________________________
121 void AliAnalysisTaskPHOSTriggerQA::UserExec(Option_t *) 
122 {
123   // Main loop, called for each event
124   // Analyze ESD/AOD  
125   
126   AliESDEvent *event = dynamic_cast<AliESDEvent*>(InputEvent());
127   
128   if (!event) {
129     Printf("ERROR: Could not retrieve event");
130     PostData(1, fOutputContainer);
131     return;
132   }
133   
134   FillHistogram("hNev",0.); // all events
135   fEventCounter++;
136   
137   AliESDCaloTrigger* trgESD = event->GetCaloTrigger("PHOS");
138   trgESD->Reset();
139   
140   if (!trgESD->GetEntries()) {
141     PostData(1, fOutputContainer);
142     return;
143   }
144   
145   FillHistogram("hNev",1.); // triggered events
146   FillHistogram("hNtr",trgESD->GetEntries());
147
148   TString trigClasses = event->GetFiredTriggerClasses();
149   printf("\nEvent %d: %d non-zero trigger digits %s\n",
150          fEventCounter,trgESD->GetEntries(),trigClasses.Data());
151
152   // Get PHOS rotation matrices from ESD and set them to the PHOS geometry
153   char key[55] ;  
154   
155   if(fEventCounter == 0) {
156     for(Int_t mod=0; mod<5; mod++) {
157       if(!event->GetPHOSMatrix(mod)) continue;
158       fPHOSGeo->SetMisalMatrix(event->GetPHOSMatrix(mod),mod) ;
159     }
160   }
161   
162   Int_t multClu = event->GetNumberOfCaloClusters();
163   AliESDCaloCells *phsCells = event->GetPHOSCells();
164  
165   Int_t inPHOS[3] = {};
166
167   //Loop over 4x4 fired regions
168   while(trgESD->Next()) {
169
170     Int_t tmod,tabsId; // "Online" module number, bottom-left 4x4 edge cell absId
171     trgESD->GetPosition(tmod,tabsId);
172     
173     Int_t trelid[4] ;
174     fPHOSGeo->AbsToRelNumbering(tabsId,trelid);
175
176     snprintf(key,55,"h4x4SM%d",trelid[0]);
177     FillHistogram(key,trelid[2]-1,trelid[3]-1);
178     
179     inPHOS[trelid[0]-1]++;
180
181     for (Int_t i=0; i<multClu; i++) {
182       
183       AliESDCaloCluster *c1 = event->GetCaloCluster(i);
184       if(!c1->IsPHOS()) continue;
185       
186       Int_t maxId, relid[4];
187       MaxEnergyCellPos(phsCells,c1,maxId);
188       
189       fPHOSGeo->AbsToRelNumbering(maxId, relid);
190       snprintf(key,55,"hPhotAllSM%d",relid[0]);
191       FillHistogram(key,c1->E());
192       
193       snprintf(key,55,"hCluSM%d",relid[0]);
194       FillHistogram(key,relid[2]-1,relid[3]-1);
195       
196       if( Matched(trelid,relid) ) {
197
198         snprintf(key,55,"hPhotTrigSM%d",relid[0]);
199         FillHistogram(key,c1->E());
200
201         snprintf(key,55,"hCluTSM%d",relid[0]);
202         FillHistogram(key,relid[2]-1,relid[3]-1);
203         
204         continue;
205       }
206       
207     }
208   } //while(trgESD->Next())
209   
210   for(Int_t sm=1; sm<4; sm++) {
211     snprintf(key,55,"hNtrSM%d",sm);
212     if(inPHOS[sm-1]) FillHistogram(key,inPHOS[sm-1]);    
213   }
214   
215   // Post output data.
216   PostData(1, fOutputContainer);
217   
218 }
219
220 //_____________________________________________________________________________
221 void AliAnalysisTaskPHOSTriggerQA::FillHistogram(const char * key,Double_t x)const{
222   //FillHistogram
223   TH1I * tmpI = dynamic_cast<TH1I*>(fOutputContainer->FindObject(key)) ;
224   if(tmpI){
225     tmpI->Fill(x) ;
226     return ;
227   }
228   TH1F * tmpF = dynamic_cast<TH1F*>(fOutputContainer->FindObject(key)) ;
229   if(tmpF){
230     tmpF->Fill(x) ;
231     return ;
232   }
233   TH1D * tmpD = dynamic_cast<TH1D*>(fOutputContainer->FindObject(key)) ;
234   if(tmpD){
235     tmpD->Fill(x) ;
236     return ;
237   }
238   AliInfo(Form("can not find histogram <%s> ",key)) ;
239 }
240 //_____________________________________________________________________________
241 void AliAnalysisTaskPHOSTriggerQA::FillHistogram(const char * key,Double_t x,Double_t y)const{
242   //FillHistogram
243   TObject * tmp = fOutputContainer->FindObject(key) ;
244   if(!tmp){
245     AliInfo(Form("can not find histogram <%s> ",key)) ;
246     return ;
247   }
248   if(tmp->IsA() == TClass::GetClass("TH1F")){
249     ((TH1F*)tmp)->Fill(x,y) ;
250     return ;
251   }
252   if(tmp->IsA() == TClass::GetClass("TH2F")){
253     ((TH2F*)tmp)->Fill(x,y) ;
254     return ;
255   }
256   AliError(Form("Calling FillHistogram with 2 parameters for histo <%s> of type %s",key,tmp->IsA()->GetName())) ;
257 }
258
259 //_____________________________________________________________________________
260 void AliAnalysisTaskPHOSTriggerQA::MaxEnergyCellPos(AliESDCaloCells *cells, AliESDCaloCluster* clu, Int_t& maxId)
261 {  
262   Double_t eMax = -111;
263   
264   for (Int_t iDig=0; iDig< clu->GetNCells(); iDig++) {
265     Int_t cellAbsId = clu->GetCellAbsId(iDig);
266     Double_t eCell = cells->GetCellAmplitude(cellAbsId)*clu->GetCellAmplitudeFraction(iDig);
267     if(eCell>eMax)  { 
268       eMax = eCell; 
269       maxId = cellAbsId;
270     }
271   }
272   
273 }
274
275 //_____________________________________________________________________________
276 Bool_t AliAnalysisTaskPHOSTriggerQA::Matched(Int_t *trig_relid, Int_t *cluster_relid)
277 {
278   //Returns kTRUE if cluster position coincides with 4x4 position.
279
280   if( trig_relid[0] != cluster_relid[0] )            return kFALSE; // different modules!
281   if( TMath::Abs(trig_relid[2]-cluster_relid[2])>3 ) return kFALSE; // X-distance too large! 
282   if( TMath::Abs(trig_relid[3]-cluster_relid[3])>3 ) return kFALSE; // Z-distance too large!
283
284   return kTRUE;
285 }