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