]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ITS/AliITSQASSDDataMakerRec.cxx
moving class
[u/mrichter/AliRoot.git] / ITS / AliITSQASSDDataMakerRec.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 2007-2009, 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/* $Id$ */
16// *************************************************************
17// Checks the quality assurance
18// by comparing with reference data
19// contained in a DB
20// -------------------------------------------------------------
21// W. Ferrarese + P. Cerello Feb 2008
22// INFN Torino
23
24// --- ROOT system ---
25#include <TH2D.h>
26#include <TTree.h>
27#include <TMath.h>
28#include <TString.h>
29#include <TSystem.h>
30
31// --- Standard library ---
32
33// --- AliRoot header files ---
34#include "AliITSQASSDDataMakerRec.h"
35#include "AliQADataMakerRec.h"
36#include "AliLog.h"
37#include "AliQA.h"
38#include "AliQAChecker.h"
39#include "AliRawReader.h"
40#include "AliRawReaderRoot.h"
41#include "AliITSRawStreamSSD.h"
42#include "AliITSgeomTGeo.h"
43#include "AliRawEventHeaderBase.h"
44#include "AliITSRecPoint.h"
45
46#include "AliCDBManager.h"
47#include "AliCDBEntry.h"
48#include "AliESDEvent.h"
49#include "AliESDtrack.h"
50
51ClassImp(AliITSQASSDDataMakerRec)
52
53AliITSQASSDDataMakerRec::AliITSQASSDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Int_t ldc) :
54TObject(),
55fAliITSQADataMakerRec(aliITSQADataMakerRec),
56fSSDEvent(0),
57fkOnline(kMode),
58fLDC(ldc),
59fSSDRawsOffset(0),
60fSSDhTask(0),
61fGenOffset(0) {
62 // Default constructor
63 //initilize the raw signal vs strip number histograms
64 Int_t fLayer = 0,fLadder = 0, fModule = 0;
65 TString fTitle; Int_t fHistCounter = 0;
66 if(fkOnline) {
67 AliCDBManager * man = AliCDBManager::Instance();
68 man->SetDefaultStorage("local://$ALICE_ROOT");
69 Int_t runNumber = atoi(gSystem->Getenv("DATE_RUN_NUMBER"));
70 if(!runNumber)
71 AliInfo("DATE_RUN_NUMBER not defined!!!\n");
72
73 man->SetRun(runNumber);
74 AliCDBEntry *geomGRP = man->Get("GRP/Geometry/Data");
75 if(!geomGRP) AliInfo("GRP geometry not found!!!\n");
76
77 for(Int_t i = 500; i < fgkSSDMODULES + 500; i++) {
78 AliITSgeomTGeo::GetModuleId(i,fLayer,fLadder,fModule);
79 fTitle = "SSD_RawSignal_Layer"; fTitle += fLayer;
80 fTitle += "_Ladder"; fTitle += fLadder;
81 fTitle += "_Module"; fTitle += fModule;
82 fHistSSDRawSignalModule[fHistCounter] = new TH1D(fTitle.Data(),fTitle.Data(),1540,0,1540);
83 fHistCounter += 1;
84 }
85 }//online flag
86 else {
87 for(Int_t i = 0; i < fgkSSDMODULES; i++) {
88 fHistSSDRawSignalModule[i]=NULL;
89 }
90 }
91}
92
93//____________________________________________________________________________
94AliITSQASSDDataMakerRec::AliITSQASSDDataMakerRec(const AliITSQASSDDataMakerRec& qadm) :
95TObject(),
96fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
97fSSDEvent(qadm.fSSDEvent),
98fkOnline(qadm.fkOnline),
99fLDC(qadm.fLDC),
100fSSDRawsOffset(qadm.fSSDRawsOffset),
101fSSDhTask(qadm.fSSDhTask),
102fGenOffset(qadm.fGenOffset) {
103 //copy ctor
104 fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ;
105 fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
106}
107
108//__________________________________________________________________
109AliITSQASSDDataMakerRec& AliITSQASSDDataMakerRec::operator = (const AliITSQASSDDataMakerRec& qac )
110{
111 // Equal operator.
112 this->~AliITSQASSDDataMakerRec();
113 new(this) AliITSQASSDDataMakerRec(qac);
114 return *this;
115}
116
117//__________________________________________________________________
118AliITSQASSDDataMakerRec::~AliITSQASSDDataMakerRec() {
119 // destructor
120 for(Int_t i = 0; i < fgkSSDMODULES; i++)
121 if(fHistSSDRawSignalModule[i]) delete fHistSSDRawSignalModule[i];
122}
123
124//____________________________________________________________________________
125void AliITSQASSDDataMakerRec::StartOfDetectorCycle()
126{
127 //Detector specific actions at start of cycle
128 AliDebug(1,"AliITSQADM::Start of SSD Cycle\n");
129}
130
131//____________________________________________________________________________
132void AliITSQASSDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
133{
134 // launch the QA checking
135 AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n");
136
137 //scaling SSD occupancy plots
138 if(fkOnline) {
139 if(fSSDEvent != 0) {
140 ((TH2F *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+28))->Scale(1./fSSDEvent);
141 ((TH2F *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+29))->Scale(1./fSSDEvent);
142 }
143 Int_t fLayer = 0, fLadder = 0;
144 Int_t fOccPosition = 0;
145 for(Int_t iLayer = 5; iLayer < 7; iLayer++) {
146 for(Int_t iLadder = 1; iLadder < AliITSgeomTGeo::GetNLadders(iLayer) + 1; iLadder++) {
147 fOccPosition = (fLayer == 5) ? fLadder : fLadder + fgkSSDLADDERSLAYER5;
148 if(fSSDEvent != 0) {
149 //P-SIDE OCCUPANCY - scaling
150 fAliITSQADataMakerRec->GetRawsData(fGenOffset + fSSDRawsOffset - 1 + 2*fOccPosition - 1)->Scale(1./fSSDEvent);
151 //N-SIDE OCCUPANCY - scaling
152 fAliITSQADataMakerRec->GetRawsData(fGenOffset + fSSDRawsOffset - 1 + 2*fOccPosition)->Scale(1./fSSDEvent);
153 }//ssd events != 0
154 }//ladder loop
155 }//layer loop
156 }//online flag for SSD
157
158 AliQAChecker::Instance()->Run( AliQA::kITS , task, list);
159}
160
161//____________________________________________________________________________
162void AliITSQASSDDataMakerRec::InitRaws() {
163
164 // Initialization for RAW data - SSD -
165 fGenOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();
166
167 if(fkOnline) {
168 AliInfo("Book Online Histograms for SSD\n");
169 }
170 else {
171 AliInfo("Book Offline Histograms for SSD\n ");
172 }
173 AliInfo(Form("Number of histograms (SPD+SDD): %d\n",fGenOffset));
174 TString fTitle = 0;
175 //book online QA histos
176 TH1F *fHistSSDEventType = new TH1F("fHistSSDEventType",
177 ";Event type;Events",
178 31,-1,30);
179 fAliITSQADataMakerRec->Add2RawsList(fHistSSDEventType,
180 fGenOffset+fSSDRawsOffset);
181 fSSDRawsOffset += 1;
182 TH1F *fHistSSDDataSize = new TH1F("fHistSSDDataSize",
183 ";log(SSD data size) [Bytes];Events",
184 100,3,8);
185 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSize,
186 fGenOffset+fSSDRawsOffset);
187 fSSDRawsOffset += 1;
188 TH1F *fHistSSDDataSizePercentage = new TH1F("fHistSSDDataSizePercentage",
189 ";SSD data size [%];Events",
190 100,0,100);
191 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePercentage,
192 fGenOffset+fSSDRawsOffset);
193 fSSDRawsOffset += 1;
194 TH1F *fHistSSDDDLId = new TH1F("fHistSSDDDLId",
195 ";DDL id;Events",20,510.5,530.5);
196 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDDLId,
197 fGenOffset+fSSDRawsOffset);
198 fSSDRawsOffset += 1;
199 TH1F *fHistSSDDataSizePerDDL = new TH1F("fHistSSDDataSizePerDDL",
200 ";DDL id;<SSD data size> [MB]",
201 20,510.5,530.5);
202 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePerDDL,
203 fGenOffset+fSSDRawsOffset);
204 fSSDRawsOffset += 1;
205 TH1F *fHistSSDDataSizeDDL[fgkNumOfDDLs];
206 for(Int_t i = 1; i < fgkNumOfDDLs+1; i++) {
207 fTitle = "fHistSSDDataSizeDDL"; fTitle += i+511;
208 fHistSSDDataSizeDDL[i-1] = new TH1F(fTitle.Data(),
209 ";log(SSD data size) [Bytes];Events",
210 100,1,8);
211 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizeDDL[i-1],
212 fGenOffset+fSSDRawsOffset);
213 fSSDRawsOffset += 1;
214 }
215
216 TH1F *fHistSSDLDCId = new TH1F("fHistSSDLDCId",";LDC id;Events",10,0.5,10.5);
217 fAliITSQADataMakerRec->Add2RawsList(fHistSSDLDCId,
218 fGenOffset+fSSDRawsOffset);
219 fSSDRawsOffset += 1;
220 TH1F *fHistSSDDataSizePerLDC = new TH1F("fHistSSDDataSizePerLDC",
221 ";LDC id;<SSD data size> [MB]",
222 100,0,20);
223 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePerLDC,
224 fGenOffset+fSSDRawsOffset);
225 fSSDRawsOffset += 1;
226 TH1F *fHistSSDDataSizeLDC[fgkNumOfLDCs];
227 for(Int_t i = 1; i < fgkNumOfLDCs+1; i++) {
228 fTitle = "fHistSSDDataSizeLDC"; fTitle += i;
229 fHistSSDDataSizeLDC[i-1] = new TH1F(fTitle.Data(),
230 ";log(SSD data size) [Bytes];Events",
231 100,1,8);
232 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizeLDC[i-1],
233 fGenOffset+fSSDRawsOffset);
234 fSSDRawsOffset += 1;
235 }
236 fSSDhTask = fSSDRawsOffset;
237 if(fkOnline) {
238 //top level occupancy plots
239 TH2D *fHistSSDOccupancyLayer5 = new TH2D("fHistSSDOccupancyLayer5",
240 ";N_{modules};N_{Ladders}",
241 fgkSSDMODULESPERLADDERLAYER5,
242 0,fgkSSDMODULESPERLADDERLAYER5,
243 3*fgkSSDLADDERSLAYER5,
244 0,fgkSSDLADDERSLAYER5);
245 Char_t fLabel[3];
246 for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER5 + 1; iBin++){
247 sprintf(fLabel,"%d",iBin);
248 fHistSSDOccupancyLayer5->GetXaxis()->SetBinLabel(iBin,fLabel);
249 }
250 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLayer5, fGenOffset+fSSDRawsOffset);
251 fSSDRawsOffset += 1;
252 TH2D *fHistSSDOccupancyLayer6 = new TH2D("fHistSSDOccupancyLayer6",
253 ";N_{modules};N_{Ladders}",
254 fgkSSDMODULESPERLADDERLAYER6,
255 0,fgkSSDMODULESPERLADDERLAYER6,
256 3*fgkSSDLADDERSLAYER6,
257 0,fgkSSDLADDERSLAYER6);
258 for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER6 + 1; iBin++){
259 sprintf(fLabel,"%d",iBin);
260 fHistSSDOccupancyLayer6->GetXaxis()->SetBinLabel(iBin,fLabel);
261 }
262 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLayer6, fGenOffset+fSSDRawsOffset);
263 fSSDRawsOffset += 1;
264
265 //Occupancy per ladder
266 fSSDhTask = fSSDRawsOffset;
267 TH1D *fHistOccupancyLadder[2*(fgkSSDLADDERSLAYER5 + fgkSSDLADDERSLAYER6)];
268 for(Int_t iLayer = 5; iLayer < 7; iLayer++) {
269 for(Int_t iLadder = 1; iLadder < AliITSgeomTGeo::GetNLadders(iLayer) + 1; iLadder++) {
270 //P-side occupancy plots
271 fTitle = "fHistSSD_Occupancy_Layer"; fTitle += iLayer;
272 fTitle += "_Ladder"; fTitle += iLadder; fTitle += "_PSide";
273 fHistOccupancyLadder[fSSDhTask] =
274 new TH1D(fTitle.Data(),
275 fTitle.Data(),
276 AliITSgeomTGeo::GetNDetectors(iLayer),
277 0.5,AliITSgeomTGeo::GetNDetectors(iLayer)+0.5);
278 fHistOccupancyLadder[fSSDhTask]->GetXaxis()->SetTitleColor(1);
279 fHistOccupancyLadder[fSSDhTask]->GetXaxis()->SetTitle("Module number");
280 fHistOccupancyLadder[fSSDhTask]->GetYaxis()->SetTitle("Occupancy [%]");
281 fAliITSQADataMakerRec->Add2RawsList(fHistOccupancyLadder[fSSDhTask],
282 fGenOffset+fSSDhTask);
283 fSSDhTask++;
284 //N-side occupancy plots
285 fTitle = "fHistSSD_Occupancy_Layer"; fTitle += iLayer;
286 fTitle += "_Ladder"; fTitle += iLadder; fTitle += "_NSide";
287 fHistOccupancyLadder[fSSDhTask] =
288 new TH1D(fTitle.Data(),
289 fTitle.Data(),
290 AliITSgeomTGeo::GetNDetectors(iLayer),
291 0.5,AliITSgeomTGeo::GetNDetectors(iLayer)+0.5);
292 fHistOccupancyLadder[fSSDhTask]->GetXaxis()->SetTitleColor(1);
293 fHistOccupancyLadder[fSSDhTask]->GetXaxis()->SetTitle("Module number");
294 fHistOccupancyLadder[fSSDhTask]->GetYaxis()->SetTitle("Occupancy [%]");
295 fAliITSQADataMakerRec->Add2RawsList(fHistOccupancyLadder[fSSDhTask],
296 fGenOffset+fSSDhTask);
297 fSSDhTask++;
298 }//ladder loop
299 }//layer loop
300 }//online flag
301 AliDebug(1,Form("%d SSD Raws histograms booked\n",fSSDhTask));
302 AliInfo(Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenOffset+fSSDhTask));
303 AliInfo(Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenOffset+fSSDRawsOffset));
304}
305
306
307//____________________________________________________________________________
308void AliITSQASSDDataMakerRec::MakeRaws(AliRawReader* rawReader) {
309 // Fill QA for RAW - SSD -
310 Int_t fStripNumber;
311 Int_t fHistPosition;
312 Int_t fLayer = 0,fLadder = 0, fModule = 0;
313
314 Double_t fSizePerDDL[fgkNumOfDDLs] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
315 Double_t fSizePerLDC[fgkNumOfLDCs] = {0.,0.,0.};
316 Double_t sumSSDDataSize = 0.0;
317 Double_t eventSize = -1.0;
318
319 if(fkOnline) {
320 //reset the signal vs strip number histograms
321 for(Int_t i = 0; i < fgkSSDMODULES; i++)
322 fHistSSDRawSignalModule[i]->Reset();
323 }
324 rawReader->Select("ITSSSD",-1,-1);
325 rawReader->Reset(); //rawReader->NextEvent();
326 (fAliITSQADataMakerRec->GetRawsData(fGenOffset))->Fill(rawReader->GetType());
327 if(rawReader->GetType() == 7)
328 fSSDEvent += 1;
329
330 AliITSRawStreamSSD fSSDStream(rawReader);
331 AliRawReaderRoot *rootReader = (AliRawReaderRoot *)rawReader;
332 if(rootReader) {
333 const AliRawEventHeaderBase *header = rootReader->GetEventHeader();
334 if(header)
335 eventSize = header->GetEventSize();
336 }
337 while (fSSDStream.Next()) {
338 if(fSSDStream.GetModuleID() < 0) continue;
339 fSizePerDDL[rawReader->GetDDLID()] = rawReader->GetDataSize();
340 fSizePerLDC[rawReader->GetLDCId()-6] = rawReader->GetDataSize();
341 AliITSgeomTGeo::GetModuleId(fSSDStream.GetModuleID(),fLayer,fLadder,fModule);
342 fStripNumber = (fSSDStream.GetSideFlag() == 0) ? fSSDStream.GetStrip() : fSSDStream.GetStrip() + fgkNumberOfPSideStrips;
343 fHistPosition = (fLayer == 5) ? ((fLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + fModule - 1) : ((fLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + fModule + fgkSSDMODULESLAYER5 - 1);
344 //AliInfo(Form("ModulePosition: %d - Layer: %d - Ladder: %d - Module: %d\n",fHistPosition,fLayer,fLadder,fModule));
345 if(fkOnline)
346 fHistSSDRawSignalModule[fHistPosition]->Fill(fStripNumber,fSSDStream.GetSignal());
347 }//streamer loop
348
349
350 for(Int_t i = 0; i < fgkNumOfDDLs; i++) {
351 sumSSDDataSize += fSizePerDDL[i];
352 if(fSizePerDDL[i] > 0) {
353 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+3))->Fill(i+512);
354 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+5+i))->Fill(TMath::Log10(fSizePerDDL[i]));
355 }
356 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+4))->Fill(i+512,fSizePerDDL[i]/1e+06);
357 }
358 for(Int_t i = 0; i < fgkNumOfLDCs; i++) {
359 if(fSizePerLDC[i] > 0) {
360 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+21))->Fill(i+6);
361 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+23+i))->Fill(TMath::Log10(fSizePerLDC[i]));
362 }
363 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+22))->Fill(i+6,fSizePerLDC[i]/1e+06);
364 }
365 if(sumSSDDataSize)
366 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+1))->Fill(TMath::Log10(sumSSDDataSize));
367 if(eventSize)
368 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+2))->Fill(100.*sumSSDDataSize/eventSize);
369 if(fkOnline) {
370 //AliInfo(Form("EVENT: %d\n",fSSDEvent));
371 Double_t occupancy = 0.0;
372 Int_t lLadderLocationY = 0;
373 Int_t fOccPosition = 0;
374 for(Int_t i = 0; i < fgkSSDMODULES; i++) {
375 AliITSgeomTGeo::GetModuleId(i+500,fLayer,fLadder,fModule);
376 fOccPosition = (fLayer == 5) ? fLadder : fLadder + fgkSSDLADDERSLAYER5;
377 //P-SIDE OCCUPANCY
378 occupancy = GetSSDOccupancyRaws(fHistSSDRawSignalModule[i],0);
379 fAliITSQADataMakerRec->GetRawsData(fGenOffset + fSSDRawsOffset - 1 + 2*fOccPosition - 1)->Fill(fModule,occupancy);
380 lLadderLocationY = 3*fLadder; // sideP=1 sideN=0
381 if(fLayer == 5)
382 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+28))->SetBinContent(fModule,lLadderLocationY,occupancy);
383 else if(fLayer == 6)
384 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+29))->SetBinContent(fModule,lLadderLocationY,occupancy);
385 //N-SIDE OCCUPANCY
386 occupancy = GetSSDOccupancyRaws(fHistSSDRawSignalModule[i],1);
387 fAliITSQADataMakerRec->GetRawsData(fGenOffset + fSSDRawsOffset - 1 + 2*fOccPosition)->Fill(fModule,occupancy);
388 if(fLayer == 5)
389 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+28))->SetBinContent(fModule,lLadderLocationY-1,occupancy);
390 else if(fLayer == 6)
391 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+29))->SetBinContent(fModule,lLadderLocationY-1,occupancy);
392 }
393 }//online flag for SSD
394}
395
396
397
398
399//____________________________________________________________________________
400Double_t AliITSQASSDDataMakerRec::GetSSDOccupancyRaws(TH1 *lHisto, Int_t stripside) {
401 // bo: TDC >0 or # of sigmas wrt noise ?
402 //stripside == 0 --> P-side
403 //stripside == 1 --> N-side
404 Int_t lNumFiredBins = 0;
405 for(Int_t iBin = 1 + stripside*fgkNumberOfPSideStrips; iBin < fgkNumberOfPSideStrips*(1 + stripside); iBin++){
406 if (lHisto->GetBinContent(iBin) > 0)
407 lNumFiredBins++;
408 }
409
410 Double_t lOccupancy = (100.*lNumFiredBins)/fgkNumberOfPSideStrips; // percentage
411 //AliInfo(Form("Fired strips: %d - Total strips: %d - Occupancy :%lf\n",lNumFiredBins,lHisto->GetNbinsX(),lOccupancy));
412
413 return lOccupancy;
414}
415
416
417//____________________________________________________________________________
418void AliITSQASSDDataMakerRec::InitRecPoints()
419{
420 // Initialization for RECPOINTS - SSD -
421 fGenOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();
422 Int_t nModuleOffset = 500;
423 Int_t nITSTotalModules = AliITSgeomTGeo::GetNModules();
424
425 TH1F *fHistModuleIdLayer5 = new TH1F("fHistModuleIdLayer5",
426 "Module Id - Layer 5;Module Id;Entries",
427 fgkSSDMODULESLAYER5,
428 nModuleOffset - 0.5,
429 nITSTotalModules-fgkSSDMODULESLAYER6+0.5);
430 fAliITSQADataMakerRec->Add2RecPointsList(fHistModuleIdLayer5,
431 fGenOffset + 0);
432 fSSDhTask += 1;
433 TH1F *fHistModuleIdLayer6 = new TH1F("fHistModuleIdLayer6",
434 "Module Id - Layer 6;Module Id;Entries",
435 fgkSSDMODULESLAYER6,
436 nModuleOffset+fgkSSDMODULESLAYER5 - 0.5,
437 nITSTotalModules + 0.5);
438 fAliITSQADataMakerRec->Add2RecPointsList(fHistModuleIdLayer6,
439 fGenOffset + 1);
440 fSSDhTask += 1;
441 TH1F *fHistLocalXLayer5 = new TH1F("fHistLocalXLayer5",
442 "Local x coord.- Layer 5;x [cm];Entries;",
443 100,-4.,4.);
444 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalXLayer5,
445 fGenOffset + 2);
446 fSSDhTask += 1;
447 TH1F *fHistLocalXLayer6 = new TH1F("fHistLocalXLayer6",
448 "Local x coord.- Layer 6;x [cm];Entries;",
449 100,-4.,4.);
450 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalXLayer6,
451 fGenOffset + 3);
452 fSSDhTask += 1;
453 TH1F *fHistLocalZLayer5 = new TH1F("fHistLocalZLayer5",
454 "Local z coord.- Layer 5;z [cm];Entries;",
455 100,-4.,4.);
456 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalZLayer5,
457 fGenOffset + 4);
458 fSSDhTask += 1;
459 TH1F *fHistLocalZLayer6 = new TH1F("fHistLocalZLayer6",
460 "Local z coord.- Layer 6;z [cm];Entries;",
461 100,-4.,4.);
462 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalZLayer6,
463 fGenOffset + 5);
464 fSSDhTask += 1;
465 TH1F *fHistGlobalXLayer5 = new TH1F("fHistGlobalXLayer5",
466 "Global x - Layer 5;x [cm];Entries;",
467 100,-40.,40.);
468 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalXLayer5,
469 fGenOffset + 6);
470 fSSDhTask += 1;
471 TH1F *fHistGlobalXLayer6 = new TH1F("fHistGlobalXLayer6",
472 "Global x - Layer 6;x [cm];Entries;",
473 100,-45.,45.);
474 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalXLayer6,
475 fGenOffset + 7);
476 fSSDhTask += 1;
477 TH1F *fHistGlobalYLayer5 = new TH1F("fHistGlobalYLayer5",
478 "Global y - Layer 5;y [cm];Entries;",
479 100,-40.,40);
480 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalYLayer5,
481 fGenOffset + 8);
482 fSSDhTask += 1;
483 TH1F *fHistGlobalYLayer6 = new TH1F("fHistGlobalYLayer6",
484 "Global y - Layer 6;y [cm];Entries;",
485 100,-45.,45.);
486 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalYLayer6,
487 fGenOffset + 9);
488 fSSDhTask += 1;
489 TH1F *fHistGlobalZLayer5 = new TH1F("fHistGlobalZLayer5",
490 "Global z - Layer 5;z [cm];Entries;",
491 100,-45.,45);
492 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalZLayer5,
493 fGenOffset + 10);
494 fSSDhTask += 1;
495 TH1F *fHistGlobalZLayer6 = new TH1F("fHistGlobalZLayer6",
496 "Global z - Layer 6;z [cm];Entries;",
497 100,-55.,55.);
498 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalZLayer6,
499 fGenOffset + 11);
500 fSSDhTask += 1;
501 TH1F *fHistPhiLayer5 = new TH1F("fHistPhiLayer5",
502 "#phi - Layer 5;#phi [rad];Entries;",
503 100,-TMath::Pi(),TMath::Pi());
504 fAliITSQADataMakerRec->Add2RecPointsList(fHistPhiLayer5,
505 fGenOffset + 12);
506 fSSDhTask += 1;
507 TH1F *fHistPhiLayer6 = new TH1F("fHistPhiLayer6",
508 "#phi - Layer 6;#phi [rad];Entries;",
509 100,-TMath::Pi(),TMath::Pi());
510 fAliITSQADataMakerRec->Add2RecPointsList(fHistPhiLayer6,
511 fGenOffset + 13);
512 fSSDhTask += 1;
513 TH1F *fHistThetaLayer5 = new TH1F("fHistThetaLayer5",
514 "#theta - Layer 5;#theta [rad];Entries;",
515 100,-TMath::Pi(),TMath::Pi());
516 fAliITSQADataMakerRec->Add2RecPointsList(fHistThetaLayer5,
517 fGenOffset + 14);
518 fSSDhTask += 1;
519 TH1F *fHistThetaLayer6 = new TH1F("fHistThetaLayer6",
520 "#theta - Layer 6;#theta [rad];Entries;",
521 100,-TMath::Pi(),TMath::Pi());
522 fAliITSQADataMakerRec->Add2RecPointsList(fHistThetaLayer6,
523 fGenOffset + 15);
524 fSSDhTask += 1;
525 TH1F *fHistRadiusLayer5 = new TH1F("fHistRadiusLayer5",
526 "r - Layer 5;r [cm];Entries;",
527 100,35.,50.);
528 fAliITSQADataMakerRec->Add2RecPointsList(fHistRadiusLayer5,
529 fGenOffset + 16);
530 fSSDhTask += 1;
531 TH1F *fHistRadiusLayer6 = new TH1F("fHistRadiusLayer6",
532 "r - Layer 6;r [cm];Entries;",
533 100,35.,50.);
534 fAliITSQADataMakerRec->Add2RecPointsList(fHistRadiusLayer6,
535 fGenOffset + 17);
536 fSSDhTask += 1;
537 TH1F *fHistClusterTypeLayer5 = new TH1F("fHistClusterTypeLayer5",
538 "CL type - Layer 5;Cluster type;Entries;",
539 150,0,150);
540 fAliITSQADataMakerRec->Add2RecPointsList(fHistClusterTypeLayer5,
541 fGenOffset + 18);
542 fSSDhTask += 1;
543 TH1F *fHistClusterTypeLayer6 = new TH1F("fHistClusterTypeLayer6",
544 "CL type - Layer 6;Cluster type;Entries;",
545 150,0,150);
546 fAliITSQADataMakerRec->Add2RecPointsList(fHistClusterTypeLayer6,
547 fGenOffset + 19);
548 fSSDhTask += 1;
549 TH1F *fHistChargeRatioLayer5 = new TH1F("fHistChargeRatioLayer5",
550 "Charge ratio - Layer 5;q_{ratio};Entries;",
551 100,-2.0,2.0);
552 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeRatioLayer5,
553 fGenOffset + 20);
554 fSSDhTask += 1;
555 TH1F *fHistChargeRatioLayer6 = new TH1F("fHistChargeRatioLayer6",
556 "Charge ratio - Layer 6;q_{ratio};Entries;",
557 100,-2.0,2.0);
558 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeRatioLayer6,
559 fGenOffset + 21);
560 fSSDhTask += 1;
561 TH1F *fHistChargekeVLayer5 = new TH1F("fHistChargekeVLayer5",
562 "Charge - Layer 5;q [keV];Entries;",
563 100,0.,300.);
564 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargekeVLayer5,
565 fGenOffset + 22);
566 fSSDhTask += 1;
567 TH1F *fHistChargekeVLayer6 = new TH1F("fHistChargekeVLayer6",
568 "Charge - Layer 6;q [keV];Entries;",
569 100,0.,300.);
570 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargekeVLayer6,
571 fGenOffset + 23);
572 fSSDhTask += 1;
573 TH1F *fHistChargeADCLayer5 = new TH1F("fHistChargeADCLayer5",
574 "Charge - Layer 5;q [ADC];Entries;",
575 100,0.,300.);
576 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeADCLayer5,
577 fGenOffset + 24);
578 fSSDhTask += 1;
579 TH1F *fHistChargeADCLayer6 = new TH1F("fHistChargeADCLayer6",
580 "Charge - Layer 6;q [ADC];Entries;",
581 100,0.,300.);
582 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeADCLayer6,
583 fGenOffset + 25);
584 fSSDhTask += 1;
585 TH2F *fHistChargeMapLayer5 = new TH2F("fHistChargeMapLayer5",
586 "Charge map;N_{modules};N_{Ladders}",
587 fgkSSDMODULESPERLADDERLAYER5,
588 -0.5,fgkSSDMODULESPERLADDERLAYER5+0.5,
589 3*fgkSSDLADDERSLAYER5,
590 -0.5,fgkSSDLADDERSLAYER5+0.5);
591 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeMapLayer5,
592 fGenOffset + 26);
593 fSSDhTask += 1;
594 TH2F *fHistChargeMapLayer6 = new TH2F("fHistChargeMapLayer6",
595 "Charge map;N_{modules};N_{Ladders}",
596 fgkSSDMODULESPERLADDERLAYER6,
597 -0.5,fgkSSDMODULESPERLADDERLAYER6+0.5,
598 3*fgkSSDLADDERSLAYER6,
599 -0.5,fgkSSDLADDERSLAYER6+0.5);
600 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeMapLayer6,
601 fGenOffset + 27);
602 fSSDhTask += 1;
603
604// custom code here
605//fSSDhTask must be incremented by one unit every time a histogram is ADDED to the QA List
606
607 AliDebug(1,Form("%d SSD Recs histograms booked\n",fSSDhTask));
608}
609
610//____________________________________________________________________________
611void AliITSQASSDDataMakerRec::MakeRecPoints(TTree *clustersTree)
612{
613 // Fill QA for recpoints - SSD -
614 Int_t fLayer = 0, fLadder = 0, fModule = 0;
615 Int_t lLadderLocationY = 0;
616 TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
617 if (!branchRecP) {
618 AliError("can't get the branch with the ITS clusters !");
619 return;
620 }
621 static TClonesArray statRecpoints("AliITSRecPoint");
622 TClonesArray *recpoints = &statRecpoints;
623 branchRecP->SetAddress(&recpoints);
624 Int_t npoints = 0;
625 Float_t cluglo[3]={0.,0.,0.};
626 for(Int_t module = 0; module < clustersTree->GetEntries(); module++){
627 branchRecP->GetEvent(module);
628 npoints += recpoints->GetEntries();
629 AliITSgeomTGeo::GetModuleId(module,fLayer,fLadder,fModule);
630 lLadderLocationY = 3*fLadder;
631
632 for(Int_t j = 0;j < recpoints->GetEntries(); j++){
633 AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);
634 Int_t layer = recp->GetLayer();
635 recp->GetGlobalXYZ(cluglo);
636 Float_t radius = TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]);
637 Float_t phi = TMath::ATan2(cluglo[1],cluglo[0]);
638 Float_t theta = TMath::ATan2(radius,cluglo[2]);
639 if(layer == 4) {
640 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 0)->Fill(module);
641 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 2)->Fill(recp->GetDetLocalX());
642 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 4)->Fill(recp->GetDetLocalZ());
643 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 6)->Fill(cluglo[0]);
644 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 8)->Fill(cluglo[1]);
645 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 10)->Fill(cluglo[2]);
646 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 12)->Fill(phi);
647 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 14)->Fill(theta);
648 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 16)->Fill(radius);
649 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 18)->Fill(recp->GetType());
650 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 20)->Fill(recp->GetChargeRatio());
651 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 22)->Fill(recp->GetQ());
652 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 26)->SetBinContent(fModule,lLadderLocationY,recp->GetQ());
653 }//layer 5 histograms
654 if(layer == 5) {
655 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 1)->Fill(module);
656 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 3)->Fill(recp->GetDetLocalX());
657 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 5)->Fill(recp->GetDetLocalZ());
658 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 7)->Fill(cluglo[0]);
659 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 9)->Fill(cluglo[1]);
660 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 11)->Fill(cluglo[2]);
661 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 13)->Fill(phi);
662 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 15)->Fill(theta);
663 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 17)->Fill(radius);
664 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 19)->Fill(recp->GetType());
665 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 21)->Fill(recp->GetChargeRatio());
666 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 23)->Fill(recp->GetQ());
667 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 27)->SetBinContent(fModule,lLadderLocationY,recp->GetQ());
668 }//layer 6 histograms
669 }//rec. points loop
670 }//module loop
671
672 statRecpoints.Clear();
673}
674
675//____________________________________________________________________________
676/*void AliITSQASSDDataMakerRec::InitESDs() {
677 // Initialization for ESD - SSD -
678 fESDsOffset = (fAliITSQADataMakerRec->fESDsQAList)->GetEntries();
679 AliDebug(1,Form("Number of ESD histograms (SPD+SDD): %d\n",fRawsOffset));
680
681 TH1F *fHistSSDTrackPt = new TH1F("fHistSSDTrackPt",
682 ";P_{T} [GeV/c];dN/dP_{T}",
683 100,0.1,10.1);
684 fAliITSQADataMakerRec->Add2ESDPointsList(fHistSSDTrackPt,
685 fESDsOffset + 0);
686 fSSDhTask += 1;
687 TH1F *fHistSSDTrackEta = new TH1F("fHistSSDTrackEta",
688 ";#eta;dN/d#eta",
689 40,-2.,2.);
690 fAliITSQADataMakerRec->Add2ESDPointsList(fHistSSDTrackEta,
691 fESDsOffset + 1);
692 fSSDhTask += 1;
693 TH1F *fHistSSDTrackPhi = new TH1F("fHistSSDTrackPhi",
694 ";#phi;dN/d#phi",
695 100,0,2.*TMath::Pi());
696 fAliITSQADataMakerRec->Add2ESDPointsList(fHistSSDTrackPhi,
697 fESDsOffset + 2);
698 fSSDhTask += 1;
699
700 AliDebug(1,Form("%d SSD ESDs histograms booked\n",fSSDhTask));
701 AliInfo(Form("Number of histograms (ITS): %d\n",fESDsOffset+fSSDhTask));
702}
703
704//____________________________________________________________________________
705void AliITSQASSDDataMakerRec::MakeESDs(AliESDEvent * esd) {
706 // make QA data from ESDs
707 Int_t clusters[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
708 for(Int_t j = 0; j < esd->GetNumberOfTracks(); j++) {
709 AliESDtrack *track = esd->GetTrack(j);
710 if ((track->GetStatus() & AliESDtrack::kITSrefit)==0) continue;
711 Int_t nClustersITS = track->GetITSclusters(clusters);
712
713 (fAliITSQADataMakerRec->GetESDsData(fESDsOffset+0))->Fill(track->Pt());
714 (fAliITSQADataMakerRec->GetESDsData(fESDsOffset+1))->Fill(track->Eta());
715 (fAliITSQADataMakerRec->GetESDsData(fESDsOffset+2))->Fill(track->Phi());
716 }//track loop
717}
718*/