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