1 /**************************************************************************
2 * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
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.
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.
44 Sumw2 for RAW histogram(s)?
45 RecPoints and ESD could have many more histograms
48 #include "AliTPCQADataMakerRec.h"
50 // --- ROOT system ---
51 #include <TClonesArray.h>
55 // --- Standard library ---
57 // --- AliRoot header files ---
58 #include "AliQAChecker.h"
59 #include "AliESDEvent.h"
60 #include "AliESDtrack.h"
62 #include "AliTPCCalPad.h"
63 #include "AliTPCCalROC.h"
64 #include "AliTPCClustersRow.h"
65 #include "AliTPCclusterMI.h"
67 ClassImp(AliTPCQADataMakerRec)
69 //____________________________________________________________________________
70 AliTPCQADataMakerRec::AliTPCQADataMakerRec() :
71 AliQADataMakerRec(AliQA::GetDetName(AliQA::kTPC),
72 "TPC Rec Quality Assurance Data Maker"),
74 fHistESDclusters(0),fHistESDratio(0), fHistESDpt(0),
75 fHistRawsOccupancy(0), fHistRawsOccupancyVsSector(0),
76 fHistRawsNClustersPerEventVsSector(0), fHistRawsQVsSector(0),
77 fHistRawsQmaxVsSector(0),
78 fHistRecPointsQmaxShort(0), fHistRecPointsQmaxMedium(0),
79 fHistRecPointsQmaxLong(0), fHistRecPointsQShort(0),
80 fHistRecPointsQMedium(0), fHistRecPointsQLong(0),
84 for(Int_t i = 0; i < 6; i++)
88 //____________________________________________________________________________
89 AliTPCQADataMakerRec::AliTPCQADataMakerRec(const AliTPCQADataMakerRec& qadm) :
92 fHistESDclusters(0),fHistESDratio(0), fHistESDpt(0),
93 fHistRawsOccupancy(0), fHistRawsOccupancyVsSector(0),
94 fHistRawsNClustersPerEventVsSector(0), fHistRawsQVsSector(0),
95 fHistRawsQmaxVsSector(0),
96 fHistRecPointsQmaxShort(0), fHistRecPointsQmaxMedium(0),
97 fHistRecPointsQmaxLong(0), fHistRecPointsQShort(0),
98 fHistRecPointsQMedium(0), fHistRecPointsQLong(0),
102 // Does not copy the calibration object, instead InitRaws have to be
104 SetName((const char*)qadm.GetName()) ;
105 SetTitle((const char*)qadm.GetTitle());
107 for(Int_t i = 0; i < 6; i++)
111 // Associate class histogram objects to the copies in the list
112 // Could also be done with the indexes
114 fHistESDclusters = (TH1F*)fESDsQAList->FindObject("hESDclusters");
115 fHistESDratio = (TH1F*)fESDsQAList->FindObject("hESDratio");
116 fHistESDpt = (TH1F*)fESDsQAList->FindObject("hESDpt");
118 fHistRawsOccupancy = (TH1F*)fRawsQAList->FindObject("hRawsOccupancy");
119 fHistRawsOccupancyVsSector =
120 (TH1F*)fRawsQAList->FindObject("hRawsOccupancyVsSector");
121 fHistRawsNClustersPerEventVsSector =
122 (TH1F*)fRawsQAList->FindObject("hRawsNClustersPerEventVsSector");
123 fHistRawsQVsSector = (TH1F*)fRawsQAList->FindObject("hRawsQVsSector");
124 fHistRawsQmaxVsSector = (TH1F*)fRawsQAList->FindObject("hRawsQmaxVsSector");
126 fHistRecPointsQmaxShort =
127 (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxShort");
128 fHistRecPointsQmaxMedium =
129 (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxMedium");
130 fHistRecPointsQmaxLong =
131 (TH1F*)fRecPointsQAList->FindObject("hRecPointsQmaxLong");
132 fHistRecPointsQShort =
133 (TH1F*)fRecPointsQAList->FindObject("hRecPointsQShort");
134 fHistRecPointsQMedium =
135 (TH1F*)fRecPointsQAList->FindObject("hRecPointsQMedium");
136 fHistRecPointsQLong =
137 (TH1F*)fRecPointsQAList->FindObject("hRecPointsQLong");
139 (TH1F*)fRecPointsQAList->FindObject("hRecPointsRow");
142 //__________________________________________________________________
143 AliTPCQADataMakerRec& AliTPCQADataMakerRec::operator = (const AliTPCQADataMakerRec& qadm )
146 this->~AliTPCQADataMakerRec();
147 new(this) AliTPCQADataMakerRec(qadm);
151 //__________________________________________________________________
152 AliTPCQADataMakerRec::~AliTPCQADataMakerRec()
157 for(Int_t i = 0; i < 6; i++)
161 //____________________________________________________________________________
162 void AliTPCQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
164 //Detector specific actions at end of cycle
166 if(fTPCdataQA) { // do the final step of the QA for Raw data
168 fTPCdataQA->Analyse(); // 31/1-08 Analyse is now protected against
169 // RAW data files with no TPC data
171 //Add2RawsList(fTPCdataQA, 0);
172 // get the histograms and add them to the output
173 // 31/8-08 Histogram is only added if the Calibration class
175 const Int_t eventCounter = fTPCdataQA->GetEventCounter();
176 if(eventCounter>0) { // some TPC data has been processed
178 // Reset histograms and refill them
179 fHistRawsOccupancy->Reset();
180 fHistRawsOccupancyVsSector->Reset();
181 fHistRawsNClustersPerEventVsSector->Reset();
182 fHistRawsQVsSector->Reset();
183 fHistRawsQmaxVsSector->Reset();
185 TH1F* hNorm72 = new TH1F("hNorm72", "histogram to normalize 72 sectors",
188 TH1F* hNorm108 = new TH1F("hNorm108", "histogram to normalize 108 sectors (medium and long pads are split up)",
192 for (Int_t iSec = 0; iSec < 72; iSec++) {
194 AliTPCCalROC* occupancyROC =
195 fTPCdataQA->GetNoThreshold()->GetCalROC(iSec);
196 AliTPCCalROC* nclusterROC =
197 fTPCdataQA->GetNLocalMaxima()->GetCalROC(iSec);
199 fTPCdataQA->GetMeanCharge()->GetCalROC(iSec);
200 AliTPCCalROC* qmaxROC =
201 fTPCdataQA->GetMaxCharge()->GetCalROC(iSec);
203 const Int_t nRows = occupancyROC->GetNrows();
204 for (Int_t iRow = 0; iRow < nRows; iRow++) {
206 Int_t helpSector = iSec;
208 helpSector += 36; // OROC (long pads)
210 const Int_t nPads = occupancyROC->GetNPads(iRow);
211 for (Int_t iPad = 0; iPad < nPads; iPad++) {
213 fHistRawsOccupancy->Fill(occupancyROC->GetValue(iRow, iPad));
215 fHistRawsOccupancyVsSector
216 ->Fill(iSec, occupancyROC->GetValue(iRow, iPad));
218 const Int_t nClusters = TMath::Nint(nclusterROC->GetValue(iRow, iPad));
222 fHistRawsNClustersPerEventVsSector->Fill(iSec, nClusters);
223 hNorm108->Fill(helpSector, nClusters);
224 fHistRawsQVsSector->Fill(helpSector,
225 nClusters*qROC->GetValue(iRow, iPad));
226 fHistRawsQmaxVsSector->Fill(helpSector,
227 nClusters*qmaxROC->GetValue(iRow, iPad));
231 } // end loop over sectors
233 // Normalize histograms
234 fHistRawsOccupancyVsSector->Divide(hNorm72);
235 fHistRawsNClustersPerEventVsSector->Scale(1.0/Float_t(eventCounter));
236 fHistRawsQVsSector->Divide(hNorm108);
237 fHistRawsQmaxVsSector->Divide(hNorm108);
244 AliQAChecker::Instance()->Run(AliQA::kTPC, task, list) ;
247 //____________________________________________________________________________
248 void AliTPCQADataMakerRec::InitESDs()
250 //create ESDs histograms in ESDs subdir
252 new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
254 fHistESDclusters->Sumw2();
255 Add2ESDsList(fHistESDclusters, 0);
258 new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
260 fHistESDratio->Sumw2();
261 Add2ESDsList(fHistESDratio, 1);
264 new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
267 Add2ESDsList(fHistESDpt, 2);
270 //____________________________________________________________________________
271 void AliTPCQADataMakerRec::InitRaws()
277 // Modified: 7/7 - 2008
278 // Laurent Aphecetche pointed out that the mapping was read from file
279 // for each event, so now we read in the map here and set if for
281 fTPCdataQA = new AliTPCdataQA();
282 LoadMaps(); // Load Altro maps
283 fTPCdataQA->SetAltroMapping(fMapping); // set Altro mapping
284 fTPCdataQA->SetRangeTime(100, 920); // set time bin interval
285 Add2RawsList(fTPCdataQA, 0); // This is used by the AMORE monitoring
288 new TH1F("hRawsOccupancy", "Occupancy (all pads); Occupancy; Counts",
290 fHistRawsOccupancy->Sumw2();
291 Add2RawsList(fHistRawsOccupancy, 1);
293 fHistRawsOccupancyVsSector =
294 new TH1F("hRawsOccupancyVsSector", "Occupancy vs sector; Sector; Occupancy",
296 fHistRawsOccupancyVsSector->Sumw2();
297 Add2RawsList(fHistRawsOccupancyVsSector, 2);
299 fHistRawsNClustersPerEventVsSector =
300 new TH1F("hRawsNClustersPerEventVsSector", "Nclusters per event vs sector; Sector; Nclusters per event",
302 fHistRawsNClustersPerEventVsSector->Sumw2();
303 Add2RawsList(fHistRawsNClustersPerEventVsSector, 3);
306 new TH1F("hRawsQVsSector", "<Q> vs sector (OROC med: 36-71, long: 72-107); Sector; <Q>",
308 fHistRawsQVsSector->Sumw2();
309 Add2RawsList(fHistRawsQVsSector, 4);
311 fHistRawsQmaxVsSector =
312 new TH1F("hRawsQmaxVsSector", "<Qmax> vs sector (OROC med: 36-71, long: 72-107); Sector; <Qmax>",
314 fHistRawsQmaxVsSector->Sumw2();
315 Add2RawsList(fHistRawsQmaxVsSector, 5);
318 //____________________________________________________________________________
319 void AliTPCQADataMakerRec::InitRecPoints()
321 fHistRecPointsQmaxShort =
322 new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
324 fHistRecPointsQmaxShort->Sumw2();
325 Add2RecPointsList(fHistRecPointsQmaxShort, 0);
327 fHistRecPointsQmaxMedium =
328 new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
330 fHistRecPointsQmaxMedium->Sumw2();
331 Add2RecPointsList(fHistRecPointsQmaxMedium, 1);
333 fHistRecPointsQmaxLong =
334 new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
336 fHistRecPointsQmaxLong->Sumw2();
337 Add2RecPointsList(fHistRecPointsQmaxLong, 2);
339 fHistRecPointsQShort =
340 new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
342 fHistRecPointsQShort->Sumw2();
343 Add2RecPointsList(fHistRecPointsQShort, 3);
345 fHistRecPointsQMedium =
346 new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
348 fHistRecPointsQMedium->Sumw2();
349 Add2RecPointsList(fHistRecPointsQMedium, 4);
351 fHistRecPointsQLong =
352 new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
354 fHistRecPointsQLong->Sumw2();
355 Add2RecPointsList(fHistRecPointsQLong, 5);
358 new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
360 fHistRecPointsRow->Sumw2();
361 Add2RecPointsList(fHistRecPointsRow, 6);
364 //____________________________________________________________________________
365 void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
367 // make QA data from ESDs
369 const Int_t nESDTracks = esd->GetNumberOfTracks();
370 Int_t nTPCtracks = 0;
371 for(Int_t i = 0; i < nESDTracks; i++) {
373 AliESDtrack * track = esd->GetTrack(i);
375 if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
380 Int_t nTPCclusters = track->GetTPCNcls();
381 Int_t nTPCclustersFindable = track->GetTPCNclsF();
382 if ( nTPCclustersFindable<=0) continue;
383 fHistESDclusters->Fill(nTPCclusters);
384 fHistESDratio->Fill(Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
385 fHistESDpt->Fill(track->Pt());
389 //____________________________________________________________________________
390 void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
393 // To make QA for the RAW data we use the TPC Calibration framework
394 // to handle the data and then in the end extract the data
397 fTPCdataQA->ProcessEvent(rawReader);
400 //____________________________________________________________________________
401 void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
403 AliTPCClustersRow *clrow = new AliTPCClustersRow();
404 clrow->SetClass("AliTPCclusterMI");
406 clrow->GetArray()->ExpandCreateFast(10000);
408 TBranch* branch = recTree->GetBranch("Segment");
409 branch->SetAddress(&clrow);
411 const Int_t nEntries = Int_t(recTree->GetEntries());
412 for (Int_t i = 0; i < nEntries; i++) {
416 const Int_t nClusters = clrow->GetArray()->GetEntriesFast();
417 for (Int_t icl=0; icl < nClusters; icl++){
419 AliTPCclusterMI* cluster =
420 (AliTPCclusterMI*)clrow->GetArray()->At(icl);
422 Float_t Qmax = cluster->GetMax();
423 Float_t Q = cluster->GetQ();
424 Int_t row = cluster->GetRow();
426 if(cluster->GetDetector()<36) { // IROC (short pads)
428 fHistRecPointsQmaxShort->Fill(Qmax);
429 fHistRecPointsQShort->Fill(Q);
430 } else { // OROC (medium and long pads)
432 if(cluster->GetRow()<64) { // medium pads
434 fHistRecPointsQmaxMedium->Fill(Qmax);
435 fHistRecPointsQMedium->Fill(Q);
436 } else { // long pads
438 fHistRecPointsQmaxLong->Fill(Qmax);
439 fHistRecPointsQLong->Fill(Q);
443 fHistRecPointsRow->Fill(row);
444 } // end loop over clusters
445 } // end loop over tree
450 //____________________________________________________________________________
451 void AliTPCQADataMakerRec::LoadMaps()
453 TString path = gSystem->Getenv("ALICE_ROOT");
454 path += "/TPC/mapping/Patch";
456 for(Int_t i = 0; i < 6; i++) {
457 TString path2 = path;
460 fMapping[i] = new AliTPCAltroMapping(path2.Data());