]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerRec.cxx
Minor issues solved. Additional printouts for debugging purposes (H. Tydesjo)
[u/mrichter/AliRoot.git] / TPC / AliTPCQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 /* $Id: $ */
18
19 /*
20   Based on AliPHOSQADataMaker
21   Produces the data needed to calculate the quality assurance. 
22   All data must be mergeable objects.
23   P. Christiansen, Lund, January 2008
24 */
25
26 /*
27   Implementation:
28
29   We have chosen to have the histograms as non-persistent meber to
30   allow better debugging. In the copy constructor we then have to
31   assign the pointers to the existing histograms in the copied
32   list. This have been implemented but not tested.
33
34   For the QA of the RAW data we use the class, AliTPCdataQA, from the
35   existing TPC Calibration framework (which is more advanced than the
36   standard QA framework) and extract the histograms at the end. This
37   has been tested with zero-suppressed data. The Analyse method of the
38   AliTPCdataQA class is called in the method, EndOfDetectorCycle, and
39   there also: 1d histogram(s) are projected and added to the QA list.
40 */
41
42 /*
43   TODO:
44   Sumw2 for RAW histogram(s)?
45   RecPoints and ESD could have many more histograms
46 */
47
48 #include "AliTPCQADataMakerRec.h"
49
50 // --- ROOT system ---
51 #include <TClonesArray.h>
52
53 // --- Standard library ---
54
55 // --- AliRoot header files ---
56 #include "AliQAChecker.h"
57 #include "AliESDEvent.h"
58 #include "AliESDtrack.h"
59 #include "AliLog.h"
60 #include "AliTPCCalPad.h"
61 #include "AliTPCClustersRow.h"
62 #include "AliTPCclusterMI.h"
63
64 ClassImp(AliTPCQADataMakerRec)
65
66 //____________________________________________________________________________ 
67 AliTPCQADataMakerRec::AliTPCQADataMakerRec() : 
68   AliQADataMakerRec(AliQA::GetDetName(AliQA::kTPC), 
69                     "TPC Rec Quality Assurance Data Maker"),
70   fTPCdataQA(0),
71   fHistESDclusters(0),fHistESDratio(0), fHistESDpt(0),
72   fHistRawsOccupancy(0),
73   fHistRecPointsQmaxShort(0), fHistRecPointsQmaxMedium(0), 
74   fHistRecPointsQmaxLong(0), fHistRecPointsQShort(0), 
75   fHistRecPointsQMedium(0), fHistRecPointsQLong(0),
76   fHistRecPointsRow(0)
77 {
78   // ctor
79 }
80
81 //____________________________________________________________________________ 
82 AliTPCQADataMakerRec::AliTPCQADataMakerRec(const AliTPCQADataMakerRec& qadm) :
83   AliQADataMakerRec()
84 {
85   //copy ctor 
86   // Does not copy the calibration object, instead InitRaws have to be
87   // called again
88   SetName((const char*)qadm.GetName()) ; 
89   SetTitle((const char*)qadm.GetTitle()); 
90
91   //
92   // Associate class histogram objects to the copies in the list
93   // Could also be done with the indexes
94   //
95   fHistESDclusters   = (TH1F*)fESDsQAList->FindObject("hESDclusters");
96   fHistESDratio      = (TH1F*)fESDsQAList->FindObject("hESDratio");
97   fHistESDpt         = (TH1F*)fESDsQAList->FindObject("hESDpt");
98
99   fHistRawsOccupancy = (TH1F*)fRawsQAList->FindObject("hRawsOccupancy");
100
101   fHistRecPointsQmaxShort  = 
102     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxShort");
103   fHistRecPointsQmaxMedium = 
104     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxMedium");
105   fHistRecPointsQmaxLong   = 
106     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxLong"); 
107   fHistRecPointsQShort     = 
108     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQShort");
109   fHistRecPointsQMedium    = 
110     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQMedium");
111   fHistRecPointsQLong      = 
112     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQLong");
113   fHistRecPointsRow      = 
114     (TH1F*)fRecPointsQAList->FindObject("hRecPointsRow");
115 }
116
117 //__________________________________________________________________
118 AliTPCQADataMakerRec& AliTPCQADataMakerRec::operator = (const AliTPCQADataMakerRec& qadm )
119 {
120   // Equal operator.
121   this->~AliTPCQADataMakerRec();
122   new(this) AliTPCQADataMakerRec(qadm);
123   return *this;
124 }
125  
126 //____________________________________________________________________________ 
127 void AliTPCQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray * list)
128 {
129   //Detector specific actions at end of cycle
130
131   if(fTPCdataQA) { // do the final step of the QA for Raw data
132
133     fTPCdataQA->Analyse(); // 31/1-08 Analyse is now protected against
134                            //         RAW data files with no TPC data
135     
136     // get the histograms and add them to the output
137     // 31/8-08 Histogram is only added if the Calibration class 
138     //         receives TPC data 
139     if(fTPCdataQA->GetNoThreshold()) { 
140       fHistRawsOccupancy = fTPCdataQA->GetNoThreshold()->MakeHisto1D(0, 1, -1);
141       Add2RawsList(fHistRawsOccupancy, 0);
142     }
143   }
144
145   AliQAChecker::Instance()->Run(AliQA::kTPC, task, list) ;  
146 }
147
148 //____________________________________________________________________________ 
149 void AliTPCQADataMakerRec::InitESDs()
150 {
151   //create ESDs histograms in ESDs subdir  
152   fHistESDclusters = 
153     new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
154              160, 0, 160);
155   fHistESDclusters->Sumw2();
156   Add2ESDsList(fHistESDclusters, 0);
157
158   fHistESDratio = 
159     new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
160              100, 0, 1);
161   fHistESDratio->Sumw2();
162   Add2ESDsList(fHistESDratio, 1);
163   
164   fHistESDpt = 
165     new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
166              50, 0, 5);
167   fHistESDpt->Sumw2();
168   Add2ESDsList(fHistESDpt, 2);
169 }
170
171 //____________________________________________________________________________ 
172 void AliTPCQADataMakerRec::InitRaws()
173 {
174   fTPCdataQA = new AliTPCdataQA();
175   fTPCdataQA->SetRangeTime(0, 999); // take all 1000 time bins
176 }
177
178 //____________________________________________________________________________ 
179 void AliTPCQADataMakerRec::InitRecPoints()
180 {
181   fHistRecPointsQmaxShort = 
182     new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
183              200, 0, 1000);
184   fHistRecPointsQmaxShort->Sumw2();
185   Add2RecPointsList(fHistRecPointsQmaxShort, 0);
186
187   fHistRecPointsQmaxMedium = 
188     new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
189              200, 0, 1000);
190   fHistRecPointsQmaxMedium->Sumw2();
191   Add2RecPointsList(fHistRecPointsQmaxMedium, 1);
192
193   fHistRecPointsQmaxLong = 
194     new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
195              200, 0, 1000);
196   fHistRecPointsQmaxLong->Sumw2();
197   Add2RecPointsList(fHistRecPointsQmaxLong, 2);
198
199   fHistRecPointsQShort = 
200     new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
201              200, 0, 5000);
202   fHistRecPointsQShort->Sumw2();
203   Add2RecPointsList(fHistRecPointsQShort, 3);
204
205   fHistRecPointsQMedium = 
206     new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
207              200, 0, 5000);
208   fHistRecPointsQMedium->Sumw2();
209   Add2RecPointsList(fHistRecPointsQMedium, 4);
210
211   fHistRecPointsQLong = 
212     new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
213              200, 0, 5000);
214   fHistRecPointsQLong->Sumw2();
215   Add2RecPointsList(fHistRecPointsQLong, 5);
216
217   fHistRecPointsRow = 
218     new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
219              159, 0, 159);
220   fHistRecPointsRow->Sumw2();
221   Add2RecPointsList(fHistRecPointsRow, 6);
222 }
223
224 //____________________________________________________________________________
225 void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
226 {
227   // make QA data from ESDs
228   
229   const Int_t nESDTracks = esd->GetNumberOfTracks();
230   Int_t nTPCtracks = 0; 
231   for(Int_t i = 0; i < nESDTracks; i++) {
232     
233     AliESDtrack * track = esd->GetTrack(i);
234     
235     if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
236       continue;
237     
238     nTPCtracks++;
239     
240     Int_t nTPCclusters         = track->GetTPCNcls();
241     Int_t nTPCclustersFindable = track->GetTPCNclsF();
242     
243     fHistESDclusters->Fill(nTPCclusters);
244     fHistESDratio->Fill(Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
245     fHistESDpt->Fill(track->Pt()); 
246   }
247 }
248
249 //____________________________________________________________________________
250 void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
251 {
252   //
253   // To make QA for the RAW data we use the TPC Calibration framework 
254   // to handle the data and then in the end extract the data
255   //
256   fTPCdataQA->ProcessEvent(rawReader);  
257 }
258
259 //____________________________________________________________________________
260 void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
261 {
262   AliTPCClustersRow *clrow = new AliTPCClustersRow();
263   clrow->SetClass("AliTPCclusterMI");
264   clrow->SetArray(0);
265   clrow->GetArray()->ExpandCreateFast(10000);
266
267   TBranch* branch = recTree->GetBranch("Segment");
268   branch->SetAddress(&clrow);
269
270   const Int_t nEntries = Int_t(recTree->GetEntries());
271   for (Int_t i = 0; i < nEntries; i++) {
272     
273     branch->GetEntry(i);
274     
275     const Int_t nClusters = clrow->GetArray()->GetEntriesFast();
276     for (Int_t icl=0; icl < nClusters; icl++){
277       
278       AliTPCclusterMI* cluster = 
279         (AliTPCclusterMI*)clrow->GetArray()->At(icl);
280       
281       Float_t Qmax = cluster->GetMax();
282       Float_t Q    = cluster->GetQ();
283       Int_t   row  = cluster->GetRow();
284
285       if(cluster->GetDetector()<36) { // IROC (short pads)
286
287         fHistRecPointsQmaxShort->Fill(Qmax);
288         fHistRecPointsQShort->Fill(Q);
289       } else { // OROC (medium and long pads)
290         row += 63;
291         if(cluster->GetRow()<64) { // medium pads
292
293           fHistRecPointsQmaxMedium->Fill(Qmax);
294           fHistRecPointsQMedium->Fill(Q);
295         } else { // long pads
296
297           fHistRecPointsQmaxLong->Fill(Qmax);
298           fHistRecPointsQLong->Fill(Q);
299         }
300       }
301       
302       fHistRecPointsRow->Fill(row);
303     } // end loop over clusters
304   } // end loop over tree
305
306   delete clrow;
307 }