]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCQADataMakerRec.cxx
Adding protection in case of missing data (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCQADataMakerRec.cxx
CommitLineData
44f32dd2 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
64ClassImp(AliTPCQADataMakerRec)
65
66//____________________________________________________________________________
67AliTPCQADataMakerRec::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//____________________________________________________________________________
82AliTPCQADataMakerRec::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//__________________________________________________________________
118AliTPCQADataMakerRec& AliTPCQADataMakerRec::operator = (const AliTPCQADataMakerRec& qadm )
119{
120 // Equal operator.
121 this->~AliTPCQADataMakerRec();
122 new(this) AliTPCQADataMakerRec(qadm);
123 return *this;
124}
125
126//____________________________________________________________________________
127void 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();
134
135 // get the histograms and add them to the output
136 fHistRawsOccupancy = fTPCdataQA->GetNoThreshold()->MakeHisto1D(0, 1, -1);
137 Add2RawsList(fHistRawsOccupancy, 0);
138 }
139
140 AliQAChecker::Instance()->Run(AliQA::kTPC, task, list) ;
141}
142
143//____________________________________________________________________________
144void AliTPCQADataMakerRec::InitESDs()
145{
146 //create ESDs histograms in ESDs subdir
147 fHistESDclusters =
148 new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
149 160, 0, 160);
150 fHistESDclusters->Sumw2();
151 Add2ESDsList(fHistESDclusters, 0);
152
153 fHistESDratio =
154 new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
155 100, 0, 1);
156 fHistESDratio->Sumw2();
157 Add2ESDsList(fHistESDratio, 1);
158
159 fHistESDpt =
160 new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
161 50, 0, 5);
162 fHistESDpt->Sumw2();
163 Add2ESDsList(fHistESDpt, 2);
164}
165
166//____________________________________________________________________________
167void AliTPCQADataMakerRec::InitRaws()
168{
169 fTPCdataQA = new AliTPCdataQA();
170 fTPCdataQA->SetRangeTime(0, 999); // take all 1000 time bins
171}
172
173//____________________________________________________________________________
174void AliTPCQADataMakerRec::InitRecPoints()
175{
176 fHistRecPointsQmaxShort =
177 new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
178 200, 0, 1000);
179 fHistRecPointsQmaxShort->Sumw2();
180 Add2RecPointsList(fHistRecPointsQmaxShort, 0);
181
182 fHistRecPointsQmaxMedium =
183 new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
184 200, 0, 1000);
185 fHistRecPointsQmaxMedium->Sumw2();
186 Add2RecPointsList(fHistRecPointsQmaxMedium, 1);
187
188 fHistRecPointsQmaxLong =
189 new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
190 200, 0, 1000);
191 fHistRecPointsQmaxLong->Sumw2();
192 Add2RecPointsList(fHistRecPointsQmaxLong, 2);
193
194 fHistRecPointsQShort =
195 new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
196 200, 0, 5000);
197 fHistRecPointsQShort->Sumw2();
198 Add2RecPointsList(fHistRecPointsQShort, 3);
199
200 fHistRecPointsQMedium =
201 new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
202 200, 0, 5000);
203 fHistRecPointsQMedium->Sumw2();
204 Add2RecPointsList(fHistRecPointsQMedium, 4);
205
206 fHistRecPointsQLong =
207 new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
208 200, 0, 5000);
209 fHistRecPointsQLong->Sumw2();
210 Add2RecPointsList(fHistRecPointsQLong, 5);
211
212 fHistRecPointsRow =
213 new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
214 159, 0, 159);
215 fHistRecPointsRow->Sumw2();
216 Add2RecPointsList(fHistRecPointsRow, 6);
217}
218
219//____________________________________________________________________________
220void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
221{
222 // make QA data from ESDs
223
224 const Int_t nESDTracks = esd->GetNumberOfTracks();
225 Int_t nTPCtracks = 0;
226 for(Int_t i = 0; i < nESDTracks; i++) {
227
228 AliESDtrack * track = esd->GetTrack(i);
229
230 if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
231 continue;
232
233 nTPCtracks++;
234
235 Int_t nTPCclusters = track->GetTPCNcls();
236 Int_t nTPCclustersFindable = track->GetTPCNclsF();
237
238 fHistESDclusters->Fill(nTPCclusters);
239 fHistESDratio->Fill(Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
240 fHistESDpt->Fill(track->Pt());
241 }
242}
243
244//____________________________________________________________________________
245void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
246{
247 //
248 // To make QA for the RAW data we use the TPC Calibration framework
249 // to handle the data and then in the end extract the data
250 //
251 fTPCdataQA->ProcessEvent(rawReader);
252}
253
254//____________________________________________________________________________
255void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
256{
257 AliTPCClustersRow *clrow = new AliTPCClustersRow();
258 clrow->SetClass("AliTPCclusterMI");
259 clrow->SetArray(0);
260 clrow->GetArray()->ExpandCreateFast(10000);
261
262 TBranch* branch = recTree->GetBranch("Segment");
263 branch->SetAddress(&clrow);
264
265 const Int_t nEntries = Int_t(recTree->GetEntries());
266 for (Int_t i = 0; i < nEntries; i++) {
267
268 branch->GetEntry(i);
269
270 const Int_t nClusters = clrow->GetArray()->GetEntriesFast();
271 for (Int_t icl=0; icl < nClusters; icl++){
272
273 AliTPCclusterMI* cluster =
274 (AliTPCclusterMI*)clrow->GetArray()->At(icl);
275
276 Float_t Qmax = cluster->GetMax();
277 Float_t Q = cluster->GetQ();
278 Int_t row = cluster->GetRow();
279
280 if(cluster->GetDetector()<36) { // IROC (short pads)
281
282 fHistRecPointsQmaxShort->Fill(Qmax);
283 fHistRecPointsQShort->Fill(Q);
284 } else { // OROC (medium and long pads)
285 row += 63;
286 if(cluster->GetRow()<64) { // medium pads
287
288 fHistRecPointsQmaxMedium->Fill(Qmax);
289 fHistRecPointsQMedium->Fill(Q);
290 } else { // long pads
291
292 fHistRecPointsQmaxLong->Fill(Qmax);
293 fHistRecPointsQLong->Fill(Q);
294 }
295 }
296
297 fHistRecPointsRow->Fill(row);
298 } // end loop over clusters
299 } // end loop over tree
300
301 delete clrow;
302}