]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSQASPDDataMakerSim.cxx
Transfer of the initialisation of the QA Data objects in the framework; clean the...
[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"
4e25ac79 36#include "AliQAv1.h"
379510c2 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),
7a0e5776 52fSPDhHTask(0),
53fSPDhSTask(0),
54fSPDhDTask(0),
3f905799 55fGenOffsetH(0),
56fGenOffsetS(0),
57fGenOffsetD(0)
379510c2 58{
59 //ctor used to discriminate OnLine-Offline analysis
60}
61
62//____________________________________________________________________________
63AliITSQASPDDataMakerSim::AliITSQASPDDataMakerSim(const AliITSQASPDDataMakerSim& qadm) :
64TObject(),
65fAliITSQADataMakerSim(qadm.fAliITSQADataMakerSim),
7a0e5776 66fSPDhHTask(qadm.fSPDhHTask),
67fSPDhSTask(qadm.fSPDhSTask),
68fSPDhDTask(qadm.fSPDhDTask),
3f905799 69fGenOffsetH(qadm.fGenOffsetH),
70fGenOffsetS(qadm.fGenOffsetS),
71fGenOffsetD(qadm.fGenOffsetD)
379510c2 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
5379c4a3 91 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SPD Cycle\n");
379510c2 92}
93
94//____________________________________________________________________________
4e25ac79 95void AliITSQASPDDataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
379510c2 96{
97 // launch the QA checking
5379c4a3 98 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n");
379510c2 99
4e25ac79 100 //AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);
379510c2 101}
102
103//____________________________________________________________________________
eca4fa66 104Int_t AliITSQASPDDataMakerSim::InitDigits()
379510c2 105{
106 // Initialization for DIGIT data - SPD -
7d297381 107 const Bool_t expert = kTRUE ;
108 const Bool_t image = kTRUE ;
eca4fa66 109 Int_t rv = 0 ;
110 //fGenOffsetD = (fAliITSQADataMakerSim->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
7a0e5776 111 //fSPDhDTask must be incremented by one unit every time a histogram is ADDED to the QA List
379510c2 112
113 Char_t name[50];
114 Char_t title[50];
115
26ee9565 116 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
379510c2 117 hlayer->GetXaxis()->SetTitle("Layer number");
118 hlayer->GetYaxis()->SetTitle("Entries");
eca4fa66 119 rv = fAliITSQADataMakerSim->Add2DigitsList(hlayer,fGenOffsetD, expert, !image);
7a0e5776 120 fSPDhDTask++;
379510c2 121
122 TH1F **hmod = new TH1F*[2];
123 for (Int_t iLay=0; iLay<2; iLay++) {
26ee9565 124 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
379510c2 125 sprintf(title,"Module map - SPD Layer %d",iLay+1);
126 hmod[iLay]=new TH1F(name,title,240,0,240);
127 hmod[iLay]->GetXaxis()->SetTitle("Module number");
128 hmod[iLay]->GetYaxis()->SetTitle("Entries");
eca4fa66 129 rv = fAliITSQADataMakerSim->Add2DigitsList(hmod[iLay],1+iLay+fGenOffsetD, !expert, image);
7a0e5776 130 fSPDhDTask++;
379510c2 131 }
132
26ee9565 133 TH1F *hcolumns = new TH1F("SPDColumns_SPD","Columns - SPD",160,0.,160.);
379510c2 134 hcolumns->GetXaxis()->SetTitle("Column number");
135 hcolumns->GetYaxis()->SetTitle("Entries");
7d297381 136 fAliITSQADataMakerSim->Add2DigitsList(hcolumns,3+fGenOffsetD, expert, !image);
7a0e5776 137 fSPDhDTask++;
379510c2 138
26ee9565 139 TH1F *hrows = new TH1F("SPDRows_SPD","Rows - SPD",256,0.,256.);
379510c2 140 hrows->GetXaxis()->SetTitle("Row number");
141 hrows->GetYaxis()->SetTitle("Entries");
eca4fa66 142 rv = fAliITSQADataMakerSim->Add2DigitsList(hrows,4+fGenOffsetD, expert, !image);
7a0e5776 143 fSPDhDTask++;
379510c2 144
145 TH1F** hMultSPDdigits = new TH1F*[2];
146 for (Int_t iLay=0; iLay<2; ++iLay) {
26ee9565 147 sprintf(name,"SPDDigitMultiplicity_SPD%d",iLay+1);
379510c2 148 sprintf(title,"Digit multiplicity - SPD Layer %d",iLay+1);
149 hMultSPDdigits[iLay]=new TH1F(name,title,200,0.,200.);
150 hMultSPDdigits[iLay]->GetXaxis()->SetTitle("Digit multiplicity");
151 hMultSPDdigits[iLay]->GetYaxis()->SetTitle("Entries");
eca4fa66 152 rv = fAliITSQADataMakerSim->Add2DigitsList(hMultSPDdigits[iLay], 5+iLay+fGenOffsetD, !expert, image);
7a0e5776 153 fSPDhDTask++;
379510c2 154 }
155
954ef57a 156 TH2F *hMultSPDdig2MultSPDdig1
26ee9565 157 = new TH2F("SPDDigitMultCorrelation_SPD","Digit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
379510c2 158 hMultSPDdig2MultSPDdig1->GetXaxis()->SetTitle("Digit multiplicity (Layer 1)");
159 hMultSPDdig2MultSPDdig1->GetYaxis()->SetTitle("Digit multiplicity (Layer 2)");
eca4fa66 160 rv = fAliITSQADataMakerSim->Add2DigitsList(hMultSPDdig2MultSPDdig1,7+fGenOffsetD, !expert, image);
7a0e5776 161 fSPDhDTask++;
379510c2 162
5379c4a3 163 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Digits histograms booked\n",fSPDhDTask));
eca4fa66 164 return rv ;
379510c2 165}
166
167//____________________________________________________________________________
eca4fa66 168Int_t AliITSQASPDDataMakerSim::MakeDigits(TTree *digits)
379510c2 169{
170 // Fill QA for DIGIT - SPD -
eca4fa66 171 Int_t rv = 0 ;
6252ceeb 172
379510c2 173 AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
174 fITS->SetTreeAddress();
175 TClonesArray *iITSdigits = fITS->DigitsAddress(0); // 0->SPD
176
177 Int_t nDigitsL1=0;
178 Int_t nDigitsL2=0;
179
180 for (Int_t imod=0; imod<240; ++imod){
181 digits->GetEvent(imod);
182 Int_t ndigits = iITSdigits->GetEntries();
183 if (imod<80) {
3f905799 184 fAliITSQADataMakerSim->GetDigitsData(0+fGenOffsetD)->Fill(0.5,ndigits);
185 fAliITSQADataMakerSim->GetDigitsData(1+fGenOffsetD)->Fill(imod,ndigits);
379510c2 186 nDigitsL1+=ndigits;
187 }
188 else {
3f905799 189 fAliITSQADataMakerSim->GetDigitsData(0+fGenOffsetD)->Fill(1,ndigits);
190 fAliITSQADataMakerSim->GetDigitsData(2+fGenOffsetD)->Fill(imod,ndigits);
379510c2 191 nDigitsL2+=ndigits;
192 }
193 for (Int_t idig=0; idig<ndigits; ++idig) {
194 AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
195 Int_t col=dig->GetCoord1(); // cell number z
196 Int_t row=dig->GetCoord2(); // cell number x
3f905799 197 fAliITSQADataMakerSim->GetDigitsData(3+fGenOffsetD)->Fill(col);
198 fAliITSQADataMakerSim->GetDigitsData(4+fGenOffsetD)->Fill(row);
379510c2 199 }
200 }
3f905799 201 fAliITSQADataMakerSim->GetDigitsData(5+fGenOffsetD)->Fill(nDigitsL1);
202 fAliITSQADataMakerSim->GetDigitsData(6+fGenOffsetD)->Fill(nDigitsL2);
203 fAliITSQADataMakerSim->GetDigitsData(7+fGenOffsetD)->Fill(nDigitsL1,nDigitsL2);
eca4fa66 204 return rv ;
379510c2 205}
206
207//____________________________________________________________________________
eca4fa66 208Int_t AliITSQASPDDataMakerSim::InitSDigits()
379510c2 209{
210 // Initialization for SDIGIT data - SPD -
7d297381 211 const Bool_t expert = kTRUE ;
212 const Bool_t image = kTRUE ;
eca4fa66 213 Int_t rv = 0 ;
214 //fGenOffsetS = (fAliITSQADataMakerSim->fSDigitsQAList[AliRecoParam::kDefault])->GetEntries();
c71529b0 215 //printf("--W-- AliITSQASPDDataMakerSim::InitSDigits() fGenOffset= %d \n",fGenOffset);
7a0e5776 216 //fSPDhSTask must be incremented by one unit every time a histogram is ADDED to the QA List
c71529b0 217
379510c2 218 Char_t name[50];
219 Char_t title[50];
220
26ee9565 221 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
379510c2 222 hlayer->GetXaxis()->SetTitle("Layer number");
223 hlayer->GetYaxis()->SetTitle("Entries");
eca4fa66 224 rv = fAliITSQADataMakerSim->Add2SDigitsList(hlayer,fGenOffsetS, expert, !image);
7a0e5776 225 fSPDhSTask++;
379510c2 226
227 TH1F **hmod = new TH1F*[2];
228 for (Int_t iLay=0; iLay<2; ++iLay) {
26ee9565 229 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
379510c2 230 sprintf(title,"Module map - SPD Layer %d",iLay+1);
231 hmod[iLay]=new TH1F(name,title,240,0,240);
232 hmod[iLay]->GetXaxis()->SetTitle("Module number");
233 hmod[iLay]->GetYaxis()->SetTitle("Entries");
eca4fa66 234 rv = fAliITSQADataMakerSim->Add2SDigitsList(hmod[iLay],1+iLay+fGenOffsetS, !expert, image);
7a0e5776 235 fSPDhSTask++;
379510c2 236 }
379510c2 237
5379c4a3 238 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD SDigits histograms booked\n",fSPDhSTask));
eca4fa66 239 return rv ;
379510c2 240}
241
242//____________________________________________________________________________
eca4fa66 243Int_t AliITSQASPDDataMakerSim::MakeSDigits(TTree *sdigits)
379510c2 244{
245 // Fill QA for SDIGIT - SPD -
eca4fa66 246 Int_t rv = 0 ;
6252ceeb 247
248 static TClonesArray * sdig ;
249 if (! sdig )
250 sdig = new TClonesArray( "AliITSpListItem",1000 );
251
379510c2 252 TBranch *brchSDigits = sdigits->GetBranch("ITS");
253 for (Int_t imod=0; imod<240; ++imod){
379510c2 254 brchSDigits->SetAddress( &sdig );
255 brchSDigits->GetEvent(imod);
256 Int_t nsdig=sdig->GetEntries();
257 if (imod<80) {
3f905799 258 fAliITSQADataMakerSim->GetSDigitsData(0+fGenOffsetS)->Fill(0.5,nsdig);
259 fAliITSQADataMakerSim->GetSDigitsData(1+fGenOffsetS)->Fill(imod,nsdig);
379510c2 260 }
261 else {
3f905799 262 fAliITSQADataMakerSim->GetSDigitsData(0+fGenOffsetS)->Fill(1,nsdig);
263 fAliITSQADataMakerSim->GetSDigitsData(2+fGenOffsetS)->Fill(imod,nsdig);
379510c2 264 }
6252ceeb 265 sdig->Clear() ;
379510c2 266 }
eca4fa66 267 return rv ;
379510c2 268}
269
270//____________________________________________________________________________
eca4fa66 271Int_t AliITSQASPDDataMakerSim::InitHits()
379510c2 272{
273 // Initialization for HITS data - SPD -
7d297381 274 const Bool_t expert = kTRUE ;
275 const Bool_t image = kTRUE ;
eca4fa66 276 Int_t rv = 0 ;
7d297381 277
eca4fa66 278 //fGenOffsetH = (fAliITSQADataMakerSim->fHitsQAList[AliRecoParam::kDefault])->GetEntries();
c71529b0 279 //printf("--W-- AliITSQASPDDataMakerSim::InitHits() fGenOffset= %d \n",fGenOffset);
7a0e5776 280 //fSPDhHTask must be incremented by one unit every time a histogram is ADDED to the QA List
379510c2 281 Char_t name[50];
282 Char_t title[50];
283
26ee9565 284 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
379510c2 285 hlayer->GetXaxis()->SetTitle("Layer number");
286 hlayer->GetYaxis()->SetTitle("Entries");
eca4fa66 287 rv = fAliITSQADataMakerSim->Add2HitsList(hlayer,fGenOffsetH, expert, !image);
7a0e5776 288 fSPDhHTask++;
379510c2 289
290 TH1F **hmod = new TH1F*[2];
291 for (Int_t iLay=0; iLay<2; ++iLay) {
26ee9565 292 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
379510c2 293 sprintf(title,"Module map - SPD Layer %d",iLay+1);
294 hmod[iLay]=new TH1F(name,title,240,0,240);
295 hmod[iLay]->GetXaxis()->SetTitle("Module number");
296 hmod[iLay]->GetYaxis()->SetTitle("Entries");
eca4fa66 297 rv = fAliITSQADataMakerSim->Add2HitsList(hmod[iLay],1+iLay+fGenOffsetH, !expert, image);
7a0e5776 298 fSPDhHTask++;
379510c2 299 }
300
db72ff3b 301 TH1F *hhitlenght = new TH1F("SPDLenght_SPD","SPD Hit lenght along y_{loc} coord",210,0.,210.);
379510c2 302 hhitlenght->GetXaxis()->SetTitle("Hit lenght [#mum]");
303 hhitlenght->GetYaxis()->SetTitle("# hits");
eca4fa66 304 rv = fAliITSQADataMakerSim->Add2HitsList(hhitlenght,3+fGenOffsetH, !expert, image);
7a0e5776 305 fSPDhHTask++;
379510c2 306
db72ff3b 307 TH1F *hEdepos = new TH1F("SPDEnergyDeposit_SPD","SPD Deposited energy distribution (y_{loc}>180 #mum)",150,0.,300.);
379510c2 308 hEdepos->GetXaxis()->SetTitle("Deposited energy [keV]");
309 hEdepos->GetYaxis()->SetTitle("# hits");
eca4fa66 310 rv = fAliITSQADataMakerSim->Add2HitsList(hEdepos,4+fGenOffsetH, !expert, image);
7a0e5776 311 fSPDhHTask++;
379510c2 312
5379c4a3 313 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Hits histograms booked\n",fSPDhHTask));
eca4fa66 314 return rv ;
379510c2 315}
316
317//____________________________________________________________________________
eca4fa66 318Int_t AliITSQASPDDataMakerSim::MakeHits(TTree *hits)
379510c2 319{
320 // Fill QA for HITS - SPD -
eca4fa66 321 Int_t rv = 0 ;
eca4fa66 322
379510c2 323 AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
324 fITS->SetTreeAddress();
325 Int_t nmodules;
326 fITS->InitModules(-1,nmodules); //-1->number of modules taken from AliITSgeom class kept in fITSgeom
327 //nmodules is set
328
329 fITS->FillModules(hits,0);
330
331 for (Int_t imod=0; imod<240; ++imod){
332 AliITSmodule *module = fITS->GetModule(imod);
333 TObjArray *arrHits = module->GetHits();
334 Int_t nhits = arrHits->GetEntriesFast();
335 if (imod<80) {
3f905799 336 fAliITSQADataMakerSim->GetHitsData(fGenOffsetH)->Fill(0.5,nhits);
337 fAliITSQADataMakerSim->GetHitsData(1+fGenOffsetH)->Fill(imod,nhits);
379510c2 338 } else {
3f905799 339 fAliITSQADataMakerSim->GetHitsData(fGenOffsetH)->Fill(1,nhits);
340 fAliITSQADataMakerSim->GetHitsData(2+fGenOffsetH)->Fill(imod,nhits);
379510c2 341 }
379510c2 342 for (Int_t iHit=0; iHit<nhits; ++iHit) {
343 AliITShit *hit = (AliITShit*) arrHits->At(iHit);
344 Double_t xl,yl,zl,xl0,yl0,zl0;
345 Double_t tof,tof0;
346 hit->GetPositionL(xl,yl,zl,tof);
347 hit->GetPositionL0(xl0,yl0,zl0,tof0);
348 Float_t dyloc=TMath::Abs(yl-yl0)*10000.;
3f905799 349 fAliITSQADataMakerSim->GetHitsData(3+fGenOffsetH)->Fill(dyloc);
379510c2 350 Float_t edep=hit->GetIonization()*1000000;
351 if(dyloc>180.){
3f905799 352 fAliITSQADataMakerSim->GetHitsData(4+fGenOffsetH)->Fill(edep);
379510c2 353 }
354 }
355 }
eca4fa66 356 return rv ;
379510c2 357}
7a0e5776 358
359
360//_______________________________________________________________
361
4e25ac79 362Int_t AliITSQASPDDataMakerSim::GetOffset(AliQAv1::TASKINDEX_t task){
7a0e5776 363 // Returns histogram offset according to the specified task
364 Int_t offset=0;
4e25ac79 365 if( task == AliQAv1::kHITS){
7a0e5776 366 offset=fGenOffsetH;
367 }
4e25ac79 368 else if( task == AliQAv1::kSDIGITS) {
7a0e5776 369 offset=fGenOffsetS;
370 }
4e25ac79 371 else if( task == AliQAv1::kDIGITS) {
7a0e5776 372 offset=fGenOffsetD;
373 }
374 else {
375 AliInfo("No task has been selected. TaskHisto set to zero.\n");
376 }
377
378 return offset;
379}
380
eca4fa66 381//____________________________________________________________________________
382void AliITSQASPDDataMakerSim::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset){
383 // Returns histogram offset according to the specified task
384 if( task == AliQAv1::kHITS){
385 fGenOffsetH = offset;
386 }
387 else if( task == AliQAv1::kSDIGITS) {
388 fGenOffsetS = offset;
389 }
390 else if( task == AliQAv1::kDIGITS) {
391 fGenOffsetD = offset;
392 }
393 else {
394 AliInfo("No task has been selected. TaskHisto set to zero.\n");
395 }
396}
7a0e5776 397
398//_______________________________________________________________
399
4e25ac79 400Int_t AliITSQASPDDataMakerSim::GetTaskHisto(AliQAv1::TASKINDEX_t task) {
7a0e5776 401 // Returns the number of booked histograms for the selected task
402 Int_t histotot=0;
4e25ac79 403 if( task == AliQAv1::kHITS) {
7a0e5776 404 histotot=fSPDhHTask ;
405 }
4e25ac79 406 else if( task == AliQAv1::kSDIGITS) {
7a0e5776 407 histotot=fSPDhSTask;
408 }
4e25ac79 409 else if( task == AliQAv1::kDIGITS) {
7a0e5776 410 histotot=fSPDhDTask ;
411 }
412 else {
413 AliInfo("No task has been selected. TaskHisto set to zero.\n");
414 }
415 return histotot;
416
417}