]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASSDDataMakerSim.cxx
Transfer of the initialisation of the QA Data objects in the framework; clean the...
[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 "AliQAv1.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(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of SSD Cycle\n");
94 }
95
96 //____________________________________________________________________________ 
97 void AliITSQASSDDataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/) {
98   // launch the QA checking
99   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); 
100   
101 //  AliQAChecker::Instance()->Run( AliQAv1::kITS , task, list);
102 }
103
104 //____________________________________________________________________________ 
105 Int_t AliITSQASSDDataMakerSim::InitDigits() { 
106   // Initialization for DIGIT data - SSD -
107   const Bool_t expert   = kTRUE ; 
108   const Bool_t image    = kTRUE ; 
109   Int_t rv = 0 ; 
110  // fGenOffsetD = (fAliITSQADataMakerSim->fDigitsQAList[AliRecoParam::kDefault])->GetEntries();
111
112   // custom code here
113   TH1F *fHistSSDModule = new TH1F("fHistSSDDigitsModule",
114                                   "SSD Digits Module;SSD Module Number;N_{DIGITS}",
115                                   1698,499.5,2197.5);  
116   rv = fAliITSQADataMakerSim->Add2DigitsList(fHistSSDModule,
117                                         fGenOffsetD + 0, !expert, image);
118   fSSDhDTask += 1;
119   TH2F *fHistSSDModuleStrip = new TH2F("fHistSSDDigitsModuleStrip",
120                                        "SSD Digits Module Strip;N_{Strip};N_{Module}",
121                                        1540,0,1540,1698,499.5,2197.5);  
122   rv = fAliITSQADataMakerSim->Add2DigitsList(fHistSSDModuleStrip,
123                                         fGenOffsetD + 1, !expert, image);
124   fSSDhDTask += 1;
125
126   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD Digits histograms booked\n",fSSDhDTask));
127   return rv ; 
128 }
129
130 //____________________________________________________________________________
131 Int_t AliITSQASSDDataMakerSim::MakeDigits(TTree *digits) { 
132   // Fill QA for DIGIT - SSD -
133   Int_t rv = 0 ; 
134
135   AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");
136   fITS->SetTreeAddress();
137   TClonesArray *iSSDdigits  = fITS->DigitsAddress(2);
138   for(Int_t iModule = 500; iModule < 2198; iModule++) {
139     iSSDdigits->Clear();
140     digits->GetEvent(iModule);    
141     Int_t ndigits = iSSDdigits->GetEntries();
142     fAliITSQADataMakerSim->GetDigitsData(fGenOffsetD + 0)->Fill(iModule,ndigits);
143     if(ndigits != 0)
144       AliDebug(AliQAv1::GetQADebugLevel(),Form("Module: %d - Digits: %d",iModule,ndigits));
145  
146     for (Int_t iDigit = 0; iDigit < ndigits; iDigit++) {
147       AliITSdigit *dig = (AliITSdigit*)iSSDdigits->UncheckedAt(iDigit);
148       Int_t fStripNumber = (dig->GetCoord1() == 0) ? dig->GetCoord2() : dig->GetCoord2() + fgkNumberOfPSideStrips;
149       ((TH2F *)fAliITSQADataMakerSim->GetDigitsData(fGenOffsetD + 1))->Fill(fStripNumber,iModule,dig->GetSignal());
150     }//digit loop
151   }//module loop
152   return rv ; 
153 }
154
155 //____________________________________________________________________________ 
156 Int_t AliITSQASSDDataMakerSim::InitSDigits() { 
157   // Initialization for SDIGIT data - SSD -
158   const Bool_t expert   = kTRUE ; 
159   const Bool_t image    = kTRUE ; 
160   Int_t rv = 0 ; 
161   //fGenOffsetS = (fAliITSQADataMakerSim->fSDigitsQAList[AliRecoParam::kDefault])->GetEntries();
162
163   // custom code here
164   TH1F *fHistSSDModule = new TH1F("fHistSSDSDigitsModule",
165                                   "SSD SDigits Module;SSD Module Number;N_{SDIGITS}",
166                                   1698,499.5,2197.5);  
167   rv = fAliITSQADataMakerSim->Add2SDigitsList(fHistSSDModule,
168                                          fGenOffsetS + 0, !expert, image);
169   fSSDhSTask += 1;  
170
171   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD SDigits histograms booked\n",fSSDhSTask));
172   return rv ; 
173 }
174
175 //____________________________________________________________________________
176 Int_t AliITSQASSDDataMakerSim::MakeSDigits(TTree *sdigits) { 
177   // Fill QA for SDIGIT - SSD -
178   Int_t rv = 0 ; 
179  
180   static TClonesArray iSSDEmpty("AliITSpListItem",10000);
181   iSSDEmpty.Clear();
182   TClonesArray *iSSDsdigits = &iSSDEmpty;
183
184   AliDebug(AliQAv1::GetQADebugLevel(), Form("Trying to access the sdigits histogram: %d\n",fGenOffsetS));
185
186   TBranch *brchSDigits = sdigits->GetBranch("ITS");
187   brchSDigits->SetAddress(&iSSDsdigits);
188   for(Int_t iModule = 500; iModule < 2198; iModule++) {
189     iSSDsdigits->Clear();
190     sdigits->GetEvent(iModule);    
191     Int_t ndigits = iSSDsdigits->GetEntries();
192     fAliITSQADataMakerSim->GetSDigitsData(fGenOffsetS + 0)->Fill(iModule,ndigits);
193     if(ndigits != 0)
194       AliDebug(AliQAv1::GetQADebugLevel(),Form("Module: %d - Digits: %d",iModule,ndigits));
195
196     for (Int_t iDigit = 0; iDigit < ndigits; iDigit++) {
197       AliITSpListItem *dig=(AliITSpListItem*)iSSDsdigits->At(iDigit);
198       dig=0;
199     }//digit loop
200   }//module loop
201   return rv ; 
202 }
203
204 //____________________________________________________________________________ 
205 Int_t AliITSQASSDDataMakerSim::InitHits() { 
206   // Initialization for HITS data - SSD -
207   const Bool_t expert   = kTRUE ; 
208   const Bool_t image    = kTRUE ; 
209   Int_t rv = 0 ; 
210
211   //fGenOffsetH = (fAliITSQADataMakerSim->fHitsQAList[fEventSpecie])->GetEntries();
212
213   // custom code here
214   TH1F *fHistSSDModule = new TH1F("fHistSSDHitsModule",
215                                   "SSD Hits Module;SDD Module Number;N_{HITS}",
216                                   1698,499.5,2197.5); 
217   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDModule,
218                                       fGenOffsetH + 0, !expert, image);
219   fSSDhHTask += 1;
220   TH1F *fHistSSDGlobalX = new TH1F("fHistSSDHitsGlobalX",
221                                    "SSD Hits Global X;x [cm];Entries",
222                                    1000,-50.,50.);
223   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalX,
224                                       fGenOffsetH + 1, !expert, image);
225   fSSDhHTask += 1;
226   TH1F *fHistSSDGlobalY = new TH1F("fHistSSDHitsGlobalY",
227                                    "SSD Hits Global Y;y [cm];Entries",
228                                    1000,-50.,50.);
229   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalY,
230                                       fGenOffsetH + 2, !expert, image);
231   fSSDhHTask += 1;
232   TH1F *fHistSSDGlobalZ = new TH1F("fHistSSDHitsGlobalZ",
233                                    "SSD Hits Global Z ;z [cm];Entries",
234                                    1000,-60.,60.);
235   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalZ,
236                                       fGenOffsetH + 3, !expert, image);
237   fSSDhHTask += 1;
238   TH1F *fHistSSDLocalX = new TH1F("fHistSSDHitsLocalX",
239                                   "SSD Hits Local X;x [cm];Entries",
240                                   1000,-4.,4.);
241   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDLocalX,
242                                       fGenOffsetH + 4, !expert, image);
243   fSSDhHTask += 1;
244   TH1F *fHistSSDLocalY = new TH1F("fHistSSDHitsLocalY",
245                                   "SSD Hits Local Y;y [cm];Entries",
246                                   1000,-0.1,0.1);
247   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDLocalY,
248                                       fGenOffsetH + 5, !expert, image);
249   fSSDhHTask += 1;
250   TH1F *fHistSSDLocalZ = new TH1F("fHistSSDHitsLocalZ",
251                                   "SSD Hits Local Z;z [cm];Entries",
252                                   1000,-4.,4.);
253   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDLocalZ,
254                                       fGenOffsetH + 6, !expert, image);
255   fSSDhHTask += 1;
256   TH1F *fHistSSDIonization = new TH1F("fHistSSDHitsIonization",
257                                       "SSD Hits Ionization;log(dE/dx) [KeV];N_{Hits}",
258                                       100,-7,-2);
259   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDIonization,
260                                       fGenOffsetH + 7, !expert, image);
261   fSSDhHTask += 1;
262   TH2F *fHistSSDGlobalXY = new TH2F("fHistSSDHitsGlobalXY",
263                                     "SSD Hits Global XY;x [cm];y [cm]",
264                                     1000,-50.,50.,
265                                     1000,-50.,50.);
266   rv = fAliITSQADataMakerSim->Add2HitsList(fHistSSDGlobalXY,
267                                       fGenOffsetH + 8, !expert, image);
268   fSSDhHTask += 1;
269  
270   AliDebug(AliQAv1::GetQADebugLevel(),Form("%d SSD Hits histograms booked\n",fSSDhHTask));
271   return rv ; 
272 }
273
274
275 //____________________________________________________________________________
276 Int_t AliITSQASSDDataMakerSim::MakeHits(TTree *hits) { 
277   // Fill QA for HITS - SSD -
278   Int_t rv = 0 ; 
279  
280   AliITS *fITS  = (AliITS*)gAlice->GetModule("ITS");
281   fITS->SetTreeAddress();
282   Int_t nmodules;
283   fITS->InitModules(-1,nmodules);
284   fITS->FillModules(hits,0);
285   for(Int_t iModule = 500; iModule < 2198; iModule++) {
286     AliITSmodule *module = fITS->GetModule(iModule);
287     TObjArray *arrHits = module->GetHits();
288     Int_t nhits = arrHits->GetEntriesFast();
289     if(nhits != 0)
290       AliDebug(AliQAv1::GetQADebugLevel(),Form("Module: %d - Hits: %d",iModule,nhits));
291     for (Int_t iHit = 0; iHit < nhits; iHit++) {
292       AliITShit *hit = (AliITShit*) arrHits->At(iHit);
293       
294       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 0)->Fill(iModule);
295       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 1)->Fill(hit->GetXG());
296       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 2)->Fill(hit->GetYG());
297       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 3)->Fill(hit->GetZG());
298       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 4)->Fill(hit->GetXL());
299       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 5)->Fill(hit->GetYL());
300       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 6)->Fill(hit->GetZL());
301       if(hit->GetIonization())
302         fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 7)->Fill(TMath::Log10(hit->GetIonization()));
303       fAliITSQADataMakerSim->GetHitsData(fGenOffsetH + 8)->Fill(hit->GetXG(),hit->GetYG());
304     }//hit loop
305   }//module loop  
306   return rv ; 
307 }
308
309 //____________________________________________________________________________ 
310 Int_t AliITSQASSDDataMakerSim::GetOffset(AliQAv1::TASKINDEX_t task){
311   // Returns histogram offset according to the specified task
312   Int_t offset=0;
313   if( task == AliQAv1::kHITS){
314     offset=fGenOffsetH;  
315   }
316   else if( task == AliQAv1::kSDIGITS) {
317     offset=fGenOffsetS;   
318   }
319   else if( task == AliQAv1::kDIGITS) {
320     offset=fGenOffsetD;   
321   }
322   else {
323     AliInfo("No task has been selected. TaskHisto set to zero.\n");
324   }
325
326   return offset;
327 }
328
329
330 //____________________________________________________________________________ 
331 void AliITSQASSDDataMakerSim::SetOffset(AliQAv1::TASKINDEX_t task, Int_t offset){
332   // Returns histogram offset according to the specified task
333   if( task == AliQAv1::kHITS){
334     fGenOffsetH = offset;  
335   }
336   else if( task == AliQAv1::kSDIGITS) {
337     fGenOffsetS = offset;   
338   }
339   else if( task == AliQAv1::kDIGITS) {
340     fGenOffsetD = offset;   
341   }
342   else {
343     AliInfo("No task has been selected. TaskHisto set to zero.\n");
344   }
345 }
346
347 //____________________________________________________________________________ 
348 Int_t AliITSQASSDDataMakerSim::GetTaskHisto(AliQAv1::TASKINDEX_t task) {
349   // Returns the number of booked histograms for the selected task
350   Int_t histotot=0;
351   if( task == AliQAv1::kHITS) {
352     histotot=fSSDhHTask ;  
353   }
354   else if( task == AliQAv1::kSDIGITS) {
355     histotot=fSSDhSTask;   
356   }
357   else if( task == AliQAv1::kDIGITS) {
358     histotot=fSSDhDTask ;   
359   }
360   else {
361     AliInfo("No task has been selected. TaskHisto set to zero.\n");
362   }
363   return histotot;
364
365 }