]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerRec.cxx
Saving some QA histograms for future corrections (Yves)
[u/mrichter/AliRoot.git] / ITS / AliITSQASDDDataMakerRec.cxx
1 /**************************************************************************\r
2  * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *\r
3  *                                                                        *\r
4  * Author: The ALICE Off-line Project.                                    *\r
5  * Contributors are mentioned in the code where appropriate.              *\r
6  *                                                                        *\r
7  * Permission to use, copy, modify and distribute this software and its   *\r
8  * documentation strictly for non-commercial purposes is hereby granted   *\r
9  * without fee, provided that the above copyright notice appears in all   *\r
10  * copies and that both the copyright notice and this permission notice   *\r
11  * appear in the supporting documentation. The authors make no claims     *\r
12  * about the suitability of this software for any purpose. It is          *\r
13  * provided "as is" without express or implied warranty.                  *\r
14  **************************************************************************/\r
15 \r
16 /* $Id$ */\r
17 \r
18 //  *************************************************************\r
19 //  Checks the quality assurance \r
20 //  by comparing with reference data\r
21 //  contained in a DB\r
22 //  -------------------------------------------------------------\r
23 //  W. Ferrarese + P. Cerello Feb 2008\r
24 //  INFN Torino\r
25 \r
26 // --- ROOT system ---\r
27 #include <TProfile2D.h>\r
28 #include <TH2D.h>\r
29 #include <TBranch.h>\r
30 #include <TTree.h>\r
31 #include <TGaxis.h>\r
32 #include <TMath.h>\r
33 // --- Standard library ---\r
34 \r
35 // --- AliRoot header files ---\r
36 #include "AliITSQASDDDataMakerRec.h"\r
37 #include "AliLog.h"\r
38 #include "AliQA.h"\r
39 #include "AliQAChecker.h"\r
40 #include "AliRawReader.h"\r
41 #include "AliITSRawStreamSDD.h"\r
42 #include "AliITSRecPoint.h"\r
43 #include "AliITSgeomTGeo.h"\r
44 \r
45 #include "AliCDBManager.h"\r
46 #include "AliCDBStorage.h"\r
47 #include "AliCDBEntry.h"\r
48 \r
49 \r
50 ClassImp(AliITSQASDDDataMakerRec)\r
51 \r
52 //____________________________________________________________________________ \r
53 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) :\r
54 TObject(),\r
55 fAliITSQADataMakerRec(aliITSQADataMakerRec),\r
56 fkOnline(kMode),\r
57 fLDC(ldc),\r
58 fSDDhTask(0),\r
59 fGenOffset(0),\r
60 fTimeBinSize(1),\r
61 fDDLModuleMap(0)\r
62 {\r
63   //ctor used to discriminate OnLine-Offline analysis\r
64   if(fLDC < 0 || fLDC > 4) {\r
65         AliError("Error: LDC number out of range; return\n");\r
66   }\r
67   //fDDLModuleMap=NULL;\r
68 }\r
69 \r
70 //____________________________________________________________________________ \r
71 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) :\r
72 TObject(),\r
73 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),\r
74 fkOnline(qadm.fkOnline),\r
75 fLDC(qadm.fLDC),\r
76 fSDDhTask(qadm.fSDDhTask),\r
77 fGenOffset(qadm.fGenOffset),\r
78 fTimeBinSize(1),\r
79 fDDLModuleMap(0)\r
80 {\r
81   //copy ctor \r
82   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; \r
83   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());\r
84   fDDLModuleMap=NULL;\r
85 }\r
86 \r
87 //____________________________________________________________________________ \r
88 AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){\r
89   // destructor\r
90   //    if(fDDLModuleMap) delete fDDLModuleMap;\r
91 }\r
92 //__________________________________________________________________\r
93 AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac )\r
94 {\r
95   // Equal operator.\r
96   this->~AliITSQASDDDataMakerRec();\r
97   new(this) AliITSQASDDDataMakerRec(qac);\r
98   return *this;\r
99 }\r
100 \r
101 //____________________________________________________________________________ \r
102 void AliITSQASDDDataMakerRec::StartOfDetectorCycle()\r
103 {\r
104   //Detector specific actions at start of cycle\r
105   AliDebug(1,"AliITSQADM::Start of SDD Cycle\n");\r
106 }\r
107 \r
108 //____________________________________________________________________________ \r
109 void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t /*task*/, TObjArray* /*list*/)\r
110 {\r
111   // launch the QA checking\r
112   AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n"); \r
113 }\r
114 \r
115 //____________________________________________________________________________ \r
116 void AliITSQASDDDataMakerRec::InitRaws()\r
117\r
118   // Initialization for RAW data - SDD -\r
119   fGenOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();\r
120 \r
121   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");\r
122   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();\r
123 \r
124   if( !ddlMapSDD){\r
125     AliError("Calibration object retrieval failed! SDD will not be processed");\r
126     fDDLModuleMap = NULL;\r
127     return;\r
128   }  \r
129   fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();\r
130   if(!cacheStatus)ddlMapSDD->SetObject(NULL);\r
131   ddlMapSDD->SetOwner(kTRUE);\r
132   if(!cacheStatus)delete ddlMapSDD;\r
133  \r
134   Int_t lay, lad, det;\r
135   Int_t LAY = -1;  //, LAD = -1;\r
136   char hname0[50];\r
137   Int_t indexlast = 0;\r
138   Int_t index1 = 0;\r
139 \r
140   if(fLDC == 1 || fLDC == 2) LAY = 2;\r
141   if(fLDC == 3 || fLDC == 4) LAY = 3;\r
142   \r
143   if(fkOnline) {\r
144     AliInfo("Book Online Histograms for SDD\n");\r
145   }\r
146   else {\r
147     AliInfo("Book Offline Histograms for SDD\n ");\r
148   }\r
149   TH1D *h0 = new TH1D("ModPattern","HW Modules pattern",fgknSDDmodules,-0.5,259.5);\r
150   h0->GetXaxis()->SetTitle("Module Number");\r
151   h0->GetYaxis()->SetTitle("Counts");\r
152   fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h0)),0+fGenOffset, kTRUE);\r
153   delete h0;\r
154   fSDDhTask++;\r
155   if(fLDC==0 || fLDC==1 || fLDC==2){\r
156     TH1D *h1 = new TH1D("LadPatternL3","Ladder pattern L3",14,0.5,14.5);  \r
157     h1->GetXaxis()->SetTitle("Ladder Number on Lay3");\r
158     h1->GetYaxis()->SetTitle("Counts");\r
159     fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h1)),1+fGenOffset, kTRUE);\r
160         delete h1;\r
161     fSDDhTask++;\r
162   }     \r
163   if(fLDC==0 || fLDC==3 || fLDC==4){\r
164     TH1D *h2 = new TH1D("LadPatternL4","Ladder pattern L4",22,0.5,22.5);  \r
165     h2->GetXaxis()->SetTitle("Ladder Number on Lay4");\r
166     h2->GetYaxis()->SetTitle("Counts");\r
167     fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h2)),2+fGenOffset, kTRUE);\r
168         delete h2;\r
169     fSDDhTask++;\r
170   }\r
171   if(fLDC==0 || fLDC==1 || fLDC==2){\r
172         for(Int_t i=1; i<=fgkLADDonLAY3; i++) {\r
173       sprintf(hname0,"ModPattern_L3_%d",i);\r
174       TH1D *h3 = new TH1D(hname0,hname0,6,0.5,6.5);\r
175       h3->GetXaxis()->SetTitle("Module Number");\r
176       h3->GetYaxis()->SetTitle("Counts");\r
177       fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h3)),i-1+3+fGenOffset, kTRUE);\r
178           delete h3;\r
179       fSDDhTask++;\r
180     }\r
181   }\r
182   if(fLDC==0 || fLDC==3 || fLDC==4){\r
183     for(Int_t i=1; i<=fgkLADDonLAY4; i++) {\r
184       sprintf(hname0,"ModPattern_L4_%d",i);\r
185           TH1D *h4 = new TH1D(hname0,hname0,8,0.5,8.5);\r
186       h4->GetXaxis()->SetTitle("Module Number");\r
187       h4->GetYaxis()->SetTitle("Counts");\r
188       fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h4)),i-1+17+fGenOffset, kTRUE);\r
189           delete h4;\r
190       fSDDhTask++;\r
191     }\r
192   }\r
193 \r
194   Int_t indexlast1 = 0;\r
195   Int_t indexlast2 = 0;\r
196 \r
197   if(fkOnline) {\r
198         fTimeBinSize = 4;\r
199     indexlast = 0;\r
200     index1 = 0;\r
201         indexlast1 = fSDDhTask;\r
202     indexlast2 = 0;\r
203     char *hname[3];\r
204     for(Int_t i=0; i<3; i++) hname[i]= new char[50];\r
205     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
206                 for(Int_t iside=0;iside<fgknSide;iside++){\r
207                         AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
208                         sprintf(hname[0],"chargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside);\r
209                         sprintf(hname[1],"ChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside);\r
210                         sprintf(hname[2],"hmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
211                         TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);\r
212                         fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin");\r
213                         fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode");\r
214                         fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMapFSE)),indexlast1 + index1 + fGenOffset);\r
215                         delete fModuleChargeMapFSE;\r
216                         \r
217                         fSDDhTask++;\r
218                         index1++;        \r
219                         indexlast2 = indexlast1 + index1;\r
220                 }\r
221         }\r
222 \r
223     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
224                 for(Int_t iside=0;iside<fgknSide;iside++){\r
225                         AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
226                         sprintf(hname[0],"chargeMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
227                         sprintf(hname[1],"ChargeMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
228                         TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);\r
229                         fModuleChargeMap->GetXaxis()->SetTitle("Time Bin");\r
230                         fModuleChargeMap->GetYaxis()->SetTitle("Anode");\r
231                         fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMap)),indexlast1 + index1 + fGenOffset);\r
232                         delete fModuleChargeMap;\r
233 \r
234                         fSDDhTask++;\r
235                         index1++;        \r
236                         indexlast2 = indexlast1 + index1;\r
237                 }\r
238         }\r
239   \r
240 }  // kONLINE\r
241 \r
242   AliDebug(1,Form("%d SDD Raws histograms booked\n",fSDDhTask));\r
243 }\r
244 \r
245 \r
246 //____________________________________________________________________________\r
247 void AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader)\r
248\r
249   // Fill QA for RAW - SDD -\r
250   if(!fDDLModuleMap){\r
251     AliError("SDD DDL module map not available - skipping SDD QA");\r
252     return;\r
253   }\r
254   if(rawReader->GetType() != 7) return;  // skips non physical triggers\r
255   AliDebug(1,"entering MakeRaws\n");                 \r
256   rawReader->SelectEquipment(17,fgkeqOffset,fgkeqOffset + fgkDDLidRange); \r
257 \r
258   rawReader->Reset();                         \r
259   AliITSRawStreamSDD s(rawReader); \r
260   s.SetDDLModuleMap(fDDLModuleMap);\r
261   Int_t lay, lad, det; \r
262 \r
263   Int_t index=0;\r
264   if(fkOnline) {\r
265     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
266         for(Int_t iside=0;iside<fgknSide;iside++) {\r
267           if(fSDDhTask > 39 + index) fAliITSQADataMakerRec->GetRawsData(39 + index +fGenOffset)->Reset();\r
268           index++;\r
269         }\r
270     }\r
271   }\r
272   Int_t cnt = 0;\r
273   Int_t ildcID = -1;\r
274   Int_t iddl = -1;\r
275   Int_t isddmod = -1;\r
276   Int_t coord1, coord2, signal, moduleSDD, ioffset, iorder, activeModule, index1;\r
277   while(s.Next()) {\r
278     ildcID = rawReader->GetLDCId();\r
279     iddl = rawReader->GetDDLID() - fgkDDLIDshift;\r
280     isddmod = s.GetModuleNumber(iddl,s.GetCarlosId());\r
281     if(isddmod==-1){\r
282       AliDebug(1,Form("Found module with iddl: %d, s.GetCarlosId: %d \n",iddl,s.GetCarlosId() ));\r
283       continue;\r
284     }\r
285     if(s.IsCompletedModule()) {\r
286       AliDebug(1,Form("IsCompletedModule == KTRUE\n"));\r
287       continue;\r
288     } \r
289     \r
290     coord1 = s.GetCoord1();\r
291     coord2 = s.GetCoord2();\r
292     signal = s.GetSignal();\r
293     \r
294     moduleSDD = isddmod - fgkmodoffset;\r
295     if(moduleSDD < 0 || moduleSDD>fgknSDDmodules) {\r
296       AliDebug(1,Form( "Module SDD = %d, resetting it to 1 \n",moduleSDD));\r
297       moduleSDD = 1;\r
298     }\r
299     fAliITSQADataMakerRec->GetRawsData(0 +fGenOffset)->Fill(moduleSDD); \r
300     \r
301     AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);\r
302     ioffset = 3;\r
303     iorder = 1;\r
304     if(lay==4) { \r
305       ioffset += 14;\r
306       iorder = 2;   \r
307     } \r
308     fAliITSQADataMakerRec->GetRawsData(iorder +fGenOffset)->Fill(lad);\r
309     fAliITSQADataMakerRec->GetRawsData(ioffset+lad-1 +fGenOffset)->Fill(det); //-1 because ladder# starts from 1    \r
310     \r
311     Short_t iside = s.GetChannel();\r
312     activeModule = moduleSDD;\r
313     index1 = activeModule * 2 + iside;\r
314     \r
315     if(index1<0){\r
316       AliDebug(1,Form("Wrong index number %d - patched to 0\n",index1));\r
317       index1 = 0;\r
318     }\r
319     \r
320     if(fkOnline) {\r
321       if(fSDDhTask > 39 + index1) {\r
322         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(39 + index1 +fGenOffset)))->Fill(coord2, coord1, signal);\r
323         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(39 + index1 + 260*2 +fGenOffset)))->Fill(coord2, coord1, signal);\r
324       }\r
325     }\r
326     cnt++;\r
327     if(!(cnt%10000)) AliDebug(1,Form(" %d raw digits read",cnt));\r
328   }\r
329   AliDebug(1,Form("Event completed, %d raw digits read",cnt)); \r
330  }\r
331 \r
332 //____________________________________________________________________________ \r
333 void AliITSQASDDDataMakerRec::InitRecPoints()\r
334 {\r
335   // Initialization for RECPOINTS - SDD -\r
336   fGenOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();\r
337   \r
338   TH1F *h0 = new TH1F("Lay3TotCh","Layer 3 total charge",1000,-0.5, 499.5);\r
339   h0->GetXaxis()->SetTitle("ADC value");\r
340   h0->GetYaxis()->SetTitle("Entries");\r
341   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h0)), 0 +fGenOffset);\r
342   delete h0;\r
343   fSDDhTask++;\r
344  \r
345   TH1F *h1 = new TH1F("Lay4TotCh","Layer 4 total charge",1000,-0.5, 499.5);\r
346   h1->GetXaxis()->SetTitle("ADC value");\r
347   h1->GetYaxis()->SetTitle("Entries");\r
348   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h1)), 1 +fGenOffset);\r
349   delete h1;\r
350   fSDDhTask++;\r
351 \r
352     \r
353   char hisnam[50];\r
354   for(Int_t i=1; i<=3; i++){\r
355     sprintf(hisnam,"Charge_L3_Strip%d",i);\r
356     TH1F *h2 = new TH1F(hisnam,hisnam,1000,-0.5, 499.5);\r
357     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h2)),i+1 +fGenOffset);\r
358         delete h2;\r
359     fSDDhTask++;\r
360   }\r
361   \r
362   for(Int_t i=1; i<=4; i++){\r
363     sprintf(hisnam,"Charge_L4_Strip%d",i);\r
364     TH1F *h3 = new TH1F(hisnam,hisnam,1000,-0.5, 499.5);\r
365     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h3)),i+4 +fGenOffset);\r
366         delete h3;\r
367     fSDDhTask++;\r
368   }\r
369   \r
370   TH1F *h4 = new TH1F("ModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); \r
371   h4->GetXaxis()->SetTitle("Module number");\r
372   h4->GetYaxis()->SetTitle("Entries");\r
373   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h4)),9 +fGenOffset);\r
374   delete h4;\r
375   fSDDhTask++;\r
376   TH1F *h5 = new TH1F("ModPatternL3 RP","Ladder pattern L3 RP",14,0.5,14.5);  \r
377   h5->GetXaxis()->SetTitle("Ladder #, Layer 3");\r
378   h5->GetYaxis()->SetTitle("Entries");\r
379   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h5)),10 +fGenOffset);\r
380   delete h5;\r
381   fSDDhTask++;\r
382   TH1F *h6 = new TH1F("ModPatternL4 RP","Ladder pattern L4 RP",22,0.5,22.5); \r
383   h6->GetXaxis()->SetTitle("Ladder #, Layer 4");\r
384   h6->GetYaxis()->SetTitle("Entries");\r
385   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h6)),11 +fGenOffset);\r
386   delete h6;\r
387   fSDDhTask++;\r
388   TH2F *h7 = new TH2F("Local Coord Distrib","Local Coord Distrib",1000,-4,4,1000,-4,4);\r
389   h7->GetXaxis()->SetTitle("X local coord, drift, cm");\r
390   h7->GetYaxis()->SetTitle("Z local coord, anode, cm");\r
391   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h7)),12 +fGenOffset);\r
392   delete h7;\r
393   fSDDhTask++;\r
394   TH2F *h8 = new TH2F("Global Coord Distrib","Global Coord Distrib",6000,-30,30,6000,-30,30);\r
395   h8->GetYaxis()->SetTitle("Y glob coord, cm");\r
396   h8->GetXaxis()->SetTitle("X glob coord, cm");\r
397   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h8)),13 +fGenOffset);\r
398   delete h8;\r
399   fSDDhTask++;\r
400    \r
401   for(Int_t iLay=0; iLay<=1; iLay++){\r
402     sprintf(hisnam,"hr_Layer%d",iLay+3);\r
403     TH1F *h9 = new TH1F(hisnam,hisnam,100,10.,30.);\r
404     h9->GetXaxis()->SetTitle("r (cm)");\r
405     h9->GetXaxis()->CenterTitle();\r
406     h9->GetYaxis()->SetTitle("Entries");\r
407     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h9)),iLay+14 +fGenOffset);\r
408         delete h9;\r
409     fSDDhTask++;\r
410   }\r
411 \r
412   for(Int_t iLay=0; iLay<=1; iLay++){\r
413     sprintf(hisnam,"hphi_Layer%d",iLay+3);\r
414     TH1F *h10 = new TH1F(hisnam,hisnam,100,-TMath::Pi(),TMath::Pi());\r
415     h10->GetXaxis()->SetTitle("#varphi (rad)");\r
416     h10->GetXaxis()->CenterTitle();\r
417     h10->GetYaxis()->SetTitle("Entries");\r
418     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h10)),iLay+16 +fGenOffset);\r
419         delete h10;\r
420     fSDDhTask++;\r
421   }\r
422 \r
423   AliDebug(1,Form("%d SDD Recs histograms booked\n",fSDDhTask));\r
424 }\r
425 \r
426 //____________________________________________________________________________ \r
427 void AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)\r
428 {\r
429   // Fill QA for RecPoints - SDD -\r
430   Int_t lay, lad, det; \r
431   TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");\r
432   if (!branchRecP) { \r
433     AliError("can't get the branch with the ITS clusters !");\r
434     return;\r
435   }\r
436   static TClonesArray statRecpoints("AliITSRecPoint") ;\r
437   TClonesArray *recpoints = &statRecpoints;\r
438   branchRecP->SetAddress(&recpoints);\r
439   Int_t npoints = 0;      \r
440   Float_t cluglo[3]={0.,0.,0.}; \r
441   for(Int_t module=0; module<clustersTree->GetEntries();module++){\r
442     branchRecP->GetEvent(module);\r
443     npoints += recpoints->GetEntries();\r
444     AliITSgeomTGeo::GetModuleId(module, lay, lad, det);\r
445     //printf("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det);\r
446     \r
447     for(Int_t j=0;j<recpoints->GetEntries();j++){\r
448       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);    \r
449       fAliITSQADataMakerRec->GetRecPointsData(9 +fGenOffset)->Fill(module);\r
450       recp->GetGlobalXYZ(cluglo);\r
451       Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); \r
452       Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);\r
453       if(recp->GetLayer() ==2) {\r
454         fAliITSQADataMakerRec->GetRecPointsData(0 +fGenOffset)->Fill(recp->GetQ()) ;\r
455         fAliITSQADataMakerRec->GetRecPointsData(10 +fGenOffset)->Fill(lad);\r
456         fAliITSQADataMakerRec->GetRecPointsData(14 +fGenOffset)->Fill(rad);\r
457         fAliITSQADataMakerRec->GetRecPointsData(16 +fGenOffset)->Fill(phi);\r
458         fAliITSQADataMakerRec->GetRecPointsData(9 +fGenOffset)->Fill(module);\r
459         fAliITSQADataMakerRec->GetRecPointsData(12 +fGenOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());\r
460         fAliITSQADataMakerRec->GetRecPointsData(13 +fGenOffset)->Fill(cluglo[0],cluglo[1]);\r
461       }\r
462       else if(recp->GetLayer() ==3) {\r
463         fAliITSQADataMakerRec->GetRecPointsData(1 +fGenOffset)->Fill(recp->GetQ()) ;\r
464         fAliITSQADataMakerRec->GetRecPointsData(11 +fGenOffset)->Fill(lad);\r
465         fAliITSQADataMakerRec->GetRecPointsData(15 +fGenOffset)->Fill(rad);\r
466         fAliITSQADataMakerRec->GetRecPointsData(17 +fGenOffset)->Fill(phi);\r
467         fAliITSQADataMakerRec->GetRecPointsData(9 +fGenOffset)->Fill(module);\r
468         fAliITSQADataMakerRec->GetRecPointsData(12 +fGenOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());\r
469         fAliITSQADataMakerRec->GetRecPointsData(13 +fGenOffset)->Fill(cluglo[0],cluglo[1]);\r
470       }\r
471     }\r
472   }\r
473   statRecpoints.Clear();\r
474 }\r
475 \r