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"
44 #include "AliITSRecPointContainer.h"
46 ClassImp(AliITSQASPDDataMakerRec)
48 //____________________________________________________________________________
49 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc, AliITSRawStreamSPDErrorLog *aliITSRawStreamSPDErrorLog) :
51 fAliITSQADataMakerRec(aliITSQADataMakerRec),
56 fSPDhRecPointsTask(0),
59 fGenRecPointsOffset(0),
60 fAdvLogger(aliITSRawStreamSPDErrorLog)
62 //ctor used to discriminate OnLine-Offline analysis
63 //AliInfo(Form("AliRecoParam::kNSpecies %d\n",AliRecoParam::kNSpecies));
64 fGenRawsOffset = new Int_t[AliRecoParam::kNSpecies];
65 fGenRecPointsOffset = new Int_t[AliRecoParam::kNSpecies];
66 fGenDigitsOffset = new Int_t[AliRecoParam::kNSpecies];
67 for(Int_t i=0; i<AliRecoParam::kNSpecies;i++) {
68 fGenRawsOffset[i] = 0;
69 fGenRecPointsOffset[i] = 0;
70 fGenDigitsOffset[i]=0;
74 //____________________________________________________________________________
75 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(const AliITSQASPDDataMakerRec& qadm) :
77 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
78 fkOnline(qadm.fkOnline),
80 fSPDhRawsTask(qadm.fSPDhRawsTask),
81 fSPDhDigitsTask(qadm.fSPDhDigitsTask),
82 fSPDhRecPointsTask(qadm.fSPDhRecPointsTask),
83 fGenRawsOffset(qadm.fGenRawsOffset),
84 fGenDigitsOffset(qadm.fGenDigitsOffset),
85 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
86 fAdvLogger(qadm.fAdvLogger)
89 fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ;
90 fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
93 //__________________________________________________________________
94 AliITSQASPDDataMakerRec::~AliITSQASPDDataMakerRec(){
98 //__________________________________________________________________
100 AliITSQASPDDataMakerRec& AliITSQASPDDataMakerRec::operator = (const AliITSQASPDDataMakerRec& qac )
103 this->~AliITSQASPDDataMakerRec();
104 new(this) AliITSQASPDDataMakerRec(qac);
108 //____________________________________________________________________________
109 void AliITSQASPDDataMakerRec::StartOfDetectorCycle()
111 //Detector specific actions at start of cycle
112 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SPD Cycle\n");
115 //____________________________________________________________________________
116 void AliITSQASPDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
118 // launch the QA checking
119 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n");
121 //AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);
124 //____________________________________________________________________________
125 Int_t AliITSQASPDDataMakerRec::InitRaws()
127 // Initialization for RAW data - SPD -
128 const Bool_t expert = kTRUE ;
129 const Bool_t saveCorr = kTRUE ;
130 const Bool_t image = kTRUE ;
132 // fGenRawsOffset = (fAliITSQADataMakerRec->fRawsQAList[AliRecoParam::kDefault])->GetEntries();
133 if(!fAdvLogger) fAdvLogger = new AliITSRawStreamSPDErrorLog();
134 AliDebug(AliQAv1::GetQADebugLevel(), "Book Offline Histograms for SPD\n ");
139 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
140 hlayer->GetXaxis()->SetTitle("Layer number");
141 hlayer->GetYaxis()->SetTitle("Entries");
142 rv = fAliITSQADataMakerRec->Add2RawsList(hlayer, 0+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
145 TH1F **hmod = new TH1F*[2];
146 TH2F **hhitMap = new TH2F*[20];
147 TH1F **herrors = new TH1F*[20];
148 for (Int_t iLay=0; iLay<2; iLay++) {
149 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
150 sprintf(title,"Module map - SPD Layer %d",iLay+1);
151 hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
152 hmod[iLay]->GetXaxis()->SetTitle("Module number");
153 hmod[iLay]->GetYaxis()->SetTitle("Entries");
154 rv = fAliITSQADataMakerRec->Add2RawsList(hmod[iLay], 1+iLay+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, !saveCorr);
157 for (Int_t iDDL=0; iDDL<20; iDDL++) {
158 sprintf(name,"SPDHitMap_SPD_DDL%d",iDDL+1);
159 sprintf(title,"Hit map - SPD DDL %d",iDDL+1);
160 hhitMap[iDDL]=new TH2F(name,title,320,0,10*32,1536,0,6*256);
161 hhitMap[iDDL]->GetXaxis()->SetTitle("Column");
162 hhitMap[iDDL]->GetYaxis()->SetTitle("Row");
163 rv = fAliITSQADataMakerRec->Add2RawsList(hhitMap[iDDL], 3+(2*iDDL)+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
165 sprintf(name,"SPDErrors_SPD_DDL%d",iDDL+1);
166 sprintf(title,"Error codes - SPD DDL %d",iDDL+1);
167 herrors[iDDL] = new TH1F (name,title,fAdvLogger->GetNrErrorCodes(),0,fAdvLogger->GetNrErrorCodes());
168 herrors[iDDL]->SetXTitle("Error Code");
169 herrors[iDDL]->SetYTitle("Nr of errors");
170 rv = fAliITSQADataMakerRec->Add2RawsList(herrors[iDDL], 4+(2*iDDL)+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
174 TH1F** hMultSPDhits = new TH1F*[2];
175 for (Int_t iLay=0; iLay<2; iLay++) {
176 sprintf(name,"SPDHitsMultiplicity_SPD%d",iLay+1);
177 sprintf(title,"Hit multiplicity - SPD Layer %d",iLay+1);
178 hMultSPDhits[iLay]=new TH1F(name,title,200,0.,200.);
179 hMultSPDhits[iLay]->GetXaxis()->SetTitle("Hit multiplicity");
180 hMultSPDhits[iLay]->GetYaxis()->SetTitle("Entries");
181 rv = fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits[iLay], 43+iLay+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image, !saveCorr);
185 TH2F *hMultSPDhits2MultSPDhits1
186 = new TH2F("SPDHitMultCorrelation_SPD","Hit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
187 hMultSPDhits2MultSPDhits1->GetXaxis()->SetTitle("Hit multiplicity (Layer 1)");
188 hMultSPDhits2MultSPDhits1->GetYaxis()->SetTitle("Hit multiplicity (Layer 2)");
189 rv = fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits2MultSPDhits1, 45+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, !saveCorr);
192 TH1F *hFastOrFiredChips = new TH1F("SPDFastOrPattern_SPD","FastOrFiredChip map - SPD",1200,0.,1200.);
193 hFastOrFiredChips->GetXaxis()->SetTitle("Chip number");
194 hFastOrFiredChips->GetYaxis()->SetTitle("Entries");
195 rv = fAliITSQADataMakerRec->Add2RawsList(hFastOrFiredChips, 46+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image, !saveCorr);
198 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Raws histograms booked\n",fSPDhRawsTask));
202 //____________________________________________________________________________
203 Int_t AliITSQASPDDataMakerRec::MakeRaws(AliRawReader* rawReader)
205 // Fill QA for RAW - SPD -
209 AliITSRawStreamSPD rawStreamSPD(rawReader);
210 rawStreamSPD.ActivateAdvancedErrorLog(kTRUE,fAdvLogger);
216 Int_t iHalfStave, iChip;
219 UInt_t module, colM, rowM;
220 while(rawStreamSPD.Next()) {
222 iEq = rawReader->GetDDLID();
223 if (iEq>=0 && iEq<20) {
224 iHalfStave = rawStreamSPD.GetHalfStaveNr();
225 iChip = rawStreamSPD.GetChipAddr();
226 col = rawStreamSPD.GetChipCol();
227 row = rawStreamSPD.GetChipRow();
229 rawStreamSPD.OnlineToOffline(iEq, iHalfStave, iChip, col, row, module, colM, rowM);
231 if (iHalfStave>=0 && iHalfStave<2) iLayer=0;
234 fAliITSQADataMakerRec->GetRawsData(0+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(iLayer);
236 fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(module);
239 fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(module);
243 fAliITSQADataMakerRec->GetRawsData((2*iEq)+3+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(colM+(module%2)*160,rowM+iHalfStave*256);
249 for (Int_t ieq=0; ieq<20; ieq++) {
250 for (UInt_t ierr=0; ierr<fAdvLogger->GetNrErrorCodes(); ierr++)
251 fAliITSQADataMakerRec->GetRawsData((2*ieq)+4+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(ierr,fAdvLogger->GetNrErrors(ierr,ieq));
253 for (Int_t ihs=0; ihs<6; ihs++) {
254 for (Int_t ichip=0; ichip<10; ichip++) {
255 if(rawStreamSPD.GetFastOrSignal(ieq,ihs,ichip)) {
256 chipKey = rawStreamSPD.GetOfflineChipKeyFromOnline(ieq,ihs,ichip);
257 fAliITSQADataMakerRec->GetRawsData(46+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(chipKey);
266 fAliITSQADataMakerRec->GetRawsData(43+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nDigitsL1);
267 fAliITSQADataMakerRec->GetRawsData(44+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nDigitsL2);
268 fAliITSQADataMakerRec->GetRawsData(45+fGenRawsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nDigitsL1,nDigitsL2);
270 AliDebug(AliQAv1::GetQADebugLevel(),Form("Event completed, %d raw digits read",nDigitsL1+nDigitsL2));
274 //____________________________________________________________________________
275 Int_t AliITSQASPDDataMakerRec::InitDigits()
277 // Initialization for DIGIT data - SPD -
278 const Bool_t expert = kTRUE ;
279 const Bool_t image = kTRUE ;
281 // fGenDigitsOffset = (fAliITSQADataMakerRec->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
282 //fSPDhDigitsTask must be incremented by one unit every time a histogram is ADDED to the QA List
287 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
288 hlayer->GetXaxis()->SetTitle("Layer number");
289 hlayer->GetYaxis()->SetTitle("Entries");
290 rv = fAliITSQADataMakerRec->Add2DigitsList(hlayer,fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
293 TH1F **hmod = new TH1F*[2];
294 for (Int_t iLay=0; iLay<2; iLay++) {
295 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
296 sprintf(title,"Module map - SPD Layer %d",iLay+1);
297 hmod[iLay]=new TH1F(name,title,240,0,240);
298 hmod[iLay]->GetXaxis()->SetTitle("Module number");
299 hmod[iLay]->GetYaxis()->SetTitle("Entries");
300 rv = fAliITSQADataMakerRec->Add2DigitsList(hmod[iLay],1+iLay+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
304 TH1F *hcolumns = new TH1F("SPDColumns_SPD","Columns - SPD",160,0.,160.);
305 hcolumns->GetXaxis()->SetTitle("Column number");
306 hcolumns->GetYaxis()->SetTitle("Entries");
307 rv = fAliITSQADataMakerRec->Add2DigitsList(hcolumns,3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
310 TH1F *hrows = new TH1F("SPDRows_SPD","Rows - SPD",256,0.,256.);
311 hrows->GetXaxis()->SetTitle("Row number");
312 hrows->GetYaxis()->SetTitle("Entries");
313 rv = fAliITSQADataMakerRec->Add2DigitsList(hrows,4+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
316 TH1F** hMultSPDdigits = new TH1F*[2];
317 for (Int_t iLay=0; iLay<2; ++iLay) {
318 sprintf(name,"SPDDigitMultiplicity_SPD%d",iLay+1);
319 sprintf(title,"Digit multiplicity - SPD Layer %d",iLay+1);
320 hMultSPDdigits[iLay]=new TH1F(name,title,200,0.,200.);
321 hMultSPDdigits[iLay]->GetXaxis()->SetTitle("Digit multiplicity");
322 hMultSPDdigits[iLay]->GetYaxis()->SetTitle("Entries");
323 rv = fAliITSQADataMakerRec->Add2DigitsList(hMultSPDdigits[iLay], 5+iLay+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
327 TH2F *hMultSPDdig2MultSPDdig1
328 = new TH2F("SPDDigitMultCorrelation_SPD","Digit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
329 hMultSPDdig2MultSPDdig1->GetXaxis()->SetTitle("Digit multiplicity (Layer 1)");
330 hMultSPDdig2MultSPDdig1->GetYaxis()->SetTitle("Digit multiplicity (Layer 2)");
331 rv = fAliITSQADataMakerRec->Add2DigitsList(hMultSPDdig2MultSPDdig1,7+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
334 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Digits histograms booked\n",fSPDhDigitsTask));
338 //____________________________________________________________________________
339 Int_t AliITSQASPDDataMakerRec::MakeDigits(TTree *digits)
341 // Fill QA for DIGIT - SPD -
344 // AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
345 // fITS->SetTreeAddress();
346 // TClonesArray *iITSdigits = fITS->DigitsAddress(0); // 0->SPD
347 TBranch *branchD = digits->GetBranch("ITSDigitsSPD");
349 AliError("can't get the branch with the SPD ITS digits !");
352 static TClonesArray statDigits("AliITSdigitSPD");
353 TClonesArray *iITSdigits = &statDigits;
354 branchD->SetAddress(&iITSdigits);
358 for (Int_t imod=0; imod<240; ++imod){
359 digits->GetEvent(imod);
360 Int_t ndigits = iITSdigits->GetEntries();
362 fAliITSQADataMakerRec->GetDigitsData(0+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(0.5,ndigits);
363 fAliITSQADataMakerRec->GetDigitsData(1+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(imod,ndigits);
367 fAliITSQADataMakerRec->GetDigitsData(0+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(1,ndigits);
368 fAliITSQADataMakerRec->GetDigitsData(2+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(imod,ndigits);
371 for (Int_t idig=0; idig<ndigits; ++idig) {
372 AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
373 Int_t col=dig->GetCoord1(); // cell number z
374 Int_t row=dig->GetCoord2(); // cell number x
375 fAliITSQADataMakerRec->GetDigitsData(3+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(col);
376 fAliITSQADataMakerRec->GetDigitsData(4+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(row);
379 fAliITSQADataMakerRec->GetDigitsData(5+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nDigitsL1);
380 fAliITSQADataMakerRec->GetDigitsData(6+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nDigitsL2);
381 fAliITSQADataMakerRec->GetDigitsData(7+fGenDigitsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nDigitsL1,nDigitsL2);
385 //____________________________________________________________________________
386 Int_t AliITSQASPDDataMakerRec::InitRecPoints()
388 // Initialization for RECPOINTS - SPD -
389 const Bool_t expert = kTRUE ;
390 const Bool_t image = kTRUE ;
392 //AliInfo(Form("fAliITSQADataMakerRec->GetEventSpecie() %d\n",fAliITSQADataMakerRec->GetEventSpecie()));
393 //AliInfo(Form("fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()] %d\n",fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]));
394 // fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
395 TH1F* hlayer= new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
396 hlayer->GetXaxis()->SetTitle("Layer number");
397 hlayer->GetYaxis()->SetTitle("Entries");
398 rv = fAliITSQADataMakerRec->Add2RecPointsList(hlayer, 0+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
399 fSPDhRecPointsTask++;
401 TH1F** hmod = new TH1F*[2];
402 TH1F** hxl = new TH1F*[2];
403 TH1F** hzl = new TH1F*[2];
404 TH1F** hxg = new TH1F*[2];
405 TH1F** hyg = new TH1F*[2];
406 TH1F** hzg = new TH1F*[2];
407 TH1F** hr = new TH1F*[2];
408 TH1F** hphi = new TH1F*[2];
409 TH1F** hMultSPDcl = new TH1F*[2];
410 TH2F** hNyNz = new TH2F*[2]; // y and z cluster length
411 TH1F** hNpixels = new TH1F*[2]; // cluster size in number of pixels
412 TH1F** hType = new TH1F*[2]; // cluster type according to conventional table
413 TH2F** hPhiZ = new TH2F*[2];
415 Float_t xlim[2]={4.5,8.};
416 Float_t zlim[2]={15.,15.};
420 for (Int_t iLay=0;iLay<2;iLay++) {
421 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
422 sprintf(title,"Module map - SPD Layer %d",iLay+1);
423 hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
424 hmod[iLay]->GetXaxis()->SetTitle("Module number");
425 hmod[iLay]->GetYaxis()->SetTitle("Entries");
426 rv = fAliITSQADataMakerRec->Add2RecPointsList(hmod[iLay], 1+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
427 fSPDhRecPointsTask++;
429 sprintf(name,"SPDxLoc_SPD%d",iLay+1);
430 sprintf(title,"Local x coordinate - SPD Layer %d",iLay+1);
431 hxl[iLay]=new TH1F(name,title,100,-4.,4.);
432 hxl[iLay]->GetXaxis()->SetTitle("Local x [cm]");
433 hxl[iLay]->GetYaxis()->SetTitle("Entries");
434 rv = fAliITSQADataMakerRec->Add2RecPointsList(hxl[iLay], 2+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
435 fSPDhRecPointsTask++;
437 sprintf(name,"SPDzLoc_SPD%d",iLay+1);
438 sprintf(title,"Local z coordinate - SPD Layer %d",iLay+1);
439 hzl[iLay]=new TH1F(name,title,100,-4.,4.);
440 hzl[iLay]->GetXaxis()->SetTitle("Local z [cm]");
441 hzl[iLay]->GetYaxis()->SetTitle("Entries");
442 rv = fAliITSQADataMakerRec->Add2RecPointsList(hzl[iLay], 3+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
443 fSPDhRecPointsTask++;
445 sprintf(name,"SPDxGlob_SPD%d",iLay+1);
446 sprintf(title,"Global x coordinate - SPD Layer %d",iLay+1);
447 hxg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
448 hxg[iLay]->GetXaxis()->SetTitle("Global x [cm]");
449 hxg[iLay]->GetYaxis()->SetTitle("Entries");
450 rv = fAliITSQADataMakerRec->Add2RecPointsList(hxg[iLay],4+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
451 fSPDhRecPointsTask++;
453 sprintf(name,"SPDyGlob_SPD%d",iLay+1);
454 sprintf(title,"Global y coordinate - SPD Layer %d",iLay+1);
455 hyg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
456 hyg[iLay]->GetXaxis()->SetTitle("Global y [cm]");
457 hyg[iLay]->GetYaxis()->SetTitle("Entries");
458 rv = fAliITSQADataMakerRec->Add2RecPointsList(hyg[iLay], 5+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
459 fSPDhRecPointsTask++;
461 sprintf(name,"SPDzGlob_SPD%d",iLay+1);
462 sprintf(title,"Global z coordinate - SPD Layer %d",iLay+1);
463 hzg[iLay]=new TH1F(name,title,150,-zlim[iLay],zlim[iLay]);
464 hzg[iLay]->GetXaxis()->SetTitle("Global z [cm]");
465 hzg[iLay]->GetYaxis()->SetTitle("Entries");
466 rv = fAliITSQADataMakerRec->Add2RecPointsList(hzg[iLay], 6+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
467 fSPDhRecPointsTask++;
469 sprintf(name,"SPDr_SPD%d",iLay+1);
470 sprintf(title,"Radius - SPD Layer %d",iLay+1);
471 hr[iLay]=new TH1F(name,title,100,0.,10.);
472 hr[iLay]->GetXaxis()->SetTitle("r [cm]");
473 hr[iLay]->GetYaxis()->SetTitle("Entries");
474 rv = fAliITSQADataMakerRec->Add2RecPointsList(hr[iLay], 7+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
475 fSPDhRecPointsTask++;
477 sprintf(name,"SPDphi_SPD%d",iLay+1);
478 sprintf(title,"#varphi - SPD Layer %d",iLay+1);
479 hphi[iLay]=new TH1F(name,title,1000,0.,2*TMath::Pi());
480 hphi[iLay]->GetXaxis()->SetTitle("#varphi [rad]");
481 hphi[iLay]->GetYaxis()->SetTitle("Entries");
482 rv = fAliITSQADataMakerRec->Add2RecPointsList(hphi[iLay], 8+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
483 fSPDhRecPointsTask++;
485 sprintf(name,"SPDSizeYvsZ_SPD%d",iLay+1);
486 sprintf(title,"Cluster dimension - SPD Layer %d",iLay+1);
487 hNyNz[iLay]=new TH2F(name,title,100,0.,100.,100,0.,100.);
488 hNyNz[iLay]->GetXaxis()->SetTitle("z length");
489 hNyNz[iLay]->GetYaxis()->SetTitle("y length");
490 rv = fAliITSQADataMakerRec->Add2RecPointsList(hNyNz[iLay], 9+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
491 fSPDhRecPointsTask++;
493 sprintf(name,"SPDSizeTot_SPD%d",iLay+1);
494 sprintf(title,"Cluster size - SPD Layer %d",iLay+1);
495 hNpixels[iLay]=new TH1F(name,title,100,0.,100.);
496 hNpixels[iLay]->GetXaxis()->SetTitle("Cluster size");
497 hNpixels[iLay]->GetYaxis()->SetTitle("Entries");
498 rv = fAliITSQADataMakerRec->Add2RecPointsList(hNpixels[iLay], 10+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
499 fSPDhRecPointsTask++;
501 sprintf(name,"SPDType_SPD%d",iLay+1);
502 sprintf(title,"Cluster type - SPD Layer %d",iLay+1);
503 hType[iLay]=new TH1F(name,title,20,0.,20.);
504 hType[iLay]->GetXaxis()->SetTitle("Cluster type");
505 hType[iLay]->GetYaxis()->SetTitle("Entries");
506 rv = fAliITSQADataMakerRec->Add2RecPointsList(hType[iLay], 11+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
507 fSPDhRecPointsTask++;
509 sprintf(name,"SPDphi_z_SPD%d",iLay+1);
510 sprintf(title,"#varphi vs z - SPD Layer %d",iLay+1);
511 hPhiZ[iLay]=new TH2F(name,title,150,-zlim[iLay],zlim[iLay],200,0.,2*TMath::Pi());
512 hPhiZ[iLay]->GetXaxis()->SetTitle("Global z [cm]");
513 hPhiZ[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
514 rv = fAliITSQADataMakerRec->Add2RecPointsList(hPhiZ[iLay], 12+(12*iLay)+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
515 fSPDhRecPointsTask++;
519 TH2F *hrPhi=new TH2F("SPDr_phi_SPD","#varphi vs r - SPD",100,0.,10.,100,0.,2*TMath::Pi());
520 hrPhi->GetXaxis()->SetTitle("r [cm]");
521 hrPhi->GetYaxis()->SetTitle("#varphi [rad]");
522 rv = fAliITSQADataMakerRec->Add2RecPointsList(hrPhi, 25+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], expert, !image);
523 fSPDhRecPointsTask++;
525 TH2F *hxy=new TH2F("SPDx_y_SPD","Global y vs x - SPD",200,-10.,10.,200,-10.,10.);
526 hxy->GetXaxis()->SetTitle("Global x [cm]");
527 hxy->GetYaxis()->SetTitle("Global y [cm]");
528 rv = fAliITSQADataMakerRec->Add2RecPointsList(hxy, 26+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
529 fSPDhRecPointsTask++;
531 for (Int_t iLay=0;iLay<2;iLay++) {
532 sprintf(name,"SPDMultiplicity_SPD%d",iLay+1);
533 sprintf(title,"Cluster multiplicity - SPD Layer %d",iLay+1);
534 hMultSPDcl[iLay]=new TH1F(name,title,200,0.,200.);
535 hMultSPDcl[iLay]->GetXaxis()->SetTitle("Cluster multiplicity");
536 hMultSPDcl[iLay]->GetYaxis()->SetTitle("Entries");
537 rv = fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl[iLay], 27+iLay+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
538 fSPDhRecPointsTask++;
541 TH2F *hMultSPDcl2MultSPDcl1 =
542 new TH2F("SPDMultCorrelation_SPD","Cluster multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
543 hMultSPDcl2MultSPDcl1->GetXaxis()->SetTitle("Clusters multiplicity (Layer 1)");
544 hMultSPDcl2MultSPDcl1->GetYaxis()->SetTitle("Clusters multiplicity (Layer 2)");
545 rv = fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl2MultSPDcl1, 29+fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()], !expert, image);
546 fSPDhRecPointsTask++;
548 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Recs histograms booked\n",fSPDhRecPointsTask));
553 //____________________________________________________________________________
554 Int_t AliITSQASPDDataMakerRec::MakeRecPoints(TTree * clusterTree)
556 // Fill QA for RecPoints - SPD -
558 AliITSRecPointContainer* rpcont=AliITSRecPointContainer::Instance();
559 TClonesArray *recpoints = rpcont->FetchClusters(0,clusterTree);
560 if(!rpcont->GetStatusOK() || !rpcont->IsSPDActive()){
561 AliError("can't get SPD clusters !");
565 //AliInfo(Form("fAliITSQADataMakerRec->GetEventSpecie() %d\n",fAliITSQADataMakerRec->GetEventSpecie()));
566 //AliInfo(Form("fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()] %d\n",fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()]));
567 Int_t nSPDmod = AliITSgeomTGeo::GetModuleIndex(3,1,1);
569 Float_t cluGlo[3] = {0.,0.,0.};
570 Int_t nClusters[2] = {0,0};
572 for (Int_t iIts=0; iIts < nSPDmod; iIts++) {
573 recpoints = rpcont->UncheckedGetClusters(iIts);
574 Int_t nCluster = recpoints->GetEntriesFast();
575 if(nCluster == 0)continue;
576 // loop over clusters
578 AliITSRecPoint* cluster =
579 (AliITSRecPoint*)recpoints->UncheckedAt(nCluster);
580 if (cluster->GetLayer()>1)continue;
581 Int_t lay=cluster->GetLayer();
582 fAliITSQADataMakerRec->GetRecPointsData(0 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(lay);
583 cluster->GetGlobalXYZ(cluGlo);
584 Float_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]);
585 Float_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
587 fAliITSQADataMakerRec->GetRecPointsData(1 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(iIts);
588 fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetDetLocalX());
589 fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetDetLocalZ());
590 fAliITSQADataMakerRec->GetRecPointsData(4 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[0]);
591 fAliITSQADataMakerRec->GetRecPointsData(5 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[1]);
592 fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[2]);
593 fAliITSQADataMakerRec->GetRecPointsData(7 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);
594 fAliITSQADataMakerRec->GetRecPointsData(8 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);
595 fAliITSQADataMakerRec->GetRecPointsData(9 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetNz(),cluster->GetNy());
596 fAliITSQADataMakerRec->GetRecPointsData(10 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetNpixels());
597 fAliITSQADataMakerRec->GetRecPointsData(11 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetSPDclusterType());
598 fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[2],phi);
600 fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(iIts);
601 fAliITSQADataMakerRec->GetRecPointsData(14 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetDetLocalX());
602 fAliITSQADataMakerRec->GetRecPointsData(15 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetDetLocalZ());
603 fAliITSQADataMakerRec->GetRecPointsData(16 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[0]);
604 fAliITSQADataMakerRec->GetRecPointsData(17 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[1]);
605 fAliITSQADataMakerRec->GetRecPointsData(18 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[2]);
606 fAliITSQADataMakerRec->GetRecPointsData(19 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad);
607 fAliITSQADataMakerRec->GetRecPointsData(20 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(phi);
608 fAliITSQADataMakerRec->GetRecPointsData(21 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetNz(),cluster->GetNy());
609 fAliITSQADataMakerRec->GetRecPointsData(22 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetNpixels());
610 fAliITSQADataMakerRec->GetRecPointsData(23 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluster->GetSPDclusterType());
611 fAliITSQADataMakerRec->GetRecPointsData(24 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[2],phi);
613 fAliITSQADataMakerRec->GetRecPointsData(25 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(rad,phi);
614 fAliITSQADataMakerRec->GetRecPointsData(26 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(cluGlo[0],cluGlo[1]);
617 } // end of cluster loop
618 } // end of its "subdetector" loop
620 for (Int_t iLay=0; iLay<2; iLay++)
621 fAliITSQADataMakerRec->GetRecPointsData(27 +iLay +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nClusters[iLay]);
623 fAliITSQADataMakerRec->GetRecPointsData(29 +fGenRecPointsOffset[fAliITSQADataMakerRec->GetEventSpecie()])->Fill(nClusters[0],nClusters[1]);
630 //_______________________________________________________________
632 Int_t AliITSQASPDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task,Int_t specie) {
633 // Returns offset number according to the specified task
635 if( task == AliQAv1::kRAWS ) {
636 offset=fGenRawsOffset[specie];
638 else if( task == AliQAv1::kDIGITSR ) {
639 offset=fGenDigitsOffset[specie];
641 else if( task == AliQAv1::kRECPOINTS ) {
642 offset=fGenRecPointsOffset[specie];
648 //_______________________________________________________________
650 void AliITSQASPDDataMakerRec::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset, Int_t specie) {
651 // Returns offset number according to the specified task
652 if( task == AliQAv1::kRAWS ) {
653 fGenRawsOffset[specie]=offset;
655 else if( task == AliQAv1::kDIGITSR ) {
656 fGenDigitsOffset[specie]=offset;
658 else if( task == AliQAv1::kRECPOINTS ) {
659 fGenRecPointsOffset[specie]=offset;
663 //_______________________________________________________________
665 Int_t AliITSQASPDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task) {
666 // Returns the number of histograms associated to the specified task
670 if( task == AliQAv1::kRAWS ) {
671 histotot=fSPDhRawsTask;
673 else if( task == AliQAv1::kDIGITSR ) {
674 histotot=fSPDhDigitsTask;
676 else if( task == AliQAv1::kRECPOINTS ){
677 histotot=fSPDhRecPointsTask;