]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerRec.cxx
Merging THbtp and HBTP in one library. Comiplation on Windows/Cygwin
[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   fTPCdataQA(0),
85   fHistESDclusters(0),  //! Clusters per ESD track
86   fHistESDratio(0),     //! Ratio of clusters to findables
87   fHistESDpt(0),        //! Pt spectrum
88   
89   fHistRawsOccupancy(0),//! Pad occupancy (1 entry per pad)
90   
91   fHistRecPointsQmaxShort(0), //! Qmax (short pads)
92   fHistRecPointsQmaxMedium(0),//! Qmax (medium pads)
93   fHistRecPointsQmaxLong(0),  //! Qmax (long pads)
94   fHistRecPointsQShort(0),    //! Q (short pads)
95   fHistRecPointsQMedium(0),   //! Q (medium pads)
96   fHistRecPointsQLong(0),     //! Q (long pads)
97   fHistRecPointsRow(0)       //! Row distribution
98 {
99   //copy ctor 
100   // Does not copy the calibration object, instead InitRaws have to be
101   // called again
102   SetName((const char*)qadm.GetName()) ; 
103   SetTitle((const char*)qadm.GetTitle()); 
104
105   //
106   // Associate class histogram objects to the copies in the list
107   // Could also be done with the indexes
108   //
109   fHistESDclusters   = (TH1F*)fESDsQAList->FindObject("hESDclusters");
110   fHistESDratio      = (TH1F*)fESDsQAList->FindObject("hESDratio");
111   fHistESDpt         = (TH1F*)fESDsQAList->FindObject("hESDpt");
112
113   fHistRawsOccupancy = (TH1F*)fRawsQAList->FindObject("hRawsOccupancy");
114
115   fHistRecPointsQmaxShort  = 
116     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxShort");
117   fHistRecPointsQmaxMedium = 
118     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxMedium");
119   fHistRecPointsQmaxLong   = 
120     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxLong"); 
121   fHistRecPointsQShort     = 
122     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQShort");
123   fHistRecPointsQMedium    = 
124     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQMedium");
125   fHistRecPointsQLong      = 
126     (TH1F*)fRecPointsQAList->FindObject("hRecPointsQLong");
127   fHistRecPointsRow      = 
128     (TH1F*)fRecPointsQAList->FindObject("hRecPointsRow");
129 }
130
131 //__________________________________________________________________
132 AliTPCQADataMakerRec& AliTPCQADataMakerRec::operator = (const AliTPCQADataMakerRec& qadm )
133 {
134   // Equal operator.
135   this->~AliTPCQADataMakerRec();
136   new(this) AliTPCQADataMakerRec(qadm);
137   return *this;
138 }
139  
140 //____________________________________________________________________________ 
141 void AliTPCQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
142 {
143   //Detector specific actions at end of cycle
144
145   if(fTPCdataQA) { // do the final step of the QA for Raw data
146
147     fTPCdataQA->Analyse(); // 31/1-08 Analyse is now protected against
148                            //         RAW data files with no TPC data
149     
150     //Add2RawsList(fTPCdataQA, 0);
151     // get the histograms and add them to the output
152     // 31/8-08 Histogram is only added if the Calibration class 
153     //         receives TPC data 
154     if(fTPCdataQA->GetNoThreshold()) { 
155       fHistRawsOccupancy = fTPCdataQA->GetNoThreshold()->MakeHisto1D(0, 1, -1);
156       //Add2RawsList(fHistRawsOccupancy, 1);
157     }
158   }
159
160   AliQAChecker::Instance()->Run(AliQA::kTPC, task, list) ;  
161 }
162
163 //____________________________________________________________________________ 
164 void AliTPCQADataMakerRec::InitESDs()
165 {
166   //create ESDs histograms in ESDs subdir  
167   fHistESDclusters = 
168     new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
169              160, 0, 160);
170   fHistESDclusters->Sumw2();
171   Add2ESDsList(fHistESDclusters, 0);
172
173   fHistESDratio = 
174     new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
175              100, 0, 1);
176   fHistESDratio->Sumw2();
177   Add2ESDsList(fHistESDratio, 1);
178   
179   fHistESDpt = 
180     new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
181              50, 0, 5);
182   fHistESDpt->Sumw2();
183   Add2ESDsList(fHistESDpt, 2);
184 }
185
186 //____________________________________________________________________________ 
187 void AliTPCQADataMakerRec::InitRaws()
188 {
189   //
190   // Adding the raw 
191   //
192   fTPCdataQA = new AliTPCdataQA();
193   fTPCdataQA->SetRangeTime(100, 920); // take all 1000 time bins 
194   Add2RawsList(fTPCdataQA, 0);
195 }
196
197 //____________________________________________________________________________ 
198 void AliTPCQADataMakerRec::InitRecPoints()
199 {
200   fHistRecPointsQmaxShort = 
201     new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
202              200, 0, 1000);
203   fHistRecPointsQmaxShort->Sumw2();
204   Add2RecPointsList(fHistRecPointsQmaxShort, 0);
205
206   fHistRecPointsQmaxMedium = 
207     new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
208              200, 0, 1000);
209   fHistRecPointsQmaxMedium->Sumw2();
210   Add2RecPointsList(fHistRecPointsQmaxMedium, 1);
211
212   fHistRecPointsQmaxLong = 
213     new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
214              200, 0, 1000);
215   fHistRecPointsQmaxLong->Sumw2();
216   Add2RecPointsList(fHistRecPointsQmaxLong, 2);
217
218   fHistRecPointsQShort = 
219     new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
220              200, 0, 5000);
221   fHistRecPointsQShort->Sumw2();
222   Add2RecPointsList(fHistRecPointsQShort, 3);
223
224   fHistRecPointsQMedium = 
225     new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
226              200, 0, 5000);
227   fHistRecPointsQMedium->Sumw2();
228   Add2RecPointsList(fHistRecPointsQMedium, 4);
229
230   fHistRecPointsQLong = 
231     new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
232              200, 0, 5000);
233   fHistRecPointsQLong->Sumw2();
234   Add2RecPointsList(fHistRecPointsQLong, 5);
235
236   fHistRecPointsRow = 
237     new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
238              159, 0, 159);
239   fHistRecPointsRow->Sumw2();
240   Add2RecPointsList(fHistRecPointsRow, 6);
241 }
242
243 //____________________________________________________________________________
244 void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
245 {
246   // make QA data from ESDs
247   
248   const Int_t nESDTracks = esd->GetNumberOfTracks();
249   Int_t nTPCtracks = 0; 
250   for(Int_t i = 0; i < nESDTracks; i++) {
251     
252     AliESDtrack * track = esd->GetTrack(i);
253     
254     if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
255       continue;
256     
257     nTPCtracks++;
258     
259     Int_t nTPCclusters         = track->GetTPCNcls();
260     Int_t nTPCclustersFindable = track->GetTPCNclsF();
261     
262     fHistESDclusters->Fill(nTPCclusters);
263     fHistESDratio->Fill(Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
264     fHistESDpt->Fill(track->Pt()); 
265   }
266 }
267
268 //____________________________________________________________________________
269 void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
270 {
271   //
272   // To make QA for the RAW data we use the TPC Calibration framework 
273   // to handle the data and then in the end extract the data
274   //
275   fTPCdataQA->ProcessEvent(rawReader);  
276 }
277
278 //____________________________________________________________________________
279 void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
280 {
281   AliTPCClustersRow *clrow = new AliTPCClustersRow();
282   clrow->SetClass("AliTPCclusterMI");
283   clrow->SetArray(0);
284   clrow->GetArray()->ExpandCreateFast(10000);
285
286   TBranch* branch = recTree->GetBranch("Segment");
287   branch->SetAddress(&clrow);
288
289   const Int_t nEntries = Int_t(recTree->GetEntries());
290   for (Int_t i = 0; i < nEntries; i++) {
291     
292     branch->GetEntry(i);
293     
294     const Int_t nClusters = clrow->GetArray()->GetEntriesFast();
295     for (Int_t icl=0; icl < nClusters; icl++){
296       
297       AliTPCclusterMI* cluster = 
298         (AliTPCclusterMI*)clrow->GetArray()->At(icl);
299       
300       Float_t Qmax = cluster->GetMax();
301       Float_t Q    = cluster->GetQ();
302       Int_t   row  = cluster->GetRow();
303
304       if(cluster->GetDetector()<36) { // IROC (short pads)
305
306         fHistRecPointsQmaxShort->Fill(Qmax);
307         fHistRecPointsQShort->Fill(Q);
308       } else { // OROC (medium and long pads)
309         row += 63;
310         if(cluster->GetRow()<64) { // medium pads
311
312           fHistRecPointsQmaxMedium->Fill(Qmax);
313           fHistRecPointsQMedium->Fill(Q);
314         } else { // long pads
315
316           fHistRecPointsQmaxLong->Fill(Qmax);
317           fHistRecPointsQLong->Fill(Q);
318         }
319       }
320       
321       fHistRecPointsRow->Fill(row);
322     } // end loop over clusters
323   } // end loop over tree
324
325   delete clrow;
326 }