]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerSim.cxx
remove obselete clean QA macros command
[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(AliQA::GetDetName(AliQA::kFMD),
50                        "FMD Quality Assurance Data Maker"),
51      fSDigitsArray("AliFMDSDigit", 1000),
52      fDigitsArray("AliFMDDigit", 1000),
53      fHitsArray("AliFMDHit", 10)
54 {
55   // ctor
56
57 }
58
59 //_____________________________________________________________________
60 AliFMDQADataMakerSim::AliFMDQADataMakerSim(const AliFMDQADataMakerSim& qadm) 
61   : AliQADataMakerSim(),
62     fSDigitsArray(qadm.fSDigitsArray),
63     fDigitsArray(qadm.fDigitsArray),
64     fHitsArray(qadm.fHitsArray)
65 {
66   //copy ctor 
67   
68   // Parameters: 
69   //    qadm    Object to copy from
70   
71 }
72 //_____________________________________________________________________
73 AliFMDQADataMakerSim& AliFMDQADataMakerSim::operator = (const AliFMDQADataMakerSim& qadm ) 
74 {
75   fSDigitsArray = qadm.fSDigitsArray;
76   fDigitsArray = qadm.fDigitsArray;
77   fHitsArray = qadm.fHitsArray;
78   
79   return *this;
80 }
81 //_____________________________________________________________________
82 AliFMDQADataMakerSim::~AliFMDQADataMakerSim()
83 {
84
85 }
86
87 //_____________________________________________________________________
88 void AliFMDQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, 
89                                               TObjArray ** list)
90 {
91   //Detector specific actions at end of cycle
92   // do the QA checking
93   AliLog::Message(5,"FMD: end of detector cycle",
94                   "AliFMDQADataMakerSim","AliFMDQADataMakerSim",
95                   "AliFMDQADataMakerSim::EndOfDetectorCycle",
96                   "AliFMDQADataMakerSim.cxx",83);
97   AliQAChecker::Instance()->Run(AliQA::kFMD, task, list) ;  
98   
99 }
100 //_____________________________________________________________________
101 void AliFMDQADataMakerSim::InitSDigits()
102 {
103   // create SDigits histograms in SDigits subdir
104   TH1I* hADCCounts = new TH1I("hADCCounts","Dist of ADC counts",1024,0,1024);
105   hADCCounts->SetXTitle("ADC counts");
106   Add2SDigitsList(hADCCounts, 0);
107 }
108
109 //____________________________________________________________________ 
110 void AliFMDQADataMakerSim::InitHits()
111 {
112   // create Digits histograms in Digits subdir
113   TH1F* hEnergyOfHits = new TH1F("hEnergyOfHits","Energy distribution",100,0,3);
114   hEnergyOfHits->SetXTitle("Edep");
115   hEnergyOfHits->SetYTitle("Counts");
116   Add2HitsList(hEnergyOfHits, 0);
117 }
118
119 //_____________________________________________________________________
120 void AliFMDQADataMakerSim::InitDigits()
121 {
122   // create Digits histograms in Digits subdir
123   TH1I* hADCCounts = new TH1I("hADCCounts","Dist of ADC counts",1024,0,1024);
124   hADCCounts->SetXTitle("ADC counts");
125   Add2DigitsList(hADCCounts, 0);
126 }
127
128 //_____________________________________________________________________
129 void AliFMDQADataMakerSim::MakeHits(TClonesArray * hits)
130 {
131   TIter next(hits);
132   AliFMDHit * hit;
133   while ((hit = static_cast<AliFMDHit *>(next()))) 
134     GetHitsData(0)->Fill(hit->Edep()/hit->Length()*0.032);
135 }
136
137 //_____________________________________________________________________
138 void AliFMDQADataMakerSim::MakeHits(TTree * hitTree)
139 {
140   // make QA data from Hit Tree
141   
142   fHitsArray.Clear();
143   
144   TBranch * branch = hitTree->GetBranch("FMD") ;
145   if (!branch) {
146     AliWarning("FMD branch in Hit Tree not found") ; 
147     return;
148   }
149   
150   TClonesArray* hitsAddress = &fHitsArray;
151   
152   branch->SetAddress(&hitsAddress) ;
153   
154   for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
155     branch->GetEntry(ientry);
156     MakeHits(hitsAddress);   //tmp); 
157   }     
158   
159 }
160
161 //_____________________________________________________________________
162 void AliFMDQADataMakerSim::MakeDigits(TClonesArray * digits)
163 {
164   // makes data from Digits
165   if(!digits) return;
166
167   for(Int_t i = 0 ; i < digits->GetEntriesFast() ; i++) {
168     //Raw ADC counts
169     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
170     GetDigitsData(0)->Fill(digit->Counts());
171   }
172 }
173
174 //_____________________________________________________________________
175 void AliFMDQADataMakerSim::MakeDigits(TTree * digitTree)
176 {
177   
178   fDigitsArray.Clear();
179   TBranch * branch = digitTree->GetBranch("FMD") ;
180   if (!branch)    {
181       AliWarning("FMD branch in Digit Tree not found") ; 
182       return;
183   } 
184   TClonesArray* digitAddress = &fDigitsArray;
185   branch->SetAddress(&digitAddress) ;
186   branch->GetEntry(0) ; 
187   MakeDigits(digitAddress) ; 
188 }
189
190 //_____________________________________________________________________
191 void AliFMDQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
192 {
193   // makes data from Digits
194   if(!sdigits) return;
195   
196   for(Int_t i = 0 ; i < sdigits->GetEntriesFast() ; i++) {
197     //Raw ADC counts
198     AliFMDSDigit* sdigit = static_cast<AliFMDSDigit*>(sdigits->At(i));
199     GetSDigitsData(0)->Fill(sdigit->Counts());
200   }
201 }
202
203 //_____________________________________________________________________
204 void AliFMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
205 {
206   
207   fSDigitsArray.Clear();
208   TBranch * branch = sdigitTree->GetBranch("FMD") ;
209   if (!branch)    {
210     AliWarning("FMD branch in SDigit Tree not found") ; 
211     return;
212   } 
213   TClonesArray* sdigitAddress = &fSDigitsArray;
214   branch->SetAddress(&sdigitAddress) ;
215   branch->GetEntry(0) ; 
216   MakeSDigits(sdigitAddress) ; 
217 }
218
219 //_____________________________________________________________________ 
220 void AliFMDQADataMakerSim::StartOfDetectorCycle()
221 {
222    
223 }
224 //_____________________________________________________________________ 
225 //
226 // EOF
227 //