]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASPDDataMakerRec.cxx
Implemented checker subdetector splitting with subdetector offsets. Checker enabled
[u/mrichter/AliRoot.git] / ITS / AliITSQASPDDataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15 /* $Id:$  */
16 //  *************************************************************
17 //  Checks the quality assurance 
18 //  by comparing with reference data
19 //  contained in a DB
20 //  -------------------------------------------------------------
21 //  W. Ferrarese + P. Cerello Feb 2008
22 //  INFN Torino
23 //  M. Nicassio D. Elia INFN Bari March 2008
24 //  maria.nicassio@ba.infn.it
25     
26
27 // --- ROOT system ---
28 #include <TTree.h>
29 #include <TH1.h>
30 #include <TH2.h>
31 #include <TMath.h>
32 // --- Standard library ---
33
34 // --- AliRoot header files ---
35 #include "AliITSQADataMakerRec.h"
36 #include "AliITSQASPDDataMakerRec.h"
37 #include "AliLog.h"
38 #include "AliQA.h"
39 #include "AliRawReader.h"
40 #include "AliITSRawStreamSPD.h"
41 #include "AliITSRawStreamSPDErrorLog.h"
42 #include "AliITSRecPoint.h"
43
44 ClassImp(AliITSQASPDDataMakerRec)
45
46 //____________________________________________________________________________
47 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc, AliITSRawStreamSPDErrorLog *aliITSRawStreamSPDErrorLog) :
48 TObject(),
49 fAliITSQADataMakerRec(aliITSQADataMakerRec),
50 fkOnline(kMode),
51 fLDC(ldc),
52 fSPDhRaws(0),
53 fSPDhRecs(0),
54 fGenOffset(0),
55 fAdvLogger(aliITSRawStreamSPDErrorLog) 
56 {
57   //ctor used to discriminate OnLine-Offline analysis  
58 }
59
60 //____________________________________________________________________________ 
61 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(const AliITSQASPDDataMakerRec& qadm) :
62 TObject(),
63 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
64 fkOnline(qadm.fkOnline),
65 fLDC(qadm.fLDC),
66 fSPDhRaws(qadm.fSPDhRaws),
67 fSPDhRecs(qadm.fSPDhRecs),
68 fGenOffset(qadm.fGenOffset),
69 fAdvLogger(qadm.fAdvLogger)
70 {
71   //copy ctor 
72   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; 
73   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
74   }
75
76 //__________________________________________________________________
77 AliITSQASPDDataMakerRec::~AliITSQASPDDataMakerRec(){
78   // destructor
79 //  delete fAdvLogger;
80 }
81 //__________________________________________________________________
82
83 AliITSQASPDDataMakerRec& AliITSQASPDDataMakerRec::operator = (const AliITSQASPDDataMakerRec& qac )
84 {
85   // Equal operator.
86   this->~AliITSQASPDDataMakerRec();
87   new(this) AliITSQASPDDataMakerRec(qac);
88   return *this;
89 }
90
91 //____________________________________________________________________________ 
92 void AliITSQASPDDataMakerRec::StartOfDetectorCycle()
93 {
94   //Detector specific actions at start of cycle
95   AliDebug(1,"AliITSQADM::Start of SPD Cycle\n");
96 }
97
98 //____________________________________________________________________________ 
99 void AliITSQASPDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t /*task*/, TObjArray* /*list*/)
100 {
101   // launch the QA checking
102   AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n"); 
103   
104   //AliQAChecker::Instance()->Run( AliQA::kITS , task, list);
105 }
106
107 //____________________________________________________________________________ 
108 void AliITSQASPDDataMakerRec::InitRaws()
109
110   // Initialization for RAW data - SPD -
111   fGenOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();
112
113   Char_t name[50];
114   Char_t title[50];
115
116   TH1F *hlayer = new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
117   hlayer->GetXaxis()->SetTitle("Layer number");
118   hlayer->GetYaxis()->SetTitle("Entries");
119   fAliITSQADataMakerRec->Add2RawsList(hlayer,fSPDhRaws+fGenOffset);
120   fSPDhRaws++;
121
122   TH1F **hmod = new TH1F*[2];
123   TH2F **hhitMap = new TH2F*[20];
124   TH1F **herrors = new TH1F*[20];
125   for (Int_t iLay=0; iLay<2; iLay++) {
126     sprintf(name,"ModPattern_SPD%d",iLay+1);
127     sprintf(title,"Module map - SPD Layer %d",iLay+1);
128     hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
129     hmod[iLay]->GetXaxis()->SetTitle("Module number");
130     hmod[iLay]->GetYaxis()->SetTitle("Entries");
131     fAliITSQADataMakerRec->Add2RawsList(hmod[iLay], fSPDhRaws +fGenOffset);
132     fSPDhRaws++;
133   }
134   fAdvLogger = new AliITSRawStreamSPDErrorLog();  
135   for (Int_t iDDL=0; iDDL<20; iDDL++) {
136     sprintf(name,"HitMap_SPD_DDL%d",iDDL+1);
137     sprintf(title,"Hit map - SPD DDL %d",iDDL+1);
138     hhitMap[iDDL]=new TH2F(name,title,320,0,10*32,1536,0,6*256);
139     hhitMap[iDDL]->GetXaxis()->SetTitle("Column");
140     hhitMap[iDDL]->GetYaxis()->SetTitle("Row");
141     fAliITSQADataMakerRec->Add2RawsList(hhitMap[iDDL], fSPDhRaws +fGenOffset);
142     fSPDhRaws++;
143     sprintf(name,"Errors_SPD_DDL%d",iDDL+1);
144     sprintf(title,"Error codes - SPD DDL %d",iDDL+1);
145     herrors[iDDL] = new TH1F (name,title,15,0,15);
146     herrors[iDDL]->SetXTitle("Error Code");
147     herrors[iDDL]->SetYTitle("Nr of errors");
148     fAliITSQADataMakerRec->Add2RawsList(herrors[iDDL], fSPDhRaws +fGenOffset);
149     fSPDhRaws++;
150   }
151
152   TH1F** hMultSPDhits = new TH1F*[2];
153   for (Int_t iLay=0; iLay<2; iLay++) {
154     sprintf(name,"HitsMultiplicity_SPD%d",iLay+1);
155     sprintf(title,"Hit multiplicity - SPD Layer %d",iLay+1);
156     hMultSPDhits[iLay]=new TH1F(name,title,200,0.,200.);
157     hMultSPDhits[iLay]->GetXaxis()->SetTitle("Hit multiplicity");
158     hMultSPDhits[iLay]->GetYaxis()->SetTitle("Entries");
159     fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits[iLay], fSPDhRaws+fGenOffset);
160     fSPDhRaws++;
161   }
162
163   TH2F *hMultSPDhits2MultSPDhits1 = new TH2F("HitMultCorrelation_SPD","Hit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
164   hMultSPDhits2MultSPDhits1->GetXaxis()->SetTitle("Hit multiplicity (Layer 1)");
165   hMultSPDhits2MultSPDhits1->GetYaxis()->SetTitle("Hit multiplicity (Layer 2)");
166   fAliITSQADataMakerRec->Add2RawsList(hMultSPDhits2MultSPDhits1, fSPDhRaws+fGenOffset);
167   fSPDhRaws++;
168  
169   AliDebug(1,Form("%d SPD Raws histograms booked\n",fSPDhRaws));
170
171 }
172
173
174 //____________________________________________________________________________
175 void AliITSQASPDDataMakerRec::MakeRaws(AliRawReader* rawReader)
176
177   // Fill QA for RAW - SPD -
178   rawReader->Reset();
179   AliITSRawStreamSPD *rawStreamSPD = new AliITSRawStreamSPD(rawReader);
180   rawStreamSPD->ActivateAdvancedErrorLog(kTRUE,fAdvLogger);
181
182   Int_t nDigitsL1 = 0;
183   Int_t nDigitsL2 = 0;
184   Int_t iEq;
185   Int_t iLayer;
186   Int_t iHalfStave, iChip;
187   Int_t col, row; 
188   UInt_t module, colM, rowM;
189   while(rawStreamSPD->Next()) {
190
191     iEq = rawReader->GetDDLID();
192     if (iEq>=0 && iEq<20) {
193       iHalfStave = rawStreamSPD->GetHalfStaveNr();
194       iChip = rawStreamSPD->GetChipAddr();
195       col  = rawStreamSPD->GetChipCol();
196       row  = rawStreamSPD->GetChipRow();
197
198       rawStreamSPD->OnlineToOffline(iEq, iHalfStave, iChip, col, row, module, colM, rowM);
199
200       if (iHalfStave>=0 && iHalfStave<2) iLayer=0;
201       else iLayer=1;
202       
203       fAliITSQADataMakerRec->GetRawsData(0+fGenOffset)->Fill(iLayer);
204       if (iLayer==0) {
205         fAliITSQADataMakerRec->GetRawsData(1+fGenOffset)->Fill(module);
206         nDigitsL1++;
207       } else {
208         fAliITSQADataMakerRec->GetRawsData(2+fGenOffset)->Fill(module);
209         nDigitsL2++;
210       }
211       
212       fAliITSQADataMakerRec->GetRawsData((2*iEq)+3+fGenOffset)->Fill(colM+(module%2)*160,rowM+iHalfStave*256);
213     }
214
215   }
216   for (Int_t ieq=0; ieq<20; ieq++)
217     for (UInt_t ierr=0; ierr<fAdvLogger->GetNrErrorCodes(); ierr++)
218       fAliITSQADataMakerRec->GetRawsData((2*ieq)+4+fGenOffset)->Fill(ierr,fAdvLogger->GetNrErrors(ierr,ieq));
219
220   fAdvLogger->Reset();
221  
222   fAliITSQADataMakerRec->GetRawsData(43+fGenOffset)->Fill(nDigitsL1);
223   fAliITSQADataMakerRec->GetRawsData(44+fGenOffset)->Fill(nDigitsL2);
224   fAliITSQADataMakerRec->GetRawsData(45+fGenOffset)->Fill(nDigitsL1,nDigitsL2);
225   
226   delete rawStreamSPD;  
227   AliDebug(1,Form("Event completed, %d raw digits read",nDigitsL1+nDigitsL2));
228 }
229
230 //____________________________________________________________________________ 
231 void AliITSQASPDDataMakerRec::InitRecPoints()
232 {
233   // Initialization for RECPOINTS - SPD -
234   fGenOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();
235
236   TH1F* hlayer= new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
237   hlayer->GetXaxis()->SetTitle("Layer number");
238   hlayer->GetYaxis()->SetTitle("Entries");
239   fAliITSQADataMakerRec->Add2RecPointsList(hlayer, fSPDhRecs+fGenOffset); 
240   fSPDhRecs++;
241
242   TH1F** hmod = new TH1F*[2];
243   TH1F** hxl  = new TH1F*[2];
244   TH1F** hzl  = new TH1F*[2];
245   TH1F** hxg  = new TH1F*[2];
246   TH1F** hyg  = new TH1F*[2];
247   TH1F** hzg  = new TH1F*[2];
248   TH1F** hr   = new TH1F*[2];
249   TH1F** hphi = new TH1F*[2];
250   TH1F** hMultSPDcl = new TH1F*[2];
251   TH2F** hNyNz = new TH2F*[2];  // y and z cluster length
252   TH2F** hPhiZ = new TH2F*[2];
253
254   Float_t xlim[2]={4.5,8.};
255   Float_t zlim[2]={15.,15.};
256
257   Char_t name[50];
258   Char_t title[50];
259   for (Int_t iLay=0;iLay<2;iLay++) {
260     sprintf(name,"ModPattern_SPD%d",iLay+1);
261     sprintf(title,"Module map - SPD Layer %d",iLay+1);
262     hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
263     hmod[iLay]->GetXaxis()->SetTitle("Module number");
264     hmod[iLay]->GetYaxis()->SetTitle("Entries");
265     fAliITSQADataMakerRec->Add2RecPointsList(hmod[iLay], fSPDhRecs +fGenOffset); 
266     fSPDhRecs++;
267
268     sprintf(name,"xLoc_SPD%d",iLay+1);
269     sprintf(title,"Local x coordinate - SPD Layer %d",iLay+1);
270     hxl[iLay]=new TH1F(name,title,100,-4.,4.);
271     hxl[iLay]->GetXaxis()->SetTitle("Local x [cm]");
272     hxl[iLay]->GetYaxis()->SetTitle("Entries");
273     fAliITSQADataMakerRec->Add2RecPointsList(hxl[iLay], fSPDhRecs +fGenOffset);
274     fSPDhRecs++;
275
276     sprintf(name,"zLoc_SPD%d",iLay+1);
277     sprintf(title,"Local z coordinate - SPD Layer %d",iLay+1);
278     hzl[iLay]=new TH1F(name,title,100,-4.,4.);
279     hzl[iLay]->GetXaxis()->SetTitle("Local z [cm]");
280     hzl[iLay]->GetYaxis()->SetTitle("Entries");
281     fAliITSQADataMakerRec->Add2RecPointsList(hzl[iLay], fSPDhRecs+fGenOffset); 
282     fSPDhRecs++;
283
284     sprintf(name,"xGlob_SPD%d",iLay+1);
285     sprintf(title,"Global x coordinate - SPD Layer %d",iLay+1);
286     hxg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
287     hxg[iLay]->GetXaxis()->SetTitle("Global x [cm]");
288     hxg[iLay]->GetYaxis()->SetTitle("Entries");
289     fAliITSQADataMakerRec->Add2RecPointsList(hxg[iLay],fSPDhRecs+fGenOffset);  
290     fSPDhRecs++;
291
292     sprintf(name,"yGlob_SPD%d",iLay+1);
293     sprintf(title,"Global y coordinate - SPD Layer %d",iLay+1);
294     hyg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
295     hyg[iLay]->GetXaxis()->SetTitle("Global y [cm]");
296     hyg[iLay]->GetYaxis()->SetTitle("Entries");
297     fAliITSQADataMakerRec->Add2RecPointsList(hyg[iLay], fSPDhRecs+fGenOffset); 
298     fSPDhRecs++;
299
300     sprintf(name,"zGlob_SPD%d",iLay+1);
301     sprintf(title,"Global z coordinate - SPD Layer %d",iLay+1);
302     hzg[iLay]=new TH1F(name,title,150,-zlim[iLay],zlim[iLay]);
303     hzg[iLay]->GetXaxis()->SetTitle("Global z [cm]");
304     hzg[iLay]->GetYaxis()->SetTitle("Entries");
305     fAliITSQADataMakerRec->Add2RecPointsList(hzg[iLay], fSPDhRecs+fGenOffset); 
306     fSPDhRecs++;
307
308     sprintf(name,"r_SPD%d",iLay+1);
309     sprintf(title,"Radius - SPD Layer %d",iLay+1);
310     hr[iLay]=new TH1F(name,title,100,0.,10.);
311     hr[iLay]->GetXaxis()->SetTitle("r [cm]");
312     hr[iLay]->GetYaxis()->SetTitle("Entries");
313     fAliITSQADataMakerRec->Add2RecPointsList(hr[iLay], fSPDhRecs+fGenOffset);  
314     fSPDhRecs++;
315
316     sprintf(name,"phi_SPD%d",iLay+1);
317     sprintf(title,"#varphi - SPD Layer %d",iLay+1);
318     hphi[iLay]=new TH1F(name,title,600,0.,2*TMath::Pi());
319     hphi[iLay]->GetXaxis()->SetTitle("#varphi [rad]");
320     hphi[iLay]->GetYaxis()->SetTitle("Entries");
321     fAliITSQADataMakerRec->Add2RecPointsList(hphi[iLay], fSPDhRecs+fGenOffset);
322     fSPDhRecs++;
323     
324     sprintf(name,"SizeYvsZ_SPD%d",iLay+1);
325     sprintf(title,"Cluster dimension - SPD Layer %d",iLay+1);
326     hNyNz[iLay]=new TH2F(name,title,100,0.,100.,100,0.,100.);
327     hNyNz[iLay]->GetXaxis()->SetTitle("z length");
328     hNyNz[iLay]->GetYaxis()->SetTitle("y length");
329     fAliITSQADataMakerRec->Add2RecPointsList(hNyNz[iLay], fSPDhRecs+fGenOffset); 
330     fSPDhRecs++;
331
332     sprintf(name,"phi_z_SPD%d",iLay+1);
333     sprintf(title,"#varphi vs z - SPD Layer %d",iLay+1);
334     hPhiZ[iLay]=new TH2F(name,title,150,-zlim[iLay],zlim[iLay],100,0.,2*TMath::Pi());
335     hPhiZ[iLay]->GetXaxis()->SetTitle("Global z [cm]");
336     hPhiZ[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
337     fAliITSQADataMakerRec->Add2RecPointsList(hPhiZ[iLay], fSPDhRecs+fGenOffset);
338     fSPDhRecs++;
339
340   }
341
342   TH2F *hrPhi=new TH2F("r_phi_SPD","#varphi vs r - SPD",100,0.,10.,100,0.,2*TMath::Pi());
343   hrPhi->GetXaxis()->SetTitle("r [cm]");
344   hrPhi->GetYaxis()->SetTitle("#varphi [rad]");
345   fAliITSQADataMakerRec->Add2RecPointsList(hrPhi, fSPDhRecs+fGenOffset);
346   fSPDhRecs++;
347
348   TH2F *hxy=new TH2F("x_y_SPD","Global y vs x - SPD",200,-10.,10.,200,-10.,10.);
349   hxy->GetXaxis()->SetTitle("Global x [cm]");
350   hxy->GetYaxis()->SetTitle("Global y [cm]");
351   fAliITSQADataMakerRec->Add2RecPointsList(hxy, fSPDhRecs+fGenOffset);
352   fSPDhRecs++;
353
354   for (Int_t iLay=0;iLay<2;iLay++) {
355     sprintf(name,"Multiplicity_SPD%d",iLay+1);
356     sprintf(title,"Cluster multiplicity - SPD Layer %d",iLay+1);
357     hMultSPDcl[iLay]=new TH1F(name,title,200,0.,200.);
358     hMultSPDcl[iLay]->GetXaxis()->SetTitle("Cluster multiplicity");
359     hMultSPDcl[iLay]->GetYaxis()->SetTitle("Entries");
360     fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl[iLay], fSPDhRecs+fGenOffset);
361     fSPDhRecs++;
362   } 
363
364   TH2F *hMultSPDcl2MultSPDcl1 =
365             new TH2F("MultCorrelation_SPD","Cluster multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
366   hMultSPDcl2MultSPDcl1->GetXaxis()->SetTitle("Clusters multiplicity (Layer 1)");
367   hMultSPDcl2MultSPDcl1->GetYaxis()->SetTitle("Clusters multiplicity (Layer 2)"); 
368   fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl2MultSPDcl1, fSPDhRecs+fGenOffset);
369   fSPDhRecs++;
370
371   AliDebug(1,Form("%d SPD Recs histograms booked\n",fSPDhRecs));
372
373 }
374
375 //____________________________________________________________________________ 
376 void AliITSQASPDDataMakerRec::MakeRecPoints(TTree * clusterTree)
377 {
378   // Fill QA for RecPoints - SPD -
379       TClonesArray* itsClusters = new TClonesArray("AliITSRecPoint");
380       TBranch* itsClusterBranch=clusterTree->GetBranch("ITSRecPoints");
381       if (!itsClusterBranch) {
382         AliError("can't get the branch with the ITS clusters !");
383         return;
384       }
385       itsClusterBranch->SetAddress(&itsClusters);
386       Int_t nItsMods = (Int_t)clusterTree->GetEntries();
387
388       Float_t cluGlo[3] = {0.,0.,0.};
389       Int_t nClusters[2] = {0,0};
390
391       for (Int_t iIts=0; iIts < nItsMods; iIts++) {
392
393       if (!clusterTree->GetEvent(iIts))    continue;
394       Int_t nCluster = itsClusters->GetEntriesFast();
395       // loop over clusters
396       while(nCluster--) {
397         AliITSRecPoint* cluster = (AliITSRecPoint*)itsClusters->UncheckedAt(nCluster);
398
399         if (cluster->GetLayer()>1)        continue;
400         Int_t lay=cluster->GetLayer();
401         fAliITSQADataMakerRec->GetRecPointsData(0 +fGenOffset)->Fill(lay);
402         cluster->GetGlobalXYZ(cluGlo);
403         Float_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]);
404         Float_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
405         if (lay==0) {
406           fAliITSQADataMakerRec->GetRecPointsData(1 +fGenOffset)->Fill(iIts);
407           fAliITSQADataMakerRec->GetRecPointsData(2 +fGenOffset)->Fill(cluster->GetDetLocalX());
408           fAliITSQADataMakerRec->GetRecPointsData(3 +fGenOffset)->Fill(cluster->GetDetLocalZ());
409           fAliITSQADataMakerRec->GetRecPointsData(4 +fGenOffset)->Fill(cluGlo[0]);
410           fAliITSQADataMakerRec->GetRecPointsData(5 +fGenOffset)->Fill(cluGlo[1]);
411           fAliITSQADataMakerRec->GetRecPointsData(6 +fGenOffset)->Fill(cluGlo[2]);
412           fAliITSQADataMakerRec->GetRecPointsData(7 +fGenOffset)->Fill(rad);
413           fAliITSQADataMakerRec->GetRecPointsData(8 +fGenOffset)->Fill(phi);
414           fAliITSQADataMakerRec->GetRecPointsData(9 +fGenOffset)->Fill(cluster->GetNz(),cluster->GetNy());
415           fAliITSQADataMakerRec->GetRecPointsData(10 +fGenOffset)->Fill(cluGlo[2],phi);
416         } else  {
417           fAliITSQADataMakerRec->GetRecPointsData(11 +fGenOffset)->Fill(iIts);
418           fAliITSQADataMakerRec->GetRecPointsData(12 +fGenOffset)->Fill(cluster->GetDetLocalX());
419           fAliITSQADataMakerRec->GetRecPointsData(13 +fGenOffset)->Fill(cluster->GetDetLocalZ());
420           fAliITSQADataMakerRec->GetRecPointsData(14 +fGenOffset)->Fill(cluGlo[0]);
421           fAliITSQADataMakerRec->GetRecPointsData(15 +fGenOffset)->Fill(cluGlo[1]);
422           fAliITSQADataMakerRec->GetRecPointsData(16 +fGenOffset)->Fill(cluGlo[2]);
423           fAliITSQADataMakerRec->GetRecPointsData(17 +fGenOffset)->Fill(rad);
424           fAliITSQADataMakerRec->GetRecPointsData(18 +fGenOffset)->Fill(phi);
425           fAliITSQADataMakerRec->GetRecPointsData(19 +fGenOffset)->Fill(cluster->GetNz(),cluster->GetNy());
426           fAliITSQADataMakerRec->GetRecPointsData(20 +fGenOffset)->Fill(cluGlo[2],phi);
427         }
428         fAliITSQADataMakerRec->GetRecPointsData(21 +fGenOffset)->Fill(rad,phi);
429         fAliITSQADataMakerRec->GetRecPointsData(22 +fGenOffset)->Fill(cluGlo[0],cluGlo[1]);
430
431         nClusters[lay]++; 
432       } // end of cluster loop
433     } // end of its "subdetector" loop
434
435     for (Int_t iLay=0; iLay<2; iLay++)
436       fAliITSQADataMakerRec->GetRecPointsData(23+iLay +fGenOffset)->Fill(nClusters[iLay]);
437
438     fAliITSQADataMakerRec->GetRecPointsData(25 +fGenOffset)->Fill(nClusters[0],nClusters[1]);
439
440     if (itsClusters) {
441       itsClusters->Delete();
442       delete itsClusters;
443       itsClusters = 0;
444     }
445
446 }