]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSQASPDDataMakerSim.cxx
updated
[u/mrichter/AliRoot.git] / ITS / AliITSQASPDDataMakerSim.cxx
CommitLineData
379510c2 1/**************************************************************************
2 * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17// *************************************************************
18// Checks the quality assurance
19// by comparing with reference data
20// contained in a DB
21// -------------------------------------------------------------
22// W. Ferrarese + P. Cerello INFN Torino Feb 2008
23// M. Nicassio D. Elia INFN Bari April 2008
24// maria.nicassio@ba.infn.it
25
26// --- ROOT system ---
27#include <TTree.h>
28#include <TH2.h>
29#include <TH1.h>
30// --- Standard library ---
31
32// --- AliRoot header files ---
33#include "AliRun.h"
34#include "AliITSQADataMakerSim.h"
35#include "AliITSQASPDDataMakerSim.h"
36#include "AliQA.h"
37#include "AliQAChecker.h"
38#include "AliITSdigit.h"
39#include "AliITSdigitSPD.h"
40#include "AliITS.h"
41#include "AliITSmodule.h"
42#include "AliITShit.h"
43#include "AliITSLoader.h"
44#include "AliRunLoader.h"
45
46ClassImp(AliITSQASPDDataMakerSim)
47
48//____________________________________________________________________________
49AliITSQASPDDataMakerSim::AliITSQASPDDataMakerSim(AliITSQADataMakerSim *aliITSQADataMakerSim) :
50TObject(),
51fAliITSQADataMakerSim(aliITSQADataMakerSim),
52fSPDhDigits(0),
53fSPDhSDigits(0),
54fSPDhHits(0),
55fDigitsOffset(0),
56fSDigitsOffset(0),
57fHitsOffset(0)
58{
59 //ctor used to discriminate OnLine-Offline analysis
60}
61
62//____________________________________________________________________________
63AliITSQASPDDataMakerSim::AliITSQASPDDataMakerSim(const AliITSQASPDDataMakerSim& qadm) :
64TObject(),
65fAliITSQADataMakerSim(qadm.fAliITSQADataMakerSim),
66fSPDhDigits(qadm.fSPDhDigits),
67fSPDhSDigits(qadm.fSPDhSDigits),
68fSPDhHits(qadm.fSPDhHits),
69fDigitsOffset(qadm.fDigitsOffset),
70fSDigitsOffset(qadm.fSDigitsOffset),
71fHitsOffset(qadm.fHitsOffset)
72{
73 //copy ctor
74 fAliITSQADataMakerSim->SetName((const char*)qadm.fAliITSQADataMakerSim->GetName()) ;
75 fAliITSQADataMakerSim->SetTitle((const char*)qadm.fAliITSQADataMakerSim->GetTitle());
76 }
77
78//__________________________________________________________________
79AliITSQASPDDataMakerSim& AliITSQASPDDataMakerSim::operator = (const AliITSQASPDDataMakerSim& qac )
80{
81 // Equal operator.
82 this->~AliITSQASPDDataMakerSim();
83 new(this) AliITSQASPDDataMakerSim(qac);
84 return *this;
85}
86
87//____________________________________________________________________________
88void AliITSQASPDDataMakerSim::StartOfDetectorCycle()
89{
90 //Detector specific actions at start of cycle
91 AliDebug(1,"AliITSQADM::Start of SPD Cycle\n");
92}
93
94//____________________________________________________________________________
95void AliITSQASPDDataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t /*task*/, TObjArray* /*list*/)
96{
97 // launch the QA checking
98 AliDebug(1,"AliITSDM instantiates checker with Run(AliQA::kITS, task, list)\n");
99
100 //AliQAChecker::Instance()->Run( AliQA::kITS , task, list);
101}
102
103//____________________________________________________________________________
104void AliITSQASPDDataMakerSim::InitDigits()
105{
106 // Initialization for DIGIT data - SPD -
107
108 fDigitsOffset = (fAliITSQADataMakerSim->fDigitsQAList)->GetEntries();
109 //fSPDhDigits must be incremented by one unit every time a histogram is ADDED to the QA List
110
111 Char_t name[50];
112 Char_t title[50];
113
114 TH1F *hlayer = new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
115 hlayer->GetXaxis()->SetTitle("Layer number");
116 hlayer->GetYaxis()->SetTitle("Entries");
117 fAliITSQADataMakerSim->Add2DigitsList(hlayer,1+fDigitsOffset);
118 fSPDhDigits++;
119
120 TH1F **hmod = new TH1F*[2];
121 for (Int_t iLay=0; iLay<2; iLay++) {
122 sprintf(name,"ModPattern_SPD%d",iLay+1);
123 sprintf(title,"Module map - SPD Layer %d",iLay+1);
124 hmod[iLay]=new TH1F(name,title,240,0,240);
125 hmod[iLay]->GetXaxis()->SetTitle("Module number");
126 hmod[iLay]->GetYaxis()->SetTitle("Entries");
127 fAliITSQADataMakerSim->Add2DigitsList(hmod[iLay],2+iLay+fDigitsOffset);
128 fSPDhDigits++;
129 }
130
131 TH1F *hcolumns = new TH1F("Columns_SPD","Columns - SPD",160,0.,160.);
132 hcolumns->GetXaxis()->SetTitle("Column number");
133 hcolumns->GetYaxis()->SetTitle("Entries");
134 fAliITSQADataMakerSim->Add2DigitsList(hcolumns,4+fDigitsOffset);
135 fSPDhDigits++;
136
137 TH1F *hrows = new TH1F("Rows_SPD","Rows - SPD",256,0.,256.);
138 hrows->GetXaxis()->SetTitle("Row number");
139 hrows->GetYaxis()->SetTitle("Entries");
140 fAliITSQADataMakerSim->Add2DigitsList(hrows,5+fDigitsOffset);
141 fSPDhDigits++;
142
143 TH1F** hMultSPDdigits = new TH1F*[2];
144 for (Int_t iLay=0; iLay<2; ++iLay) {
145 sprintf(name,"DigitMultiplicity_SPD%d",iLay+1);
146 sprintf(title,"Digit multiplicity - SPD Layer %d",iLay+1);
147 hMultSPDdigits[iLay]=new TH1F(name,title,200,0.,200.);
148 hMultSPDdigits[iLay]->GetXaxis()->SetTitle("Digit multiplicity");
149 hMultSPDdigits[iLay]->GetYaxis()->SetTitle("Entries");
150 fAliITSQADataMakerSim->Add2DigitsList(hMultSPDdigits[iLay], 6+iLay+fDigitsOffset);
151 fSPDhDigits++;
152 }
153
154 TH2F *hMultSPDdig2MultSPDdig1 = new TH2F("DigitMultCorrelation_SPD","Digit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
155 hMultSPDdig2MultSPDdig1->GetXaxis()->SetTitle("Digit multiplicity (Layer 1)");
156 hMultSPDdig2MultSPDdig1->GetYaxis()->SetTitle("Digit multiplicity (Layer 2)");
157 fAliITSQADataMakerSim->Add2DigitsList(hMultSPDdig2MultSPDdig1,8+fDigitsOffset);
158 fSPDhDigits++;
159
160 AliDebug(1,Form("%d SPD Digits histograms booked\n",fSPDhDigits));
161
162}
163
164//____________________________________________________________________________
165void AliITSQASPDDataMakerSim::MakeDigits(TTree *digits)
166{
167 // Fill QA for DIGIT - SPD -
168 AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
169 fITS->SetTreeAddress();
170 TClonesArray *iITSdigits = fITS->DigitsAddress(0); // 0->SPD
171
172 Int_t nDigitsL1=0;
173 Int_t nDigitsL2=0;
174
175 for (Int_t imod=0; imod<240; ++imod){
176 digits->GetEvent(imod);
177 Int_t ndigits = iITSdigits->GetEntries();
178 if (imod<80) {
179 fAliITSQADataMakerSim->GetDigitsData(1+fDigitsOffset)->Fill(0.5,ndigits);
180 fAliITSQADataMakerSim->GetDigitsData(2+fDigitsOffset)->Fill(imod,ndigits);
181 nDigitsL1+=ndigits;
182 }
183 else {
184 fAliITSQADataMakerSim->GetDigitsData(1+fDigitsOffset)->Fill(1,ndigits);
185 fAliITSQADataMakerSim->GetDigitsData(3+fDigitsOffset)->Fill(imod,ndigits);
186 nDigitsL2+=ndigits;
187 }
188 for (Int_t idig=0; idig<ndigits; ++idig) {
189 AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
190 Int_t col=dig->GetCoord1(); // cell number z
191 Int_t row=dig->GetCoord2(); // cell number x
192 fAliITSQADataMakerSim->GetDigitsData(4+fDigitsOffset)->Fill(col);
193 fAliITSQADataMakerSim->GetDigitsData(5+fDigitsOffset)->Fill(row);
194 }
195 }
196 fAliITSQADataMakerSim->GetDigitsData(6+fDigitsOffset)->Fill(nDigitsL1);
197 fAliITSQADataMakerSim->GetDigitsData(7+fDigitsOffset)->Fill(nDigitsL2);
198 fAliITSQADataMakerSim->GetDigitsData(8+fDigitsOffset)->Fill(nDigitsL1,nDigitsL2);
199}
200
201//____________________________________________________________________________
202void AliITSQASPDDataMakerSim::InitSDigits()
203{
204 // Initialization for SDIGIT data - SPD -
205 fSDigitsOffset = (fAliITSQADataMakerSim->fSDigitsQAList)->GetEntries();
206
207 //fSPDhSDigits must be incremented by one unit every time a histogram is ADDED to the QA List
208 Char_t name[50];
209 Char_t title[50];
210
211 TH1F *hlayer = new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
212 hlayer->GetXaxis()->SetTitle("Layer number");
213 hlayer->GetYaxis()->SetTitle("Entries");
214 fAliITSQADataMakerSim->Add2SDigitsList(hlayer,1+fSDigitsOffset);
215 fSPDhSDigits++;
216
217 TH1F **hmod = new TH1F*[2];
218 for (Int_t iLay=0; iLay<2; ++iLay) {
219 sprintf(name,"ModPattern_SPD%d",iLay+1);
220 sprintf(title,"Module map - SPD Layer %d",iLay+1);
221 hmod[iLay]=new TH1F(name,title,240,0,240);
222 hmod[iLay]->GetXaxis()->SetTitle("Module number");
223 hmod[iLay]->GetYaxis()->SetTitle("Entries");
224 fAliITSQADataMakerSim->Add2SDigitsList(hmod[iLay],2+iLay+fSDigitsOffset);
225 fSPDhSDigits++;
226 }
227
228
229 AliDebug(1,Form("%d SPD SDigits histograms booked\n",fSPDhSDigits));
230
231}
232
233//____________________________________________________________________________
234void AliITSQASPDDataMakerSim::MakeSDigits(TTree *sdigits)
235{
236 // Fill QA for SDIGIT - SPD -
237 TBranch *brchSDigits = sdigits->GetBranch("ITS");
238 for (Int_t imod=0; imod<240; ++imod){
239 TClonesArray * sdig = new TClonesArray( "AliITSpListItem",1000 );
240 brchSDigits->SetAddress( &sdig );
241 brchSDigits->GetEvent(imod);
242 Int_t nsdig=sdig->GetEntries();
243 if (imod<80) {
244 fAliITSQADataMakerSim->GetSDigitsData(1+fSDigitsOffset)->Fill(0.5,nsdig);
245 fAliITSQADataMakerSim->GetSDigitsData(2+fSDigitsOffset)->Fill(imod,nsdig);
246 }
247 else {
248 fAliITSQADataMakerSim->GetSDigitsData(1+fSDigitsOffset)->Fill(1,nsdig);
249 fAliITSQADataMakerSim->GetSDigitsData(3+fSDigitsOffset)->Fill(imod,nsdig);
250 }
251 delete sdig;
252 }
253
254}
255
256//____________________________________________________________________________
257void AliITSQASPDDataMakerSim::InitHits()
258{
259 // Initialization for HITS data - SPD -
260 fHitsOffset = (fAliITSQADataMakerSim->fHitsQAList)->GetEntries();
261
262 //fSPDhHits must be incremented by one unit every time a histogram is ADDED to the QA List
263 Char_t name[50];
264 Char_t title[50];
265
266 TH1F *hlayer = new TH1F("LayPattern_SPD","Layer map - SPD",6,0.,6.);
267 hlayer->GetXaxis()->SetTitle("Layer number");
268 hlayer->GetYaxis()->SetTitle("Entries");
269 fAliITSQADataMakerSim->Add2HitsList(hlayer,1+fHitsOffset);
270 fSPDhHits++;
271
272 TH1F **hmod = new TH1F*[2];
273 for (Int_t iLay=0; iLay<2; ++iLay) {
274 sprintf(name,"ModPattern_SPD%d",iLay+1);
275 sprintf(title,"Module map - SPD Layer %d",iLay+1);
276 hmod[iLay]=new TH1F(name,title,240,0,240);
277 hmod[iLay]->GetXaxis()->SetTitle("Module number");
278 hmod[iLay]->GetYaxis()->SetTitle("Entries");
279 fAliITSQADataMakerSim->Add2HitsList(hmod[iLay],2+iLay+fHitsOffset);
280 fSPDhHits++;
281 }
282
283 TH1F *hhitlenght = new TH1F("Lenght_SPD","Hit lenght along y_{loc} coord",210,0.,210.);
284 hhitlenght->GetXaxis()->SetTitle("Hit lenght [#mum]");
285 hhitlenght->GetYaxis()->SetTitle("# hits");
286 fAliITSQADataMakerSim->Add2HitsList(hhitlenght,4+fHitsOffset);
287 fSPDhHits++;
288
289 TH1F *hEdepos = new TH1F("EnergyDeposit_SPD","Deposited energy distribution (y_{loc}>180 #mum)",150,0.,300.);
290 hEdepos->GetXaxis()->SetTitle("Deposited energy [keV]");
291 hEdepos->GetYaxis()->SetTitle("# hits");
292 fAliITSQADataMakerSim->Add2HitsList(hEdepos,5+fHitsOffset);
293 fSPDhHits++;
294
295 AliDebug(1,Form("%d SPD Hits histograms booked\n",fSPDhHits));
296
297}
298
299//____________________________________________________________________________
300void AliITSQASPDDataMakerSim::MakeHits(TTree *hits)
301{
302 // Fill QA for HITS - SPD -
303 AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
304 fITS->SetTreeAddress();
305 Int_t nmodules;
306 fITS->InitModules(-1,nmodules); //-1->number of modules taken from AliITSgeom class kept in fITSgeom
307 //nmodules is set
308
309 fITS->FillModules(hits,0);
310
311 for (Int_t imod=0; imod<240; ++imod){
312 AliITSmodule *module = fITS->GetModule(imod);
313 TObjArray *arrHits = module->GetHits();
314 Int_t nhits = arrHits->GetEntriesFast();
315 if (imod<80) {
316 fAliITSQADataMakerSim->GetHitsData(1+fHitsOffset)->Fill(0.5,nhits);
317 fAliITSQADataMakerSim->GetHitsData(2+fHitsOffset)->Fill(imod,nhits);
318 } else {
319 fAliITSQADataMakerSim->GetHitsData(1+fHitsOffset)->Fill(1,nhits);
320 fAliITSQADataMakerSim->GetHitsData(3+fHitsOffset)->Fill(imod,nhits);
321 }
322 //printf("--w--AliITSQASPDDataMakerSim::MakeHits nhits = %d\n",nhits);
323 for (Int_t iHit=0; iHit<nhits; ++iHit) {
324 AliITShit *hit = (AliITShit*) arrHits->At(iHit);
325 Double_t xl,yl,zl,xl0,yl0,zl0;
326 Double_t tof,tof0;
327 hit->GetPositionL(xl,yl,zl,tof);
328 hit->GetPositionL0(xl0,yl0,zl0,tof0);
329 Float_t dyloc=TMath::Abs(yl-yl0)*10000.;
330 fAliITSQADataMakerSim->GetHitsData(4+fHitsOffset)->Fill(dyloc);
331 Float_t edep=hit->GetIonization()*1000000;
332 if(dyloc>180.){
333 fAliITSQADataMakerSim->GetHitsData(5+fHitsOffset)->Fill(edep);
334 }
335 }
336 }
337}