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