1 /**************************************************************************
2 * Copyright(c) 2007-2009, 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 **************************************************************************/
16 // *************************************************************
17 // Checks the quality assurance
18 // by comparing with reference data
20 // -------------------------------------------------------------
21 // W. Ferrarese + P. Cerello Feb 2008
24 // --- ROOT system ---
31 // --- Standard library ---
33 // --- AliRoot header files ---
34 #include "AliITSQASSDDataMakerRec.h"
35 #include "AliQADataMakerRec.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 #include "AliITSBadChannelsSSDv2.h"
47 #include "AliCDBManager.h"
48 #include "AliCDBEntry.h"
50 ClassImp(AliITSQASSDDataMakerRec)
52 AliITSQASSDDataMakerRec::AliITSQASSDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Int_t ldc) :
54 fAliITSQADataMakerRec(aliITSQADataMakerRec),
59 fSSDRawsOffset(0), fSSDRawsDAOffset(0),
60 fSSDRawsCommonLevelOffset(0),
64 // Default constructor
65 //initilize the raw signal vs strip number histograms
67 fCDBManager = AliCDBManager::Instance();
68 //fCDBManager->SetDefaultStorage("local://$ALICE_ROOT");
69 fCDBManager->SetDefaultStorage(gSystem->Getenv("AMORE_CDB_URI"));
70 Int_t runNumber = atoi(gSystem->Getenv("DATE_RUN_NUMBER"));
72 AliInfo("DATE_RUN_NUMBER not defined!!!\n");
74 fCDBManager->SetRun(runNumber);
75 AliCDBEntry *geomGRP = fCDBManager->Get("GRP/Geometry/Data");
76 if(!geomGRP) AliInfo("GRP geometry not found!!!\n");
78 Int_t gLayer = 0,gLadder = 0, gModule = 0;
79 Int_t gHistCounter = 0;
82 for(Int_t iModule = 500; iModule < fgkSSDMODULES + 500; iModule++) {
83 AliITSgeomTGeo::GetModuleId(iModule,gLayer,gLadder,gModule);
84 gTitle = "SSD_RawSignal_Layer"; gTitle += gLayer;
85 gTitle += "_Ladder"; gTitle += gLadder;
86 gTitle += "_Module"; gTitle += gModule;
87 fHistSSDRawSignalModule[gHistCounter] = new TH1D(gTitle.Data(),gTitle.Data(),
88 2*fgkNumberOfPSideStrips,0,2*fgkNumberOfPSideStrips);
91 for(Int_t iStrip = 0; iStrip < 2*fgkNumberOfPSideStrips; iStrip++)
92 fOccupancyMatrix[iModule-500][iStrip] = 0;
96 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++)
97 fHistSSDRawSignalModule[iModule]=NULL;
102 //____________________________________________________________________________
103 AliITSQASSDDataMakerRec::AliITSQASSDDataMakerRec(const AliITSQASSDDataMakerRec& qadm) :
105 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
106 fSSDEvent(qadm.fSSDEvent),
107 fSSDEventPerCycle(qadm.fSSDEventPerCycle),
108 fkOnline(qadm.fkOnline),
110 fSSDRawsOffset(qadm.fSSDRawsOffset), fSSDRawsDAOffset(qadm.fSSDRawsDAOffset),
111 fSSDRawsCommonLevelOffset(qadm.fSSDRawsCommonLevelOffset),
112 fSSDhTask(qadm.fSSDhTask),
113 fGenOffset(qadm.fGenOffset),
114 fCDBManager(qadm.fCDBManager) {
116 fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ;
117 fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
120 //__________________________________________________________________
121 AliITSQASSDDataMakerRec& AliITSQASSDDataMakerRec::operator = (const AliITSQASSDDataMakerRec& qac )
124 this->~AliITSQASSDDataMakerRec();
125 new(this) AliITSQASSDDataMakerRec(qac);
129 //__________________________________________________________________
130 AliITSQASSDDataMakerRec::~AliITSQASSDDataMakerRec() {
132 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++)
133 if(fHistSSDRawSignalModule[iModule]) delete fHistSSDRawSignalModule[iModule];
134 if(fCDBManager) delete fCDBManager;
137 //____________________________________________________________________________
138 void AliITSQASSDDataMakerRec::StartOfDetectorCycle()
140 //Detector specific actions at start of cycle
141 AliDebug(1,"AliITSQADM::Start of SSD Cycle\n");
144 for(Int_t iModule = 500; iModule < fgkSSDMODULES + 500; iModule++) {
145 for(Int_t iStrip = 0; iStrip < 2*fgkNumberOfPSideStrips; iStrip++)
146 fOccupancyMatrix[iModule-500][iStrip] = 0;
149 Int_t gHistPositionOccupancyPerLadder = 0;
150 Int_t gLayer = 0, gLadder = 0, gModule = 0;
151 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
152 AliITSgeomTGeo::GetModuleId(iModule+500,gLayer,gLadder,gModule);
154 gHistPositionOccupancyPerLadder = (gLayer == 5) ? 2*(gLadder - 1) : 2*(gLadder - 1 + fgkSSDLADDERSLAYER5);
157 fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder)->Reset();
159 fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder+1)->Reset();
161 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6))->Reset();
162 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+1))->Reset();
167 //____________________________________________________________________________
168 void AliITSQASSDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
170 // launch the QA checking
171 AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n");
176 MonitorOCDBObjects();
178 Int_t gHistPositionOccupancyPerModule = 0;
179 Int_t gLayer = 0, gLadder = 0, gModule = 0;
180 //occupancy per module
181 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
182 AliITSgeomTGeo::GetModuleId(iModule+500,gLayer,gLadder,gModule);
184 gHistPositionOccupancyPerModule = (gLayer == 5) ? ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + gModule - 1) : ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + gModule + fgkSSDMODULESLAYER5 - 1);
185 for(Int_t iBins = 1; iBins < fHistSSDRawSignalModule[iModule]->GetXaxis()->GetNbins(); iBins++)
186 fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule)->SetBinContent(iBins,fOccupancyMatrix[iModule][iBins-1]);
188 if(fSSDEventPerCycle != 0)
189 ((TH1D *)(fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule)))->Scale(100./fSSDEventPerCycle);
192 //occupancy per ladder
193 Int_t gHistPositionOccupancyPerLadder = 0;
194 Int_t lLadderLocationY = 0;
195 Double_t occupancy = 0.0;
196 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
197 AliITSgeomTGeo::GetModuleId(iModule+500,gLayer,gLadder,gModule);
199 gHistPositionOccupancyPerModule = (gLayer == 5) ? ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + gModule - 1) : ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + gModule + fgkSSDMODULESLAYER5 - 1);
200 gHistPositionOccupancyPerLadder = (gLayer == 5) ? 2*(gLadder - 1) : 2*(gLadder - 1 + fgkSSDLADDERSLAYER5);
203 occupancy = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),0);
204 fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder)->Fill(gModule,occupancy);
205 lLadderLocationY = 3*gLadder; // sideP=1 sideN=0
207 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6))->SetBinContent(gModule,lLadderLocationY,occupancy);
209 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+1))->SetBinContent(gModule,lLadderLocationY,occupancy);
212 occupancy = GetOccupancyModule((TH1 *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+gHistPositionOccupancyPerModule),1);
213 fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+gHistPositionOccupancyPerLadder+1)->Fill(gModule,occupancy);
215 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6))->SetBinContent(gModule,lLadderLocationY-1,occupancy);
217 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsCommonLevelOffset+fgkSSDMODULES+2*fgkSSDLADDERSLAYER5+2*fgkSSDLADDERSLAYER6+1))->SetBinContent(gModule,lLadderLocationY-1,occupancy);
219 }//online flag for SSD
221 fSSDEventPerCycle = 0;
223 AliQAChecker::Instance()->Run( AliQA::kITS , task, list);
226 //____________________________________________________________________________
227 void AliITSQASSDDataMakerRec::InitRaws() {
228 // Initialization for RAW data - SSD -
229 fGenOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();
232 AliInfo("Book Online Histograms for SSD\n");
235 AliInfo("Book Offline Histograms for SSD\n ");
237 AliInfo(Form("Number of histograms (SPD+SDD): %d\n",fGenOffset));
239 //book online-offline QA histos
240 TH1F *fHistSSDEventType = new TH1F("fHistSSDEventType",
241 ";Event type;Events",
243 fAliITSQADataMakerRec->Add2RawsList(fHistSSDEventType,
244 fGenOffset+fSSDRawsOffset);
246 TH1F *fHistSSDDataSize = new TH1F("fHistSSDDataSize",
247 ";log(SSD data size) [Bytes];Events",
249 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSize,
250 fGenOffset+fSSDRawsOffset);
252 TH1F *fHistSSDDataSizePercentage = new TH1F("fHistSSDDataSizePercentage",
253 ";SSD data size [%];Events",
255 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePercentage,
256 fGenOffset+fSSDRawsOffset);
258 TH1F *fHistSSDDDLId = new TH1F("fHistSSDDDLId",
259 ";DDL id;Events",20,510.5,530.5);
260 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDDLId,
261 fGenOffset+fSSDRawsOffset);
263 TH1F *fHistSSDDataSizePerDDL = new TH1F("fHistSSDDataSizePerDDL",
264 ";DDL id;<SSD data size> [MB]",
266 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePerDDL,
267 fGenOffset+fSSDRawsOffset);
269 TH1F *fHistSSDDataSizeDDL[fgkNumOfDDLs];
270 for(Int_t i = 1; i < fgkNumOfDDLs+1; i++) {
271 gTitle = "fHistSSDDataSizeDDL"; gTitle += i+511;
272 fHistSSDDataSizeDDL[i-1] = new TH1F(gTitle.Data(),
273 ";log(SSD data size) [Bytes];Events",
275 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizeDDL[i-1],
276 fGenOffset+fSSDRawsOffset);
280 TH1F *fHistSSDLDCId = new TH1F("fHistSSDLDCId",";LDC id;Events",10,0.5,10.5);
281 fAliITSQADataMakerRec->Add2RawsList(fHistSSDLDCId,
282 fGenOffset+fSSDRawsOffset);
284 TH1F *fHistSSDDataSizePerLDC = new TH1F("fHistSSDDataSizePerLDC",
285 ";LDC id;<SSD data size> [MB]",
287 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizePerLDC,
288 fGenOffset+fSSDRawsOffset);
290 TH1F *fHistSSDDataSizeLDC[fgkNumOfLDCs];
291 for(Int_t i = 1; i < fgkNumOfLDCs+1; i++) {
292 gTitle = "fHistSSDDataSizeLDC"; gTitle += i;
293 fHistSSDDataSizeLDC[i-1] = new TH1F(gTitle.Data(),
294 ";log(SSD data size) [Bytes];Events",
296 fAliITSQADataMakerRec->Add2RawsList(fHistSSDDataSizeLDC[i-1],
297 fGenOffset+fSSDRawsOffset);
300 fSSDRawsCommonLevelOffset = fSSDRawsOffset;
303 Int_t gLayer = 0, gLadder = 0, gModule = 0;
304 //occupancy per SSD module
305 TH1D *fHistSSDOccupancyModule[fgkSSDMODULES];
306 for(Int_t i = 500; i < fgkSSDMODULES + 500; i++) {
307 AliITSgeomTGeo::GetModuleId(i,gLayer,gLadder,gModule);
308 gTitle = "fHistSSD_Occupancy_Layer"; gTitle += gLayer;
311 gTitle += 499+gLadder;
313 gTitle += 599+gLadder;
314 gTitle += "_Module"; gTitle += gModule;
315 fHistSSDOccupancyModule[i-500] = new TH1D(gTitle.Data(),gTitle.Data(),
316 2*fgkNumberOfPSideStrips,0,2*fgkNumberOfPSideStrips);
317 fHistSSDOccupancyModule[i-500]->GetXaxis()->SetTitleColor(1);
318 fHistSSDOccupancyModule[i-500]->GetXaxis()->SetTitle("N_{strip}");
319 fHistSSDOccupancyModule[i-500]->GetYaxis()->SetTitle("Occupancy [%]");
320 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyModule[i-500], fGenOffset+fSSDRawsOffset);
324 //Occupancy per SSD ladder
325 Int_t occupancyCounter = 0;
326 TH1D *fHistSSDOccupancyLadder[2*(fgkSSDLADDERSLAYER5 + fgkSSDLADDERSLAYER6)];
327 for(Int_t iLayer = 5; iLayer < 7; iLayer++) {
328 for(Int_t iLadder = 1; iLadder < AliITSgeomTGeo::GetNLadders(iLayer) + 1; iLadder++) {
329 //P-side occupancy plots
330 gTitle = "fHistSSD_Occupancy_Layer"; gTitle += iLayer;
333 gTitle += 499+iLadder;
335 gTitle += 599+iLadder;
337 fHistSSDOccupancyLadder[occupancyCounter] = new TH1D(gTitle.Data(),
339 AliITSgeomTGeo::GetNDetectors(iLayer),
340 0.5,AliITSgeomTGeo::GetNDetectors(iLayer)+0.5);
341 fHistSSDOccupancyLadder[occupancyCounter]->GetXaxis()->SetTitleColor(1);
342 fHistSSDOccupancyLadder[occupancyCounter]->GetXaxis()->SetTitle("Module number");
343 fHistSSDOccupancyLadder[occupancyCounter]->GetYaxis()->SetTitle("Occupancy [%]");
344 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLadder[occupancyCounter],
345 fGenOffset+fSSDRawsOffset);
346 occupancyCounter += 1; fSSDRawsOffset += 1;
347 //N-side occupancy plots
348 gTitle = "fHistSSD_Occupancy_Layer"; gTitle += iLayer;
351 gTitle += 499+iLadder;
353 gTitle += 599+iLadder;
355 fHistSSDOccupancyLadder[occupancyCounter] = new TH1D(gTitle.Data(),
357 AliITSgeomTGeo::GetNDetectors(iLayer),
358 0.5,AliITSgeomTGeo::GetNDetectors(iLayer)+0.5);
359 fHistSSDOccupancyLadder[occupancyCounter]->GetXaxis()->SetTitleColor(1);
360 fHistSSDOccupancyLadder[occupancyCounter]->GetXaxis()->SetTitle("Module number");
361 fHistSSDOccupancyLadder[occupancyCounter]->GetYaxis()->SetTitle("Occupancy [%]");
362 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLadder[occupancyCounter],
363 fGenOffset+fSSDRawsOffset);
364 occupancyCounter += 1; fSSDRawsOffset += 1;
368 //top level occupancy plots
369 TH2D *fHistSSDOccupancyLayer5 = new TH2D("fHistSSDOccupancyLayer5",
370 ";N_{modules};N_{Ladders}",
371 fgkSSDMODULESPERLADDERLAYER5,
372 0,fgkSSDMODULESPERLADDERLAYER5,
373 3*fgkSSDLADDERSLAYER5,
374 0,fgkSSDLADDERSLAYER5);
376 for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER5 + 1; iBin++){
377 sprintf(fLabel,"%d",iBin);
378 fHistSSDOccupancyLayer5->GetXaxis()->SetBinLabel(iBin,fLabel);
380 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLayer5, fGenOffset+fSSDRawsOffset);
382 TH2D *fHistSSDOccupancyLayer6 = new TH2D("fHistSSDOccupancyLayer6",
383 ";N_{modules};N_{Ladders}",
384 fgkSSDMODULESPERLADDERLAYER6,
385 0,fgkSSDMODULESPERLADDERLAYER6,
386 3*fgkSSDLADDERSLAYER6,
387 0,fgkSSDLADDERSLAYER6);
388 for(Int_t iBin = 1; iBin < fgkSSDMODULESPERLADDERLAYER6 + 1; iBin++){
389 sprintf(fLabel,"%d",iBin);
390 fHistSSDOccupancyLayer6->GetXaxis()->SetBinLabel(iBin,fLabel);
392 fAliITSQADataMakerRec->Add2RawsList(fHistSSDOccupancyLayer6, fGenOffset+fSSDRawsOffset);
396 TH2D *fHistPSideBadChannelMapLayer5 = new TH2D("fHistPSideBadChannelMapLayer5",
397 "Layer 5;N_{module};N_{ladder}",
400 fHistPSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
401 fHistPSideBadChannelMapLayer5->SetStats(kFALSE);
402 fHistPSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
403 fHistPSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
404 fHistPSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
405 fHistPSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
406 fHistPSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
407 fHistPSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
408 fHistPSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
409 fAliITSQADataMakerRec->Add2RawsList(fHistPSideBadChannelMapLayer5, fGenOffset+fSSDRawsOffset);
410 fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
412 TH2D *fHistNSideBadChannelMapLayer5 = new TH2D("fHistNSideBadChannelMapLayer5",
413 "Layer 5;N_{module};N_{ladder}",
416 fHistNSideBadChannelMapLayer5->GetXaxis()->SetTitleColor(1);
417 fHistNSideBadChannelMapLayer5->SetStats(kFALSE);
418 fHistNSideBadChannelMapLayer5->GetYaxis()->SetTitleOffset(1.8);
419 fHistNSideBadChannelMapLayer5->GetXaxis()->SetNdivisions(22);
420 fHistNSideBadChannelMapLayer5->GetYaxis()->SetNdivisions(34);
421 fHistNSideBadChannelMapLayer5->GetXaxis()->SetLabelSize(0.03);
422 fHistNSideBadChannelMapLayer5->GetYaxis()->SetLabelSize(0.03);
423 fHistNSideBadChannelMapLayer5->GetZaxis()->SetTitleOffset(1.6);
424 fHistNSideBadChannelMapLayer5->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
425 fAliITSQADataMakerRec->Add2RawsList(fHistNSideBadChannelMapLayer5, fGenOffset+fSSDRawsOffset);
426 fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
428 TH2D *fHistPSideBadChannelMapLayer6 = new TH2D("fHistPSideBadChannelMapLayer6",
429 "Layer 6;N_{module};N_{ladder}",
432 fHistPSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
433 fHistPSideBadChannelMapLayer6->SetStats(kFALSE);
434 fHistPSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
435 fHistPSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
436 fHistPSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
437 fHistPSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
438 fHistPSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
439 fHistPSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
440 fHistPSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (p-side)[%]");
441 fAliITSQADataMakerRec->Add2RawsList(fHistPSideBadChannelMapLayer6, fGenOffset+fSSDRawsOffset);
442 fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
444 TH2D *fHistNSideBadChannelMapLayer6 = new TH2D("fHistNSideBadChannelMapLayer6",
445 "Layer 6;N_{module};N_{ladder}",
448 fHistNSideBadChannelMapLayer6->GetXaxis()->SetTitleColor(1);
449 fHistNSideBadChannelMapLayer6->SetStats(kFALSE);
450 fHistNSideBadChannelMapLayer6->GetYaxis()->SetTitleOffset(1.8);
451 fHistNSideBadChannelMapLayer6->GetXaxis()->SetNdivisions(25);
452 fHistNSideBadChannelMapLayer6->GetYaxis()->SetNdivisions(38);
453 fHistNSideBadChannelMapLayer6->GetXaxis()->SetLabelSize(0.03);
454 fHistNSideBadChannelMapLayer6->GetYaxis()->SetLabelSize(0.03);
455 fHistNSideBadChannelMapLayer6->GetZaxis()->SetTitleOffset(1.6);
456 fHistNSideBadChannelMapLayer6->GetZaxis()->SetTitle("Bad channels (n-side)[%]");
457 fAliITSQADataMakerRec->Add2RawsList(fHistNSideBadChannelMapLayer6, fGenOffset+fSSDRawsOffset);
458 fSSDRawsOffset += 1; fSSDRawsDAOffset += 1;
460 fSSDhTask = fSSDRawsOffset;
461 AliDebug(1,Form("%d SSD Raws histograms booked\n",fSSDhTask));
462 AliInfo(Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenOffset+fSSDhTask));
463 AliDebug(1,Form("Number of histograms (SPD+SDD+SSD): %d\n",fGenOffset+fSSDRawsOffset));
466 //____________________________________________________________________________
467 void AliITSQASSDDataMakerRec::MakeRaws(AliRawReader* rawReader) {
468 // Fill QA for RAW - SSD -
471 Int_t gLayer = 0,gLadder = 0, gModule = 0;
473 Double_t gSizePerDDL[fgkNumOfDDLs] = {0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.};
474 Double_t gSizePerLDC[fgkNumOfLDCs] = {0.,0.,0.};
475 Double_t sumSSDDataSize = 0.0;
476 Double_t eventSize = -1.0;
479 //reset the signal vs strip number histograms
480 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++)
481 fHistSSDRawSignalModule[iModule]->Reset();
484 rawReader->Select("ITSSSD",-1,-1);
485 rawReader->Reset(); //rawReader->NextEvent();
486 (fAliITSQADataMakerRec->GetRawsData(fGenOffset))->Fill(rawReader->GetType());
487 if(rawReader->GetType() == 7) {
489 fSSDEventPerCycle += 1;
492 AliITSRawStreamSSD gSSDStream(rawReader);
493 AliRawReaderRoot *rootReader = (AliRawReaderRoot *)rawReader;
495 const AliRawEventHeaderBase *header = rootReader->GetEventHeader();
497 eventSize = header->GetEventSize();
499 while (gSSDStream.Next()) {
500 if(gSSDStream.GetModuleID() < 0) continue;
501 /*cout<<"DDL: "<<rawReader->GetDDLID()<<
502 " - LDC: "<<rawReader->GetLDCId()<<
503 " - Size: "<<rawReader->GetDataSize()<<
504 " - Equipment size: "<<rawReader->GetEquipmentSize()<<endl;*/
505 gSizePerDDL[rawReader->GetDDLID()] = rawReader->GetDataSize();
506 gSizePerLDC[rawReader->GetLDCId()-6] = rawReader->GetDataSize();
507 AliITSgeomTGeo::GetModuleId(gSSDStream.GetModuleID(),gLayer,gLadder,gModule);
508 if(gSSDStream.GetStrip() < 0) continue;
509 gStripNumber = (gSSDStream.GetSideFlag() == 0) ? gSSDStream.GetStrip() : -gSSDStream.GetStrip() + 2*fgkNumberOfPSideStrips;
510 gHistPosition = (gLayer == 5) ? ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER5 + gModule - 1) : ((gLadder - 1)*fgkSSDMODULESPERLADDERLAYER6 + gModule + fgkSSDMODULESLAYER5 - 1);
511 //AliInfo(Form("ModulePosition: %d - Layer: %d - Ladder: %d - Module: %d\n",gHistPosition,gLayer,gLadder,gModule));
513 fHistSSDRawSignalModule[gHistPosition]->Fill(gStripNumber,gSSDStream.GetSignal());
514 //fAliITSQADataMakerRec->GetRawsData(fGenOffset+gHistPosition+fSSDRawsCommonLevelOffset)->Fill(gStripNumber,gSSDStream.GetSignal());
517 //event size calculation and filling info
518 for(Int_t i = 0; i < fgkNumOfDDLs; i++) {
519 sumSSDDataSize += gSizePerDDL[i];
520 if(gSizePerDDL[i] > 0) {
521 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+3))->Fill(i+512);
522 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+5+i))->Fill(TMath::Log10(gSizePerDDL[i]));
524 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+4))->Fill(i+512,gSizePerDDL[i]/1e+06);
526 for(Int_t i = 0; i < fgkNumOfLDCs; i++) {
527 if(gSizePerLDC[i] > 0) {
528 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+21))->Fill(i+6);
529 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+23+i))->Fill(TMath::Log10(gSizePerLDC[i]));
530 //cout<<"Event: "<<fSSDEventPerCycle<<" - LDC: "<<i+6<<
531 //" - Data size: "<<gSizePerLDC[i]<<endl;
533 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+22))->Fill(i+6,gSizePerLDC[i]/1e+06);
536 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+1))->Fill(TMath::Log10(sumSSDDataSize));
538 (fAliITSQADataMakerRec->GetRawsData(fGenOffset+2))->Fill(100.*sumSSDDataSize/eventSize);
540 //Occupancy calculation
542 for(Int_t iModule = 0; iModule < fgkSSDMODULES; iModule++) {
543 GetOccupancyStrip(fHistSSDRawSignalModule[iModule],fOccupancyMatrix[iModule]);
544 //if(iModule == 156) {
545 //cout<<"========================================"<<endl;
546 //AliITSgeomTGeo::GetModuleId(656,gLayer,gLadder,gModule);
547 /*for(Int_t iBin = 1; iBin < fHistSSDRawSignalModule[iModule]->GetXaxis()->GetNbins(); iBin++) {
548 if((iBin >= 750)&&(iBin <= 780))
549 cout<<"Event: "<<fSSDEventPerCycle<<" - Ladder: "<<gLadder+499<<
550 " - Module: "<<gModule<<" - Strip: "<<iBin<<
551 " - Signal: "<<fHistSSDRawSignalModule[iModule]->GetBinContent(iBin)<<endl;
552 }*///strip loop --> to be removed
553 //}//module cut --> to be removed
555 }//online flag for SSD
558 //____________________________________________________________________________ //
559 void AliITSQASSDDataMakerRec::GetOccupancyStrip(TH1 *lHisto, Int_t *occupancyMatrix) {
560 //Increments the entries in the occupancy matrix based
561 //on whether the signal for each strip is larger than the cutoff
562 //Currently the cutoff is at 0 which means that if ZS
563 //works, whatever comes from the FEROM is considered as "signal"
564 // Double_t cutoff = 0.0;
565 TString histname = lHisto->GetName();
566 //cout<<histname.Data()<<endl;
567 //if(histname.Contains("Layer5_Ladder8_Module3")) {
568 for(Int_t iBin = 1; iBin < lHisto->GetXaxis()->GetNbins(); iBin++) {
569 Double_t y = lHisto->GetBinContent(iBin);
571 occupancyMatrix[iBin-1] += 1;
573 //if((iBin >= 750)&&(iBin <= 780))
574 //cout<<"Strip: "<<iBin<<
576 //" - Occupancy: "<<occupancyMatrix[iBin-1]<<endl;
582 //____________________________________________________________________________
583 Double_t AliITSQASSDDataMakerRec::GetOccupancyModule(TH1 *lHisto, Int_t stripside) {
584 // bo: TDC >0 or # of sigmas wrt noise ?
585 //stripside == 0 --> P-side
586 //stripside == 1 --> N-side
587 Int_t lNumFiredBins = 0;
588 TString histname = lHisto->GetName();
589 for(Int_t iBin = 1 + stripside*fgkNumberOfPSideStrips; iBin < fgkNumberOfPSideStrips*(1 + stripside); iBin++){
590 /*if(histname.Contains("Layer5_Ladder507_Module3")) {
591 cout<<lHisto->GetName()<<
593 " - Bin content: "<<lHisto->GetBinContent(iBin)<<endl;
595 if (lHisto->GetBinContent(iBin) > 0)
599 Double_t lOccupancy = (100.*lNumFiredBins)/fgkNumberOfPSideStrips; // percentage
601 /*if(histname.Contains("Layer5_Ladder507_Module3"))
602 cout<<"Fired strips: "<<lNumFiredBins<<
603 " - Occupancy: "<<lOccupancy<<endl;*/
604 //AliInfo(Form("Fired strips: %d - Total strips: %d - Occupancy :%lf\n",lNumFiredBins,lHisto->GetNbinsX(),lOccupancy));
609 //____________________________________________________________________________
610 void AliITSQASSDDataMakerRec::MonitorOCDBObjects() {
611 //Monitor in AMORE the output of the DA
612 //Currently only the bad channel list is monitored
613 //Todo: Noise - Pedestal
614 AliCDBEntry *entryBadChannelsSSD = fCDBManager->Get("ITS/Calib/BadChannelsSSD");
615 if(!entryBadChannelsSSD)
616 AliError("OCDB entry for the bad channel list is not valid!");
617 AliITSBadChannelsSSDv2 *badchannelsSSD = (AliITSBadChannelsSSDv2 *)entryBadChannelsSSD->GetObject();
619 AliError("Bad channel list object is not a valid AliITSBadChannelsSSD object!");
621 //_____________________________________________________________________________//
622 Int_t nBadPSideChannels = 0, nBadNSideChannels = 0;
623 Int_t layer = 0, ladder = 0, module = 0;
624 Int_t nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
625 Int_t nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
626 //_____________________________________________________________________________//
628 for(Int_t i = 0; i < fgkSSDMODULES; i++) {
629 AliITSgeomTGeo::GetModuleId(i+500,layer,ladder,module);
630 nBadPSideChannels = 0, nBadNSideChannels = 0;
631 nPSideChannelsLayer5 = 0, nNSideChannelsLayer5 = 0;
632 nPSideChannelsLayer6 = 0, nNSideChannelsLayer6 = 0;
634 Int_t badChannel = 0;
635 for(Int_t j = 0; j < fgkNumberOfPSideStrips; j++) {
636 badChannel = (Int_t)(badchannelsSSD->GetBadChannelP(i,j));
637 if(badChannel != 0) {
639 nPSideChannelsLayer5 += 1;
641 nPSideChannelsLayer6 += 1;
642 nBadPSideChannels += 1;
643 }//badchannel flag != 0
644 badChannel = (Int_t)(badchannelsSSD->GetBadChannelN(i,j));
645 if(badChannel != 0) {
647 nNSideChannelsLayer5 += 1;
649 nNSideChannelsLayer6 += 1;
650 nBadNSideChannels += 1;
651 }//badchannel flag != 0
654 if(nPSideChannelsLayer5 > 0)
655 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset))->Fill(module,499+ladder,
656 100.*nPSideChannelsLayer5/fgkNumberOfPSideStrips);
657 else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset))->Fill(module,499+ladder,0.0001);
658 if(nNSideChannelsLayer5 > 0)
659 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset+1))->Fill(module,499+ladder,
660 100.*nNSideChannelsLayer5/fgkNumberOfPSideStrips);
661 else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset+1))->Fill(module,499+ladder,0.0001);
664 if(nPSideChannelsLayer6 > 0)
665 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset+2))->Fill(module,599+ladder,
666 100.*nPSideChannelsLayer6/fgkNumberOfPSideStrips);
667 else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset+2))->Fill(module,599+ladder,0.0001);
668 if(nNSideChannelsLayer6 > 0)
669 ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset+3))->Fill(module,599+ladder,
670 100.*nNSideChannelsLayer6/fgkNumberOfPSideStrips);
671 else ((TH2D *)fAliITSQADataMakerRec->GetRawsData(fGenOffset+fSSDRawsOffset-fSSDRawsDAOffset+3))->Fill(module,599+ladder,0.0001);
676 //____________________________________________________________________________
677 void AliITSQASSDDataMakerRec::InitRecPoints()
679 // Initialization for RECPOINTS - SSD -
681 fGenOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();
682 Int_t nModuleOffset = 500;
683 Int_t nITSTotalModules = AliITSgeomTGeo::GetNModules();
685 TH1F *fHistModuleIdLayer5 = new TH1F("fHistModuleIdLayer5",
686 "Module Id - Layer 5;Module Id;Entries",
689 nITSTotalModules-fgkSSDMODULESLAYER6+0.5);
690 fAliITSQADataMakerRec->Add2RecPointsList(fHistModuleIdLayer5,
693 TH1F *fHistModuleIdLayer6 = new TH1F("fHistModuleIdLayer6",
694 "Module Id - Layer 6;Module Id;Entries",
696 nModuleOffset+fgkSSDMODULESLAYER5 - 0.5,
697 nITSTotalModules + 0.5);
698 fAliITSQADataMakerRec->Add2RecPointsList(fHistModuleIdLayer6,
701 TH1F *fHistClusterPerEventLayer5 = new TH1F("fHistClusterPerEventLayer5",
702 "N_{clusters} - Layer 5;N_{clusters};Entries;",
704 fAliITSQADataMakerRec->Add2RecPointsList(fHistClusterPerEventLayer5,
707 TH1F *fHistClusterPerEventLayer6 = new TH1F("fHistClusterPerEventLayer6",
708 "N_{clusters} - Layer 6;N_{clusters};Entries;",
710 fAliITSQADataMakerRec->Add2RecPointsList(fHistClusterPerEventLayer6,
713 TH1F *fHistLocalXLayer5 = new TH1F("fHistLocalXLayer5",
714 "Local x coord.- Layer 5;x [cm];Entries;",
716 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalXLayer5,
719 TH1F *fHistLocalXLayer6 = new TH1F("fHistLocalXLayer6",
720 "Local x coord.- Layer 6;x [cm];Entries;",
722 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalXLayer6,
725 TH1F *fHistLocalZLayer5 = new TH1F("fHistLocalZLayer5",
726 "Local z coord.- Layer 5;z [cm];Entries;",
728 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalZLayer5,
731 TH1F *fHistLocalZLayer6 = new TH1F("fHistLocalZLayer6",
732 "Local z coord.- Layer 6;z [cm];Entries;",
734 fAliITSQADataMakerRec->Add2RecPointsList(fHistLocalZLayer6,
737 TH1F *fHistGlobalXLayer5 = new TH1F("fHistGlobalXLayer5",
738 "Global x - Layer 5;x [cm];Entries;",
740 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalXLayer5,
743 TH1F *fHistGlobalXLayer6 = new TH1F("fHistGlobalXLayer6",
744 "Global x - Layer 6;x [cm];Entries;",
746 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalXLayer6,
749 TH1F *fHistGlobalYLayer5 = new TH1F("fHistGlobalYLayer5",
750 "Global y - Layer 5;y [cm];Entries;",
752 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalYLayer5,
755 TH1F *fHistGlobalYLayer6 = new TH1F("fHistGlobalYLayer6",
756 "Global y - Layer 6;y [cm];Entries;",
758 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalYLayer6,
761 TH1F *fHistGlobalZLayer5 = new TH1F("fHistGlobalZLayer5",
762 "Global z - Layer 5;z [cm];Entries;",
764 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalZLayer5,
767 TH1F *fHistGlobalZLayer6 = new TH1F("fHistGlobalZLayer6",
768 "Global z - Layer 6;z [cm];Entries;",
770 fAliITSQADataMakerRec->Add2RecPointsList(fHistGlobalZLayer6,
773 TH1F *fHistPhiLayer5 = new TH1F("fHistPhiLayer5",
774 "#phi - Layer 5;#phi [rad];Entries;",
775 100,-TMath::Pi(),TMath::Pi());
776 fAliITSQADataMakerRec->Add2RecPointsList(fHistPhiLayer5,
779 TH1F *fHistPhiLayer6 = new TH1F("fHistPhiLayer6",
780 "#phi - Layer 6;#phi [rad];Entries;",
781 100,-TMath::Pi(),TMath::Pi());
782 fAliITSQADataMakerRec->Add2RecPointsList(fHistPhiLayer6,
785 TH1F *fHistThetaLayer5 = new TH1F("fHistThetaLayer5",
786 "#theta - Layer 5;#theta [rad];Entries;",
787 100,-TMath::Pi(),TMath::Pi());
788 fAliITSQADataMakerRec->Add2RecPointsList(fHistThetaLayer5,
791 TH1F *fHistThetaLayer6 = new TH1F("fHistThetaLayer6",
792 "#theta - Layer 6;#theta [rad];Entries;",
793 100,-TMath::Pi(),TMath::Pi());
794 fAliITSQADataMakerRec->Add2RecPointsList(fHistThetaLayer6,
797 TH1F *fHistRadiusLayer5 = new TH1F("fHistRadiusLayer5",
798 "r - Layer 5;r [cm];Entries;",
800 fAliITSQADataMakerRec->Add2RecPointsList(fHistRadiusLayer5,
803 TH1F *fHistRadiusLayer6 = new TH1F("fHistRadiusLayer6",
804 "r - Layer 6;r [cm];Entries;",
806 fAliITSQADataMakerRec->Add2RecPointsList(fHistRadiusLayer6,
809 TH1F *fHistClusterTypeLayer5 = new TH1F("fHistClusterTypeLayer5",
810 "CL type - Layer 5;Cluster type;Entries;",
812 fAliITSQADataMakerRec->Add2RecPointsList(fHistClusterTypeLayer5,
815 TH1F *fHistClusterTypeLayer6 = new TH1F("fHistClusterTypeLayer6",
816 "CL type - Layer 6;Cluster type;Entries;",
818 fAliITSQADataMakerRec->Add2RecPointsList(fHistClusterTypeLayer6,
821 TH1F *fHistChargeRatioLayer5 = new TH1F("fHistChargeRatioLayer5",
822 "Charge ratio - Layer 5;q_{ratio};Entries;",
824 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeRatioLayer5,
827 TH1F *fHistChargeRatioLayer6 = new TH1F("fHistChargeRatioLayer6",
828 "Charge ratio - Layer 6;q_{ratio};Entries;",
830 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeRatioLayer6,
833 TH1F *fHistChargekeVLayer5 = new TH1F("fHistChargekeVLayer5",
834 "Charge - Layer 5;q [keV];Entries;",
836 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargekeVLayer5,
839 TH1F *fHistChargekeVLayer6 = new TH1F("fHistChargekeVLayer6",
840 "Charge - Layer 6;q [keV];Entries;",
842 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargekeVLayer6,
845 TH1F *fHistChargePSideLayer5 = new TH1F("fHistChargePSideLayer5",
846 "Charge P- Layer 5;q_{P} [keV];Entries;",
848 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargePSideLayer5,
851 TH1F *fHistChargePSideLayer6 = new TH1F("fHistChargePSideLayer6",
852 "Charge P- Layer 6;q_{P} [keV];Entries;",
854 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargePSideLayer6,
857 TH1F *fHistChargeNSideLayer5 = new TH1F("fHistChargeNSideLayer5",
858 "Charge N- Layer 5;q_{N} [keV];Entries;",
860 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeNSideLayer5,
863 TH1F *fHistChargeNSideLayer6 = new TH1F("fHistChargeNSideLayer6",
864 "Charge N- Layer 6;q_{N} [keV];Entries;",
866 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeNSideLayer6,
869 TH1F *fHistChargeRatio2Layer5 = new TH1F("fHistChargeRatio2Layer5",
870 "Charge Ratio qN/qP - Layer 5;q_{N}/q_{P};Entries;",
872 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeRatio2Layer5,
875 TH1F *fHistChargeRatio2Layer6 = new TH1F("fHistChargeRatio2Layer6",
876 "Charge Ratio qN/qP - Layer 6;q_{N}/q_{P};Entries;",
878 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeRatio2Layer6,
881 TH2F *fHistChargePNSideLayer5 = new TH2F("fHistChargePNSideLayer5",
882 "Charge correlation - Layer 5;q_{P} [keV];q_{N} [keV]",
885 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargePNSideLayer5,
888 TH2F *fHistChargePNSideLayer6 = new TH2F("fHistChargePNSideLayer6",
889 "Charge correlation - Layer 6;q_{P} [keV];q_{N} [keV]",
892 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargePNSideLayer6,
895 TH2F *fHistChargeMapLayer5 = new TH2F("fHistChargeMapLayer5",
896 "Charge map;N_{modules};N_{Ladders}",
897 fgkSSDMODULESPERLADDERLAYER5,
898 -0.5,fgkSSDMODULESPERLADDERLAYER5+0.5,
899 3*fgkSSDLADDERSLAYER5,
900 -0.5,fgkSSDLADDERSLAYER5+0.5);
901 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeMapLayer5,
904 TH2F *fHistChargeMapLayer6 = new TH2F("fHistChargeMapLayer6",
905 "Charge map;N_{modules};N_{Ladders}",
906 fgkSSDMODULESPERLADDERLAYER6,
907 -0.5,fgkSSDMODULESPERLADDERLAYER6+0.5,
908 3*fgkSSDLADDERSLAYER6,
909 -0.5,fgkSSDLADDERSLAYER6+0.5);
910 fAliITSQADataMakerRec->Add2RecPointsList(fHistChargeMapLayer6,
914 AliDebug(1,Form("%d SSD Recs histograms booked\n",fSSDhTask));
917 //____________________________________________________________________________
918 void AliITSQASSDDataMakerRec::MakeRecPoints(TTree *clustersTree)
920 // Fill QA for recpoints - SSD -
922 Int_t gLayer = 0, gLadder = 0, gModule = 0;
923 Int_t lLadderLocationY = 0;
924 TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
926 AliError("can't get the branch with the ITS clusters !");
929 static TClonesArray statRecpoints("AliITSRecPoint");
930 TClonesArray *recpoints = &statRecpoints;
931 branchRecP->SetAddress(&recpoints);
932 Int_t nClustersLayer5 = 0, nClustersLayer6 = 0;
934 Float_t cluglo[3]={0.,0.,0.};
935 for(Int_t module = 0; module < clustersTree->GetEntries(); module++){
936 branchRecP->GetEvent(module);
937 npoints += recpoints->GetEntries();
938 AliITSgeomTGeo::GetModuleId(module,gLayer,gLadder,gModule);
939 lLadderLocationY = 3*gLadder;
941 for(Int_t j = 0;j < recpoints->GetEntries(); j++){
942 AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);
943 Int_t layer = recp->GetLayer();
944 recp->GetGlobalXYZ(cluglo);
945 Float_t radius = TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]);
946 Float_t phi = TMath::ATan2(cluglo[1],cluglo[0]);
947 Float_t theta = TMath::ATan2(radius,cluglo[2]);
948 Double_t chargeRatio = recp->GetChargeRatio();
949 Double_t clusterCharge = recp->GetQ();
950 Double_t chargePSide = clusterCharge*(1. + chargeRatio);
951 Double_t chargeNSide = clusterCharge*(1. - chargeRatio);
953 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 0)->Fill(module);
954 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 4)->Fill(recp->GetDetLocalX());
955 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 6)->Fill(recp->GetDetLocalZ());
956 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 8)->Fill(cluglo[0]);
957 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 10)->Fill(cluglo[1]);
958 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 12)->Fill(cluglo[2]);
959 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 14)->Fill(phi);
960 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 16)->Fill(theta);
961 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 18)->Fill(radius);
962 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 20)->Fill(recp->GetType());
963 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 22)->Fill(recp->GetChargeRatio());
964 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 24)->Fill(recp->GetQ());
965 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 26)->Fill(chargePSide);
966 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 28)->Fill(chargeNSide);
967 if(chargePSide != 0.) fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 30)->Fill(chargeNSide/chargePSide);
968 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 32)->Fill(chargePSide,chargeNSide);
969 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 34)->SetBinContent(gModule,lLadderLocationY,recp->GetQ());
970 nClustersLayer5 += 1;
971 }//layer 5 histograms
973 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 1)->Fill(module);
974 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 5)->Fill(recp->GetDetLocalX());
975 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 7)->Fill(recp->GetDetLocalZ());
976 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 9)->Fill(cluglo[0]);
977 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 11)->Fill(cluglo[1]);
978 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 13)->Fill(cluglo[2]);
979 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 15)->Fill(phi);
980 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 17)->Fill(theta);
981 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 19)->Fill(radius);
982 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 21)->Fill(recp->GetType());
983 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 23)->Fill(recp->GetChargeRatio());
984 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 25)->Fill(recp->GetQ());
985 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 27)->Fill(chargePSide);
986 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 29)->Fill(chargeNSide);
987 if(chargePSide != 0.) fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 31)->Fill(chargeNSide/chargePSide);
988 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 33)->Fill(chargePSide,chargeNSide);
989 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 35)->SetBinContent(gModule,lLadderLocationY,recp->GetQ());
990 nClustersLayer6 += 1;
991 }//layer 6 histograms
995 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 2)->Fill(nClustersLayer5);
996 fAliITSQADataMakerRec->GetRecPointsData(fGenOffset + 3)->Fill(nClustersLayer6);
998 statRecpoints.Clear();