1 /**************************************************************************
2 * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 // *************************************************************
18 // Checks the quality assurance
19 // by comparing with reference data
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
26 // --- ROOT system ---
30 // --- Standard library ---
32 // --- AliRoot header files ---
34 #include "AliITSQADataMakerSim.h"
35 #include "AliITSQASPDDataMakerSim.h"
37 #include "AliQAChecker.h"
38 #include "AliITSdigit.h"
39 #include "AliITSdigitSPD.h"
41 #include "AliITSmodule.h"
42 #include "AliITShit.h"
43 #include "AliITSLoader.h"
44 #include "AliRunLoader.h"
46 ClassImp(AliITSQASPDDataMakerSim)
48 //____________________________________________________________________________
49 AliITSQASPDDataMakerSim::AliITSQASPDDataMakerSim(AliITSQADataMakerSim *aliITSQADataMakerSim) :
51 fAliITSQADataMakerSim(aliITSQADataMakerSim),
59 //ctor used to discriminate OnLine-Offline analysis
60 fGenOffsetH= new Int_t[AliRecoParam::kNSpecies];
61 fGenOffsetS= new Int_t[AliRecoParam::kNSpecies];
62 fGenOffsetD= new Int_t[AliRecoParam::kNSpecies];
63 for(Int_t i=0; i<AliRecoParam::kNSpecies; i++) {
70 //____________________________________________________________________________
71 AliITSQASPDDataMakerSim::AliITSQASPDDataMakerSim(const AliITSQASPDDataMakerSim& qadm) :
73 fAliITSQADataMakerSim(qadm.fAliITSQADataMakerSim),
74 fSPDhHTask(qadm.fSPDhHTask),
75 fSPDhSTask(qadm.fSPDhSTask),
76 fSPDhDTask(qadm.fSPDhDTask),
77 fGenOffsetH(qadm.fGenOffsetH),
78 fGenOffsetS(qadm.fGenOffsetS),
79 fGenOffsetD(qadm.fGenOffsetD)
82 fAliITSQADataMakerSim->SetName((const char*)qadm.fAliITSQADataMakerSim->GetName()) ;
83 fAliITSQADataMakerSim->SetTitle((const char*)qadm.fAliITSQADataMakerSim->GetTitle());
87 //__________________________________________________________________
88 AliITSQASPDDataMakerSim& AliITSQASPDDataMakerSim::operator = (const AliITSQASPDDataMakerSim& qac )
91 this->~AliITSQASPDDataMakerSim();
92 new(this) AliITSQASPDDataMakerSim(qac);
96 //____________________________________________________________________________
97 void AliITSQASPDDataMakerSim::StartOfDetectorCycle()
99 //Detector specific actions at start of cycle
100 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SPD Cycle\n");
103 //____________________________________________________________________________
104 void AliITSQASPDDataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
106 // launch the QA checking
107 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n");
109 //AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);
112 //____________________________________________________________________________
113 Int_t AliITSQASPDDataMakerSim::InitDigits()
115 // Initialization for DIGIT data - SPD -
116 const Bool_t expert = kTRUE ;
117 const Bool_t image = kTRUE ;
119 //fGenOffsetD = (fAliITSQADataMakerSim->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
120 //fSPDhDTask must be incremented by one unit every time a histogram is ADDED to the QA List
125 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
126 hlayer->GetXaxis()->SetTitle("Layer number");
127 hlayer->GetYaxis()->SetTitle("Entries");
128 rv = fAliITSQADataMakerSim->Add2DigitsList(hlayer,fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], expert, !image);
131 TH1F **hmod = new TH1F*[2];
132 for (Int_t iLay=0; iLay<2; iLay++) {
133 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
134 sprintf(title,"Module map - SPD Layer %d",iLay+1);
135 hmod[iLay]=new TH1F(name,title,240,0,240);
136 hmod[iLay]->GetXaxis()->SetTitle("Module number");
137 hmod[iLay]->GetYaxis()->SetTitle("Entries");
138 rv = fAliITSQADataMakerSim->Add2DigitsList(hmod[iLay],1+iLay+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
142 TH1F *hcolumns = new TH1F("SPDColumns_SPD","Columns - SPD",160,0.,160.);
143 hcolumns->GetXaxis()->SetTitle("Column number");
144 hcolumns->GetYaxis()->SetTitle("Entries");
145 fAliITSQADataMakerSim->Add2DigitsList(hcolumns,3+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], expert, !image);
148 TH1F *hrows = new TH1F("SPDRows_SPD","Rows - SPD",256,0.,256.);
149 hrows->GetXaxis()->SetTitle("Row number");
150 hrows->GetYaxis()->SetTitle("Entries");
151 rv = fAliITSQADataMakerSim->Add2DigitsList(hrows,4+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], expert, !image);
154 TH1F** hMultSPDdigits = new TH1F*[2];
155 for (Int_t iLay=0; iLay<2; ++iLay) {
156 sprintf(name,"SPDDigitMultiplicity_SPD%d",iLay+1);
157 sprintf(title,"Digit multiplicity - SPD Layer %d",iLay+1);
158 hMultSPDdigits[iLay]=new TH1F(name,title,200,0.,200.);
159 hMultSPDdigits[iLay]->GetXaxis()->SetTitle("Digit multiplicity");
160 hMultSPDdigits[iLay]->GetYaxis()->SetTitle("Entries");
161 rv = fAliITSQADataMakerSim->Add2DigitsList(hMultSPDdigits[iLay], 5+iLay+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
165 TH2F *hMultSPDdig2MultSPDdig1
166 = new TH2F("SPDDigitMultCorrelation_SPD","Digit multiplicity correlation - SPD",200,0.,200.,200,0.,200.);
167 hMultSPDdig2MultSPDdig1->GetXaxis()->SetTitle("Digit multiplicity (Layer 1)");
168 hMultSPDdig2MultSPDdig1->GetYaxis()->SetTitle("Digit multiplicity (Layer 2)");
169 rv = fAliITSQADataMakerSim->Add2DigitsList(hMultSPDdig2MultSPDdig1,7+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
172 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Digits histograms booked\n",fSPDhDTask));
176 //____________________________________________________________________________
177 Int_t AliITSQASPDDataMakerSim::MakeDigits(TTree *digits)
179 // Fill QA for DIGIT - SPD -
182 AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
183 fITS->SetTreeAddress();
184 TClonesArray *iITSdigits = fITS->DigitsAddress(0); // 0->SPD
189 for (Int_t imod=0; imod<240; ++imod){
190 digits->GetEvent(imod);
191 Int_t ndigits = iITSdigits->GetEntries();
193 fAliITSQADataMakerSim->GetDigitsData(0+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(0.5,ndigits);
194 fAliITSQADataMakerSim->GetDigitsData(1+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(imod,ndigits);
198 fAliITSQADataMakerSim->GetDigitsData(0+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(1,ndigits);
199 fAliITSQADataMakerSim->GetDigitsData(2+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(imod,ndigits);
202 for (Int_t idig=0; idig<ndigits; ++idig) {
203 AliITSdigit *dig=(AliITSdigit*)iITSdigits->UncheckedAt(idig);
204 Int_t col=dig->GetCoord1(); // cell number z
205 Int_t row=dig->GetCoord2(); // cell number x
206 fAliITSQADataMakerSim->GetDigitsData(3+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(col);
207 fAliITSQADataMakerSim->GetDigitsData(4+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(row);
210 fAliITSQADataMakerSim->GetDigitsData(5+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(nDigitsL1);
211 fAliITSQADataMakerSim->GetDigitsData(6+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(nDigitsL2);
212 fAliITSQADataMakerSim->GetDigitsData(7+fGenOffsetD[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(nDigitsL1,nDigitsL2);
216 //____________________________________________________________________________
217 Int_t AliITSQASPDDataMakerSim::InitSDigits()
219 // Initialization for SDIGIT data - SPD -
220 const Bool_t expert = kTRUE ;
221 const Bool_t image = kTRUE ;
223 //fGenOffsetS = (fAliITSQADataMakerSim->fSDigitsQAList[AliRecoParam::kDefault])->GetEntries();
224 //printf("--W-- AliITSQASPDDataMakerSim::InitSDigits() fGenOffset= %d \n",fGenOffset);
225 //fSPDhSTask must be incremented by one unit every time a histogram is ADDED to the QA List
230 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
231 hlayer->GetXaxis()->SetTitle("Layer number");
232 hlayer->GetYaxis()->SetTitle("Entries");
233 rv = fAliITSQADataMakerSim->Add2SDigitsList(hlayer,fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()], expert, !image);
236 TH1F **hmod = new TH1F*[2];
237 for (Int_t iLay=0; iLay<2; ++iLay) {
238 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
239 sprintf(title,"Module map - SPD Layer %d",iLay+1);
240 hmod[iLay]=new TH1F(name,title,240,0,240);
241 hmod[iLay]->GetXaxis()->SetTitle("Module number");
242 hmod[iLay]->GetYaxis()->SetTitle("Entries");
243 rv = fAliITSQADataMakerSim->Add2SDigitsList(hmod[iLay],1+iLay+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
247 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD SDigits histograms booked\n",fSPDhSTask));
251 //____________________________________________________________________________
252 Int_t AliITSQASPDDataMakerSim::MakeSDigits(TTree *sdigits)
254 // Fill QA for SDIGIT - SPD -
257 static TClonesArray * sdig ;
259 sdig = new TClonesArray( "AliITSpListItem",1000 );
261 TBranch *brchSDigits = sdigits->GetBranch("ITS");
262 for (Int_t imod=0; imod<240; ++imod){
263 brchSDigits->SetAddress( &sdig );
264 brchSDigits->GetEvent(imod);
265 Int_t nsdig=sdig->GetEntries();
267 fAliITSQADataMakerSim->GetSDigitsData(0+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(0.5,nsdig);
268 fAliITSQADataMakerSim->GetSDigitsData(1+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(imod,nsdig);
271 fAliITSQADataMakerSim->GetSDigitsData(0+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(1,nsdig);
272 fAliITSQADataMakerSim->GetSDigitsData(2+fGenOffsetS[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(imod,nsdig);
279 //____________________________________________________________________________
280 Int_t AliITSQASPDDataMakerSim::InitHits()
282 // Initialization for HITS data - SPD -
283 const Bool_t expert = kTRUE ;
284 const Bool_t image = kTRUE ;
287 //fGenOffsetH = (fAliITSQADataMakerSim->fHitsQAList[AliRecoParam::kDefault])->GetEntries();
288 //printf("--W-- AliITSQASPDDataMakerSim::InitHits() fGenOffset= %d \n",fGenOffset);
289 //fSPDhHTask must be incremented by one unit every time a histogram is ADDED to the QA List
293 TH1F *hlayer = new TH1F("SPDLayPattern_SPD","Layer map - SPD",6,0.,6.);
294 hlayer->GetXaxis()->SetTitle("Layer number");
295 hlayer->GetYaxis()->SetTitle("Entries");
296 rv = fAliITSQADataMakerSim->Add2HitsList(hlayer,fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], expert, !image);
299 TH1F **hmod = new TH1F*[2];
300 for (Int_t iLay=0; iLay<2; ++iLay) {
301 sprintf(name,"SPDModPattern_SPD%d",iLay+1);
302 sprintf(title,"Module map - SPD Layer %d",iLay+1);
303 hmod[iLay]=new TH1F(name,title,240,0,240);
304 hmod[iLay]->GetXaxis()->SetTitle("Module number");
305 hmod[iLay]->GetYaxis()->SetTitle("Entries");
306 rv = fAliITSQADataMakerSim->Add2HitsList(hmod[iLay],1+iLay+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
310 TH1F *hhitlenght = new TH1F("SPDLenght_SPD","SPD Hit lenght along y_{loc} coord",210,0.,210.);
311 hhitlenght->GetXaxis()->SetTitle("Hit lenght [#mum]");
312 hhitlenght->GetYaxis()->SetTitle("# hits");
313 rv = fAliITSQADataMakerSim->Add2HitsList(hhitlenght,3+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
316 TH1F *hEdepos = new TH1F("SPDEnergyDeposit_SPD","SPD Deposited energy distribution (y_{loc}>180 #mum)",150,0.,300.);
317 hEdepos->GetXaxis()->SetTitle("Deposited energy [keV]");
318 hEdepos->GetYaxis()->SetTitle("# hits");
319 rv = fAliITSQADataMakerSim->Add2HitsList(hEdepos,4+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()], !expert, image);
322 AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SPD Hits histograms booked\n",fSPDhHTask));
326 //____________________________________________________________________________
327 Int_t AliITSQASPDDataMakerSim::MakeHits(TTree *hits)
329 // Fill QA for HITS - SPD -
332 AliITS *fITS = (AliITS*)gAlice->GetModule("ITS");
333 fITS->SetTreeAddress();
335 fITS->InitModules(-1,nmodules); //-1->number of modules taken from AliITSgeom class kept in fITSgeom
338 fITS->FillModules(hits,0);
340 for (Int_t imod=0; imod<240; ++imod){
341 AliITSmodule *module = fITS->GetModule(imod);
342 TObjArray *arrHits = module->GetHits();
343 Int_t nhits = arrHits->GetEntriesFast();
345 fAliITSQADataMakerSim->GetHitsData(fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(0.5,nhits);
346 fAliITSQADataMakerSim->GetHitsData(1+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(imod,nhits);
348 fAliITSQADataMakerSim->GetHitsData(fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(1,nhits);
349 fAliITSQADataMakerSim->GetHitsData(2+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(imod,nhits);
351 for (Int_t iHit=0; iHit<nhits; ++iHit) {
352 AliITShit *hit = (AliITShit*) arrHits->At(iHit);
353 Double_t xl,yl,zl,xl0,yl0,zl0;
355 hit->GetPositionL(xl,yl,zl,tof);
356 hit->GetPositionL0(xl0,yl0,zl0,tof0);
357 Float_t dyloc=TMath::Abs(yl-yl0)*10000.;
358 fAliITSQADataMakerSim->GetHitsData(3+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(dyloc);
359 Float_t edep=hit->GetIonization()*1000000;
361 fAliITSQADataMakerSim->GetHitsData(4+fGenOffsetH[fAliITSQADataMakerSim->GetEventSpecie()])->Fill(edep);
369 //_______________________________________________________________
371 Int_t AliITSQASPDDataMakerSim::GetOffset(AliQAv1::TASKINDEX_t task,Int_t specie){
372 // Returns histogram offset according to the specified task
374 if( task == AliQAv1::kHITS){
375 offset=fGenOffsetH[specie];
377 else if( task == AliQAv1::kSDIGITS) {
378 offset=fGenOffsetS[specie];
380 else if( task == AliQAv1::kDIGITS) {
381 offset=fGenOffsetD[specie];
384 AliInfo("No task has been selected. TaskHisto set to zero.\n");
390 //____________________________________________________________________________
391 void AliITSQASPDDataMakerSim::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset,Int_t specie ){
392 // Returns histogram offset according to the specified task
393 if( task == AliQAv1::kHITS){
394 fGenOffsetH[specie] = offset;
396 else if( task == AliQAv1::kSDIGITS) {
397 fGenOffsetS[specie] = offset;
399 else if( task == AliQAv1::kDIGITS) {
400 fGenOffsetD[specie] = offset;
403 AliInfo("No task has been selected. TaskHisto set to zero.\n");
407 //_______________________________________________________________
409 Int_t AliITSQASPDDataMakerSim::GetTaskHisto(AliQAv1::TASKINDEX_t task) {
410 // Returns the number of booked histograms for the selected task
412 if( task == AliQAv1::kHITS) {
413 histotot=fSPDhHTask ;
415 else if( task == AliQAv1::kSDIGITS) {
418 else if( task == AliQAv1::kDIGITS) {
419 histotot=fSPDhDTask ;
422 AliInfo("No task has been selected. TaskHisto set to zero.\n");