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
23 // M. Nicassio D. Elia INFN Bari March 2008
24 // maria.nicassio@ba.infn.it
27 // --- ROOT system ---
32 // --- Standard library ---
34 // --- AliRoot header files ---
35 #include "AliITSQADataMakerRec.h"
36 #include "AliITSQASPDDataMakerRec.h"
39 #include "AliRawReader.h"
40 #include "AliITSRawStreamSPD.h"
41 #include "AliITSRawStreamSPDErrorLog.h"
42 #include "AliITSdigitSPD.h"
43 #include "AliITSRecPoint.h"
45 ClassImp(AliITSQASPDDataMakerRec)
47 //____________________________________________________________________________
48 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc, AliITSRawStreamSPDErrorLog *aliITSRawStreamSPDErrorLog) :
50 fAliITSQADataMakerRec(aliITSQADataMakerRec),
55 fSPDhRecPointsTask(0),
58 fGenRecPointsOffset(0),
59 fAdvLogger(aliITSRawStreamSPDErrorLog)
61 //ctor used to discriminate OnLine-Offline analysis
64 //____________________________________________________________________________
65 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(const AliITSQASPDDataMakerRec& qadm) :
67 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
68 fkOnline(qadm.fkOnline),
70 fSPDhRawsTask(qadm.fSPDhRawsTask),
71 fSPDhDigitsTask(qadm.fSPDhDigitsTask),
72 fSPDhRecPointsTask(qadm.fSPDhRecPointsTask),
73 fGenRawsOffset(qadm.fGenRawsOffset),
74 fGenDigitsOffset(qadm.fGenDigitsOffset),
75 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
76 fAdvLogger(qadm.fAdvLogger)
79 fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ;
80 fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
83 //__________________________________________________________________
84 AliITSQASPDDataMakerRec::~AliITSQASPDDataMakerRec(){
88 //__________________________________________________________________
90 AliITSQASPDDataMakerRec& AliITSQASPDDataMakerRec::operator = (const AliITSQASPDDataMakerRec& qac )
93 this->~AliITSQASPDDataMakerRec();
94 new(this) AliITSQASPDDataMakerRec(qac);
98 //____________________________________________________________________________
99 void AliITSQASPDDataMakerRec::StartOfDetectorCycle()
101 //Detector specific actions at start of cycle
102 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SPD Cycle\n");
105 //____________________________________________________________________________
106 void AliITSQASPDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
108 // launch the QA checking
109 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n");
111 //AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);
114 //____________________________________________________________________________
115 Int_t AliITSQASPDDataMakerRec::InitRaws()
117 // Initialization for RAW data - SPD -
118 const Bool_t expert = kTRUE ;
119 const Bool_t saveCorr = kTRUE ;
120 const Bool_t image = kTRUE ;
122 // fGenRawsOffset = (fAliITSQADataMakerRec->fRawsQAList[AliRecoParam::kDefault])->GetEntries();
123 fAdvLogger = new AliITSRawStreamSPDErrorLog();
124 AliDebug(AliQAv1::GetQADebugLevel(), "Book Offline Histograms for SPD\n ");
129 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
130 hlayer->GetXaxis()->SetTitle("Layer number");
131 hlayer->GetYaxis()->SetTitle("Entries");
132 rv = fAliITSQADataMakerRec->Add2RawsList(hlayer, 0+fGenRawsOffset, expert, !image, !saveCorr);
135 TH1F **hmod = new TH1F*[2];
136 TH2F **hhitMap = new TH2F*[20];
137 TH1F **herrors = new TH1F*[20];
138 for (Int_t iLay=0; iLay<2; iLay++) {
139 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
140 sprintf(title,"Module map - SPD Layer %d",iLay+1);
141 hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
142 hmod[iLay]->GetXaxis()->SetTitle("Module number");
143 hmod[iLay]->GetYaxis()->SetTitle("Entries");
144 rv = fAliITSQADataMakerRec->Add2RawsList(hmod[iLay], 1+iLay+fGenRawsOffset, !expert, image, !saveCorr);
147 for (Int_t iDDL=0; iDDL<20; iDDL++) {
148 sprintf(name,"SPDHitMap_SPD_DDL%d",iDDL+1);
149 sprintf(title,"Hit map - SPD DDL %d",iDDL+1);
150 hhitMap[iDDL]=new TH2F(name,title,320,0,10*32,1536,0,6*256);
151 hhitMap[iDDL]->GetXaxis()->SetTitle("Column");
152 hhitMap[iDDL]->GetYaxis()->SetTitle("Row");
153 rv = fAliITSQADataMakerRec->Add2RawsList(hhitMap[iDDL], 3+(2*iDDL)+fGenRawsOffset, expert, !image, !saveCorr);
155 sprintf(name,"SPDErrors_SPD_DDL%d",iDDL+1);
156 sprintf(title,"Error codes - SPD DDL %d",iDDL+1);
157 herrors[iDDL] = new TH1F (name,title,fAdvLogger->GetNrErrorCodes(),0,fAdvLogger->GetNrErrorCodes());
158 herrors[iDDL]->SetXTitle("Error Code");
159 herrors[iDDL]->SetYTitle("Nr of errors");
160 rv = fAliITSQADataMakerRec->Add2RawsList(herrors[iDDL], 4+(2*iDDL)+fGenRawsOffset, expert, !image, !saveCorr);
164 TH1F** hMultSPDhits = new TH1F*[2];
165 for (Int_t iLay=0; iLay<2; iLay++) {
166 sprintf(name,"SPDHitsMultiplicity_SPD%d",iLay+1);
167 sprintf(title,"Hit multiplicity - SPD Layer %d",iLay+1);
168 hMultSPDhits[iLay]=new TH1F(name,title,200,0.,200.);
169 hMultSPDhits[iLay]->GetXaxis()->SetTitle("Hit multiplicity");
170 hMultSPDhits[iLay]->GetYaxis()->SetTitle("Entries");
171 rv = fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits[iLay], 43+iLay+fGenRawsOffset, expert, !image, !saveCorr);
175 TH2F *hMultSPDhits2MultSPDhits1
176 = new TH2F("SPDHitMultCorrelation_SPD","Hit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
177 hMultSPDhits2MultSPDhits1->GetXaxis()->SetTitle("Hit multiplicity (Layer 1)");
178 hMultSPDhits2MultSPDhits1->GetYaxis()->SetTitle("Hit multiplicity (Layer 2)");
179 rv = fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits2MultSPDhits1, 45+fGenRawsOffset, !expert, image, !saveCorr);
182 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Raws histograms booked\n",fSPDhRawsTask));
186 //____________________________________________________________________________
187 Int_t AliITSQASPDDataMakerRec::MakeRaws(AliRawReader* rawReader)
189 // Fill QA for RAW - SPD -
191 // Check id histograms already created for this Event Specie
192 if ( ! fAliITSQADataMakerRec->GetRawsData(fGenRawsOffset) )
196 AliITSRawStreamSPD *rawStreamSPD = new AliITSRawStreamSPD(rawReader);
197 rawStreamSPD->ActivateAdvancedErrorLog(kTRUE,fAdvLogger);
203 Int_t iHalfStave, iChip;
205 UInt_t module, colM, rowM;
206 while(rawStreamSPD->Next()) {
208 iEq = rawReader->GetDDLID();
209 if (iEq>=0 && iEq<20) {
210 iHalfStave = rawStreamSPD->GetHalfStaveNr();
211 iChip = rawStreamSPD->GetChipAddr();
212 col = rawStreamSPD->GetChipCol();
213 row = rawStreamSPD->GetChipRow();
215 rawStreamSPD->OnlineToOffline(iEq, iHalfStave, iChip, col, row, module, colM, rowM);
217 if (iHalfStave>=0 && iHalfStave<2) iLayer=0;
220 fAliITSQADataMakerRec->GetRawsData(0+fGenRawsOffset)->Fill(iLayer);
222 fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset)->Fill(module);
225 fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset)->Fill(module);
229 fAliITSQADataMakerRec->GetRawsData((2*iEq)+3+fGenRawsOffset)->Fill(colM+(module%2)*160,rowM+iHalfStave*256);
233 for (Int_t ieq=0; ieq<20; ieq++)
234 for (UInt_t ierr=0; ierr<fAdvLogger->GetNrErrorCodes(); ierr++)
235 fAliITSQADataMakerRec->GetRawsData((2*ieq)+4+fGenRawsOffset)->Fill(ierr,fAdvLogger->GetNrErrors(ierr,ieq));
239 fAliITSQADataMakerRec->GetRawsData(43+fGenRawsOffset)->Fill(nDigitsL1);
240 fAliITSQADataMakerRec->GetRawsData(44+fGenRawsOffset)->Fill(nDigitsL2);
241 fAliITSQADataMakerRec->GetRawsData(45+fGenRawsOffset)->Fill(nDigitsL1,nDigitsL2);
244 AliDebug(AliQAv1::GetQADebugLevel(),Form("Event completed, %d raw digits read",nDigitsL1+nDigitsL2));
248 //____________________________________________________________________________
249 Int_t AliITSQASPDDataMakerRec::InitDigits()
251 // Initialization for DIGIT data - SPD -
252 const Bool_t expert = kTRUE ;
253 const Bool_t image = kTRUE ;
255 // fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
256 //fSPDhDigitsTask must be incremented by one unit every time a histogram is ADDED to the QA List
261 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
262 hlayer->GetXaxis()->SetTitle("Layer number");
263 hlayer->GetYaxis()->SetTitle("Entries");
264 rv = fAliITSQADataMakerRec->Add2DigitsList(hlayer,fGenDigitsOffset, expert, !image);
267 TH1F **hmod = new TH1F*[2];
268 for (Int_t iLay=0; iLay<2; iLay++) {
269 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
270 sprintf(title,"Module map - SPD Layer %d",iLay+1);
271 hmod[iLay]=new TH1F(name,title,240,0,240);
272 hmod[iLay]->GetXaxis()->SetTitle("Module number");
273 hmod[iLay]->GetYaxis()->SetTitle("Entries");
274 rv = fAliITSQADataMakerRec->Add2DigitsList(hmod[iLay],1+iLay+fGenDigitsOffset, !expert, image);
278 TH1F *hcolumns = new TH1F("SPDColumns_SPD","Columns - SPD",160,0.,160.);
279 hcolumns->GetXaxis()->SetTitle("Column number");
280 hcolumns->GetYaxis()->SetTitle("Entries");
281 rv = fAliITSQADataMakerRec->Add2DigitsList(hcolumns,3+fGenDigitsOffset, expert, !image);
284 TH1F *hrows = new TH1F("SPDRows_SPD","Rows - SPD",256,0.,256.);
285 hrows->GetXaxis()->SetTitle("Row number");
286 hrows->GetYaxis()->SetTitle("Entries");
287 rv = fAliITSQADataMakerRec->Add2DigitsList(hrows,4+fGenDigitsOffset, expert, !image);
290 TH1F** hMultSPDdigits = new TH1F*[2];
291 for (Int_t iLay=0; iLay<2; ++iLay) {
292 sprintf(name,"SPDDigitMultiplicity_SPD%d",iLay+1);
293 sprintf(title,"Digit multiplicity - SPD Layer %d",iLay+1);
294 hMultSPDdigits[iLay]=new TH1F(name,title,200,0.,200.);
295 hMultSPDdigits[iLay]->GetXaxis()->SetTitle("Digit multiplicity");
296 hMultSPDdigits[iLay]->GetYaxis()->SetTitle("Entries");
297 rv = fAliITSQADataMakerRec->Add2DigitsList(hMultSPDdigits[iLay], 5+iLay+fGenDigitsOffset, !expert, image);
301 TH2F *hMultSPDdig2MultSPDdig1
302 = new TH2F("SPDDigitMultCorrelation_SPD","Digit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
303 hMultSPDdig2MultSPDdig1->GetXaxis()->SetTitle("Digit multiplicity (Layer 1)");
304 hMultSPDdig2MultSPDdig1->GetYaxis()->SetTitle("Digit multiplicity (Layer 2)");
305 rv = fAliITSQADataMakerRec->Add2DigitsList(hMultSPDdig2MultSPDdig1,7+fGenDigitsOffset, !expert, image);
308 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Digits histograms booked\n",fSPDhDigitsTask));
312 //____________________________________________________________________________
313 Int_t AliITSQASPDDataMakerRec::MakeDigits(TTree *digits)
315 // Fill QA for DIGIT - SPD -
317 // Check id histograms already created for this Event Specie
318 if ( ! fAliITSQADataMakerRec->GetDigitsData(fGenDigitsOffset) )
321 // AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
322 // fITS->SetTreeAddress();
323 // TClonesArray *iITSdigits = fITS->DigitsAddress(0); // 0->SPD
324 TBranch *branchD = digits->GetBranch("ITSDigitsSPD");
326 AliError("can't get the branch with the SPD ITS digits !");
329 static TClonesArray statDigits("AliITSdigitSPD");
330 TClonesArray *iITSdigits = &statDigits;
331 branchD->SetAddress(&iITSdigits);
335 for (Int_t imod=0; imod<240; ++imod){
336 digits->GetEvent(imod);
337 Int_t ndigits = iITSdigits->GetEntries();
339 fAliITSQADataMakerRec->GetDigitsData(0+fGenDigitsOffset)->Fill(0.5,ndigits);
340 fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset)->Fill(imod,ndigits);
344 fAliITSQADataMakerRec->GetDigitsData(0+fGenDigitsOffset)->Fill(1,ndigits);
345 fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset)->Fill(imod,ndigits);
348 for (Int_t idig=0; idig<ndigits; ++idig) {
349 AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
350 Int_t col=dig->GetCoord1(); // cell number z
351 Int_t row=dig->GetCoord2(); // cell number x
352 fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset)->Fill(col);
353 fAliITSQADataMakerRec->GetDigitsData(4+fGenDigitsOffset)->Fill(row);
356 fAliITSQADataMakerRec->GetDigitsData(5+fGenDigitsOffset)->Fill(nDigitsL1);
357 fAliITSQADataMakerRec->GetDigitsData(6+fGenDigitsOffset)->Fill(nDigitsL2);
358 fAliITSQADataMakerRec->GetDigitsData(7+fGenDigitsOffset)->Fill(nDigitsL1,nDigitsL2);
362 //____________________________________________________________________________
363 Int_t AliITSQASPDDataMakerRec::InitRecPoints()
365 // Initialization for RECPOINTS - SPD -
366 const Bool_t expert = kTRUE ;
367 const Bool_t image = kTRUE ;
369 // fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
370 TH1F* hlayer= new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
371 hlayer->GetXaxis()->SetTitle("Layer number");
372 hlayer->GetYaxis()->SetTitle("Entries");
373 rv = fAliITSQADataMakerRec->Add2RecPointsList(hlayer, 0+fGenRecPointsOffset, expert, !image);
374 fSPDhRecPointsTask++;
376 TH1F** hmod = new TH1F*[2];
377 TH1F** hxl = new TH1F*[2];
378 TH1F** hzl = new TH1F*[2];
379 TH1F** hxg = new TH1F*[2];
380 TH1F** hyg = new TH1F*[2];
381 TH1F** hzg = new TH1F*[2];
382 TH1F** hr = new TH1F*[2];
383 TH1F** hphi = new TH1F*[2];
384 TH1F** hMultSPDcl = new TH1F*[2];
385 TH2F** hNyNz = new TH2F*[2]; // y and z cluster length
386 TH2F** hPhiZ = new TH2F*[2];
388 Float_t xlim[2]={4.5,8.};
389 Float_t zlim[2]={15.,15.};
393 for (Int_t iLay=0;iLay<2;iLay++) {
394 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
395 sprintf(title,"Module map - SPD Layer %d",iLay+1);
396 hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
397 hmod[iLay]->GetXaxis()->SetTitle("Module number");
398 hmod[iLay]->GetYaxis()->SetTitle("Entries");
399 rv = fAliITSQADataMakerRec->Add2RecPointsList(hmod[iLay], 1+(10*iLay)+fGenRecPointsOffset, expert, !image);
400 fSPDhRecPointsTask++;
402 sprintf(name,"SPDxLoc_SPD%d",iLay+1);
403 sprintf(title,"Local x coordinate - SPD Layer %d",iLay+1);
404 hxl[iLay]=new TH1F(name,title,100,-4.,4.);
405 hxl[iLay]->GetXaxis()->SetTitle("Local x [cm]");
406 hxl[iLay]->GetYaxis()->SetTitle("Entries");
407 rv = fAliITSQADataMakerRec->Add2RecPointsList(hxl[iLay], 2+(10*iLay)+fGenRecPointsOffset, expert, !image);
408 fSPDhRecPointsTask++;
410 sprintf(name,"SPDzLoc_SPD%d",iLay+1);
411 sprintf(title,"Local z coordinate - SPD Layer %d",iLay+1);
412 hzl[iLay]=new TH1F(name,title,100,-4.,4.);
413 hzl[iLay]->GetXaxis()->SetTitle("Local z [cm]");
414 hzl[iLay]->GetYaxis()->SetTitle("Entries");
415 rv = fAliITSQADataMakerRec->Add2RecPointsList(hzl[iLay], 3+(10*iLay)+fGenRecPointsOffset, expert, !image);
416 fSPDhRecPointsTask++;
418 sprintf(name,"SPDxGlob_SPD%d",iLay+1);
419 sprintf(title,"Global x coordinate - SPD Layer %d",iLay+1);
420 hxg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
421 hxg[iLay]->GetXaxis()->SetTitle("Global x [cm]");
422 hxg[iLay]->GetYaxis()->SetTitle("Entries");
423 rv = fAliITSQADataMakerRec->Add2RecPointsList(hxg[iLay],4+(10*iLay)+fGenRecPointsOffset, expert, !image);
424 fSPDhRecPointsTask++;
426 sprintf(name,"SPDyGlob_SPD%d",iLay+1);
427 sprintf(title,"Global y coordinate - SPD Layer %d",iLay+1);
428 hyg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
429 hyg[iLay]->GetXaxis()->SetTitle("Global y [cm]");
430 hyg[iLay]->GetYaxis()->SetTitle("Entries");
431 rv = fAliITSQADataMakerRec->Add2RecPointsList(hyg[iLay], 5+(10*iLay)+fGenRecPointsOffset, expert, !image);
432 fSPDhRecPointsTask++;
434 sprintf(name,"SPDzGlob_SPD%d",iLay+1);
435 sprintf(title,"Global z coordinate - SPD Layer %d",iLay+1);
436 hzg[iLay]=new TH1F(name,title,150,-zlim[iLay],zlim[iLay]);
437 hzg[iLay]->GetXaxis()->SetTitle("Global z [cm]");
438 hzg[iLay]->GetYaxis()->SetTitle("Entries");
439 rv = fAliITSQADataMakerRec->Add2RecPointsList(hzg[iLay], 6+(10*iLay)+fGenRecPointsOffset, expert, !image);
440 fSPDhRecPointsTask++;
442 sprintf(name,"SPDr_SPD%d",iLay+1);
443 sprintf(title,"Radius - SPD Layer %d",iLay+1);
444 hr[iLay]=new TH1F(name,title,100,0.,10.);
445 hr[iLay]->GetXaxis()->SetTitle("r [cm]");
446 hr[iLay]->GetYaxis()->SetTitle("Entries");
447 rv = fAliITSQADataMakerRec->Add2RecPointsList(hr[iLay], 7+(10*iLay)+fGenRecPointsOffset, expert, !image);
448 fSPDhRecPointsTask++;
450 sprintf(name,"SPDphi_SPD%d",iLay+1);
451 sprintf(title,"#varphi - SPD Layer %d",iLay+1);
452 hphi[iLay]=new TH1F(name,title,1000,0.,2*TMath::Pi());
453 hphi[iLay]->GetXaxis()->SetTitle("#varphi [rad]");
454 hphi[iLay]->GetYaxis()->SetTitle("Entries");
455 rv = fAliITSQADataMakerRec->Add2RecPointsList(hphi[iLay], 8+(10*iLay)+fGenRecPointsOffset, expert, !image);
456 fSPDhRecPointsTask++;
458 sprintf(name,"SPDSizeYvsZ_SPD%d",iLay+1);
459 sprintf(title,"Cluster dimension - SPD Layer %d",iLay+1);
460 hNyNz[iLay]=new TH2F(name,title,100,0.,100.,100,0.,100.);
461 hNyNz[iLay]->GetXaxis()->SetTitle("z length");
462 hNyNz[iLay]->GetYaxis()->SetTitle("y length");
463 rv = fAliITSQADataMakerRec->Add2RecPointsList(hNyNz[iLay], 9+(10*iLay)+fGenRecPointsOffset, expert, !image);
464 fSPDhRecPointsTask++;
466 sprintf(name,"SPDphi_z_SPD%d",iLay+1);
467 sprintf(title,"#varphi vs z - SPD Layer %d",iLay+1);
468 hPhiZ[iLay]=new TH2F(name,title,150,-zlim[iLay],zlim[iLay],200,0.,2*TMath::Pi());
469 hPhiZ[iLay]->GetXaxis()->SetTitle("Global z [cm]");
470 hPhiZ[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
471 rv = fAliITSQADataMakerRec->Add2RecPointsList(hPhiZ[iLay], 10+(10*iLay)+fGenRecPointsOffset, !expert, image);
472 fSPDhRecPointsTask++;
476 TH2F *hrPhi=new TH2F("SPDr_phi_SPD","#varphi vs r - SPD",100,0.,10.,100,0.,2*TMath::Pi());
477 hrPhi->GetXaxis()->SetTitle("r [cm]");
478 hrPhi->GetYaxis()->SetTitle("#varphi [rad]");
479 rv = fAliITSQADataMakerRec->Add2RecPointsList(hrPhi, 21+fGenRecPointsOffset, expert, !image);
480 fSPDhRecPointsTask++;
482 TH2F *hxy=new TH2F("SPDx_y_SPD","Global y vs x - SPD",200,-10.,10.,200,-10.,10.);
483 hxy->GetXaxis()->SetTitle("Global x [cm]");
484 hxy->GetYaxis()->SetTitle("Global y [cm]");
485 rv = fAliITSQADataMakerRec->Add2RecPointsList(hxy, 22+fGenRecPointsOffset, !expert, image);
486 fSPDhRecPointsTask++;
488 for (Int_t iLay=0;iLay<2;iLay++) {
489 sprintf(name,"SPDMultiplicity_SPD%d",iLay+1);
490 sprintf(title,"Cluster multiplicity - SPD Layer %d",iLay+1);
491 hMultSPDcl[iLay]=new TH1F(name,title,200,0.,200.);
492 hMultSPDcl[iLay]->GetXaxis()->SetTitle("Cluster multiplicity");
493 hMultSPDcl[iLay]->GetYaxis()->SetTitle("Entries");
494 rv = fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl[iLay], 23+iLay+fGenRecPointsOffset, !expert, image);
495 fSPDhRecPointsTask++;
498 TH2F *hMultSPDcl2MultSPDcl1 =
499 new TH2F("SPDMultCorrelation_SPD","Cluster multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
500 hMultSPDcl2MultSPDcl1->GetXaxis()->SetTitle("Clusters multiplicity (Layer 1)");
501 hMultSPDcl2MultSPDcl1->GetYaxis()->SetTitle("Clusters multiplicity (Layer 2)");
502 rv = fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl2MultSPDcl1, 25+fGenRecPointsOffset, !expert, image);
503 fSPDhRecPointsTask++;
505 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Recs histograms booked\n",fSPDhRecPointsTask));
510 //____________________________________________________________________________
511 Int_t AliITSQASPDDataMakerRec::MakeRecPoints(TTree * clusterTree)
513 // Fill QA for RecPoints - SPD -
515 static TClonesArray statITSCluster("AliITSRecPoint");
516 TClonesArray *ITSCluster = &statITSCluster;
517 TBranch* itsClusterBranch=clusterTree->GetBranch("ITSRecPoints");
518 if (!itsClusterBranch) {
519 AliError("can't get the branch with the ITS clusters !");
522 // Check id histograms already created for this Event Specie
523 if ( ! fAliITSQADataMakerRec->GetRecPointsData(fGenRecPointsOffset) )
524 rv = InitRecPoints() ;
526 itsClusterBranch->SetAddress(&ITSCluster);
527 Int_t nItsMods = (Int_t)clusterTree->GetEntries();
529 Float_t cluGlo[3] = {0.,0.,0.};
530 Int_t nClusters[2] = {0,0};
532 for (Int_t iIts=0; iIts < nItsMods; iIts++) {
534 if (!clusterTree->GetEvent(iIts)) continue;
535 Int_t nCluster = ITSCluster->GetEntriesFast();
536 // loop over clusters
538 AliITSRecPoint* cluster = (AliITSRecPoint*)ITSCluster->UncheckedAt(nCluster);
540 if (cluster->GetLayer()>1) continue;
541 Int_t lay=cluster->GetLayer();
542 fAliITSQADataMakerRec->GetRecPointsData(0 +fGenRecPointsOffset)->Fill(lay);
543 cluster->GetGlobalXYZ(cluGlo);
544 Float_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]);
545 Float_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
547 fAliITSQADataMakerRec->GetRecPointsData(1 +fGenRecPointsOffset)->Fill(iIts);
548 fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset)->Fill(cluster->GetDetLocalX());
549 fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset)->Fill(cluster->GetDetLocalZ());
550 fAliITSQADataMakerRec->GetRecPointsData(4 +fGenRecPointsOffset)->Fill(cluGlo[0]);
551 fAliITSQADataMakerRec->GetRecPointsData(5 +fGenRecPointsOffset)->Fill(cluGlo[1]);
552 fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset)->Fill(cluGlo[2]);
553 fAliITSQADataMakerRec->GetRecPointsData(7 +fGenRecPointsOffset)->Fill(rad);
554 fAliITSQADataMakerRec->GetRecPointsData(8 +fGenRecPointsOffset)->Fill(phi);
555 fAliITSQADataMakerRec->GetRecPointsData(9 +fGenRecPointsOffset)->Fill(cluster->GetNz(),cluster->GetNy());
556 fAliITSQADataMakerRec->GetRecPointsData(10 +fGenRecPointsOffset)->Fill(cluGlo[2],phi);
558 fAliITSQADataMakerRec->GetRecPointsData(11 +fGenRecPointsOffset)->Fill(iIts);
559 fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset)->Fill(cluster->GetDetLocalX());
560 fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset)->Fill(cluster->GetDetLocalZ());
561 fAliITSQADataMakerRec->GetRecPointsData(14 +fGenRecPointsOffset)->Fill(cluGlo[0]);
562 fAliITSQADataMakerRec->GetRecPointsData(15 +fGenRecPointsOffset)->Fill(cluGlo[1]);
563 fAliITSQADataMakerRec->GetRecPointsData(16 +fGenRecPointsOffset)->Fill(cluGlo[2]);
564 fAliITSQADataMakerRec->GetRecPointsData(17 +fGenRecPointsOffset)->Fill(rad);
565 fAliITSQADataMakerRec->GetRecPointsData(18 +fGenRecPointsOffset)->Fill(phi);
566 fAliITSQADataMakerRec->GetRecPointsData(19 +fGenRecPointsOffset)->Fill(cluster->GetNz(),cluster->GetNy());
567 fAliITSQADataMakerRec->GetRecPointsData(20 +fGenRecPointsOffset)->Fill(cluGlo[2],phi);
569 fAliITSQADataMakerRec->GetRecPointsData(21 +fGenRecPointsOffset)->Fill(rad,phi);
570 fAliITSQADataMakerRec->GetRecPointsData(22 +fGenRecPointsOffset)->Fill(cluGlo[0],cluGlo[1]);
573 } // end of cluster loop
574 } // end of its "subdetector" loop
576 for (Int_t iLay=0; iLay<2; iLay++)
577 fAliITSQADataMakerRec->GetRecPointsData(23+iLay +fGenRecPointsOffset)->Fill(nClusters[iLay]);
579 fAliITSQADataMakerRec->GetRecPointsData(25 +fGenRecPointsOffset)->Fill(nClusters[0],nClusters[1]);
581 statITSCluster.Clear();
587 //_______________________________________________________________
589 Int_t AliITSQASPDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task) {
590 // Returns offset number according to the specified task
592 if( task == AliQAv1::kRAWS ) {
593 offset=fGenRawsOffset;
595 else if( task == AliQAv1::kDIGITSR ) {
596 offset=fGenDigitsOffset;
598 else if( task == AliQAv1::kRECPOINTS ) {
599 offset=fGenRecPointsOffset;
601 else if( task == AliQAv1::kDIGITS ) {
602 offset=fGenDigitsOffset;
605 AliInfo("No task has been selected. Offset set to zero.\n");
611 //_______________________________________________________________
613 void AliITSQASPDDataMakerRec::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset) {
614 // Returns offset number according to the specified task
615 if( task == AliQAv1::kRAWS ) {
616 fGenRawsOffset=offset;
618 else if( task == AliQAv1::kDIGITSR ) {
619 fGenDigitsOffset=offset;
621 else if( task == AliQAv1::kRECPOINTS ) {
622 fGenRecPointsOffset=offset;
625 AliInfo("No task has been selected. Offset set to zero.\n");
629 //_______________________________________________________________
631 Int_t AliITSQASPDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task) {
632 // Returns the number of histograms associated to the specified task
636 if( task == AliQAv1::kRAWS ) {
637 histotot=fSPDhRawsTask;
639 else if( task == AliQAv1::kDIGITSR ) {
640 histotot=fSPDhDigitsTask;
642 else if( task == AliQAv1::kRECPOINTS ){
643 histotot=fSPDhRecPointsTask;
645 else if( task == AliQAv1::kDIGITS ){
646 histotot=fSPDhDigitsTask;
649 AliInfo("No task has been selected. TaskHisto set to zero.\n");