]>
Commit | Line | Data |
---|---|---|
61899827 | 1 | /* $Id: AliPhysicsSelection.cxx 35782 2009-10-22 11:54:31Z jgrosseo $ */ |
2 | ||
3 | /************************************************************************** | |
4 | * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. * | |
5 | * * | |
6 | * Author: The ALICE Off-line Project. * | |
7 | * Contributors are mentioned in the code where appropriate. * | |
8 | * * | |
9 | * Permission to use, copy, modify and distribute this software and its * | |
10 | * documentation strictly for non-commercial purposes is hereby granted * | |
11 | * without fee, provided that the above copyright notice appears in all * | |
12 | * copies and that both the copyright notice and this permission notice * | |
13 | * appear in the supporting documentation. The authors make no claims * | |
14 | * about the suitability of this software for any purpose. It is * | |
15 | * provided "as is" without express or implied warranty. * | |
16 | **************************************************************************/ | |
17 | ||
18 | //------------------------------------------------------------------------- | |
19 | // Implementation of Class AliPhysicsSelection | |
528640ed | 20 | // This class selects collision candidates from data runs, applying selection cuts on triggers |
21 | // and background rejection based on the content of the ESD | |
22 | // | |
23 | // Usage: | |
24 | // | |
25 | // Create the object and initialize it with the correct run number: | |
26 | // fPhysicsSelection = new AliPhysicsSelection; | |
27 | // fPhysicsSelection->Initialize(104160); | |
28 | // | |
29 | // To check if an event is a collision candidate, use: | |
30 | // fPhysicsSelection->IsCollisionCandidate(fESD) | |
31 | // | |
32 | // After processing save the resulting histograms to a file with (a folder physics_selection | |
33 | // will be created that contains the histograms): | |
34 | // fPhysicsSelection->SaveHistograms("physics_selection") | |
35 | // | |
36 | // To print statistics after processing use: | |
37 | // fPhysicsSelection->Print(); | |
38 | // | |
61899827 | 39 | // Origin: Jan Fiete Grosse-Oetringhaus, CERN |
40 | //------------------------------------------------------------------------- | |
41 | ||
42 | #include <Riostream.h> | |
43 | #include <TH1F.h> | |
44 | #include <TH2F.h> | |
45 | #include <TList.h> | |
46 | #include <TIterator.h> | |
47 | #include <TDirectory.h> | |
48 | ||
49 | #include <AliPhysicsSelection.h> | |
50 | ||
51 | #include <AliTriggerAnalysis.h> | |
52 | #include <AliLog.h> | |
53 | ||
54 | #include <AliESDEvent.h> | |
55 | ||
56 | ClassImp(AliPhysicsSelection) | |
57 | ||
58 | AliPhysicsSelection::AliPhysicsSelection() : | |
59 | fTriggerAnalysis(0), | |
60 | fHistStatistics(0), | |
61 | fHistBunchCrossing(0) | |
62 | { | |
63 | // constructor | |
64 | ||
65 | AliLog::SetClassDebugLevel("AliPhysicsSelection", AliLog::kWarning); | |
66 | } | |
67 | ||
68 | AliPhysicsSelection::~AliPhysicsSelection() | |
69 | { | |
70 | // destructor | |
71 | ||
72 | if (fTriggerAnalysis) | |
73 | { | |
74 | delete fTriggerAnalysis; | |
75 | fTriggerAnalysis = 0; | |
76 | } | |
77 | ||
78 | if (fHistStatistics) | |
79 | { | |
80 | delete fHistStatistics; | |
81 | fHistStatistics = 0; | |
82 | } | |
83 | ||
84 | if (fHistBunchCrossing) | |
85 | { | |
86 | delete fHistBunchCrossing; | |
87 | fHistBunchCrossing = 0; | |
88 | } | |
89 | } | |
90 | ||
91 | Bool_t AliPhysicsSelection::IsCollisionCandidate(const AliESDEvent* aEsd) | |
92 | { | |
93 | // checks if the given event is a collision candidate | |
94 | ||
95 | if (!fTriggerAnalysis) | |
96 | AliFatal("Not initialized!"); | |
97 | ||
98 | const AliESDHeader* esdHeader = aEsd->GetHeader(); | |
99 | if (!esdHeader) | |
100 | { | |
101 | AliError("ESD Header could not be retrieved"); | |
102 | return kFALSE; | |
103 | } | |
104 | ||
105 | // check event type (should be PHYSICS = 7) | |
106 | if (esdHeader->GetEventType() != 7) | |
107 | return kFALSE; | |
108 | ||
109 | fHistStatistics->Fill(1); | |
110 | ||
111 | fTriggerAnalysis->FillTriggerClasses(aEsd); | |
112 | ||
113 | for (Int_t i=0; i < fRequTrigClasses.GetEntries(); i++) | |
114 | { | |
115 | const char* triggerClass = ((TObjString*) fRequTrigClasses.At(i))->String(); | |
116 | if (!aEsd->IsTriggerClassFired(triggerClass)) | |
117 | { | |
118 | AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is not present", triggerClass)); | |
119 | return kFALSE; | |
120 | } | |
121 | } | |
122 | ||
123 | for (Int_t i=0; i < fRejTrigClasses.GetEntries(); i++) | |
124 | { | |
125 | const char* triggerClass = ((TObjString*) fRejTrigClasses.At(i))->String(); | |
126 | if (aEsd->IsTriggerClassFired(triggerClass)) | |
127 | { | |
128 | AliDebug(AliLog::kDebug, Form("Rejecting event because trigger class %s is present", triggerClass)); | |
129 | return kFALSE; | |
130 | } | |
131 | } | |
132 | ||
133 | fTriggerAnalysis->FillHistograms(aEsd); | |
134 | ||
135 | fHistStatistics->Fill(2); | |
136 | ||
137 | Bool_t fastOR = fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kSPDGFO); | |
138 | Bool_t v0BB = fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0A) || fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0C); | |
139 | Bool_t v0BG = fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0ABG) || fTriggerAnalysis->IsOfflineTriggerFired(aEsd, AliTriggerAnalysis::kV0CBG); | |
140 | ||
141 | if (fastOR) | |
142 | fHistStatistics->Fill(3); | |
143 | if (v0BB) | |
144 | fHistStatistics->Fill(4); | |
145 | if (v0BG) | |
146 | fHistStatistics->Fill(5); | |
147 | ||
148 | if (fastOR || v0BB) | |
149 | fHistStatistics->Fill(6); | |
150 | ||
151 | if (!fastOR && !v0BB) | |
152 | { | |
153 | AliDebug(AliLog::kDebug, "Rejecting event because neither FO nor V0 has triggered"); | |
154 | return kFALSE; | |
155 | } | |
156 | ||
157 | if (v0BG) | |
158 | { | |
159 | AliDebug(AliLog::kDebug, "Rejecting event because of V0 BG flag"); | |
160 | return kFALSE; | |
161 | } | |
162 | ||
163 | fHistStatistics->Fill(7); | |
164 | ||
165 | // TODO additional background identification | |
166 | ||
167 | fHistStatistics->Fill(9); | |
168 | ||
169 | fHistBunchCrossing->Fill(aEsd->GetBunchCrossNumber()); | |
170 | ||
171 | AliDebug(AliLog::kDebug, "Accepted event"); | |
172 | ||
173 | return kTRUE; | |
174 | } | |
175 | ||
176 | Bool_t AliPhysicsSelection::Initialize(UInt_t runNumber) | |
177 | { | |
178 | // initializes the object for the given run | |
179 | // TODO having the run number here and parameters hardcoded is clearly temporary, a way needs to be found to have a CDB-like configuration also for analysis | |
180 | ||
181 | AliInfo(Form("Initializing for run %d", runNumber)); | |
182 | ||
183 | fRequTrigClasses.Clear(); | |
184 | fRejTrigClasses.Clear(); | |
185 | ||
186 | fRequTrigClasses.Add(new TObjString("CINT1B-ABCE-NOPF-ALL")); | |
187 | ||
188 | if (!fTriggerAnalysis) | |
189 | { | |
190 | fTriggerAnalysis = new AliTriggerAnalysis; | |
191 | fTriggerAnalysis->EnableHistograms(); | |
192 | } | |
193 | ||
194 | fTriggerAnalysis->SetSPDGFOThreshhold(1); | |
195 | fTriggerAnalysis->SetV0TimeOffset(0); | |
196 | ||
197 | if (runNumber == 104321) | |
198 | fTriggerAnalysis->SetV0TimeOffset(7.5); | |
199 | ||
200 | if (fHistStatistics) | |
201 | { | |
202 | fHistStatistics->Reset(); | |
203 | } | |
204 | else | |
205 | { | |
206 | fHistStatistics = new TH1F("fHistStatistics", "fHistStatistics;;event count", 10, 0.5, 10.5); | |
207 | ||
208 | Int_t n = 1; | |
209 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "Total"); | |
210 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "Correct trigger class(es)"); | |
211 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO"); | |
212 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0 BB"); | |
213 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "V0 BG"); | |
214 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "FO | V0 BB"); | |
215 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "(FO | V0 BB) & !V0 BG"); | |
216 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "Background identification"); | |
217 | fHistStatistics->GetXaxis()->SetBinLabel(n++, "Accepted"); | |
218 | } | |
219 | ||
220 | if (fHistBunchCrossing) | |
221 | { | |
222 | fHistBunchCrossing->Reset(); | |
223 | } | |
224 | else | |
225 | fHistBunchCrossing = new TH1F("fHistBunchCrossing", "fHistBunchCrossing;bunch crossing number;accepted events", 4000, -0.5, 3999.5); | |
226 | ||
227 | return kTRUE; | |
228 | } | |
229 | ||
230 | void AliPhysicsSelection::Print(Option_t* /* option */) const | |
231 | { | |
232 | // print the configuration | |
233 | ||
528640ed | 234 | AliInfo("Configuration:"); |
61899827 | 235 | |
236 | TString str("Required trigger classes: "); | |
237 | for (Int_t i=0; i < fRequTrigClasses.GetEntries(); i++) | |
238 | str += ((TObjString*) fRequTrigClasses.At(i))->String() + " "; | |
239 | AliInfo(str); | |
240 | ||
241 | str = "Rejected trigger classes: "; | |
242 | for (Int_t i=0; i < fRejTrigClasses.GetEntries(); i++) | |
243 | str += ((TObjString*) fRejTrigClasses.At(i))->String() + " "; | |
244 | AliInfo(str); | |
245 | ||
246 | AliInfo(Form("Requiring %d FO chips (offline) or V0A or V0C and no V0 BG flag", fTriggerAnalysis->GetSPDGFOThreshhold())); | |
247 | ||
248 | if (fTriggerAnalysis->GetV0TimeOffset() > 0) | |
249 | AliInfo(Form("V0 time offset active: %.2f ns", fTriggerAnalysis->GetV0TimeOffset())); | |
250 | ||
528640ed | 251 | AliInfo(""); |
252 | ||
253 | AliInfo("Selection statistics:"); | |
254 | ||
255 | AliInfo("Total available events:"); | |
61899827 | 256 | fTriggerAnalysis->PrintTriggerClasses(); |
528640ed | 257 | |
258 | AliInfo(Form("Total events: %d", (Int_t) fHistStatistics->GetBinContent(1))); | |
259 | AliInfo(Form("Correct trigger class: %d", (Int_t) fHistStatistics->GetBinContent(2))); | |
260 | AliInfo(Form("Selected collision candidates: %d", (Int_t) fHistStatistics->GetBinContent(9))); | |
61899827 | 261 | } |
262 | ||
263 | Long64_t AliPhysicsSelection::Merge(TCollection* list) | |
264 | { | |
265 | // Merge a list of AliMultiplicityCorrection objects with this (needed for | |
266 | // PROOF). | |
267 | // Returns the number of merged objects (including this). | |
268 | ||
269 | if (!list) | |
270 | return 0; | |
271 | ||
272 | if (list->IsEmpty()) | |
273 | return 1; | |
274 | ||
275 | TIterator* iter = list->MakeIterator(); | |
276 | TObject* obj; | |
277 | ||
278 | // collections of all histograms | |
279 | const Int_t nHists = 9; | |
280 | TList collections[nHists]; | |
281 | ||
282 | Int_t count = 0; | |
283 | while ((obj = iter->Next())) { | |
284 | ||
285 | AliPhysicsSelection* entry = dynamic_cast<AliPhysicsSelection*> (obj); | |
286 | if (entry == 0) | |
287 | continue; | |
288 | ||
289 | Int_t n = 0; | |
290 | collections[n++].Add(entry->fTriggerAnalysis); | |
291 | collections[n++].Add(entry->fHistStatistics); | |
292 | collections[n++].Add(entry->fHistBunchCrossing); | |
293 | ||
294 | count++; | |
295 | } | |
296 | ||
297 | Int_t n = 0; | |
298 | fTriggerAnalysis->Merge(&collections[n++]); | |
299 | fHistStatistics->Merge(&collections[n++]); | |
300 | fHistBunchCrossing->Merge(&collections[n++]); | |
301 | ||
302 | delete iter; | |
303 | ||
304 | return count+1; | |
305 | } | |
306 | ||
307 | void AliPhysicsSelection::SaveHistograms(const char* folder) const | |
308 | { | |
309 | // write histograms to current directory | |
310 | ||
311 | if (!fHistStatistics) | |
312 | return; | |
313 | ||
314 | if (folder) | |
315 | { | |
316 | gDirectory->mkdir(folder); | |
317 | gDirectory->cd(folder); | |
318 | } | |
319 | ||
320 | fHistStatistics->Write(); | |
321 | fHistBunchCrossing->Write(); | |
322 | ||
323 | gDirectory->mkdir("trigger_histograms"); | |
324 | gDirectory->cd("trigger_histograms"); | |
325 | ||
326 | fTriggerAnalysis->SaveHistograms(); | |
327 | ||
328 | gDirectory->cd(".."); | |
329 | ||
330 | if (folder) | |
331 | gDirectory->cd(".."); | |
332 | } |