]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASSDDataMakerSim.cxx
Updated macro
[u/mrichter/AliRoot.git] / ITS / AliITSQASSDDataMakerSim.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 //  SSD QA part: P. Christakoglou - NIKHEF/UU\r
26 \r
27 // --- ROOT system ---\r
28 #include <TTree.h>\r
29 #include <TH1.h>\r
30 #include <TH2.h>\r
31 #include <TMath.h>\r
32 // --- Standard library ---\r
33 \r
34 // --- AliRoot header files ---\r
35 #include "AliITS.h"\r
36 #include "AliITSmodule.h"\r
37 #include "AliITShit.h"\r
38 #include "AliITSdigit.h"\r
39 #include "AliITSpListItem.h"\r
40 #include "AliRun.h"\r
41 #include "AliITSQADataMakerSim.h"\r
42 #include "AliITSQASSDDataMakerSim.h"\r
43 #include "AliLog.h"\r
44 #include "AliQA.h"\r
45 #include "AliQAChecker.h"\r
46 #include "AliRawReader.h"\r
47 \r
48 ClassImp(AliITSQASSDDataMakerSim)\r
49 \r
50 //____________________________________________________________________________ \r
51 AliITSQASSDDataMakerSim::AliITSQASSDDataMakerSim(AliITSQADataMakerSim *aliITSQADataMakerSim) :\r
52 TObject(),\r
53 fAliITSQADataMakerSim(aliITSQADataMakerSim),\r
54 fSSDhTask(0),\r
55 fGenOffset(0) {\r
56   //ctor used to discriminate OnLine-Offline analysis   \r
57 }\r
58 \r
59 //____________________________________________________________________________ \r
60 AliITSQASSDDataMakerSim::AliITSQASSDDataMakerSim(const AliITSQASSDDataMakerSim& qadm) :\r
61 TObject(),\r
62 fAliITSQADataMakerSim(qadm.fAliITSQADataMakerSim),\r
63 fSSDhTask(qadm.fSSDhTask),\r
64 fGenOffset(qadm.fGenOffset) {\r
65   //copy ctor \r
66   fAliITSQADataMakerSim->SetName((const char*)qadm.fAliITSQADataMakerSim->GetName()) ; \r
67   fAliITSQADataMakerSim->SetTitle((const char*)qadm.fAliITSQADataMakerSim->GetTitle());\r
68   }\r
69 \r
70 //__________________________________________________________________\r
71 AliITSQASSDDataMakerSim& AliITSQASSDDataMakerSim::operator = (const AliITSQASSDDataMakerSim& qac ) {\r
72   // Equal operator.\r
73   this->~AliITSQASSDDataMakerSim();\r
74   new(this) AliITSQASSDDataMakerSim(qac);\r
75   return *this;\r
76 }\r
77 \r
78 //____________________________________________________________________________ \r
79 void AliITSQASSDDataMakerSim::StartOfDetectorCycle() {\r
80   //Detector specific actions at start of cycle\r
81   AliDebug(1,"AliITSQADM::Start of SSD Cycle\n");\r
82 }\r
83 \r
84 //____________________________________________________________________________ \r
85 void AliITSQASSDDataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list) {\r
86   // launch the QA checking\r
87   AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n"); \r
88   \r
89   AliQAChecker::Instance()->Run( AliQA::kITS , task, list);\r
90 }\r
91 \r
92 //____________________________________________________________________________ \r
93 void AliITSQASSDDataMakerSim::InitDigits() { \r
94   // Initialization for DIGIT data - SSD -\r
95   fGenOffset = (fAliITSQADataMakerSim->fDigitsQAList)->GetEntries();\r
96 \r
97   // custom code here\r
98   TH1F *fHistSSDModule = new TH1F("fHistSSDDigitsModule",\r
99                                   ";SSD Module Number;N_{DIGITS}",\r
100                                   1698,499.5,2197.5);  \r
101   fAliITSQADataMakerSim->Add2DigitsList(fHistSSDModule,\r
102                                         fGenOffset + 0);\r
103   fSSDhTask += 1;\r
104   TH2F *fHistSSDModuleStrip = new TH2F("fHistSSDDigitsModuleStrip",\r
105                                        ";N_{Strip};N_{Module}",\r
106                                        1540,0,1540,1698,499.5,2197.5);  \r
107   fAliITSQADataMakerSim->Add2DigitsList(fHistSSDModuleStrip,\r
108                                         fGenOffset + 1);\r
109   fSSDhTask += 1;\r
110 \r
111   AliDebug(1,Form("%d SSD Digits histograms booked\n",fSSDhTask));\r
112 \r
113 }\r
114 \r
115 //____________________________________________________________________________\r
116 void AliITSQASSDDataMakerSim::MakeDigits(TTree *digits) { \r
117   // Fill QA for DIGIT - SSD -\r
118   AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");\r
119   fITS->SetTreeAddress();\r
120   TClonesArray *iSSDdigits  = fITS->DigitsAddress(2);\r
121   for(Int_t iModule = 500; iModule < 2198; iModule++) {\r
122     iSSDdigits->Clear();\r
123     digits->GetEvent(iModule);    \r
124     Int_t ndigits = iSSDdigits->GetEntries();\r
125     fAliITSQADataMakerSim->GetDigitsData(fGenOffset + 0)->Fill(iModule,ndigits);\r
126     if(ndigits != 0)\r
127       AliDebug(1,Form("Module: %d - Digits: %d",iModule,ndigits));\r
128  \r
129     for (Int_t iDigit = 0; iDigit < ndigits; iDigit++) {\r
130       AliITSdigit *dig = (AliITSdigit*)iSSDdigits->UncheckedAt(iDigit);\r
131       Int_t fStripNumber = (dig->GetCoord1() == 0) ? dig->GetCoord2() : dig->GetCoord2() + fgkNumberOfPSideStrips;\r
132       ((TH2F *)fAliITSQADataMakerSim->GetDigitsData(fGenOffset + 1))->Fill(fStripNumber,iModule,dig->GetSignal());\r
133     }//digit loop\r
134   }//module loop\r
135 }\r
136 \r
137 //____________________________________________________________________________ \r
138 void AliITSQASSDDataMakerSim::InitSDigits() { \r
139   // Initialization for SDIGIT data - SSD -\r
140   fGenOffset = (fAliITSQADataMakerSim->fSDigitsQAList)->GetEntries();\r
141 \r
142   // custom code here\r
143   TH1F *fHistSSDModule = new TH1F("fHistSSDSDigitsModule",\r
144                                   ";SSD Module Number;N_{SDIGITS}",\r
145                                   1698,499.5,2197.5);  \r
146   fAliITSQADataMakerSim->Add2SDigitsList(fHistSSDModule,\r
147                                          fGenOffset + 0);\r
148   fSSDhTask += 1;  \r
149 \r
150   AliDebug(1,Form("%d SSD SDigits histograms booked\n",fSSDhTask));\r
151 }\r
152 \r
153 //____________________________________________________________________________\r
154 void AliITSQASSDDataMakerSim::MakeSDigits(TTree *sdigits) { \r
155   // Fill QA for SDIGIT - SSD -\r
156   static TClonesArray iSSDEmpty("AliITSpListItem",10000);\r
157   iSSDEmpty.Clear();\r
158   TClonesArray *iSSDsdigits = &iSSDEmpty;\r
159 \r
160   AliInfo(Form("Trying to access the sdigits histogram: %d\n",fGenOffset));\r
161 \r
162   TBranch *brchSDigits = sdigits->GetBranch("ITS");\r
163   brchSDigits->SetAddress(&iSSDsdigits);\r
164   for(Int_t iModule = 500; iModule < 2198; iModule++) {\r
165     iSSDsdigits->Clear();\r
166     sdigits->GetEvent(iModule);    \r
167     Int_t ndigits = iSSDsdigits->GetEntries();\r
168     fAliITSQADataMakerSim->GetSDigitsData(fGenOffset + 0)->Fill(iModule,ndigits);\r
169     if(ndigits != 0)\r
170       AliDebug(1,Form("Module: %d - Digits: %d",iModule,ndigits));\r
171 \r
172     for (Int_t iDigit = 0; iDigit < ndigits; iDigit++) {\r
173       AliITSpListItem *dig=(AliITSpListItem*)iSSDsdigits->At(iDigit);\r
174       dig=0;\r
175     }//digit loop\r
176   }//module loop\r
177 }\r
178 \r
179 //____________________________________________________________________________ \r
180 void AliITSQASSDDataMakerSim::InitHits() { \r
181   // Initialization for HITS data - SSD -\r
182   fGenOffset = (fAliITSQADataMakerSim->fHitsQAList)->GetEntries();\r
183 \r
184   // custom code here\r
185   TH1F *fHistSSDModule = new TH1F("fHistSSDHitsModule",\r
186                                   ";SDD Module Number;N_{HITS}",\r
187                                   1698,499.5,2197.5); \r
188   fAliITSQADataMakerSim->Add2HitsList(fHistSSDModule,\r
189                                       fGenOffset + 0);\r
190   fSSDhTask += 1;\r
191   TH1F *fHistSSDGlobalX = new TH1F("fHistSSDHitsGlobalX",\r
192                                    ";x [cm];Entries",\r
193                                    1000,-50.,50.);\r
194   fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalX,\r
195                                       fGenOffset + 1);\r
196   fSSDhTask += 1;\r
197   TH1F *fHistSSDGlobalY = new TH1F("fHistSSDHitsGlobalY",\r
198                                    ";y [cm];Entries",\r
199                                    1000,-50.,50.);\r
200   fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalY,\r
201                                       fGenOffset + 2);\r
202   fSSDhTask += 1;\r
203   TH1F *fHistSSDGlobalZ = new TH1F("fHistSSDHitsGlobalZ",\r
204                                    ";z [cm];Entries",\r
205                                    1000,-60.,60.);\r
206   fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalZ,\r
207                                       fGenOffset + 3);\r
208   fSSDhTask += 1;\r
209   TH1F *fHistSSDLocalX = new TH1F("fHistSSDHitsLocalX",\r
210                                   ";x [cm];Entries",\r
211                                   1000,-4.,4.);\r
212   fAliITSQADataMakerSim->Add2HitsList(fHistSSDLocalX,\r
213                                       fGenOffset + 4);\r
214   fSSDhTask += 1;\r
215   TH1F *fHistSSDLocalY = new TH1F("fHistSSDHitsLocalY",\r
216                                   ";y [cm];Entries",\r
217                                   1000,-0.1,0.1);\r
218   fAliITSQADataMakerSim->Add2HitsList(fHistSSDLocalY,\r
219                                       fGenOffset + 5);\r
220   fSSDhTask += 1;\r
221   TH1F *fHistSSDLocalZ = new TH1F("fHistSSDHitsLocalZ",\r
222                                   ";z [cm];Entries",\r
223                                   1000,-4.,4.);\r
224   fAliITSQADataMakerSim->Add2HitsList(fHistSSDLocalZ,\r
225                                       fGenOffset + 6);\r
226   fSSDhTask += 1;\r
227   TH1F *fHistSSDIonization = new TH1F("fHistSSDHitsIonization",\r
228                                       ";log(dE/dx) [KeV];N_{Hits}",\r
229                                       100,-7,-2);\r
230   fAliITSQADataMakerSim->Add2HitsList(fHistSSDIonization,\r
231                                       fGenOffset + 7 );\r
232   fSSDhTask += 1;\r
233   TH2F *fHistSSDGlobalXY = new TH2F("fHistSSDHitsGlobalXY",\r
234                                     ";x [cm];y [cm]",\r
235                                     1000,-50.,50.,\r
236                                     1000,-50.,50.);\r
237   fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalXY,\r
238                                       fGenOffset + 8 );\r
239   fSSDhTask += 1;\r
240  \r
241   AliDebug(1,Form("%d SSD Hits histograms booked\n",fSSDhTask));\r
242 }\r
243 \r
244 \r
245 //____________________________________________________________________________\r
246 void AliITSQASSDDataMakerSim::MakeHits(TTree *hits) { \r
247   // Fill QA for HITS - SSD -\r
248   AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");\r
249   fITS->SetTreeAddress();\r
250   Int_t nmodules;\r
251   fITS->InitModules(-1,nmodules);\r
252   fITS->FillModules(hits,0);\r
253   for(Int_t iModule = 500; iModule < 2198; iModule++) {\r
254     AliITSmodule *module = fITS->GetModule(iModule);\r
255     TObjArray *arrHits = module->GetHits();\r
256     Int_t nhits = arrHits->GetEntriesFast();\r
257     if(nhits != 0)\r
258       AliDebug(1,Form("Module: %d - Hits: %d",iModule,nhits));\r
259     for (Int_t iHit = 0; iHit < nhits; iHit++) {\r
260       AliITShit *hit = (AliITShit*) arrHits->At(iHit);\r
261       \r
262       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 0)->Fill(iModule);\r
263       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 1)->Fill(hit->GetXG());\r
264       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 2)->Fill(hit->GetYG());\r
265       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 3)->Fill(hit->GetZG());\r
266       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 4)->Fill(hit->GetXL());\r
267       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 5)->Fill(hit->GetYL());\r
268       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 6)->Fill(hit->GetZL());\r
269       if(hit->GetIonization())\r
270         fAliITSQADataMakerSim->GetHitsData(fGenOffset + 7)->Fill(TMath::Log10(hit->GetIonization()));\r
271       fAliITSQADataMakerSim->GetHitsData(fGenOffset + 8)->Fill(hit->GetXG(),hit->GetYG());\r
272     }//hit loop\r
273   }//module loop  \r
274 }\r