]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliEMCALQADataMaker.cxx
comment out signal->Reset(), at least for now
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALQADataMaker.cxx
CommitLineData
94594e5d 1/**************************************************************************
2 * Copyright(c) 1998-1999, 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/*
17 Produces the data needed to calculate the quality assurance.
18 All data must be mergeable objects.
19
20 Based on PHOS code written by
21 Y. Schutz CERN July 2007
22*/
23
24// --- ROOT system ---
25#include <TClonesArray.h>
26#include <TFile.h>
27#include <TH1F.h>
28#include <TH1I.h>
29#include <TH2F.h>
30
31// --- Standard library ---
32
33// --- AliRoot header files ---
34#include "AliESDCaloCluster.h"
35#include "AliESDEvent.h"
36#include "AliLog.h"
37#include "AliEMCALDigit.h"
38#include "AliEMCALHit.h"
39#include "AliEMCALQADataMaker.h"
40#include "AliQAChecker.h"
41#include "AliEMCALRecPoint.h"
42#include "AliEMCALRawUtils.h"
43#include "AliEMCALReconstructor.h"
44#include "AliEMCALRecParam.h"
45
46ClassImp(AliEMCALQADataMaker)
47
48//____________________________________________________________________________
49 AliEMCALQADataMaker::AliEMCALQADataMaker() :
50 AliQADataMaker(AliQA::GetDetName(AliQA::kEMCAL), "EMCAL Quality Assurance Data Maker")
51{
52 // ctor
53}
54
55//____________________________________________________________________________
56AliEMCALQADataMaker::AliEMCALQADataMaker(const AliEMCALQADataMaker& qadm) :
57 AliQADataMaker()
58{
59 //copy ctor
60 SetName((const char*)qadm.GetName()) ;
61 SetTitle((const char*)qadm.GetTitle());
62}
63
64//__________________________________________________________________
65AliEMCALQADataMaker& AliEMCALQADataMaker::operator = (const AliEMCALQADataMaker& qadm )
66{
67 // Equal operator.
68 this->~AliEMCALQADataMaker();
69 new(this) AliEMCALQADataMaker(qadm);
70 return *this;
71}
72
73//____________________________________________________________________________
74void AliEMCALQADataMaker::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray * list)
75{
76 //Detector specific actions at end of cycle
77 // do the QA checking
78 AliQAChecker::Instance()->Run(AliQA::kEMCAL, task, list) ;
79}
80
81//____________________________________________________________________________
82void AliEMCALQADataMaker::InitESDs()
83{
84 //create ESDs histograms in ESDs subdir
85
86 TH1F * h0 = new TH1F("hEmcalESDs", "ESDs energy distribution in EMCAL", 100, 0., 100.) ;
87 h0->Sumw2() ;
88 Add2ESDsList(h0, 0) ;
89 TH1I * h1 = new TH1I("hEmcalESDsMul", "ESDs multiplicity distribution in EMCAL", 100, 0., 100) ;
90 h1->Sumw2() ;
91 Add2ESDsList(h1, 1) ;
92}
93
94//____________________________________________________________________________
95void AliEMCALQADataMaker::InitHits()
96{
97 // create Hits histograms in Hits subdir
98 TH1F * h0 = new TH1F("hEmcalHits", "Hits energy distribution in EMCAL", 100, 0., 100.) ;
99 h0->Sumw2() ;
100 Add2HitsList(h0, 0) ;
101 TH1I * h1 = new TH1I("hEmcalHitsMul", "Hits multiplicity distribution in EMCAL", 500, 0., 10000) ;
102 h1->Sumw2() ;
103 Add2HitsList(h1, 1) ;
104}
105
106//____________________________________________________________________________
107void AliEMCALQADataMaker::InitDigits()
108{
109 // create Digits histograms in Digits subdir
110 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL", 500, 0, 5000) ;
111 h0->Sumw2() ;
112 Add2DigitsList(h0, 0) ;
113 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL", 500, 0, 1000) ;
114 h1->Sumw2() ;
115 Add2DigitsList(h1, 1) ;
116}
117
118//____________________________________________________________________________
119void AliEMCALQADataMaker::InitRecPoints()
120{
121 // create Reconstructed Points histograms in RecPoints subdir
122 TH1F * h0 = new TH1F("hEmcalRecPoints", "RecPoints energy distribution in EMCAL", 100, 0., 100.) ;
123 h0->Sumw2() ;
124 Add2RecPointsList(h0, 0) ;
125 TH1I * h1 = new TH1I("hEmcalRecPointsMul", "RecPoints multiplicity distribution in EMCAL", 100, 0, 100) ;
126 h1->Sumw2() ;
127 Add2RecPointsList(h1, 1) ;
128
129}
130
131//____________________________________________________________________________
132void AliEMCALQADataMaker::InitRaws()
133{
134 AliInfo(Form("Raw QA infor for EMCAL not yet implemented"));
135}
136
137//____________________________________________________________________________
138void AliEMCALQADataMaker::InitSDigits()
139{
140 // create SDigits histograms in SDigits subdir
141 TH1F * h0 = new TH1F("hEmcalSDigits", "SDigits energy distribution in EMCAL", 100, 0., 100.) ;
142 h0->Sumw2() ;
143 Add2SDigitsList(h0, 0) ;
144 TH1I * h1 = new TH1I("hEmcalSDigitsMul", "SDigits multiplicity distribution in EMCAL", 500, 0, 10000) ;
145 h1->Sumw2() ;
146 Add2SDigitsList(h1, 1) ;
147}
148
149//____________________________________________________________________________
150void AliEMCALQADataMaker::MakeESDs(AliESDEvent * esd)
151{
152 // make QA data from ESDs
153
154 Int_t count = 0 ;
155 for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
156 AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
157 if ( clu->IsEMCAL() ) {
158 GetESDsData(0)->Fill(clu->E()) ;
159 count++ ;
160 }
161 }
162 GetESDsData(1)->Fill(count) ;
163}
164
165//____________________________________________________________________________
166void AliEMCALQADataMaker::MakeHits(TClonesArray * hits)
167{
168 //make QA data from Hits
169
170 GetHitsData(1)->Fill(hits->GetEntriesFast()) ;
171 TIter next(hits) ;
172 AliEMCALHit * hit ;
173 while ( (hit = dynamic_cast<AliEMCALHit *>(next())) ) {
174 GetHitsData(0)->Fill( hit->GetEnergy()) ;
175 }
176}
177
178//____________________________________________________________________________
179void AliEMCALQADataMaker::MakeHits(TTree * hitTree)
180{
181 // make QA data from Hit Tree
182
183 TClonesArray * hits = new TClonesArray("AliEMCALHit", 1000);
184
185 TBranch * branch = hitTree->GetBranch("EMCAL") ;
186 if ( ! branch ) {
187 AliWarning("EMCAL branch in Hit Tree not found") ;
188 } else {
189 TClonesArray * tmp = new TClonesArray("AliEMCALHit", 1000) ;
190 branch->SetAddress(&tmp) ;
191 Int_t index = 0 ;
192 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
193 branch->GetEntry(ientry) ;
194 for (Int_t ihit = 0 ; ihit < tmp->GetEntries() ; ihit++) {
195 AliEMCALHit * hit = dynamic_cast<AliEMCALHit *> (tmp->At(ihit)) ;
196 new((*hits)[index]) AliEMCALHit(*hit) ;
197 index++ ;
198 }
199 }
200 tmp->Delete() ;
201 delete tmp ;
202 MakeHits(hits) ;
203 }
204}
205
206//____________________________________________________________________________
207void AliEMCALQADataMaker::MakeDigits(TClonesArray * digits)
208{
209 // makes data from Digits
210
211 GetDigitsData(1)->Fill(digits->GetEntriesFast()) ;
212 TIter next(digits) ;
213 AliEMCALDigit * digit ;
214 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
215 GetDigitsData(0)->Fill( digit->GetAmp()) ;
216 }
217}
218
219//____________________________________________________________________________
220void AliEMCALQADataMaker::MakeDigits(TTree * digitTree)
221{
222 // makes data from Digit Tree
223 TClonesArray * digits = new TClonesArray("AliEMCALDigit", 1000) ;
224
225 TBranch * branch = digitTree->GetBranch("EMCAL") ;
226 if ( ! branch ) {
227 AliWarning("EMCAL branch in Digit Tree not found") ;
228 } else {
229 branch->SetAddress(&digits) ;
230 branch->GetEntry(0) ;
231 MakeDigits(digits) ;
232 }
233}
234
235//____________________________________________________________________________
236void AliEMCALQADataMaker::MakeRaws(AliRawReader* rawReader)
237{
238 //Raw QA info not yet implemented for EMCAL
239
240}
241
242//____________________________________________________________________________
243void AliEMCALQADataMaker::MakeRecPoints(TTree * clustersTree)
244{
245 // makes data from RecPoints
246 TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
247 if (!emcbranch) {
248 AliError("can't get the branch with the EMCAL clusters !");
249 return;
250 }
251 TObjArray * emcrecpoints = new TObjArray(100) ;
252 emcbranch->SetAddress(&emcrecpoints);
253 emcbranch->GetEntry(0);
254
255 GetRecPointsData(1)->Fill(emcrecpoints->GetEntriesFast()) ;
256 TIter next(emcrecpoints) ;
257 AliEMCALRecPoint * rp ;
258 while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
259 GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
260 }
261 emcrecpoints->Delete();
262 delete emcrecpoints;
263
264}
265
266//____________________________________________________________________________
267void AliEMCALQADataMaker::MakeSDigits(TClonesArray * sdigits)
268{
269 // makes data from SDigits
270
271 GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ;
272 TIter next(sdigits) ;
273 AliEMCALDigit * sdigit ;
274 while ( (sdigit = dynamic_cast<AliEMCALDigit *>(next())) ) {
275 GetSDigitsData(0)->Fill( sdigit->GetAmp()) ;
276 }
277}
278
279//____________________________________________________________________________
280void AliEMCALQADataMaker::MakeSDigits(TTree * sdigitTree)
281{
282 // makes data from SDigit Tree
283 TClonesArray * sdigits = new TClonesArray("AliEMCALDigit", 1000) ;
284
285 TBranch * branch = sdigitTree->GetBranch("EMCAL") ;
286 if ( ! branch ) {
287 AliWarning("EMCAL branch in SDigit Tree not found") ;
288 } else {
289 branch->SetAddress(&sdigits) ;
290 branch->GetEntry(0) ;
291 MakeSDigits(sdigits) ;
292 }
293}
294
295//____________________________________________________________________________
296void AliEMCALQADataMaker::StartOfDetectorCycle()
297{
298 //Detector specific actions at start of cycle
299
300}