]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerSim.cxx
Fixed Coverity issues
[u/mrichter/AliRoot.git] / FMD / AliFMDQADataMakerSim.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, 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 // --- ROOT system ---
16 #include <iostream>
17 #include <TClonesArray.h>
18 #include <TFile.h> 
19 #include <TH1F.h> 
20 #include <TH1I.h> 
21
22 // --- AliRoot header files ---
23 #include "AliESDEvent.h"
24 #include "AliLog.h"
25 #include "AliFMDQADataMakerSim.h"
26 #include "AliFMDDigit.h"
27 #include "AliFMDHit.h"
28 #include "AliQAChecker.h"
29 #include "AliFMDParameters.h"
30 #include "AliFMDSDigit.h"
31
32 //_____________________________________________________________________
33 // This is the class that collects the QA data for the FMD during simulation.
34 // The following data types are picked up:
35 // - hits
36 // - digits
37 // The following data types are not supported (yet):
38 // - raws
39 // - sdigits
40 // Author : Hans Hjersing Dalsgaard, Niels Bohr Institute, hans.dalsgaard@cern.ch
41 //_____________________________________________________________________
42
43 ClassImp(AliFMDQADataMakerSim)
44 #if 0
45 ; // This line is for Emacs - do not delete!
46 #endif
47 //_____________________________________________________________________
48 AliFMDQADataMakerSim::AliFMDQADataMakerSim() 
49   :  AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kFMD),
50                        "FMD Quality Assurance Data Maker")
51 {
52   // ctor
53
54 }
55
56 //_____________________________________________________________________
57 AliFMDQADataMakerSim::AliFMDQADataMakerSim(const AliFMDQADataMakerSim& /*qadm*/) 
58   : AliQADataMakerSim()
59 {
60   // copy ctor 
61   // 
62   // Parameters: 
63   //    qadm    Object to copy from
64   
65 }
66 //_____________________________________________________________________
67 AliFMDQADataMakerSim& 
68 AliFMDQADataMakerSim::operator = (const AliFMDQADataMakerSim& ) 
69 {
70   
71   return *this;
72 }
73 //_____________________________________________________________________
74 AliFMDQADataMakerSim::~AliFMDQADataMakerSim()
75 {
76
77 }
78
79 //_____________________________________________________________________
80 void AliFMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, 
81                                               TObjArray ** list)
82 {
83   //Detector specific actions at end of cycle
84   // do the QA checking
85   ResetEventTrigClasses(); // reset triggers list to select all histos
86   AliLog::Message(5,"FMD: end of detector cycle",
87                   "AliFMDQADataMakerSim","AliFMDQADataMakerSim",
88                   "AliFMDQADataMakerSim::EndOfDetectorCycle",
89                   "AliFMDQADataMakerSim.cxx",83);
90   AliQAChecker::Instance()->Run(AliQAv1::kFMD, task, list) ;  
91   
92 }
93 //_____________________________________________________________________
94 void AliFMDQADataMakerSim::InitSDigits()
95 {
96   // create SDigits histograms in SDigits subdir
97   const Bool_t expert   = kTRUE ; 
98   const Bool_t image    = kTRUE ; 
99   
100   TH1I* hADCCounts = new TH1I("hADCCounts",
101                               "Dist of ADC counts;ADC counts;Entries",
102                               1024,0,1024);
103   hADCCounts->SetXTitle("ADC counts");
104   Add2SDigitsList(hADCCounts, 0, !expert, image);
105   //
106   ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
107 }
108
109 //____________________________________________________________________ 
110 void AliFMDQADataMakerSim::InitHits()
111 {
112   // create Digits histograms in Digits subdir
113   const Bool_t expert   = kTRUE ; 
114   const Bool_t image    = kTRUE ; 
115   
116   TH1F* hEnergyOfHits = new TH1F("hEnergyOfHits",
117                                  "Energy distribution;Energy [MeV];Counts",
118                                  100,0,3);
119   hEnergyOfHits->SetXTitle("Edep");
120   hEnergyOfHits->SetYTitle("Counts");
121   Add2HitsList(hEnergyOfHits, 0, !expert, image);
122   //
123   ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
124 }
125
126 //_____________________________________________________________________
127 void AliFMDQADataMakerSim::InitDigits()
128 {
129   // create Digits histograms in Digits subdir
130   const Bool_t expert   = kTRUE ; 
131   const Bool_t image    = kTRUE ; 
132   
133   TH1I* hADCCounts = new TH1I("hADCCounts",
134                               "Dist of ADC counts; ADC counts;Entries",
135                               1024,0,1024);
136   hADCCounts->SetXTitle("ADC counts");
137   Add2DigitsList(hADCCounts, 0, !expert, image);
138   //
139   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
140 }
141
142 //_____________________________________________________________________
143 void AliFMDQADataMakerSim::MakeHits()
144 {
145   // Check id histograms already created for this Event Specie
146   if ( ! GetHitsData(0) )
147     InitHits() ;
148
149   TIter next(fHitsArray);
150   AliFMDHit * hit;
151   while ((hit = static_cast<AliFMDHit *>(next()))) FillHitsData(0,hit->Edep()/hit->Length()*0.032);
152   //
153 }
154
155 //_____________________________________________________________________
156 void AliFMDQADataMakerSim::MakeHits(TTree * hitTree)
157 {
158   // make QA data from Hit Tree
159   // 
160   // Parameters: 
161   //   hitTree    Hits container 
162   //
163   if (!fHitsArray) 
164     fHitsArray = new TClonesArray("AliFMDHit", 1000) ; 
165   fHitsArray->Clear() ; 
166   
167   TBranch * branch = hitTree->GetBranch("FMD") ;
168   if (!branch) {
169     AliWarning("FMD branch in Hit Tree not found") ; 
170     return;
171   }
172     
173   branch->SetAddress(&fHitsArray) ;
174   //  
175   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
176     branch->GetEntry(ientry);
177     MakeHits();   //tmp); 
178     fHitsArray->Clear() ; 
179   } 
180   //
181   IncEvCountCycleHits();
182   IncEvCountTotalHits();
183   //    
184 }
185
186 //_____________________________________________________________________
187 void AliFMDQADataMakerSim::MakeDigits()
188 {
189   // makes data from Digits
190   // 
191   // Parameters: 
192   //    none
193   if(!fDigitsArray) return;
194   
195   for(Int_t i = 0 ; i < fDigitsArray->GetEntriesFast() ; i++) {
196     //Raw ADC counts
197     AliFMDDigit* digit = static_cast<AliFMDDigit*>(fDigitsArray->At(i));
198     FillDigitsData(0,digit->Counts());
199   }
200   //
201 }
202
203 //_____________________________________________________________________
204 void AliFMDQADataMakerSim::MakeDigits(TTree * digitTree)
205 {
206   // Make data from digits. 
207   // 
208   // Parameters: 
209   //    digitTree    Tree holding digits. 
210   // RS: counters are incremented in MakeDigits()
211   
212   if (!fDigitsArray) 
213     fDigitsArray = new TClonesArray("AliFMDDigit", 1000) ; 
214   fDigitsArray->Clear();
215   
216   TBranch * branch = digitTree->GetBranch("FMD") ;
217   if (!branch)    {
218       AliWarning("FMD branch in Digit Tree not found") ; 
219       return;
220   } 
221   branch->SetAddress(&fDigitsArray) ;
222
223   if (fDigitsArray) fDigitsArray->Clear();
224
225   branch->GetEntry(0) ; 
226   MakeDigits() ; 
227   //
228   IncEvCountCycleDigits();
229   IncEvCountTotalDigits();
230   //
231 }
232
233 //_____________________________________________________________________
234 void AliFMDQADataMakerSim::MakeSDigits()
235 {
236   // makes data from Digits
237   // 
238   // Parameters: 
239   //   none 
240   if(!fSDigitsArray) return;
241   
242   for(Int_t i = 0 ; i < fSDigitsArray->GetEntriesFast() ; i++) {
243     //Raw ADC counts
244     AliFMDSDigit* sdigit = static_cast<AliFMDSDigit*>(fSDigitsArray->At(i));
245     FillSDigitsData(0,sdigit->Counts());
246   }
247   //
248 }
249
250 //_____________________________________________________________________
251 void AliFMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
252 {
253   // Make data from digits. 
254   // 
255   // Parameters: 
256   //    digitTree    Tree holding digits. 
257   //
258   if (!fSDigitsArray) 
259     fSDigitsArray = new TClonesArray("AliFMDSDigit", 1000) ; 
260   fSDigitsArray->Clear() ;
261
262   TBranch * branch = sdigitTree->GetBranch("FMD") ;
263   if (!branch)    {
264     AliWarning("FMD branch in SDigit Tree not found") ; 
265     return;
266   } 
267   branch->SetAddress(&fSDigitsArray) ;
268   branch->GetEntry(0) ; 
269   MakeSDigits() ; 
270   //
271   IncEvCountCycleSDigits();
272   IncEvCountTotalSDigits();
273   //
274 }
275
276 //_____________________________________________________________________ 
277 void AliFMDQADataMakerSim::StartOfDetectorCycle()
278 {
279   // Does 
280   // not 
281   // do 
282   // anything 
283 }
284 //_____________________________________________________________________ 
285 //
286 // EOF
287 //