]>
Commit | Line | Data |
---|---|---|
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 | #include <TTree.h> | |
31 | ||
32 | // --- Standard library --- | |
33 | ||
34 | // --- AliRoot header files --- | |
35 | #include "AliESDCaloCluster.h" | |
36 | #include "AliLog.h" | |
37 | #include "AliEMCALDigit.h" | |
38 | #include "AliEMCALHit.h" | |
39 | #include "AliEMCALQADataMakerSim.h" | |
40 | #include "AliQAChecker.h" | |
41 | #include "AliEMCALSDigitizer.h" | |
42 | ||
43 | ClassImp(AliEMCALQADataMakerSim) | |
44 | ||
45 | //____________________________________________________________________________ | |
46 | AliEMCALQADataMakerSim::AliEMCALQADataMakerSim() : | |
47 | AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kEMCAL), "EMCAL Quality Assurance Data Maker") | |
48 | { | |
49 | // ctor | |
50 | } | |
51 | ||
52 | //____________________________________________________________________________ | |
53 | AliEMCALQADataMakerSim::AliEMCALQADataMakerSim(const AliEMCALQADataMakerSim& qadm) : | |
54 | AliQADataMakerSim() | |
55 | { | |
56 | //copy ctor | |
57 | SetName((const char*)qadm.GetName()) ; | |
58 | SetTitle((const char*)qadm.GetTitle()); | |
59 | } | |
60 | ||
61 | //__________________________________________________________________ | |
62 | AliEMCALQADataMakerSim& AliEMCALQADataMakerSim::operator = (const AliEMCALQADataMakerSim& qadm ) | |
63 | { | |
64 | // Assign operator. | |
65 | this->~AliEMCALQADataMakerSim(); | |
66 | new(this) AliEMCALQADataMakerSim(qadm); | |
67 | return *this; | |
68 | } | |
69 | ||
70 | //____________________________________________________________________________ | |
71 | void AliEMCALQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list) | |
72 | { | |
73 | //Detector specific actions at end of cycle | |
74 | // do the QA checking | |
75 | AliQAChecker::Instance()->Run(AliQAv1::kEMCAL, task, list) ; | |
76 | } | |
77 | ||
78 | //____________________________________________________________________________ | |
79 | void AliEMCALQADataMakerSim::InitHits() | |
80 | { | |
81 | // create Hits histograms in Hits subdir | |
82 | const Bool_t expert = kTRUE ; | |
83 | const Bool_t image = kTRUE ; | |
84 | ||
85 | TH1F * h0 = new TH1F("hEmcalHits", "Hits energy distribution in EMCAL;Energy [MeV];Counts", 200, 0., 2.) ; //GeV | |
86 | h0->Sumw2() ; | |
87 | Add2HitsList(h0, 0, !expert, image) ; | |
88 | TH1I * h1 = new TH1I("hEmcalHitsMul", "Hits multiplicity distribution in EMCAL;# of Hits;Entries", 1000, 0, 10000) ; | |
89 | h1->Sumw2() ; | |
90 | Add2HitsList(h1, 1, !expert, image) ; | |
91 | } | |
92 | ||
93 | //____________________________________________________________________________ | |
94 | void AliEMCALQADataMakerSim::InitDigits() | |
95 | { | |
96 | // create Digits histograms in Digits subdir | |
97 | const Bool_t expert = kTRUE ; | |
98 | const Bool_t image = kTRUE ; | |
99 | ||
100 | TH1I * h0 = new TH1I("hEmcalDigits", "Digits amplitude distribution in EMCAL;Amplitude [ADC counts];Counts", 500, 0, 500) ; | |
101 | h0->Sumw2() ; | |
102 | Add2DigitsList(h0, 0, !expert, image) ; | |
103 | TH1I * h1 = new TH1I("hEmcalDigitsMul", "Digits multiplicity distribution in EMCAL;# of Digits;Entries", 200, 0, 2000) ; | |
104 | h1->Sumw2() ; | |
105 | Add2DigitsList(h1, 1, !expert, image) ; | |
106 | } | |
107 | ||
108 | //____________________________________________________________________________ | |
109 | void AliEMCALQADataMakerSim::InitSDigits() | |
110 | { | |
111 | // create SDigits histograms in SDigits subdir | |
112 | const Bool_t expert = kTRUE ; | |
113 | const Bool_t image = kTRUE ; | |
114 | ||
115 | TH1F * h0 = new TH1F("hEmcalSDigits", "SDigits energy distribution in EMCAL;Energy [MeV];Counts", 200, 0., 20.) ; | |
116 | h0->Sumw2() ; | |
117 | Add2SDigitsList(h0, 0, !expert, image) ; | |
118 | TH1I * h1 = new TH1I("hEmcalSDigitsMul", "SDigits multiplicity distribution in EMCAL;# of SDigits;Entries", 500, 0, 5000) ; | |
119 | h1->Sumw2() ; | |
120 | Add2SDigitsList(h1, 1, !expert, image) ; | |
121 | } | |
122 | ||
123 | //____________________________________________________________________________ | |
124 | void AliEMCALQADataMakerSim::MakeHits() | |
125 | { | |
126 | //make QA data from Hits | |
127 | ||
128 | GetHitsData(1)->Fill(fHitsArray->GetEntriesFast()) ; | |
129 | TIter next(fHitsArray) ; | |
130 | AliEMCALHit * hit ; | |
131 | while ( (hit = dynamic_cast<AliEMCALHit *>(next())) ) { | |
132 | GetHitsData(0)->Fill( hit->GetEnergy()) ; | |
133 | } | |
134 | ||
135 | } | |
136 | ||
137 | //____________________________________________________________________________ | |
138 | void AliEMCALQADataMakerSim::MakeHits(TTree * hitTree) | |
139 | { | |
140 | // make QA data from Hit Tree | |
141 | ||
142 | if (fHitsArray) | |
143 | fHitsArray->Clear() ; | |
144 | else | |
145 | fHitsArray = new TClonesArray("AliEMCALHit", 1000); | |
146 | ||
147 | TBranch * branch = hitTree->GetBranch("EMCAL") ; | |
148 | if ( ! branch ) { | |
149 | AliWarning("EMCAL branch in Hit Tree not found") ; | |
150 | } else { | |
151 | branch->SetAddress(&fHitsArray) ; | |
152 | for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) { | |
153 | branch->GetEntry(ientry) ; | |
154 | MakeHits() ; | |
155 | fHitsArray->Clear() ; | |
156 | } | |
157 | } | |
158 | } | |
159 | ||
160 | //____________________________________________________________________________ | |
161 | void AliEMCALQADataMakerSim::MakeDigits() | |
162 | { | |
163 | // makes data from Digits | |
164 | ||
165 | GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ; | |
166 | TIter next(fDigitsArray) ; | |
167 | AliEMCALDigit * digit ; | |
168 | while ( (digit = dynamic_cast<AliEMCALDigit *>(next())) ) { | |
169 | GetDigitsData(0)->Fill( digit->GetAmp()) ; | |
170 | } | |
171 | ||
172 | } | |
173 | ||
174 | //____________________________________________________________________________ | |
175 | void AliEMCALQADataMakerSim::MakeDigits(TTree * digitTree) | |
176 | { | |
177 | // makes data from Digit Tree | |
178 | if (fDigitsArray) | |
179 | fDigitsArray->Clear("C") ; | |
180 | else | |
181 | fDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ; | |
182 | ||
183 | TBranch * branch = digitTree->GetBranch("EMCAL") ; | |
184 | if ( ! branch ) { | |
185 | AliWarning("EMCAL branch in Digit Tree not found") ; | |
186 | } else { | |
187 | branch->SetAddress(&fDigitsArray) ; | |
188 | branch->GetEntry(0) ; | |
189 | MakeDigits() ; | |
190 | } | |
191 | ||
192 | } | |
193 | ||
194 | //____________________________________________________________________________ | |
195 | void AliEMCALQADataMakerSim::MakeSDigits() | |
196 | { | |
197 | // makes data from SDigits | |
198 | //Need a copy of the SDigitizer to calibrate the sdigit amplitude to | |
199 | //energy in GeV | |
200 | ||
201 | AliEMCALSDigitizer* sDigitizer = new AliEMCALSDigitizer(); | |
202 | ||
203 | GetSDigitsData(1)->Fill(fSDigitsArray->GetEntriesFast()) ; | |
204 | TIter next(fSDigitsArray) ; | |
205 | AliEMCALDigit * sdigit ; | |
206 | while ( (sdigit = dynamic_cast<AliEMCALDigit *>(next())) ) { | |
207 | GetSDigitsData(0)->Fill( sDigitizer->Calibrate(sdigit->GetAmp())) ; | |
208 | } | |
209 | ||
210 | delete sDigitizer; | |
211 | } | |
212 | ||
213 | //____________________________________________________________________________ | |
214 | void AliEMCALQADataMakerSim::MakeSDigits(TTree * sdigitTree) | |
215 | { | |
216 | // makes data from SDigit Tree | |
217 | if (fSDigitsArray) | |
218 | fSDigitsArray->Clear("C") ; | |
219 | else | |
220 | fSDigitsArray = new TClonesArray("AliEMCALDigit", 1000) ; | |
221 | ||
222 | TBranch * branch = sdigitTree->GetBranch("EMCAL") ; | |
223 | if ( ! branch ) { | |
224 | AliWarning("EMCAL branch in SDigit Tree not found") ; | |
225 | } else { | |
226 | branch->SetAddress(&fSDigitsArray) ; | |
227 | branch->GetEntry(0) ; | |
228 | MakeSDigits() ; | |
229 | } | |
230 | ||
231 | } | |
232 | ||
233 | //____________________________________________________________________________ | |
234 | void AliEMCALQADataMakerSim::StartOfDetectorCycle() | |
235 | { | |
236 | //Detector specific actions at start of cycle | |
237 | ||
238 | } |