1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 Produces the data needed to calculate the quality assurance.
18 All data must be mergeable objects.
20 Based on PHOS code written by
21 Y. Schutz CERN July 2007
24 // --- ROOT system ---
25 #include <TClonesArray.h>
32 // --- Standard library ---
34 // --- AliRoot header files ---
35 #include "AliESDCaloCluster.h"
37 #include "AliEMCALDigit.h"
38 #include "AliEMCALHit.h"
39 #include "AliEMCALQADataMakerSim.h"
40 #include "AliQAChecker.h"
41 #include "AliEMCALSDigitizer.h"
43 ClassImp(AliEMCALQADataMakerSim)
45 //____________________________________________________________________________
46 AliEMCALQADataMakerSim::AliEMCALQADataMakerSim() :
47 AliQADataMakerSim(AliQA::GetDetName(AliQA::kEMCAL), "EMCAL Quality Assurance Data Maker")
52 //____________________________________________________________________________
53 AliEMCALQADataMakerSim::AliEMCALQADataMakerSim(const AliEMCALQADataMakerSim& qadm) :
57 SetName((const char*)qadm.GetName()) ;
58 SetTitle((const char*)qadm.GetTitle());
61 //__________________________________________________________________
62 AliEMCALQADataMakerSim& AliEMCALQADataMakerSim::operator = (const AliEMCALQADataMakerSim& qadm )
65 this->~AliEMCALQADataMakerSim();
66 new(this) AliEMCALQADataMakerSim(qadm);
70 //____________________________________________________________________________
71 void AliEMCALQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
73 //Detector specific actions at end of cycle
75 AliQAChecker::Instance()->Run(AliQA::kEMCAL, task, list) ;
78 //____________________________________________________________________________
79 void AliEMCALQADataMakerSim::InitHits()
81 // create Hits histograms in Hits subdir
82 TH1F * h0 = new TH1F("hEmcalHits", "Hits energy distribution in EMCAL", 200, 0., 2.) ; //GeV
85 TH1I * h1 = new TH1I("hEmcalHitsMul", "Hits multiplicity distribution in EMCAL", 1000, 0, 10000) ;
90 //____________________________________________________________________________
91 void AliEMCALQADataMakerSim::InitDigits()
93 // create Digits histograms in Digits subdir
94 TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL", 500, 0, 500) ;
96 Add2DigitsList(h0, 0) ;
97 TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL", 200, 0, 2000) ;
99 Add2DigitsList(h1, 1) ;
102 //____________________________________________________________________________
103 void AliEMCALQADataMakerSim::InitSDigits()
105 // create SDigits histograms in SDigits subdir
106 TH1F * h0 = new TH1F("hEmcalSDigits", "SDigits energy distribution in EMCAL", 200, 0., 20.) ;
108 Add2SDigitsList(h0, 0) ;
109 TH1I * h1 = new TH1I("hEmcalSDigitsMul", "SDigits multiplicity distribution in EMCAL", 500, 0, 5000) ;
111 Add2SDigitsList(h1, 1) ;
114 //____________________________________________________________________________
115 void AliEMCALQADataMakerSim::MakeHits(TClonesArray * hits)
117 //make QA data from Hits
119 GetHitsData(1)->Fill(hits->GetEntriesFast()) ;
122 while ( (hit = dynamic_cast<AliEMCALHit *>(next())) ) {
123 GetHitsData(0)->Fill( hit->GetEnergy()) ;
128 //____________________________________________________________________________
129 void AliEMCALQADataMakerSim::MakeHits(TTree * hitTree)
131 // make QA data from Hit Tree
133 TClonesArray * hits = new TClonesArray("AliEMCALHit", 1000);
135 TBranch * branch = hitTree->GetBranch("EMCAL") ;
137 AliWarning("EMCAL branch in Hit Tree not found") ;
139 TClonesArray * tmp = new TClonesArray("AliEMCALHit", 1000) ;
140 branch->SetAddress(&tmp) ;
142 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
143 branch->GetEntry(ientry) ;
144 for (Int_t ihit = 0 ; ihit < tmp->GetEntries() ; ihit++) {
145 AliEMCALHit * hit = dynamic_cast<AliEMCALHit *> (tmp->At(ihit)) ;
146 new((*hits)[index]) AliEMCALHit(*hit) ;
156 //____________________________________________________________________________
157 void AliEMCALQADataMakerSim::MakeDigits(TClonesArray * digits)
159 // makes data from Digits
161 GetDigitsData(1)->Fill(digits->GetEntriesFast()) ;
163 AliEMCALDigit * digit ;
164 while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) {
165 GetDigitsData(0)->Fill( digit->GetAmp()) ;
170 //____________________________________________________________________________
171 void AliEMCALQADataMakerSim::MakeDigits(TTree * digitTree)
173 // makes data from Digit Tree
174 TClonesArray * digits = new TClonesArray("AliEMCALDigit", 1000) ;
176 TBranch * branch = digitTree->GetBranch("EMCAL") ;
178 AliWarning("EMCAL branch in Digit Tree not found") ;
180 branch->SetAddress(&digits) ;
181 branch->GetEntry(0) ;
187 //____________________________________________________________________________
188 void AliEMCALQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
190 // makes data from SDigits
191 //Need a copy of the SDigitizer to calibrate the sdigit amplitude to
193 AliEMCALSDigitizer* sDigitizer = new AliEMCALSDigitizer();
195 GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ;
196 TIter next(sdigits) ;
197 AliEMCALDigit * sdigit ;
198 while ( (sdigit = dynamic_cast<AliEMCALDigit *>(next())) ) {
199 GetSDigitsData(0)->Fill( sDigitizer->Calibrate(sdigit->GetAmp())) ;
205 //____________________________________________________________________________
206 void AliEMCALQADataMakerSim::MakeSDigits(TTree * sdigitTree)
208 // makes data from SDigit Tree
209 TClonesArray * sdigits = new TClonesArray("AliEMCALDigit", 1000) ;
211 TBranch * branch = sdigitTree->GetBranch("EMCAL") ;
213 AliWarning("EMCAL branch in SDigit Tree not found") ;
215 branch->SetAddress(&sdigits) ;
216 branch->GetEntry(0) ;
217 MakeSDigits(sdigits) ;
222 //____________________________________________________________________________
223 void AliEMCALQADataMakerSim::StartOfDetectorCycle()
225 //Detector specific actions at start of cycle