]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASPDDataMakerRec.cxx
Added method SetTaskOffset in ITS checkers. Updated data members containing the numbe...
[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 fSPDhTask(0),
53 fGenOffset(0),
54 fAdvLogger(aliITSRawStreamSPDErrorLog) 
55 {
56   //ctor used to discriminate OnLine-Offline analysis  
57 }
58
59 //____________________________________________________________________________ 
60 AliITSQASPDDataMakerRec::AliITSQASPDDataMakerRec(const AliITSQASPDDataMakerRec& qadm) :
61 TObject(),
62 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
63 fkOnline(qadm.fkOnline),
64 fLDC(qadm.fLDC),
65 fSPDhTask(qadm.fSPDhTask),
66 fGenOffset(qadm.fGenOffset),
67 fAdvLogger(qadm.fAdvLogger)
68 {
69   //copy ctor 
70   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; 
71   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
72   }
73
74 //__________________________________________________________________
75 AliITSQASPDDataMakerRec::~AliITSQASPDDataMakerRec(){
76   // destructor
77 //  delete fAdvLogger;
78 }
79 //__________________________________________________________________
80
81 AliITSQASPDDataMakerRec& AliITSQASPDDataMakerRec::operator = (const AliITSQASPDDataMakerRec& qac )
82 {
83   // Equal operator.
84   this->~AliITSQASPDDataMakerRec();
85   new(this) AliITSQASPDDataMakerRec(qac);
86   return *this;
87 }
88
89 //____________________________________________________________________________ 
90 void AliITSQASPDDataMakerRec::StartOfDetectorCycle()
91 {
92   //Detector specific actions at start of cycle
93   AliDebug(1,"AliITSQADM::Start of SPD Cycle\n");
94 }
95
96 //____________________________________________________________________________ 
97 void AliITSQASPDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t /*task*/, TObjArray* /*list*/)
98 {
99   // launch the QA checking
100   AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n"); 
101   
102   //AliQAChecker::Instance()->Run( AliQA::kITS , task, list);
103 }
104
105 //____________________________________________________________________________ 
106 void AliITSQASPDDataMakerRec::InitRaws()
107
108   // Initialization for RAW data - SPD -
109   fGenOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();
110
111   Char_t name[50];
112   Char_t title[50];
113
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);
118   fSPDhTask++;
119
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);
130     fSPDhTask++;
131   }
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);
140     fSPDhTask++;
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);
147     fSPDhTask++;
148   }
149
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);
158     fSPDhTask++;
159   }
160
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);
165   fSPDhTask++;
166  
167   AliDebug(1,Form("%d SPD Raws histograms booked\n",fSPDhTask));
168 }
169
170
171 //____________________________________________________________________________
172 void AliITSQASPDDataMakerRec::MakeRaws(AliRawReader* rawReader)
173
174   // Fill QA for RAW - SPD -
175   rawReader->Reset();
176   AliITSRawStreamSPD *rawStreamSPD = new AliITSRawStreamSPD(rawReader);
177   rawStreamSPD->ActivateAdvancedErrorLog(kTRUE,fAdvLogger);
178
179   Int_t nDigitsL1 = 0;
180   Int_t nDigitsL2 = 0;
181   Int_t iEq;
182   Int_t iLayer;
183   Int_t iHalfStave, iChip;
184   Int_t col, row; 
185   UInt_t module, colM, rowM;
186   while(rawStreamSPD->Next()) {
187
188     iEq = rawReader->GetDDLID();
189     if (iEq>=0 && iEq<20) {
190       iHalfStave = rawStreamSPD->GetHalfStaveNr();
191       iChip = rawStreamSPD->GetChipAddr();
192       col  = rawStreamSPD->GetChipCol();
193       row  = rawStreamSPD->GetChipRow();
194
195       rawStreamSPD->OnlineToOffline(iEq, iHalfStave, iChip, col, row, module, colM, rowM);
196
197       if (iHalfStave>=0 && iHalfStave<2) iLayer=0;
198       else iLayer=1;
199       
200       fAliITSQADataMakerRec->GetRawsData(0+fGenOffset)->Fill(iLayer);
201       if (iLayer==0) {
202         fAliITSQADataMakerRec->GetRawsData(1+fGenOffset)->Fill(module);
203         nDigitsL1++;
204       } else {
205         fAliITSQADataMakerRec->GetRawsData(2+fGenOffset)->Fill(module);
206         nDigitsL2++;
207       }
208       
209       fAliITSQADataMakerRec->GetRawsData((2*iEq)+3+fGenOffset)->Fill(colM+(module%2)*160,rowM+iHalfStave*256);
210     }
211
212   }
213   for (Int_t ieq=0; ieq<20; ieq++)
214     for (UInt_t ierr=0; ierr<fAdvLogger->GetNrErrorCodes(); ierr++)
215       fAliITSQADataMakerRec->GetRawsData((2*ieq)+4+fGenOffset)->Fill(ierr,fAdvLogger->GetNrErrors(ierr,ieq));
216
217   fAdvLogger->Reset();
218  
219   fAliITSQADataMakerRec->GetRawsData(43+fGenOffset)->Fill(nDigitsL1);
220   fAliITSQADataMakerRec->GetRawsData(44+fGenOffset)->Fill(nDigitsL2);
221   fAliITSQADataMakerRec->GetRawsData(45+fGenOffset)->Fill(nDigitsL1,nDigitsL2);
222   
223   delete rawStreamSPD;  
224   AliDebug(1,Form("Event completed, %d raw digits read",nDigitsL1+nDigitsL2));
225 }
226
227 //____________________________________________________________________________ 
228 void AliITSQASPDDataMakerRec::InitRecPoints()
229 {
230   // Initialization for RECPOINTS - SPD -
231   fGenOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();
232
233   TH1F* hlayer= new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
234   hlayer->GetXaxis()->SetTitle("Layer number");
235   hlayer->GetYaxis()->SetTitle("Entries");
236   fAliITSQADataMakerRec->Add2RecPointsList(hlayer, fSPDhTask+fGenOffset); 
237   fSPDhTask++;
238
239   TH1F** hmod = new TH1F*[2];
240   TH1F** hxl  = new TH1F*[2];
241   TH1F** hzl  = new TH1F*[2];
242   TH1F** hxg  = new TH1F*[2];
243   TH1F** hyg  = new TH1F*[2];
244   TH1F** hzg  = new TH1F*[2];
245   TH1F** hr   = new TH1F*[2];
246   TH1F** hphi = new TH1F*[2];
247   TH1F** hMultSPDcl = new TH1F*[2];
248   TH2F** hNyNz = new TH2F*[2];  // y and z cluster length
249   TH2F** hPhiZ = new TH2F*[2];
250
251   Float_t xlim[2]={4.5,8.};
252   Float_t zlim[2]={15.,15.};
253
254   Char_t name[50];
255   Char_t title[50];
256   for (Int_t iLay=0;iLay<2;iLay++) {
257     sprintf(name,"ModPattern_SPD%d",iLay+1);
258     sprintf(title,"Module map - SPD Layer %d",iLay+1);
259     hmod[iLay]=new TH1F(name,title,fgknSPDmodules,0,fgknSPDmodules);
260     hmod[iLay]->GetXaxis()->SetTitle("Module number");
261     hmod[iLay]->GetYaxis()->SetTitle("Entries");
262     fAliITSQADataMakerRec->Add2RecPointsList(hmod[iLay], fSPDhTask +fGenOffset); 
263     fSPDhTask++;
264
265     sprintf(name,"xLoc_SPD%d",iLay+1);
266     sprintf(title,"Local x coordinate - SPD Layer %d",iLay+1);
267     hxl[iLay]=new TH1F(name,title,100,-4.,4.);
268     hxl[iLay]->GetXaxis()->SetTitle("Local x [cm]");
269     hxl[iLay]->GetYaxis()->SetTitle("Entries");
270     fAliITSQADataMakerRec->Add2RecPointsList(hxl[iLay], fSPDhTask +fGenOffset);
271     fSPDhTask++;
272
273     sprintf(name,"zLoc_SPD%d",iLay+1);
274     sprintf(title,"Local z coordinate - SPD Layer %d",iLay+1);
275     hzl[iLay]=new TH1F(name,title,100,-4.,4.);
276     hzl[iLay]->GetXaxis()->SetTitle("Local z [cm]");
277     hzl[iLay]->GetYaxis()->SetTitle("Entries");
278     fAliITSQADataMakerRec->Add2RecPointsList(hzl[iLay], fSPDhTask+fGenOffset); 
279     fSPDhTask++;
280
281     sprintf(name,"xGlob_SPD%d",iLay+1);
282     sprintf(title,"Global x coordinate - SPD Layer %d",iLay+1);
283     hxg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
284     hxg[iLay]->GetXaxis()->SetTitle("Global x [cm]");
285     hxg[iLay]->GetYaxis()->SetTitle("Entries");
286     fAliITSQADataMakerRec->Add2RecPointsList(hxg[iLay],fSPDhTask+fGenOffset);  
287     fSPDhTask++;
288
289     sprintf(name,"yGlob_SPD%d",iLay+1);
290     sprintf(title,"Global y coordinate - SPD Layer %d",iLay+1);
291     hyg[iLay]=new TH1F(name,title,100,-xlim[iLay],xlim[iLay]);
292     hyg[iLay]->GetXaxis()->SetTitle("Global y [cm]");
293     hyg[iLay]->GetYaxis()->SetTitle("Entries");
294     fAliITSQADataMakerRec->Add2RecPointsList(hyg[iLay], fSPDhTask+fGenOffset); 
295     fSPDhTask++;
296
297     sprintf(name,"zGlob_SPD%d",iLay+1);
298     sprintf(title,"Global z coordinate - SPD Layer %d",iLay+1);
299     hzg[iLay]=new TH1F(name,title,150,-zlim[iLay],zlim[iLay]);
300     hzg[iLay]->GetXaxis()->SetTitle("Global z [cm]");
301     hzg[iLay]->GetYaxis()->SetTitle("Entries");
302     fAliITSQADataMakerRec->Add2RecPointsList(hzg[iLay], fSPDhTask+fGenOffset); 
303     fSPDhTask++;
304
305     sprintf(name,"r_SPD%d",iLay+1);
306     sprintf(title,"Radius - SPD Layer %d",iLay+1);
307     hr[iLay]=new TH1F(name,title,100,0.,10.);
308     hr[iLay]->GetXaxis()->SetTitle("r [cm]");
309     hr[iLay]->GetYaxis()->SetTitle("Entries");
310     fAliITSQADataMakerRec->Add2RecPointsList(hr[iLay], fSPDhTask+fGenOffset);  
311     fSPDhTask++;
312
313     sprintf(name,"phi_SPD%d",iLay+1);
314     sprintf(title,"#varphi - SPD Layer %d",iLay+1);
315     hphi[iLay]=new TH1F(name,title,600,0.,2*TMath::Pi());
316     hphi[iLay]->GetXaxis()->SetTitle("#varphi [rad]");
317     hphi[iLay]->GetYaxis()->SetTitle("Entries");
318     fAliITSQADataMakerRec->Add2RecPointsList(hphi[iLay], fSPDhTask+fGenOffset);
319     fSPDhTask++;
320     
321     sprintf(name,"SizeYvsZ_SPD%d",iLay+1);
322     sprintf(title,"Cluster dimension - SPD Layer %d",iLay+1);
323     hNyNz[iLay]=new TH2F(name,title,100,0.,100.,100,0.,100.);
324     hNyNz[iLay]->GetXaxis()->SetTitle("z length");
325     hNyNz[iLay]->GetYaxis()->SetTitle("y length");
326     fAliITSQADataMakerRec->Add2RecPointsList(hNyNz[iLay], fSPDhTask+fGenOffset); 
327     fSPDhTask++;
328
329     sprintf(name,"phi_z_SPD%d",iLay+1);
330     sprintf(title,"#varphi vs z - SPD Layer %d",iLay+1);
331     hPhiZ[iLay]=new TH2F(name,title,150,-zlim[iLay],zlim[iLay],100,0.,2*TMath::Pi());
332     hPhiZ[iLay]->GetXaxis()->SetTitle("Global z [cm]");
333     hPhiZ[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
334     fAliITSQADataMakerRec->Add2RecPointsList(hPhiZ[iLay], fSPDhTask+fGenOffset);
335     fSPDhTask++;
336
337   }
338
339   TH2F *hrPhi=new TH2F("r_phi_SPD","#varphi vs r - SPD",100,0.,10.,100,0.,2*TMath::Pi());
340   hrPhi->GetXaxis()->SetTitle("r [cm]");
341   hrPhi->GetYaxis()->SetTitle("#varphi [rad]");
342   fAliITSQADataMakerRec->Add2RecPointsList(hrPhi, fSPDhTask+fGenOffset);
343   fSPDhTask++;
344
345   TH2F *hxy=new TH2F("x_y_SPD","Global y vs x - SPD",200,-10.,10.,200,-10.,10.);
346   hxy->GetXaxis()->SetTitle("Global x [cm]");
347   hxy->GetYaxis()->SetTitle("Global y [cm]");
348   fAliITSQADataMakerRec->Add2RecPointsList(hxy, fSPDhTask+fGenOffset);
349   fSPDhTask++;
350
351   for (Int_t iLay=0;iLay<2;iLay++) {
352     sprintf(name,"Multiplicity_SPD%d",iLay+1);
353     sprintf(title,"Cluster multiplicity - SPD Layer %d",iLay+1);
354     hMultSPDcl[iLay]=new TH1F(name,title,200,0.,200.);
355     hMultSPDcl[iLay]->GetXaxis()->SetTitle("Cluster multiplicity");
356     hMultSPDcl[iLay]->GetYaxis()->SetTitle("Entries");
357     fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl[iLay], fSPDhTask+fGenOffset);
358     fSPDhTask++;
359   } 
360
361   TH2F *hMultSPDcl2MultSPDcl1 =
362             new TH2F("MultCorrelation_SPD","Cluster multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
363   hMultSPDcl2MultSPDcl1->GetXaxis()->SetTitle("Clusters multiplicity (Layer 1)");
364   hMultSPDcl2MultSPDcl1->GetYaxis()->SetTitle("Clusters multiplicity (Layer 2)"); 
365   fAliITSQADataMakerRec->Add2RecPointsList(hMultSPDcl2MultSPDcl1, fSPDhTask+fGenOffset);
366   fSPDhTask++;
367
368   AliDebug(1,Form("%d SPD Recs histograms booked\n",fSPDhTask));
369
370 }
371
372 //____________________________________________________________________________ 
373 void AliITSQASPDDataMakerRec::MakeRecPoints(TTree * clusterTree)
374 {
375   // Fill QA for RecPoints - SPD -
376   static TClonesArray statITSCluster("AliITSRecPoint");
377   TClonesArray *ITSCluster = &statITSCluster;
378   TBranch* itsClusterBranch=clusterTree->GetBranch("ITSRecPoints");
379   if (!itsClusterBranch) {
380     AliError("can't get the branch with the ITS clusters !");
381     return;
382   }
383   itsClusterBranch->SetAddress(&ITSCluster);
384   Int_t nItsMods = (Int_t)clusterTree->GetEntries();
385   
386   Float_t cluGlo[3] = {0.,0.,0.};
387   Int_t nClusters[2] = {0,0};
388   
389   for (Int_t iIts=0; iIts < nItsMods; iIts++) {
390     
391     if (!clusterTree->GetEvent(iIts))    continue;
392     Int_t nCluster = ITSCluster->GetEntriesFast();
393     // loop over clusters
394     while(nCluster--) {
395       AliITSRecPoint* cluster = (AliITSRecPoint*)ITSCluster->UncheckedAt(nCluster);
396       
397       if (cluster->GetLayer()>1)        continue;
398       Int_t lay=cluster->GetLayer();
399       fAliITSQADataMakerRec->GetRecPointsData(0 +fGenOffset)->Fill(lay);
400       cluster->GetGlobalXYZ(cluGlo);
401       Float_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]);
402         Float_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
403         if (lay==0) {
404           fAliITSQADataMakerRec->GetRecPointsData(1 +fGenOffset)->Fill(iIts);
405           fAliITSQADataMakerRec->GetRecPointsData(2 +fGenOffset)->Fill(cluster->GetDetLocalX());
406           fAliITSQADataMakerRec->GetRecPointsData(3 +fGenOffset)->Fill(cluster->GetDetLocalZ());
407           fAliITSQADataMakerRec->GetRecPointsData(4 +fGenOffset)->Fill(cluGlo[0]);
408           fAliITSQADataMakerRec->GetRecPointsData(5 +fGenOffset)->Fill(cluGlo[1]);
409           fAliITSQADataMakerRec->GetRecPointsData(6 +fGenOffset)->Fill(cluGlo[2]);
410           fAliITSQADataMakerRec->GetRecPointsData(7 +fGenOffset)->Fill(rad);
411           fAliITSQADataMakerRec->GetRecPointsData(8 +fGenOffset)->Fill(phi);
412           fAliITSQADataMakerRec->GetRecPointsData(9 +fGenOffset)->Fill(cluster->GetNz(),cluster->GetNy());
413           fAliITSQADataMakerRec->GetRecPointsData(10 +fGenOffset)->Fill(cluGlo[2],phi);
414         } else  {
415           fAliITSQADataMakerRec->GetRecPointsData(11 +fGenOffset)->Fill(iIts);
416           fAliITSQADataMakerRec->GetRecPointsData(12 +fGenOffset)->Fill(cluster->GetDetLocalX());
417           fAliITSQADataMakerRec->GetRecPointsData(13 +fGenOffset)->Fill(cluster->GetDetLocalZ());
418           fAliITSQADataMakerRec->GetRecPointsData(14 +fGenOffset)->Fill(cluGlo[0]);
419           fAliITSQADataMakerRec->GetRecPointsData(15 +fGenOffset)->Fill(cluGlo[1]);
420           fAliITSQADataMakerRec->GetRecPointsData(16 +fGenOffset)->Fill(cluGlo[2]);
421           fAliITSQADataMakerRec->GetRecPointsData(17 +fGenOffset)->Fill(rad);
422           fAliITSQADataMakerRec->GetRecPointsData(18 +fGenOffset)->Fill(phi);
423           fAliITSQADataMakerRec->GetRecPointsData(19 +fGenOffset)->Fill(cluster->GetNz(),cluster->GetNy());
424           fAliITSQADataMakerRec->GetRecPointsData(20 +fGenOffset)->Fill(cluGlo[2],phi);
425         }
426         fAliITSQADataMakerRec->GetRecPointsData(21 +fGenOffset)->Fill(rad,phi);
427         fAliITSQADataMakerRec->GetRecPointsData(22 +fGenOffset)->Fill(cluGlo[0],cluGlo[1]);
428         
429         nClusters[lay]++; 
430     } // end of cluster loop
431   } // end of its "subdetector" loop
432   
433   for (Int_t iLay=0; iLay<2; iLay++)
434     fAliITSQADataMakerRec->GetRecPointsData(23+iLay +fGenOffset)->Fill(nClusters[iLay]);
435   
436   fAliITSQADataMakerRec->GetRecPointsData(25 +fGenOffset)->Fill(nClusters[0],nClusters[1]);
437   
438   statITSCluster.Clear();
439 }