]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerRec.cxx
Switch off the QA
[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 fSDDhRaws(0),\r
59 fSDDhRecs(0),\r
60 fRawsOffset(0),\r
61 fRecsOffset(0),\r
62 fDDLModuleMap(0)\r
63 {\r
64   //ctor used to discriminate OnLine-Offline analysis\r
65   if(fLDC < 0 || fLDC > 4) {\r
66         AliError("Error: LDC number out of range; return\n");\r
67   }\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 fSDDhRaws(qadm.fSDDhRaws),\r
77 fSDDhRecs(qadm.fSDDhRecs),\r
78 fRawsOffset(qadm.fRawsOffset),\r
79 fRecsOffset(qadm.fRecsOffset),\r
80 fDDLModuleMap(0)\r
81 {\r
82   //copy ctor \r
83   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; \r
84   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());\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   fRawsOffset = (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+fRawsOffset);\r
153   delete h0;\r
154   fSDDhRaws++;\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+fRawsOffset);\r
160         delete h1;\r
161     fSDDhRaws++;\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+fRawsOffset);\r
168         delete h2;\r
169     fSDDhRaws++;\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+fRawsOffset);\r
178           delete h3;\r
179       fSDDhRaws++;\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+fRawsOffset);\r
189           delete h4;\r
190       fSDDhRaws++;\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     indexlast = 0;\r
199     index1 = 0;\r
200         indexlast1 = fSDDhRaws;\r
201     indexlast2 = 0;\r
202     char *hname[3];\r
203     for(Int_t i=0; i<3; i++) hname[i]= new char[50];\r
204     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
205                 for(Int_t iside=0;iside<fgknSide;iside++){\r
206                         AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
207                         sprintf(hname[0],"chargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside);\r
208                         sprintf(hname[1],"ChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside);\r
209                         sprintf(hname[2],"hmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
210                         TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],128,-0.5,255.5,256,-0.5,255.5);\r
211                         fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin");\r
212                         fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode");\r
213                         fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMapFSE)),indexlast1 + index1 + fRawsOffset);\r
214                         delete fModuleChargeMapFSE;\r
215                         \r
216                         fSDDhRaws++;\r
217                         index1++;        \r
218                         indexlast2 = indexlast1 + index1;\r
219                 }\r
220         }\r
221 \r
222     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
223                 for(Int_t iside=0;iside<fgknSide;iside++){\r
224                         AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
225                         sprintf(hname[0],"chargeMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
226                         sprintf(hname[1],"ChargeMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
227                         TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],128,-0.5,255.5,256,-0.5,255.5);\r
228                         fModuleChargeMap->GetXaxis()->SetTitle("Time Bin");\r
229                         fModuleChargeMap->GetYaxis()->SetTitle("Anode");\r
230                         fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMap)),indexlast1 + index1 + fRawsOffset);\r
231                         delete fModuleChargeMap;\r
232 \r
233                         fSDDhRaws++;\r
234                         index1++;        \r
235                         indexlast2 = indexlast1 + index1;\r
236                 }\r
237         }\r
238   \r
239 }  // kONLINE\r
240 \r
241 \r
242   AliDebug(1,Form("%d SDD Raws histograms booked\n",fSDDhRaws));\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   /*\r
259   if(rawReader->GetEventId()!=fEvtId){\r
260     TFile *DAoutput = new TFile::Open(filename);\r
261     TH1F *BLhisto;\r
262     for(Int_t imod=0; imod<nSDDmodules; imod++){ \r
263       BLhisto = (TH1F*)DAoutput->Get("BLhistoname[imod]");\r
264       mean = BLhisto->take mean;\r
265       fAliITSQADataMakerRec->GetRawsData(i+1887)->Fill(mean, imod);\r
266       Noisehisto....;\r
267       Vdrifthisto...;\r
268       <Q>histo....;\r
269     }\r
270   }\r
271   fEvtId==rawReader->GetEventId();\r
272   */\r
273 \r
274   rawReader->Reset();                         \r
275   AliITSRawStreamSDD s(rawReader); \r
276   s.SetDDLModuleMap(fDDLModuleMap);\r
277   Int_t lay, lad, det; \r
278 \r
279   Int_t index=0;\r
280   if(fkOnline) {\r
281     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
282         for(Int_t iside=0;iside<fgknSide;iside++) {\r
283           if(fSDDhRaws > 39 + index) fAliITSQADataMakerRec->GetRawsData(39 + index +fRawsOffset)->Reset();\r
284           index++;\r
285         }\r
286     }\r
287   }\r
288   Int_t cnt = 0;\r
289   Int_t ildcID = -1;\r
290   Int_t iddl = -1;\r
291   Int_t isddmod = -1;\r
292   Int_t coord1, coord2, signal, moduleSDD, ioffset, iorder, activeModule, index1;\r
293   while(s.Next()) {\r
294     ildcID = rawReader->GetLDCId();\r
295     iddl = rawReader->GetDDLID() - fgkDDLIDshift;\r
296     isddmod = s.GetModuleNumber(iddl,s.GetCarlosId());\r
297     if(isddmod==-1){\r
298       AliDebug(1,Form("Found module with iddl: %d, s.GetCarlosId: %d \n",iddl,s.GetCarlosId() ));\r
299       continue;\r
300     }\r
301     if(s.IsCompletedModule()) {\r
302       AliDebug(1,Form("IsCompletedModule == KTRUE\n"));\r
303       continue;\r
304     } \r
305     \r
306     coord1 = s.GetCoord1();\r
307     coord2 = s.GetCoord2();\r
308     signal = s.GetSignal();\r
309     \r
310     moduleSDD = isddmod - fgkmodoffset;\r
311     if(moduleSDD < 0 || moduleSDD>fgknSDDmodules) {\r
312       AliDebug(1,Form( "Module SDD = %d, resetting it to 1 \n",moduleSDD));\r
313       moduleSDD = 1;\r
314     }\r
315     fAliITSQADataMakerRec->GetRawsData(0 +fRawsOffset)->Fill(moduleSDD); \r
316     \r
317     AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);\r
318     ioffset = 3;\r
319     iorder = 1;\r
320     if(lay==4) { \r
321       ioffset += 14;\r
322       iorder = 2;   \r
323     } \r
324     fAliITSQADataMakerRec->GetRawsData(iorder +fRawsOffset)->Fill(lad);\r
325     fAliITSQADataMakerRec->GetRawsData(ioffset+lad-1 +fRawsOffset)->Fill(det); //-1 because ladder# starts from 1    \r
326     \r
327     Short_t iside = s.GetChannel();\r
328     activeModule = moduleSDD;\r
329     index1 = activeModule * 2 + iside;\r
330     \r
331     if(index1<0){\r
332       AliDebug(1,Form("Wrong index number %d - patched to 0\n",index1));\r
333       index1 = 0;\r
334     }\r
335     \r
336     if(fkOnline) {\r
337       if(fSDDhRaws > 39 + index1) {\r
338         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(39 + index1 +fRawsOffset)))->Fill(coord2, coord1, signal);\r
339         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(39 + index1 + 260*2 +fRawsOffset)))->Fill(coord2, coord1, signal);\r
340       }\r
341     }\r
342     cnt++;\r
343     if(!(cnt%10000)) AliDebug(1,Form(" %d raw digits read",cnt));\r
344   }\r
345   AliDebug(1,Form("Event completed, %d raw digits read",cnt));  \r
346 \r
347 }\r
348 \r
349 //____________________________________________________________________________ \r
350 void AliITSQASDDDataMakerRec::InitRecPoints()\r
351 {\r
352   // Initialization for RECPOINTS - SDD -\r
353   fRecsOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();\r
354   \r
355   TH1F *h0 = new TH1F("Lay3TotCh","Layer 3 total charge",1000,-0.5, 499.5);\r
356   h0->GetXaxis()->SetTitle("ADC value");\r
357   h0->GetYaxis()->SetTitle("Entries");\r
358   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h0)), 0 +fRecsOffset);\r
359   delete h0;\r
360   fSDDhRecs++;\r
361  \r
362   TH1F *h1 = new TH1F("Lay4TotCh","Layer 4 total charge",1000,-0.5, 499.5);\r
363   h1->GetXaxis()->SetTitle("ADC value");\r
364   h1->GetYaxis()->SetTitle("Entries");\r
365   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h1)), 1 +fRecsOffset);\r
366   delete h1;\r
367   fSDDhRecs++;\r
368 \r
369     \r
370   char hisnam[50];\r
371   for(Int_t i=1; i<=3; i++){\r
372     sprintf(hisnam,"Charge_L3_Strip%d",i);\r
373     TH1F *h2 = new TH1F(hisnam,hisnam,1000,-0.5, 499.5);\r
374     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h2)),i+1 +fRecsOffset);\r
375         delete h2;\r
376     fSDDhRecs++;\r
377   }\r
378   \r
379   for(Int_t i=1; i<=4; i++){\r
380     sprintf(hisnam,"Charge_L4_Strip%d",i);\r
381     TH1F *h3 = new TH1F(hisnam,hisnam,1000,-0.5, 499.5);\r
382     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h3)),i+4 +fRecsOffset);\r
383         delete h3;\r
384     fSDDhRecs++;\r
385   }\r
386   \r
387   TH1F *h4 = new TH1F("ModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); \r
388   h4->GetXaxis()->SetTitle("Module number");\r
389   h4->GetYaxis()->SetTitle("Entries");\r
390   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h4)),9 +fRecsOffset);\r
391   delete h4;\r
392   fSDDhRecs++;\r
393   TH1F *h5 = new TH1F("ModPatternL3 RP","Ladder pattern L3 RP",14,0.5,14.5);  \r
394   h5->GetXaxis()->SetTitle("Ladder #, Layer 3");\r
395   h5->GetYaxis()->SetTitle("Entries");\r
396   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h5)),10 +fRecsOffset);\r
397   delete h5;\r
398   fSDDhRecs++;\r
399   TH1F *h6 = new TH1F("ModPatternL4 RP","Ladder pattern L4 RP",22,0.5,22.5); \r
400   h6->GetXaxis()->SetTitle("Ladder #, Layer 4");\r
401   h6->GetYaxis()->SetTitle("Entries");\r
402   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h6)),11 +fRecsOffset);\r
403   delete h6;\r
404   fSDDhRecs++;\r
405   TH2F *h7 = new TH2F("Local Coord Distrib","Local Coord Distrib",1000,-4,4,1000,-4,4);\r
406   h7->GetXaxis()->SetTitle("X local coord, drift, cm");\r
407   h7->GetYaxis()->SetTitle("Z local coord, anode, cm");\r
408   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h7)),12 +fRecsOffset);\r
409   delete h7;\r
410   fSDDhRecs++;\r
411   TH2F *h8 = new TH2F("Global Coord Distrib","Global Coord Distrib",6000,-30,30,6000,-30,30);\r
412   h8->GetYaxis()->SetTitle("Y glob coord, cm");\r
413   h8->GetXaxis()->SetTitle("X glob coord, cm");\r
414   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h8)),13 +fRecsOffset);\r
415   delete h8;\r
416   fSDDhRecs++;\r
417    \r
418   for(Int_t iLay=0; iLay<=1; iLay++){\r
419     sprintf(hisnam,"hr_Layer%d",iLay+3);\r
420     TH1F *h9 = new TH1F(hisnam,hisnam,100,10.,30.);\r
421     h9->GetXaxis()->SetTitle("r (cm)");\r
422     h9->GetXaxis()->CenterTitle();\r
423     h9->GetYaxis()->SetTitle("Entries");\r
424     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h9)),iLay+14 +fRecsOffset);\r
425         delete h9;\r
426     fSDDhRecs++;\r
427   }\r
428 \r
429   for(Int_t iLay=0; iLay<=1; iLay++){\r
430     sprintf(hisnam,"hphi_Layer%d",iLay+3);\r
431     TH1F *h10 = new TH1F(hisnam,hisnam,100,-TMath::Pi(),TMath::Pi());\r
432     h10->GetXaxis()->SetTitle("#varphi (rad)");\r
433     h10->GetXaxis()->CenterTitle();\r
434     h10->GetYaxis()->SetTitle("Entries");\r
435     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h10)),iLay+16 +fRecsOffset);\r
436         delete h10;\r
437     fSDDhRecs++;\r
438   }\r
439 \r
440   AliDebug(1,Form("%d SDD Recs histograms booked\n",fSDDhRecs));\r
441 }\r
442 \r
443 //____________________________________________________________________________ \r
444 void AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)\r
445 {\r
446   // Fill QA for RecPoints - SDD -\r
447   Int_t lay, lad, det; \r
448   TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");\r
449   if (!branchRecP) { \r
450     AliError("can't get the branch with the ITS clusters !");\r
451     return;\r
452   }\r
453   TClonesArray * recpoints = new TClonesArray("AliITSRecPoint") ;\r
454   branchRecP->SetAddress(&recpoints);\r
455   Int_t npoints = 0;      \r
456   Float_t cluglo[3]={0.,0.,0.}; \r
457   for(Int_t module=0; module<clustersTree->GetEntries();module++){\r
458     branchRecP->GetEvent(module);\r
459     npoints += recpoints->GetEntries();\r
460     AliITSgeomTGeo::GetModuleId(module, lay, lad, det);\r
461     //printf("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det);\r
462     \r
463     for(Int_t j=0;j<recpoints->GetEntries();j++){\r
464       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);    \r
465       fAliITSQADataMakerRec->GetRecPointsData(9 +fRecsOffset)->Fill(module);\r
466       recp->GetGlobalXYZ(cluglo);\r
467       Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); \r
468       Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);\r
469       if(recp->GetLayer() ==2) {\r
470         fAliITSQADataMakerRec->GetRecPointsData(0 +fRecsOffset)->Fill(recp->GetQ()) ;\r
471         fAliITSQADataMakerRec->GetRecPointsData(10 +fRecsOffset)->Fill(lad);\r
472         fAliITSQADataMakerRec->GetRecPointsData(14 +fRecsOffset)->Fill(rad);\r
473         fAliITSQADataMakerRec->GetRecPointsData(16 +fRecsOffset)->Fill(phi);\r
474         fAliITSQADataMakerRec->GetRecPointsData(9 +fRecsOffset)->Fill(module);\r
475         fAliITSQADataMakerRec->GetRecPointsData(12 +fRecsOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());\r
476         fAliITSQADataMakerRec->GetRecPointsData(13 +fRecsOffset)->Fill(cluglo[0],cluglo[1]);\r
477       }\r
478       else if(recp->GetLayer() ==3) {\r
479         fAliITSQADataMakerRec->GetRecPointsData(1 +fRecsOffset)->Fill(recp->GetQ()) ;\r
480         fAliITSQADataMakerRec->GetRecPointsData(11 +fRecsOffset)->Fill(lad);\r
481         fAliITSQADataMakerRec->GetRecPointsData(15 +fRecsOffset)->Fill(rad);\r
482         fAliITSQADataMakerRec->GetRecPointsData(17 +fRecsOffset)->Fill(phi);\r
483         fAliITSQADataMakerRec->GetRecPointsData(9 +fRecsOffset)->Fill(module);\r
484         fAliITSQADataMakerRec->GetRecPointsData(12 +fRecsOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());\r
485         fAliITSQADataMakerRec->GetRecPointsData(13 +fRecsOffset)->Fill(cluglo[0],cluglo[1]);\r
486       }\r
487     }\r
488   }\r
489   recpoints->Delete();\r
490   delete recpoints;\r
491 \r
492 }\r
493 \r