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