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 "AliITSRecPoint.h"
44 ClassImp(AliITSQASPDDataMakerRec)
46 //____________________________________________________________________________
47 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc, AliITSRawStreamSPDErrorLog *aliITSRawStreamSPDErrorLog) :
49 fAliITSQADataMakerRec(aliITSQADataMakerRec),
54 fAdvLogger(aliITSRawStreamSPDErrorLog)
56 //ctor used to discriminate OnLine-Offline analysis
59 //____________________________________________________________________________
60 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(const AliITSQASPDDataMakerRec& qadm) :
62 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
63 fkOnline(qadm.fkOnline),
65 fSPDhTask(qadm.fSPDhTask),
66 fGenOffset(qadm.fGenOffset),
67 fAdvLogger(qadm.fAdvLogger)
70 fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ;
71 fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
74 //__________________________________________________________________
75 AliITSQASPDDataMakerRec::~AliITSQASPDDataMakerRec(){
79 //__________________________________________________________________
81 AliITSQASPDDataMakerRec& AliITSQASPDDataMakerRec::operator = (const AliITSQASPDDataMakerRec& qac )
84 this->~AliITSQASPDDataMakerRec();
85 new(this) AliITSQASPDDataMakerRec(qac);
89 //____________________________________________________________________________
90 void AliITSQASPDDataMakerRec::StartOfDetectorCycle()
92 //Detector specific actions at start of cycle
93 AliDebug(1,"AliITSQADM::Start of SPD Cycle\n");
96 //____________________________________________________________________________
97 void AliITSQASPDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t /*task*/, TObjArray* /*list*/)
99 // launch the QA checking
100 AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n");
102 //AliQAChecker::Instance()->Run( AliQA::kITS , task, list);
105 //____________________________________________________________________________
106 void AliITSQASPDDataMakerRec::InitRaws()
108 // Initialization for RAW data - SPD -
109 fGenOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();
114 TH1F *hlayer = new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
115 hlayer->GetXaxis()->SetTitle("Layer number");
116 hlayer->GetYaxis()->SetTitle("Entries");
117 fAliITSQADataMakerRec->Add2RawsList(hlayer,fSPDhTask+fGenOffset);
120 TH1F **hmod = new TH1F*[2];
121 TH2F **hhitMap = new TH2F*[20];
122 TH1F **herrors = new TH1F*[20];
123 for (Int_t iLay=0; iLay<2; iLay++) {
124 sprintf(name,"ModPattern_SPD%d",iLay+1);
125 sprintf(title,"Module map - SPD Layer %d",iLay+1);
126 hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
127 hmod[iLay]->GetXaxis()->SetTitle("Module number");
128 hmod[iLay]->GetYaxis()->SetTitle("Entries");
129 fAliITSQADataMakerRec->Add2RawsList(hmod[iLay], fSPDhTask +fGenOffset);
132 fAdvLogger = new AliITSRawStreamSPDErrorLog();
133 for (Int_t iDDL=0; iDDL<20; iDDL++) {
134 sprintf(name,"HitMap_SPD_DDL%d",iDDL+1);
135 sprintf(title,"Hit map - SPD DDL %d",iDDL+1);
136 hhitMap[iDDL]=new TH2F(name,title,320,0,10*32,1536,0,6*256);
137 hhitMap[iDDL]->GetXaxis()->SetTitle("Column");
138 hhitMap[iDDL]->GetYaxis()->SetTitle("Row");
139 fAliITSQADataMakerRec->Add2RawsList(hhitMap[iDDL], fSPDhTask +fGenOffset);
141 sprintf(name,"Errors_SPD_DDL%d",iDDL+1);
142 sprintf(title,"Error codes - SPD DDL %d",iDDL+1);
143 herrors[iDDL] = new TH1F (name,title,15,0,15);
144 herrors[iDDL]->SetXTitle("Error Code");
145 herrors[iDDL]->SetYTitle("Nr of errors");
146 fAliITSQADataMakerRec->Add2RawsList(herrors[iDDL], fSPDhTask +fGenOffset);
150 TH1F** hMultSPDhits = new TH1F*[2];
151 for (Int_t iLay=0; iLay<2; iLay++) {
152 sprintf(name,"HitsMultiplicity_SPD%d",iLay+1);
153 sprintf(title,"Hit multiplicity - SPD Layer %d",iLay+1);
154 hMultSPDhits[iLay]=new TH1F(name,title,200,0.,200.);
155 hMultSPDhits[iLay]->GetXaxis()->SetTitle("Hit multiplicity");
156 hMultSPDhits[iLay]->GetYaxis()->SetTitle("Entries");
157 fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits[iLay], fSPDhTask+fGenOffset);
161 TH2F *hMultSPDhits2MultSPDhits1 = new TH2F("HitMultCorrelation_SPD","Hit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
162 hMultSPDhits2MultSPDhits1->GetXaxis()->SetTitle("Hit multiplicity (Layer 1)");
163 hMultSPDhits2MultSPDhits1->GetYaxis()->SetTitle("Hit multiplicity (Layer 2)");
164 fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits2MultSPDhits1, fSPDhTask+fGenOffset);
167 AliDebug(1,Form("%d SPD Raws histograms booked\n",fSPDhTask));
172 //____________________________________________________________________________
173 void AliITSQASPDDataMakerRec::MakeRaws(AliRawReader* rawReader)
175 // Fill QA for RAW - SPD -
177 AliITSRawStreamSPD *rawStreamSPD = new AliITSRawStreamSPD(rawReader);
178 rawStreamSPD->ActivateAdvancedErrorLog(kTRUE,fAdvLogger);
184 Int_t iHalfStave, iChip;
186 UInt_t module, colM, rowM;
187 while(rawStreamSPD->Next()) {
189 iEq = rawReader->GetDDLID();
190 if (iEq>=0 && iEq<20) {
191 iHalfStave = rawStreamSPD->GetHalfStaveNr();
192 iChip = rawStreamSPD->GetChipAddr();
193 col = rawStreamSPD->GetChipCol();
194 row = rawStreamSPD->GetChipRow();
196 rawStreamSPD->OnlineToOffline(iEq, iHalfStave, iChip, col, row, module, colM, rowM);
198 if (iHalfStave>=0 && iHalfStave<2) iLayer=0;
201 fAliITSQADataMakerRec->GetRawsData(0+fGenOffset)->Fill(iLayer);
203 fAliITSQADataMakerRec->GetRawsData(1+fGenOffset)->Fill(module);
206 fAliITSQADataMakerRec->GetRawsData(2+fGenOffset)->Fill(module);
210 fAliITSQADataMakerRec->GetRawsData((2*iEq)+3+fGenOffset)->Fill(colM+(module%2)*160,rowM+iHalfStave*256);
214 for (Int_t ieq=0; ieq<20; ieq++)
215 for (UInt_t ierr=0; ierr<fAdvLogger->GetNrErrorCodes(); ierr++)
216 fAliITSQADataMakerRec->GetRawsData((2*ieq)+4+fGenOffset)->Fill(ierr,fAdvLogger->GetNrErrors(ierr,ieq));
220 fAliITSQADataMakerRec->GetRawsData(43+fGenOffset)->Fill(nDigitsL1);
221 fAliITSQADataMakerRec->GetRawsData(44+fGenOffset)->Fill(nDigitsL2);
222 fAliITSQADataMakerRec->GetRawsData(45+fGenOffset)->Fill(nDigitsL1,nDigitsL2);
225 AliDebug(1,Form("Event completed, %d raw digits read",nDigitsL1+nDigitsL2));
228 //____________________________________________________________________________
229 void AliITSQASPDDataMakerRec::InitRecPoints()
231 // Initialization for RECPOINTS - SPD -
232 fGenOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();
234 TH1F* hlayer= new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
235 hlayer->GetXaxis()->SetTitle("Layer number");
236 hlayer->GetYaxis()->SetTitle("Entries");
237 fAliITSQADataMakerRec->Add2RecPointsList(hlayer, fSPDhTask+fGenOffset);
240 TH1F** hmod = new TH1F*[2];
241 TH1F** hxl = new TH1F*[2];
242 TH1F** hzl = new TH1F*[2];
243 TH1F** hxg = new TH1F*[2];
244 TH1F** hyg = new TH1F*[2];
245 TH1F** hzg = new TH1F*[2];
246 TH1F** hr = new TH1F*[2];
247 TH1F** hphi = new TH1F*[2];
248 TH1F** hMultSPDcl = new TH1F*[2];
249 TH2F** hNyNz = new TH2F*[2]; // y and z cluster length
250 TH2F** hPhiZ = new TH2F*[2];
252 Float_t xlim[2]={4.5,8.};
253 Float_t zlim[2]={15.,15.};
257 for (Int_t iLay=0;iLay<2;iLay++) {
258 sprintf(name,"ModPattern_SPD%d",iLay+1);
259 sprintf(title,"Module map - SPD Layer %d",iLay+1);
260 hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
261 hmod[iLay]->GetXaxis()->SetTitle("Module number");
262 hmod[iLay]->GetYaxis()->SetTitle("Entries");
263 fAliITSQADataMakerRec->Add2RecPointsList(hmod[iLay], fSPDhTask +fGenOffset);
266 sprintf(name,"xLoc_SPD%d",iLay+1);
267 sprintf(title,"Local x coordinate - SPD Layer %d",iLay+1);
268 hxl[iLay]=new TH1F(name,title,100,-4.,4.);
269 hxl[iLay]->GetXaxis()->SetTitle("Local x [cm]");
270 hxl[iLay]->GetYaxis()->SetTitle("Entries");
271 fAliITSQADataMakerRec->Add2RecPointsList(hxl[iLay], fSPDhTask +fGenOffset);
274 sprintf(name,"zLoc_SPD%d",iLay+1);
275 sprintf(title,"Local z coordinate - SPD Layer %d",iLay+1);
276 hzl[iLay]=new TH1F(name,title,100,-4.,4.);
277 hzl[iLay]->GetXaxis()->SetTitle("Local z [cm]");
278 hzl[iLay]->GetYaxis()->SetTitle("Entries");
279 fAliITSQADataMakerRec->Add2RecPointsList(hzl[iLay], fSPDhTask+fGenOffset);
282 sprintf(name,"xGlob_SPD%d",iLay+1);
283 sprintf(title,"Global x coordinate - SPD Layer %d",iLay+1);
284 hxg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
285 hxg[iLay]->GetXaxis()->SetTitle("Global x [cm]");
286 hxg[iLay]->GetYaxis()->SetTitle("Entries");
287 fAliITSQADataMakerRec->Add2RecPointsList(hxg[iLay],fSPDhTask+fGenOffset);
290 sprintf(name,"yGlob_SPD%d",iLay+1);
291 sprintf(title,"Global y coordinate - SPD Layer %d",iLay+1);
292 hyg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
293 hyg[iLay]->GetXaxis()->SetTitle("Global y [cm]");
294 hyg[iLay]->GetYaxis()->SetTitle("Entries");
295 fAliITSQADataMakerRec->Add2RecPointsList(hyg[iLay], fSPDhTask+fGenOffset);
298 sprintf(name,"zGlob_SPD%d",iLay+1);
299 sprintf(title,"Global z coordinate - SPD Layer %d",iLay+1);
300 hzg[iLay]=new TH1F(name,title,150,-zlim[iLay],zlim[iLay]);
301 hzg[iLay]->GetXaxis()->SetTitle("Global z [cm]");
302 hzg[iLay]->GetYaxis()->SetTitle("Entries");
303 fAliITSQADataMakerRec->Add2RecPointsList(hzg[iLay], fSPDhTask+fGenOffset);
306 sprintf(name,"r_SPD%d",iLay+1);
307 sprintf(title,"Radius - SPD Layer %d",iLay+1);
308 hr[iLay]=new TH1F(name,title,100,0.,10.);
309 hr[iLay]->GetXaxis()->SetTitle("r [cm]");
310 hr[iLay]->GetYaxis()->SetTitle("Entries");
311 fAliITSQADataMakerRec->Add2RecPointsList(hr[iLay], fSPDhTask+fGenOffset);
314 sprintf(name,"phi_SPD%d",iLay+1);
315 sprintf(title,"#varphi - SPD Layer %d",iLay+1);
316 hphi[iLay]=new TH1F(name,title,600,0.,2*TMath::Pi());
317 hphi[iLay]->GetXaxis()->SetTitle("#varphi [rad]");
318 hphi[iLay]->GetYaxis()->SetTitle("Entries");
319 fAliITSQADataMakerRec->Add2RecPointsList(hphi[iLay], fSPDhTask+fGenOffset);
322 sprintf(name,"SizeYvsZ_SPD%d",iLay+1);
323 sprintf(title,"Cluster dimension - SPD Layer %d",iLay+1);
324 hNyNz[iLay]=new TH2F(name,title,100,0.,100.,100,0.,100.);
325 hNyNz[iLay]->GetXaxis()->SetTitle("z length");
326 hNyNz[iLay]->GetYaxis()->SetTitle("y length");
327 fAliITSQADataMakerRec->Add2RecPointsList(hNyNz[iLay], fSPDhTask+fGenOffset);
330 sprintf(name,"phi_z_SPD%d",iLay+1);
331 sprintf(title,"#varphi vs z - SPD Layer %d",iLay+1);
332 hPhiZ[iLay]=new TH2F(name,title,150,-zlim[iLay],zlim[iLay],100,0.,2*TMath::Pi());
333 hPhiZ[iLay]->GetXaxis()->SetTitle("Global z [cm]");
334 hPhiZ[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
335 fAliITSQADataMakerRec->Add2RecPointsList(hPhiZ[iLay], fSPDhTask+fGenOffset);
340 TH2F *hrPhi=new TH2F("r_phi_SPD","#varphi vs r - SPD",100,0.,10.,100,0.,2*TMath::Pi());
341 hrPhi->GetXaxis()->SetTitle("r [cm]");
342 hrPhi->GetYaxis()->SetTitle("#varphi [rad]");
343 fAliITSQADataMakerRec->Add2RecPointsList(hrPhi, fSPDhTask+fGenOffset);
346 TH2F *hxy=new TH2F("x_y_SPD","Global y vs x - SPD",200,-10.,10.,200,-10.,10.);
347 hxy->GetXaxis()->SetTitle("Global x [cm]");
348 hxy->GetYaxis()->SetTitle("Global y [cm]");
349 fAliITSQADataMakerRec->Add2RecPointsList(hxy, fSPDhTask+fGenOffset);
352 for (Int_t iLay=0;iLay<2;iLay++) {
353 sprintf(name,"Multiplicity_SPD%d",iLay+1);
354 sprintf(title,"Cluster multiplicity - SPD Layer %d",iLay+1);
355 hMultSPDcl[iLay]=new TH1F(name,title,200,0.,200.);
356 hMultSPDcl[iLay]->GetXaxis()->SetTitle("Cluster multiplicity");
357 hMultSPDcl[iLay]->GetYaxis()->SetTitle("Entries");
358 fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl[iLay], fSPDhTask+fGenOffset);
362 TH2F *hMultSPDcl2MultSPDcl1 =
363 new TH2F("MultCorrelation_SPD","Cluster multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
364 hMultSPDcl2MultSPDcl1->GetXaxis()->SetTitle("Clusters multiplicity (Layer 1)");
365 hMultSPDcl2MultSPDcl1->GetYaxis()->SetTitle("Clusters multiplicity (Layer 2)");
366 fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl2MultSPDcl1, fSPDhTask+fGenOffset);
369 AliDebug(1,Form("%d SPD Recs histograms booked\n",fSPDhTask));
375 //____________________________________________________________________________
376 void AliITSQASPDDataMakerRec::MakeRecPoints(TTree * clusterTree)
378 // Fill QA for RecPoints - SPD -
379 static TClonesArray statITSCluster("AliITSRecPoint");
380 TClonesArray *ITSCluster = &statITSCluster;
381 TBranch* itsClusterBranch=clusterTree->GetBranch("ITSRecPoints");
382 if (!itsClusterBranch) {
383 AliError("can't get the branch with the ITS clusters !");
386 itsClusterBranch->SetAddress(&ITSCluster);
387 Int_t nItsMods = (Int_t)clusterTree->GetEntries();
389 Float_t cluGlo[3] = {0.,0.,0.};
390 Int_t nClusters[2] = {0,0};
392 for (Int_t iIts=0; iIts < nItsMods; iIts++) {
394 if (!clusterTree->GetEvent(iIts)) continue;
395 Int_t nCluster = ITSCluster->GetEntriesFast();
396 // loop over clusters
398 AliITSRecPoint* cluster = (AliITSRecPoint*)ITSCluster->UncheckedAt(nCluster);
400 if (cluster->GetLayer()>1) continue;
401 Int_t lay=cluster->GetLayer();
402 fAliITSQADataMakerRec->GetRecPointsData(0 +fGenOffset)->Fill(lay);
403 cluster->GetGlobalXYZ(cluGlo);
404 Float_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]);
405 Float_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
407 fAliITSQADataMakerRec->GetRecPointsData(1 +fGenOffset)->Fill(iIts);
408 fAliITSQADataMakerRec->GetRecPointsData(2 +fGenOffset)->Fill(cluster->GetDetLocalX());
409 fAliITSQADataMakerRec->GetRecPointsData(3 +fGenOffset)->Fill(cluster->GetDetLocalZ());
410 fAliITSQADataMakerRec->GetRecPointsData(4 +fGenOffset)->Fill(cluGlo[0]);
411 fAliITSQADataMakerRec->GetRecPointsData(5 +fGenOffset)->Fill(cluGlo[1]);
412 fAliITSQADataMakerRec->GetRecPointsData(6 +fGenOffset)->Fill(cluGlo[2]);
413 fAliITSQADataMakerRec->GetRecPointsData(7 +fGenOffset)->Fill(rad);
414 fAliITSQADataMakerRec->GetRecPointsData(8 +fGenOffset)->Fill(phi);
415 fAliITSQADataMakerRec->GetRecPointsData(9 +fGenOffset)->Fill(cluster->GetNz(),cluster->GetNy());
416 fAliITSQADataMakerRec->GetRecPointsData(10 +fGenOffset)->Fill(cluGlo[2],phi);
418 fAliITSQADataMakerRec->GetRecPointsData(11 +fGenOffset)->Fill(iIts);
419 fAliITSQADataMakerRec->GetRecPointsData(12 +fGenOffset)->Fill(cluster->GetDetLocalX());
420 fAliITSQADataMakerRec->GetRecPointsData(13 +fGenOffset)->Fill(cluster->GetDetLocalZ());
421 fAliITSQADataMakerRec->GetRecPointsData(14 +fGenOffset)->Fill(cluGlo[0]);
422 fAliITSQADataMakerRec->GetRecPointsData(15 +fGenOffset)->Fill(cluGlo[1]);
423 fAliITSQADataMakerRec->GetRecPointsData(16 +fGenOffset)->Fill(cluGlo[2]);
424 fAliITSQADataMakerRec->GetRecPointsData(17 +fGenOffset)->Fill(rad);
425 fAliITSQADataMakerRec->GetRecPointsData(18 +fGenOffset)->Fill(phi);
426 fAliITSQADataMakerRec->GetRecPointsData(19 +fGenOffset)->Fill(cluster->GetNz(),cluster->GetNy());
427 fAliITSQADataMakerRec->GetRecPointsData(20 +fGenOffset)->Fill(cluGlo[2],phi);
429 fAliITSQADataMakerRec->GetRecPointsData(21 +fGenOffset)->Fill(rad,phi);
430 fAliITSQADataMakerRec->GetRecPointsData(22 +fGenOffset)->Fill(cluGlo[0],cluGlo[1]);
433 } // end of cluster loop
434 } // end of its "subdetector" loop
436 for (Int_t iLay=0; iLay<2; iLay++)
437 fAliITSQADataMakerRec->GetRecPointsData(23+iLay +fGenOffset)->Fill(nClusters[iLay]);
439 fAliITSQADataMakerRec->GetRecPointsData(25 +fGenOffset)->Fill(nClusters[0],nClusters[1]);
441 statITSCluster.Clear();