]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerRec.cxx
Possibility to reconstruct tracks with 1 point + vertex, and possibility to reuse...
[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 <TH2D.h>\r
28 #include <TBranch.h>\r
29 #include <TTree.h>\r
30 #include <TGaxis.h>\r
31 #include <TMath.h>\r
32 // --- Standard library ---\r
33 \r
34 // --- AliRoot header files ---\r
35 #include "AliITSQASDDDataMakerRec.h"\r
36 #include "AliLog.h"\r
37 #include "AliQA.h"\r
38 #include "AliQAChecker.h"\r
39 #include "AliRawReader.h"\r
40 #include "AliITSRawStreamSDD.h"\r
41 #include "AliITSRecPoint.h"\r
42 #include "AliITSgeomTGeo.h"\r
43 \r
44 #include "AliCDBManager.h"\r
45 #include "AliCDBStorage.h"\r
46 #include "AliCDBEntry.h"\r
47 \r
48 \r
49 ClassImp(AliITSQASDDDataMakerRec)\r
50 \r
51 //____________________________________________________________________________ \r
52 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) :\r
53 TObject(),\r
54 fAliITSQADataMakerRec(aliITSQADataMakerRec),\r
55 fkOnline(kMode),\r
56 fLDC(ldc),\r
57 fSDDhRaws(0),\r
58 fSDDhRecs(0),\r
59 fRawsOffset(0),\r
60 fRecsOffset(0),\r
61 fSDDDDLModuleMap(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   for(Int_t i=0;i<2*fgknSDDmodules;i++){\r
68     fModuleChargeMap[i] = NULL;\r
69   }\r
70 }\r
71 \r
72 //____________________________________________________________________________ \r
73 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) :\r
74 TObject(),\r
75 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),\r
76 fkOnline(qadm.fkOnline),\r
77 fLDC(qadm.fLDC),\r
78 fSDDhRaws(qadm.fSDDhRaws),\r
79 fSDDhRecs(qadm.fSDDhRecs),\r
80 fRawsOffset(qadm.fRawsOffset),\r
81 fRecsOffset(qadm.fRecsOffset),\r
82 fSDDDDLModuleMap(0)\r
83 {\r
84   //copy ctor \r
85   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; \r
86   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());\r
87 }\r
88 \r
89 //____________________________________________________________________________ \r
90 AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){\r
91   // destructor\r
92 \r
93   for(Int_t i=0;i<2*fgknSDDmodules;i++){\r
94     if(fModuleChargeMap[i]) delete fModuleChargeMap[i];\r
95   }\r
96   \r
97 }\r
98 //__________________________________________________________________\r
99 AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac )\r
100 {\r
101   // Equal operator.\r
102   this->~AliITSQASDDDataMakerRec();\r
103   new(this) AliITSQASDDDataMakerRec(qac);\r
104   return *this;\r
105 }\r
106 \r
107 //____________________________________________________________________________ \r
108 void AliITSQASDDDataMakerRec::StartOfDetectorCycle()\r
109 {\r
110   //Detector specific actions at start of cycle\r
111   AliDebug(1,"AliITSQADM::Start of SDD Cycle\n");\r
112 }\r
113 \r
114 //____________________________________________________________________________ \r
115 void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t /*task*/, TObjArray* /*list*/)\r
116 {\r
117   // launch the QA checking\r
118   AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n"); \r
119   \r
120   //AliQAChecker::Instance()->Run( AliQA::kITS , task, list);\r
121 }\r
122 \r
123 //____________________________________________________________________________ \r
124 void AliITSQASDDDataMakerRec::InitRaws()\r
125\r
126   // Initialization for RAW data - SDD -\r
127   fRawsOffset = (fAliITSQADataMakerRec->fRawsQAList)->GetEntries();\r
128 \r
129   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");\r
130   if( !ddlMapSDD){\r
131     AliFatal("Calibration object retrieval failed! ");\r
132   }  \r
133   AliITSDDLModuleMapSDD *ddlsdd=(AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();\r
134   ddlMapSDD->SetOwner(kTRUE);\r
135   fSDDDDLModuleMap = ddlsdd;\r
136  \r
137   Int_t lay, lad, det;\r
138   Int_t LAY = -1;  //, LAD = -1;\r
139   char hname0[50];\r
140   Int_t indexlast = 0;\r
141   Int_t index1 = 0;\r
142 \r
143   if(fLDC == 1 || fLDC == 2) LAY = 2;\r
144   if(fLDC == 3 || fLDC == 4) LAY = 3;\r
145   \r
146   if(fkOnline) {\r
147     AliInfo("Book Online Histograms for SDD\n");\r
148   }\r
149   else {\r
150     AliInfo("Book Offline Histograms for SDD\n ");\r
151   }\r
152   TH1F *h0 = new TH1F("ModPattern","Modules pattern",fgknSDDmodules,-0.5,259.5);\r
153   h0->GetXaxis()->SetTitle("Module Number");\r
154   h0->GetYaxis()->SetTitle("Counts");\r
155   fAliITSQADataMakerRec->Add2RawsList(h0,0+fRawsOffset);\r
156   fSDDhRaws++;\r
157   if(fLDC==0 || fLDC==1 || fLDC==2){\r
158     TH1F *h1 = new TH1F("LadPatternL3","Ladder pattern L3",14,0.5,14.5);  \r
159     h1->GetXaxis()->SetTitle("Ladder Number on Lay3");\r
160     h1->GetYaxis()->SetTitle("Counts");\r
161     fAliITSQADataMakerRec->Add2RawsList(h1,1+fRawsOffset);\r
162     fSDDhRaws++;\r
163   }     \r
164   if(fLDC==0 || fLDC==3 || fLDC==4){\r
165     TH1F *h2 = new TH1F("LadPatternL4","Ladder pattern L4",22,0.5,22.5);  \r
166     h2->GetXaxis()->SetTitle("Ladder Number on Lay4");\r
167     h2->GetYaxis()->SetTitle("Counts");\r
168     fAliITSQADataMakerRec->Add2RawsList(h2,2+fRawsOffset);\r
169     fSDDhRaws++;\r
170   }\r
171   if(fLDC==0 || fLDC==1 || fLDC==2){\r
172     TH1D *h3[fgkLADDonLAY3] ; \r
173     for(Int_t i=1; i<=fgkLADDonLAY3; i++) {\r
174       sprintf(hname0,"ModPattern_L3_%d",i);\r
175       h3[i-1] = new TH1D(hname0,hname0,6,0.5,6.5);\r
176       h3[i-1]->GetXaxis()->SetTitle("Channel Number");\r
177       h3[i-1]->GetYaxis()->SetTitle("Counts");\r
178       fAliITSQADataMakerRec->Add2RawsList(h3[i-1],i-1+3+fRawsOffset);\r
179       fSDDhRaws++;\r
180     }\r
181   }\r
182   if(fLDC==0 || fLDC==3 || fLDC==4){\r
183     TH1D *h4[fgkLADDonLAY4] ; \r
184     for(Int_t i=1; i<=fgkLADDonLAY4; i++) {\r
185       sprintf(hname0,"ModPattern_L4_%d",i);\r
186       h4[i-1] = new TH1D(hname0,hname0,8,0.5,8.5);\r
187       h4[i-1]->GetXaxis()->SetTitle("Channel Number");\r
188       h4[i-1]->GetYaxis()->SetTitle("Counts");\r
189       fAliITSQADataMakerRec->Add2RawsList(h4[i-1],i-1+17+fRawsOffset);\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     TH1D *h5[fgknSDDmodules*2] ; \r
200     index1 = 0;\r
201 \r
202 /*\r
203     if(fLDC == 0 || fLDC == 1 || fLDC == 4){\r
204       \r
205         for(Int_t lad=0; lad<fgkLADDonLAY3; lad++){\r
206         for(Int_t ichan=0; det<3; det++){\r
207           for(Int_t iside=0; iside<fgknSide; iside++){\r
208             sprintf(hname0,"ProjYMap_L%d_%d_%d_%d",lay,lad,det+sideflag*3,iside);  //need to define "lay"\r
209             h5[index1] = new TH1D(hname0,hname0,256,-0.5,255.5);\r
210             h5[index1]->GetXaxis()->SetTitle("Anode Number");\r
211             h5[index1]->GetYaxis()->SetTitle("Counts");\r
212             fAliITSQADataMakerRec->Add2RawsList(h5[index1],index1+39+fRawsOffset);\r
213                 fSDDhRaws++;\r
214             index1++;\r
215             indexlast = index1+39;\r
216           }\r
217         }       \r
218       }\r
219     }\r
220 */  \r
221     for(Int_t moduleSDD=0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
222       if((moduleSDD >= 0 && moduleSDD < 36) || (moduleSDD >= 84 && moduleSDD < 180)) {\r
223         for(Int_t iside=0; iside<fgknSide; iside++){\r
224           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
225           //if(fLDC == 0 || ((fLDC == 1 || fLDC == 2) && lay == 2) || ((fLDC == 3 || fLDC == 4) && lay == 3)\r
226           sprintf(hname0,"ProjYMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
227           h5[index1] = new TH1D(hname0,hname0,256,-0.5,255.5);\r
228           h5[index1]->GetXaxis()->SetTitle("Anode Number");\r
229           h5[index1]->GetYaxis()->SetTitle("Counts");\r
230           fAliITSQADataMakerRec->Add2RawsList(h5[index1],index1+39+fRawsOffset);\r
231           fSDDhRaws++;\r
232           index1++;\r
233           indexlast = index1+39;\r
234         }\r
235       }\r
236     }\r
237     \r
238     TH1D *h6[fgknSDDmodules*8] ; \r
239     Int_t indextot = 0;\r
240     indexlast1 = 0;\r
241     for(Int_t htype=0; htype<4; htype++){\r
242       for(Int_t moduleSDD=0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
243         if((moduleSDD >= 0 && moduleSDD < 36) || (moduleSDD >= 84 && moduleSDD < 180)) {         \r
244           for(Int_t iside=0; iside<fgknSide; iside++){\r
245             AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
246             if(htype == 0) sprintf(hname0,"CountsVSAnode_L%d_%d_%d_%d",lay,lad,det,iside);\r
247             if(htype == 1) sprintf(hname0,"ChargeVSAnode_L%d_%d_%d_%d",lay,lad,det,iside);\r
248             if(htype == 2) sprintf(hname0,"CountsVSTbin_L%d_%d_%d_%d",lay,lad,det,iside);\r
249             if(htype == 3) sprintf(hname0,"ChargeVSTbin_L%d_%d_%d_%d",lay,lad,det,iside);\r
250             h6[indextot] = new TH1D(hname0,hname0,256,-0.5,255.5);\r
251             if(htype == 0){\r
252               h6[indextot]->GetXaxis()->SetTitle("Anode Number");\r
253               h6[indextot]->GetYaxis()->SetTitle("Counts");\r
254             }\r
255             if(htype == 1){\r
256               h6[indextot]->GetXaxis()->SetTitle("Anode Number");\r
257               h6[indextot]->GetYaxis()->SetTitle("Charge (ADC)");\r
258             }\r
259             if(htype == 2){\r
260               h6[indextot]->GetXaxis()->SetTitle("Time bin");\r
261               h6[indextot]->GetYaxis()->SetTitle("Counts");\r
262             }\r
263             if(htype == 3){\r
264               h6[indextot]->GetXaxis()->SetTitle("Time bin");\r
265               h6[indextot]->GetYaxis()->SetTitle("Charge (ADC)");\r
266             }\r
267             fAliITSQADataMakerRec->Add2RawsList(h6[indextot],indexlast + indextot + fRawsOffset);\r
268             fSDDhRaws++;\r
269             indextot++;\r
270             indexlast1 = indexlast + indextot;\r
271           }\r
272         }\r
273       }\r
274     }\r
275 \r
276 \r
277     indexlast2 = 0;\r
278     index1 = 0;\r
279     char *hname[3];\r
280     for(Int_t i=0; i<3; i++) hname[i]= new char[50];\r
281     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
282       if((moduleSDD >= 0 && moduleSDD < 36) || (moduleSDD >= 84 && moduleSDD < 180)) {\r
283         for(Int_t iside=0;iside<fgknSide;iside++){\r
284         AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
285         sprintf(hname[0],"chargeMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
286         sprintf(hname[1],"TotalCharge_L%d_%d_%d_%d",lay,lad,det,iside);\r
287         sprintf(hname[2],"hmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);\r
288         fModuleChargeMap[index1] = new TH2D(hname[0],hname[1],256,-0.5,255.5,256,-0.5,255.5);\r
289         fModuleChargeMap[index1]->GetXaxis()->SetTitle("Time Bin");\r
290         fModuleChargeMap[index1]->GetYaxis()->SetTitle("Anode");\r
291         fAliITSQADataMakerRec->Add2RawsList((new TH2D(*fModuleChargeMap[index1])),indexlast1 + index1 + fRawsOffset);\r
292         fSDDhRaws++;\r
293         index1++;        \r
294         indexlast2 = indexlast1 + index1;\r
295       }\r
296     }\r
297   }\r
298   \r
299   \r
300   TH2D *h7[fgknSDDmodules*2] ;\r
301   index1 = 0;\r
302   for(Int_t moduleSDD=0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
303     if((moduleSDD >= 0 && moduleSDD < 36) || (moduleSDD >= 84 && moduleSDD < 180)) {\r
304       for(Int_t iside=0; iside<fgknSide; iside++){\r
305         AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);\r
306         sprintf(hname0,"Anode_vs_Charge_L%d_%d_%d_%d",lay,lad,det,iside);\r
307         h7[index1] = new TH2D(hname0,hname0,fgknSDDmodules*2,-0.5,-0.5+fgknSDDmodules*2,256,0.5,256.5);\r
308         h7[index1]->GetXaxis()->SetTitle("Charge Value (ADC)");\r
309         h7[index1]->GetYaxis()->SetTitle("Anode");\r
310         fAliITSQADataMakerRec->Add2RawsList(h7[index1],indexlast2 + index1 + fRawsOffset);\r
311         fSDDhRaws++;\r
312         index1++;\r
313       }\r
314     }\r
315   }\r
316   /*\r
317   TH2F *h8[3];\r
318   for(Int_t htype =0; htype<3; htype++){\r
319     if(htype == 0) sprintf(hname0,"BLcompare");\r
320     if(htype == 1) sprintf(hname0,"NoiseCompare");\r
321     if(htype == 2) sprintf(hname0,"VdriftCompare");\r
322     h8[htype]= new TH2F(hname0,hname0,520,0,260,100,-0.5,99.5);\r
323     if(htype == 0){\r
324       h8[htype]->GetXaxis()->SetTitle("Module Number");\r
325       h8[htype]->GetYaxis()->SetTitle("Time #");\r
326     }\r
327     if(htype == 1){\r
328       h8[htype]->GetXaxis()->SetTitle("Module Number");\r
329       h8[htype]->GetYaxis()->SetTitle("Time #)");\r
330     }\r
331     if(htype == 2){\r
332       h8[htype]->GetXaxis()->SetTitle("Module Number");\r
333       h8[htype]->GetYaxis()->SetTitle("Time #");\r
334     }\r
335 \r
336     fAliITSQADataMakerRec->Add2RawsList(h8[htype],indexlast2 + index1 + htype+1 + fRawsOffset);\r
337     fSDDhRaws++;\r
338     printf("added histos %c, at position %d \n", hname0, indexlast2 + index1 + htype+1 + fRawsOffset);\r
339   }\r
340   */\r
341 }  // kONLINE\r
342 \r
343 \r
344   AliDebug(1,Form("%d SDD Raws histograms booked\n",fSDDhRaws));\r
345 }\r
346 \r
347 \r
348 //____________________________________________________________________________\r
349 void AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader)\r
350\r
351   // Fill QA for RAW - SDD -\r
352   if(rawReader->GetType() != 7) return;  // skips non physical triggers\r
353   AliDebug(1,"entering MakeRaws\n");                 \r
354   rawReader->SelectEquipment(17,fgkeqOffset,fgkeqOffset + fgkDDLidRange); \r
355 \r
356   /*\r
357   if(rawReader->GetEventId()!=fEvtId){\r
358     TFile *DAoutput = new TFile::Open(filename);\r
359     TH1F *BLhisto;\r
360     for(Int_t imod=0; imod<nSDDmodules; imod++){ \r
361       BLhisto = (TH1F*)DAoutput->Get("BLhistoname[imod]");\r
362       mean = BLhisto->take mean;\r
363       fAliITSQADataMakerRec->GetRawsData(i+1887)->Fill(mean, imod);\r
364       Noisehisto....;\r
365       Vdrifthisto...;\r
366       <Q>histo....;\r
367     }\r
368   }\r
369   fEvtId==rawReader->GetEventId();\r
370   */\r
371 \r
372   rawReader->Reset();                         \r
373   AliITSRawStreamSDD s(rawReader); \r
374   s.SetDDLModuleMap(fSDDDDLModuleMap);\r
375   Int_t lay, lad, det; \r
376 \r
377   Int_t index=0;\r
378   if(fkOnline) {\r
379     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
380       if((moduleSDD >= 0 && moduleSDD < 36) || (moduleSDD >= 84 && moduleSDD < 180)) {\r
381         for(Int_t iside=0;iside<fgknSide;iside++) {\r
382           if(fSDDhRaws > 39+12 * 132 + index) fAliITSQADataMakerRec->GetRawsData(39+12 * 132 + index +fRawsOffset)->Reset();\r
383           index++;\r
384         }\r
385       }\r
386     }\r
387   }\r
388   Int_t cnt = 0;\r
389   Int_t ildcID = -1;\r
390   Int_t iddl = -1;\r
391   Int_t isddmod = -1;\r
392   Int_t coord1, coord2, signal, moduleSDD, ioffset, iorder, activeModule, index1 ,nBins , iside, bin ;\r
393   while(s.Next()) {\r
394     ildcID = rawReader->GetLDCId();\r
395     iddl = rawReader->GetDDLID() - fgkDDLIDshift;\r
396     isddmod = s.GetModuleNumber(iddl,s.GetCarlosId());\r
397     if(isddmod==-1){\r
398       AliDebug(1,Form("Found module with iddl: %d, s.GetCarlosId: %d \n",iddl,s.GetCarlosId() ));\r
399       continue;\r
400     }\r
401     if(s.IsCompletedModule()) {\r
402       AliDebug(1,Form("IsCompletedModule == KTRUE\n"));\r
403       continue;\r
404     } \r
405     \r
406     coord1 = s.GetCoord1();\r
407     coord2 = s.GetCoord2();\r
408     signal = s.GetSignal();\r
409     \r
410     moduleSDD = isddmod - fgkmodoffset;\r
411     if(moduleSDD < 0 || moduleSDD>fgknSDDmodules) {\r
412       AliDebug(1,Form( "Module SDD = %d, resetting it to 1 \n",moduleSDD));\r
413       moduleSDD = 1;\r
414     }\r
415     fAliITSQADataMakerRec->GetRawsData(0 +fRawsOffset)->Fill(moduleSDD); \r
416     \r
417     AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);\r
418     ioffset = 3;\r
419     iorder = 1;\r
420     if(lay==4) { \r
421       ioffset += 14;\r
422       iorder = 2;   \r
423     } \r
424     fAliITSQADataMakerRec->GetRawsData(iorder +fRawsOffset)->Fill(lad);\r
425     fAliITSQADataMakerRec->GetRawsData(ioffset+lad-1 +fRawsOffset)->Fill(det); //-1 because ladder# starts from 1    \r
426     \r
427     Short_t iside = s.GetChannel();\r
428     activeModule = moduleSDD;\r
429     if(moduleSDD > 35) activeModule -= 48;\r
430     index1 = activeModule * 2 + iside;\r
431     \r
432     if(index1<0){\r
433       AliDebug(1,Form("Wrong index number %d - patched to 0\n",index1));\r
434       index1 = 0;\r
435     }\r
436     \r
437     if(fkOnline) {\r
438       if(fSDDhRaws > 39+12 * 132 + index1) {\r
439         \r
440         fAliITSQADataMakerRec->GetRawsData(39+ 2 * 132 + index1 +fRawsOffset)->Fill(coord1); \r
441         fAliITSQADataMakerRec->GetRawsData(39+ 4 * 132 + index1 +fRawsOffset)->Fill(coord1,signal); \r
442         fAliITSQADataMakerRec->GetRawsData(39+ 6 * 132 + index1 +fRawsOffset)->Fill(coord2);    \r
443         fAliITSQADataMakerRec->GetRawsData(39+ 8 * 132 + index1 +fRawsOffset)->Fill(coord2,signal); \r
444         fAliITSQADataMakerRec->GetRawsData(39+12 * 132 + index1 +fRawsOffset)->Fill(signal,coord1);\r
445         ((TH2D *)(fAliITSQADataMakerRec->GetRawsData(39+10 * 132 + index1 +fRawsOffset)))->Fill(coord2, coord1, signal);\r
446 \r
447       }\r
448     }\r
449     cnt++;\r
450     if(!(cnt%10000)) AliDebug(1,Form(" %d raw digits read",cnt));\r
451   }\r
452   AliDebug(1,Form("Event completed, %d raw digits read",cnt));  \r
453 \r
454   if(fkOnline) {\r
455     TH1D *ptr = NULL;\r
456     nBins = 256;\r
457     for(moduleSDD=0; moduleSDD<fgknSDDmodules; moduleSDD++){\r
458       if((moduleSDD >= 0 && moduleSDD < 36) || (moduleSDD >= 84 && moduleSDD < 180)) {\r
459         for(iside=0; iside<fgknSide; iside++){\r
460           activeModule = moduleSDD;\r
461           if(moduleSDD > 35) activeModule -= 48;\r
462           index1 = activeModule * 2 + iside;\r
463           if(fSDDhRaws > 39 + 2 * 132 + index1) {\r
464             ptr = ((TH2D *) (fAliITSQADataMakerRec->GetRawsData(39+10 * 132 + index1 +fRawsOffset)))->ProjectionY();\r
465             for(bin=0; bin<nBins; bin++) fAliITSQADataMakerRec->GetRawsData(index1+39 +fRawsOffset)->Fill(bin,ptr->GetBinContent(bin+1) );\r
466           }  \r
467         }\r
468       }\r
469     }\r
470     for(Int_t i=0; i<fSDDhRaws; i++){\r
471       Int_t entries = static_cast<Int_t>(fAliITSQADataMakerRec->GetRawsData(i +fRawsOffset)->GetEntries());\r
472       if(entries != 0)\r
473         AliDebug(1,Form("histo %d, name %s , entries %d ",i,fAliITSQADataMakerRec->GetRawsData(i +fRawsOffset)->GetName(),entries));\r
474     }\r
475   }\r
476 }\r
477 \r
478 //____________________________________________________________________________ \r
479 void AliITSQASDDDataMakerRec::InitRecPoints()\r
480 {\r
481   // Initialization for RECPOINTS - SDD -\r
482   fRecsOffset = (fAliITSQADataMakerRec->fRecPointsQAList)->GetEntries();\r
483   \r
484   TH1F * h0 = new TH1F("Lay3TotCh","Layer 3 total charge",1000,-0.5, 499.5);\r
485   h0->GetXaxis()->SetTitle("ADC value");\r
486   h0->GetYaxis()->SetTitle("Entries");\r
487   fAliITSQADataMakerRec->Add2RecPointsList(h0, 0 +fRecsOffset);\r
488   fSDDhRecs++;\r
489  \r
490   TH1F * h1 = new TH1F("Lay4TotCh","Layer 4 total charge",1000,-0.5, 499.5);\r
491   h1->GetXaxis()->SetTitle("ADC value");\r
492   h1->GetYaxis()->SetTitle("Entries");\r
493   fAliITSQADataMakerRec->Add2RecPointsList(h1, 1 +fRecsOffset);\r
494   fSDDhRecs++;\r
495 \r
496     \r
497   char hisnam[50];\r
498   TH1F *h2[3]; \r
499   for(Int_t i=1; i<=3; i++){\r
500     sprintf(hisnam,"Charge_L3_Strip%d",i);\r
501     h2[i] = new TH1F(hisnam,hisnam,1000,-0.5, 499.5);\r
502     fAliITSQADataMakerRec->Add2RecPointsList(h2[i],i+1 +fRecsOffset);\r
503     fSDDhRecs++;\r
504   }\r
505   \r
506   TH1F *h3[4]; \r
507   for(Int_t i=1; i<=4; i++){\r
508     sprintf(hisnam,"Charge_L4_Strip%d",i);\r
509     h3[i] = new TH1F(hisnam,hisnam,1000,-0.5, 499.5);\r
510     fAliITSQADataMakerRec->Add2RecPointsList(h3[i],i+4 +fRecsOffset);\r
511     fSDDhRecs++;\r
512   }\r
513   \r
514   TH1F *h4 = new TH1F("ModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); \r
515   h4->GetXaxis()->SetTitle("Module number");\r
516   h4->GetYaxis()->SetTitle("Entries");\r
517   fAliITSQADataMakerRec->Add2RecPointsList(h4,9 +fRecsOffset);\r
518   fSDDhRecs++;\r
519   TH1F *h5 = new TH1F("ModPatternL3 RP","Ladder pattern L3 RP",14,0.5,14.5);  \r
520   h5->GetXaxis()->SetTitle("Ladder #, Layer 3");\r
521   h5->GetYaxis()->SetTitle("Entries");\r
522   fAliITSQADataMakerRec->Add2RecPointsList(h5,10 +fRecsOffset);\r
523   fSDDhRecs++;\r
524   TH1F *h6 = new TH1F("ModPatternL4 RP","Ladder pattern L4 RP",22,0.5,22.5); \r
525   h6->GetXaxis()->SetTitle("Ladder #, Layer 4");\r
526   h6->GetYaxis()->SetTitle("Entries");\r
527   fAliITSQADataMakerRec->Add2RecPointsList(h6,11 +fRecsOffset);\r
528   fSDDhRecs++;\r
529   TH2F *h7 = new TH2F("Local Coord Distrib","Local Coord Distrib",1000,-4,4,1000,-4,4);\r
530   h7->GetXaxis()->SetTitle("X local coord, drift, cm");\r
531   h7->GetYaxis()->SetTitle("Z local coord, anode, cm");\r
532   fAliITSQADataMakerRec->Add2RecPointsList(h7,12 +fRecsOffset);\r
533   fSDDhRecs++;\r
534   TH2F *h8 = new TH2F("Global Coord Distrib","Global Coord Distrib",6000,-30,30,6000,-30,30);\r
535   h8->GetYaxis()->SetTitle("Y glob coord, cm");\r
536   h8->GetXaxis()->SetTitle("X glob coord, cm");\r
537   fAliITSQADataMakerRec->Add2RecPointsList(h8,13 +fRecsOffset);\r
538   fSDDhRecs++;\r
539   \r
540   TH1F *h9[2]; \r
541   for(Int_t iLay=0; iLay<=1; iLay++){\r
542     sprintf(hisnam,"hr_Layer%d",iLay+3);\r
543     h9[iLay] = new TH1F(hisnam,hisnam,100,10.,30.);\r
544     h9[iLay]->GetXaxis()->SetTitle("r (cm)");\r
545     h9[iLay]->GetXaxis()->CenterTitle();\r
546     h9[iLay]->GetYaxis()->SetTitle("Entries");\r
547     fAliITSQADataMakerRec->Add2RecPointsList(h9[iLay],iLay+14 +fRecsOffset);\r
548     fSDDhRecs++;\r
549   }\r
550   TH1F *h10[2]; \r
551   for(Int_t iLay=0; iLay<=1; iLay++){\r
552     sprintf(hisnam,"hphi_Layer%d",iLay+3);\r
553     h10[iLay] = new TH1F(hisnam,hisnam,100,-TMath::Pi(),TMath::Pi());\r
554     h10[iLay]->GetXaxis()->SetTitle("#varphi (rad)");\r
555     h10[iLay]->GetXaxis()->CenterTitle();\r
556     h10[iLay]->GetYaxis()->SetTitle("Entries");\r
557     fAliITSQADataMakerRec->Add2RecPointsList(h10[iLay],iLay+16 +fRecsOffset);\r
558     fSDDhRecs++;\r
559   }\r
560 \r
561   AliDebug(1,Form("%d SDD Recs histograms booked\n",fSDDhRecs));\r
562 }\r
563 \r
564 //____________________________________________________________________________ \r
565 void AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)\r
566 {\r
567   // Fill QA for RecPoints - SDD -\r
568   Int_t lay, lad, det; \r
569   TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");\r
570   if (!branchRecP) { \r
571     AliError("can't get the branch with the ITS clusters !");\r
572     return;\r
573   }\r
574   TClonesArray * recpoints = new TClonesArray("AliITSRecPoint") ;\r
575   branchRecP->SetAddress(&recpoints);\r
576   Int_t npoints = 0;      \r
577   Float_t cluglo[3]={0.,0.,0.}; \r
578   for(Int_t module=0; module<clustersTree->GetEntries();module++){\r
579     branchRecP->GetEvent(module);\r
580     npoints += recpoints->GetEntries();\r
581     AliITSgeomTGeo::GetModuleId(module, lay, lad, det);\r
582     //printf("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det);\r
583     \r
584     for(Int_t j=0;j<recpoints->GetEntries();j++){\r
585       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);    \r
586       fAliITSQADataMakerRec->GetRecPointsData(9 +fRecsOffset)->Fill(module);\r
587       recp->GetGlobalXYZ(cluglo);\r
588       Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); \r
589       Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);\r
590       if(recp->GetLayer() ==2) {\r
591         fAliITSQADataMakerRec->GetRecPointsData(0 +fRecsOffset)->Fill(recp->GetQ()) ;\r
592         fAliITSQADataMakerRec->GetRecPointsData(10 +fRecsOffset)->Fill(lad);\r
593         fAliITSQADataMakerRec->GetRecPointsData(14 +fRecsOffset)->Fill(rad);\r
594         fAliITSQADataMakerRec->GetRecPointsData(16 +fRecsOffset)->Fill(phi);\r
595         fAliITSQADataMakerRec->GetRecPointsData(9 +fRecsOffset)->Fill(module);\r
596         fAliITSQADataMakerRec->GetRecPointsData(12 +fRecsOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());\r
597         fAliITSQADataMakerRec->GetRecPointsData(13 +fRecsOffset)->Fill(cluglo[0],cluglo[1]);\r
598       }\r
599       else if(recp->GetLayer() ==3) {\r
600         fAliITSQADataMakerRec->GetRecPointsData(1 +fRecsOffset)->Fill(recp->GetQ()) ;\r
601         fAliITSQADataMakerRec->GetRecPointsData(11 +fRecsOffset)->Fill(lad);\r
602         fAliITSQADataMakerRec->GetRecPointsData(15 +fRecsOffset)->Fill(rad);\r
603         fAliITSQADataMakerRec->GetRecPointsData(17 +fRecsOffset)->Fill(phi);\r
604         fAliITSQADataMakerRec->GetRecPointsData(9 +fRecsOffset)->Fill(module);\r
605         fAliITSQADataMakerRec->GetRecPointsData(12 +fRecsOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());\r
606         fAliITSQADataMakerRec->GetRecPointsData(13 +fRecsOffset)->Fill(cluglo[0],cluglo[1]);\r
607       }\r
608     }\r
609   }\r
610   recpoints->Delete();\r
611   delete recpoints;\r
612 \r
613 }\r
614 \r