]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerSim.cxx
clusterizer,reconstructor + many fixes (Magnus,Stefan)
[u/mrichter/AliRoot.git] / ITS / AliITSQASDDDataMakerSim.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 <TTree.h>\r
28 // --- Standard library ---\r
29 \r
30 // --- AliRoot header files ---\r
31 #include "AliITSQASDDDataMakerSim.h"\r
32 #include "AliLog.h"\r
33 #include "AliQAv1.h"\r
34 #include "AliQAChecker.h"\r
35 #include "AliQADataMakerSim.h"\r
36 #include "AliITSQADataMakerSim.h"\r
37 #include "AliRawReader.h"\r
38 #include "AliITSdigit.h"\r
39 #include "AliITS.h"\r
40 #include "AliITSmodule.h"\r
41 #include "AliITShit.h"\r
42 #include "AliITSLoader.h"\r
43 #include "AliRunLoader.h"\r
44 #include "AliRun.h"\r
45 #include "AliITSsegmentationSDD.h"\r
46 #include "AliITSpList.h"\r
47 \r
48 ClassImp(AliITSQASDDDataMakerSim)\r
49 \r
50 //____________________________________________________________________________ \r
51 AliITSQASDDDataMakerSim::AliITSQASDDDataMakerSim(AliITSQADataMakerSim *aliITSQADataMakerSim) :\r
52 TObject(),\r
53 fAliITSQADataMakerSim(aliITSQADataMakerSim),\r
54 fSDDhHTask(0),\r
55 fSDDhSTask(0),\r
56 fSDDhDTask(0),\r
57 fGenOffsetH(0),\r
58 fGenOffsetS(0),\r
59 fGenOffsetD(0)\r
60 {\r
61   //ctor used to discriminate OnLine-Offline analysis   \r
62   fGenOffsetH=  new Int_t[AliRecoParam::kNSpecies];                       \r
63   fGenOffsetS=  new Int_t[AliRecoParam::kNSpecies];                           \r
64   fGenOffsetD=  new Int_t[AliRecoParam::kNSpecies];\r
65   for(Int_t i=0; i<AliRecoParam::kNSpecies; i++) \r
66     {\r
67       fGenOffsetH[i]= 0;\r
68       fGenOffsetS[i]= 0;\r
69       fGenOffsetD[i]= 0;\r
70     } \r
71 }\r
72 \r
73 //____________________________________________________________________________ \r
74 AliITSQASDDDataMakerSim::AliITSQASDDDataMakerSim(const AliITSQASDDDataMakerSim& qadm) :\r
75 TObject(),\r
76 fAliITSQADataMakerSim(qadm.fAliITSQADataMakerSim),\r
77 fSDDhHTask(qadm.fSDDhHTask),\r
78 fSDDhSTask(qadm.fSDDhSTask),\r
79 fSDDhDTask(qadm.fSDDhDTask),\r
80 fGenOffsetH(qadm.fGenOffsetH),\r
81 fGenOffsetS(qadm.fGenOffsetS),\r
82 fGenOffsetD(qadm.fGenOffsetD)\r
83 {\r
84   //copy ctor \r
85   fAliITSQADataMakerSim->SetName((const char*)qadm.fAliITSQADataMakerSim->GetName()) ; \r
86   fAliITSQADataMakerSim->SetTitle((const char*)qadm.fAliITSQADataMakerSim->GetTitle());\r
87 \r
88 }\r
89 \r
90 //__________________________________________________________________\r
91 AliITSQASDDDataMakerSim& AliITSQASDDDataMakerSim::operator = (const AliITSQASDDDataMakerSim& qac )\r
92 {\r
93   // Equal operator.\r
94   this->~AliITSQASDDDataMakerSim();\r
95   new(this) AliITSQASDDDataMakerSim(qac);\r
96   return *this;\r
97 }\r
98 \r
99 //____________________________________________________________________________ \r
100 void AliITSQASDDDataMakerSim::StartOfDetectorCycle()\r
101 {\r
102   //Detector specific actions at start of cycle\r
103   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SDD Cycle\n");\r
104 }\r
105 \r
106 //____________________________________________________________________________ \r
107 void AliITSQASDDDataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray** /*list*/)\r
108 {\r
109   // launch the QA checking\r
110   //\r
111   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); \r
112   //AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);\r
113 }\r
114 \r
115 //____________________________________________________________________________ \r
116 Int_t AliITSQASDDDataMakerSim::InitDigits()\r
117\r
118   // Initialization for DIGIT data - SDD -  \r
119   const Bool_t expert   = kTRUE ; \r
120   const Bool_t image    = kTRUE ;\r
121   Int_t rv = 0 ; \r
122   //fGenOffsetD = (fAliITSQADataMakerSim->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();\r
123   //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List\r
124   TH1F* h0=new TH1F("SDD DIGITS Module Pattern","SDD DIGITS Module Pattern",260,239.5,499.5);       //hmod\r
125   h0->GetXaxis()->SetTitle("SDD Module Number");\r
126   h0->GetYaxis()->SetTitle("# DIGITS");\r
127   rv = fAliITSQADataMakerSim->Add2DigitsList(h0,fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
128   fSDDhDTask ++;\r
129   TH1F* h1=new TH1F("SDD Anode Distribution","SDD DIGITS Anode Distribution",512,-0.5,511.5);      //hanocc\r
130   h1->GetXaxis()->SetTitle("Anode Number");\r
131   h1->GetYaxis()->SetTitle("# DIGITS");\r
132   rv = fAliITSQADataMakerSim->Add2DigitsList(h1,1+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
133   fSDDhDTask ++;\r
134   TH1F* h2=new TH1F("SDD Tbin Distribution","SDD DIGITS Tbin Distribution",256,-0.5,255.5);      //htbocc\r
135   h2->GetXaxis()->SetTitle("Tbin Number");\r
136   h2->GetYaxis()->SetTitle("# DIGITS");\r
137   rv = fAliITSQADataMakerSim->Add2DigitsList(h2,2+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
138   fSDDhDTask ++;\r
139   TH1F* h3=new TH1F("SDD ADC Counts Distribution","SDD DIGITS ADC Counts Distribution",200,0.,1024.);          //hsig\r
140   h3->GetXaxis()->SetTitle("ADC Value");\r
141   h3->GetYaxis()->SetTitle("# DIGITS");\r
142   rv = fAliITSQADataMakerSim->Add2DigitsList(h3,3+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
143   fSDDhDTask ++;\r
144   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Digits histograms booked\n",fSDDhDTask));\r
145   return rv ; \r
146 }\r
147 \r
148 //____________________________________________________________________________\r
149 Int_t AliITSQASDDDataMakerSim::MakeDigits(TTree * digits)\r
150\r
151 \r
152   // Fill QA for DIGIT - SDD -\r
153   Int_t rv = 0 ; \r
154 \r
155   AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");\r
156   fITS->SetTreeAddress();\r
157   TClonesArray *iITSdigits  = fITS->DigitsAddress(1);\r
158   for(Int_t i=0; i<260; i++){\r
159     Int_t nmod=i+240;\r
160     digits->GetEvent(nmod);\r
161     Int_t ndigits = iITSdigits->GetEntries();\r
162     fAliITSQADataMakerSim->FillDigitsData(fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()],nmod,ndigits);\r
163     for (Int_t idig=0; idig<ndigits; idig++) {\r
164       AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);\r
165       Int_t iz=dig->GetCoord1();  // cell number z\r
166       Int_t ix=dig->GetCoord2();  // cell number x\r
167       Int_t sig=dig->GetSignal();\r
168       fAliITSQADataMakerSim->FillDigitsData(1+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()],iz);\r
169       fAliITSQADataMakerSim->FillDigitsData(2+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()],ix);\r
170       fAliITSQADataMakerSim->FillDigitsData(3+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()],sig);\r
171     }\r
172   }\r
173   //\r
174   return rv ; \r
175 }\r
176 \r
177 //____________________________________________________________________________ \r
178 Int_t AliITSQASDDDataMakerSim::InitSDigits()\r
179\r
180   // Initialization for SDIGIT data - SDD -\r
181   const Bool_t expert   = kTRUE ; \r
182   const Bool_t image    = kTRUE ;\r
183   Int_t rv = 0 ; \r
184   //fGenOffsetS = (fAliITSQADataMakerSim->fSDigitsQAList[AliRecoParam::kDefault])->GetEntries();\r
185   //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List\r
186   TH1F* h0=new TH1F("SDD SDIGITS Module Pattern","SDIGITS SDD Module Pattern",260,239.5,499.5);       //hmod\r
187   h0->GetXaxis()->SetTitle("SDD Module Number");\r
188   h0->GetYaxis()->SetTitle("# SDIGITS");\r
189   rv = fAliITSQADataMakerSim->Add2SDigitsList(h0,fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
190   fSDDhSTask ++;\r
191   TH1F* h1=new TH1F("SDD Anode Distribution","SDIGITS SDD Anode Distribution",512,-0.5,511.5);      //hanocc\r
192   h1->GetXaxis()->SetTitle("Anode Number");\r
193   h1->GetYaxis()->SetTitle("# SDIGITS");\r
194   rv = fAliITSQADataMakerSim->Add2SDigitsList(h1,1+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
195   fSDDhSTask ++;\r
196   TH1F* h2=new TH1F("SDD Tbin Distribution","SDIGITS SDD Tbin Distribution",256,-0.5,255.5);      //htbocc\r
197   h2->GetXaxis()->SetTitle("Tbin Number");\r
198   h2->GetYaxis()->SetTitle("# SDIGITS");\r
199   rv = fAliITSQADataMakerSim->Add2SDigitsList(h2,2+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()]);\r
200   fSDDhSTask ++;\r
201   TH1F* h3=new TH1F("SDD ADC Counts Distribution","SDIGITS SDD ADC Counts Distribution",200,0.,1024.);          //hsig\r
202   h3->GetXaxis()->SetTitle("ADC Value");\r
203   h3->GetYaxis()->SetTitle("# SDIGITS");\r
204   rv = fAliITSQADataMakerSim->Add2SDigitsList(h3,3+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
205   fSDDhSTask ++;\r
206 \r
207   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD SDigits histograms booked\n",fSDDhSTask));\r
208   return rv ; \r
209 }\r
210 \r
211 //____________________________________________________________________________\r
212 Int_t AliITSQASDDDataMakerSim::MakeSDigits(TTree * sdigits)\r
213\r
214   // Fill QA for SDIGIT - SDD -\r
215   Int_t rv = 0 ; \r
216   \r
217    AliITSsegmentationSDD* seg = new AliITSsegmentationSDD();\r
218   Int_t nan=seg->Npz();\r
219   Int_t ntb=seg->Npx();\r
220   Int_t scaleSize=4;\r
221   AliITSpList* list=new AliITSpList(nan,ntb*scaleSize);\r
222 \r
223   //AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");\r
224   //fITS->SetTreeAddress();\r
225   //TClonesArray *ITSdigits  = fITS->DigitsAddress(1);\r
226   //TFile *sper = new TFile("sper.root","CREATE"); //agginto a mano x prova\r
227   //digits->Write();\r
228   //sper->Close();\r
229 \r
230 \r
231   TBranch *brchSDigits = sdigits->GetBranch("ITS");\r
232   static TClonesArray * sdig ; \r
233   if (! sdig )\r
234     sdig = new TClonesArray( "AliITSpListItem",1000 );\r
235   for(Int_t id=0; id<260; id++){\r
236     Int_t nmod=id+240;\r
237     brchSDigits->SetAddress( &sdig );\r
238     brchSDigits->GetEvent(nmod);\r
239     Int_t nsdig=sdig->GetEntries();\r
240     fAliITSQADataMakerSim->FillSDigitsData(fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()],nmod,nsdig);\r
241     for(Int_t i=0;i<nsdig;i++){\r
242       AliITSpListItem *cell=(AliITSpListItem*)sdig->At(i);\r
243       Float_t sig=cell->GetSignal();\r
244       Int_t idx=cell->GetIndex();\r
245       Int_t ia,it;\r
246       list->GetCell(idx,ia,it);\r
247       fAliITSQADataMakerSim->FillSDigitsData(1+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()],ia);\r
248       fAliITSQADataMakerSim->FillSDigitsData(2+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()],it);\r
249       fAliITSQADataMakerSim->FillSDigitsData(3+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()],sig);\r
250     }\r
251     sdig->Clear();\r
252   }\r
253   //\r
254   return rv ; \r
255 }\r
256 \r
257 //____________________________________________________________________________ \r
258 Int_t AliITSQASDDDataMakerSim::InitHits()\r
259\r
260 \r
261   // Initialization for HITS data - SDD -\r
262   const Bool_t expert   = kTRUE ; \r
263   const Bool_t image    = kTRUE ;\r
264   Int_t rv = 0 ; \r
265 \r
266   //fGenOffsetH = (fAliITSQADataMakerSim->fHitsQAList[AliRecoParam::kDefault])->GetEntries();\r
267   //fSDDhTask must be incremented by one unit every time a histogram is ADDED to the QA List\r
268   //printf("AliITSQASDDDataMakerSim::InitHits called \n");\r
269   TH1F *h0=new TH1F("SDD HITS Module Pattern","SDD HITS Module Pattern",260,239.5,499.5);  \r
270   h0->GetXaxis()->SetTitle("SDD Module Number");\r
271   h0->GetYaxis()->SetTitle("# HITS");\r
272   rv = fAliITSQADataMakerSim->Add2HitsList(h0,fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
273   fSDDhHTask ++;\r
274   TH1F *h1=new TH1F("SDD HIT lenght along local Y Coord","SDD HIT lenght along local Y Coord",200,0.,350.);\r
275   h1->GetXaxis()->SetTitle("HIT lenght (um)");\r
276   h1->GetYaxis()->SetTitle("# HITS");\r
277   rv = fAliITSQADataMakerSim->Add2HitsList(h1,1+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
278   fSDDhHTask ++;\r
279   TH1F *h2=new TH1F("SDD HIT lenght along local Y Coord - Zoom","SDD HIT lenght along local Y Coord - Zoom",200,250.,350.);\r
280   h2->GetXaxis()->SetTitle("HIT lenght (um)");\r
281   h2->GetYaxis()->SetTitle("# HITS");\r
282   rv = fAliITSQADataMakerSim->Add2HitsList(h2,2+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
283   fSDDhHTask ++;\r
284   TH1F *h3=new TH1F("SDD Deposited Energy Distribution (loc Y > 200um)","SDD HITS Deposited Energy Distribution (loc Y > 200um)",200,0.,350.);\r
285   h3->GetXaxis()->SetTitle("ADC counts ");\r
286   h3->GetYaxis()->SetTitle("# HITS");\r
287   rv = fAliITSQADataMakerSim->Add2HitsList(h3,3+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);\r
288   fSDDhHTask ++;\r
289   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SDD Hits histograms booked\n",fSDDhHTask));\r
290   return rv ; \r
291 }\r
292 \r
293 //____________________________________________________________________________\r
294 Int_t AliITSQASDDDataMakerSim::MakeHits(TTree * hits)\r
295\r
296   // Fill QA for HITS - SDD -\r
297   Int_t rv = 0 ; \r
298 \r
299   AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");\r
300   fITS->SetTreeAddress();\r
301   Int_t nmodules;\r
302   if(!(fITS->InitModules(-1,nmodules))){\r
303     AliError("ITS geometry not available - nothing done");\r
304     return rv;\r
305   }\r
306  \r
307   fITS->FillModules(hits,0);\r
308 \r
309   for(Int_t i=0; i<260; i++){\r
310     Int_t nmod=i+240;\r
311     AliITSmodule *modu = fITS->GetModule(nmod);\r
312     TObjArray *arrHits = modu->GetHits();\r
313     Int_t nhits = arrHits->GetEntriesFast();\r
314     ////printf("--w--AliITSQASDDDataMakerSim::MakeHits  nhits = %d\n",nhits);\r
315     for (Int_t iHit=0;iHit<nhits;iHit++) {\r
316       AliITShit *hit = (AliITShit*) arrHits->At(iHit);\r
317       fAliITSQADataMakerSim->FillHitsData(fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()],nmod);\r
318       Double_t xl,yl,zl,xl0,yl0,zl0;\r
319       Double_t tof,tof0;\r
320       hit->GetPositionL(xl,yl,zl,tof);\r
321       hit->GetPositionL0(xl0,yl0,zl0,tof0);\r
322       Float_t dyloc=TMath::Abs(yl-yl0)*10000.;\r
323       fAliITSQADataMakerSim->FillHitsData(1+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()],dyloc);\r
324       Float_t edep=hit->GetIonization()*1000000;\r
325       if(dyloc>200.){ \r
326         fAliITSQADataMakerSim->FillHitsData(2+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()],edep);\r
327         fAliITSQADataMakerSim->FillHitsData(3+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()],dyloc);\r
328       }\r
329     }\r
330   }\r
331   //\r
332   return rv ; \r
333 }\r
334 \r
335 //_______________________________________________________________\r
336 \r
337 Int_t AliITSQASDDDataMakerSim::GetOffset(AliQAv1::TASKINDEX_t task, Int_t specie){\r
338   // Returns histogram offset according to the specified task\r
339   Int_t offset=0;\r
340   if( task == AliQAv1::kHITS){\r
341     offset=fGenOffsetH[specie];  \r
342   }\r
343   else if( task == AliQAv1::kSDIGITS) {\r
344     offset=fGenOffsetS[specie];   \r
345   }\r
346   else if( task == AliQAv1::kDIGITS) {\r
347     offset=fGenOffsetD[specie];   \r
348   }\r
349   else {\r
350     AliInfo("No task has been selected. TaskHisto set to zero.\n");\r
351   }\r
352 \r
353   return offset;\r
354 }\r
355 \r
356 //____________________________________________________________________________ \r
357 void AliITSQASDDDataMakerSim::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset,Int_t specie ){\r
358   // Returns histogram offset according to the specified task\r
359   if( task == AliQAv1::kHITS){\r
360     fGenOffsetH[specie] = offset;  \r
361   }\r
362   else if( task == AliQAv1::kSDIGITS) {\r
363     fGenOffsetS[specie] = offset;   \r
364   }\r
365   else if( task == AliQAv1::kDIGITS) {\r
366     fGenOffsetD[specie] = offset;   \r
367   }\r
368   else {\r
369     AliInfo("No task has been selected. TaskHisto set to zero.\n");\r
370   }\r
371 }\r
372 \r
373 //_______________________________________________________________\r
374 \r
375 Int_t AliITSQASDDDataMakerSim::GetTaskHisto(AliQAv1::TASKINDEX_t task) {\r
376   // Returns the number of booked histograms for the selected task\r
377   Int_t histotot=0;\r
378   if( task == AliQAv1::kHITS) {\r
379     histotot=fSDDhHTask ;  \r
380   }\r
381   else if( task == AliQAv1::kSDIGITS) {\r
382     histotot=fSDDhSTask;   \r
383   }\r
384   else if( task == AliQAv1::kDIGITS) {\r
385     histotot=fSDDhDTask ;   \r
386   }\r
387   else {\r
388     AliInfo("No task has been selected. TaskHisto set to zero.\n");\r
389   }\r
390   return histotot;\r
391 \r
392 }\r