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