]>
Commit | Line | Data |
---|---|---|
c52c2132 | 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 | /* $Id$ */ | |
17 | // Author: Andrei Gheata, 31/05/2006 | |
18 | ||
19 | //============================================================================== | |
20 | // AliAnalysisSelector - A transparent selector to be created by | |
21 | // AliAnalysisManager to handle analysis. | |
22 | //============================================================================== | |
23 | ||
c2e40c93 | 24 | #include <Riostream.h> |
25 | #include <TProcessID.h> | |
c52c2132 | 26 | |
27 | #include "AliAnalysisManager.h" | |
28 | #include "AliAnalysisTask.h" | |
29 | #include "AliAnalysisDataContainer.h" | |
30 | #include "AliAnalysisSelector.h" | |
31 | ||
32 | ClassImp(AliAnalysisSelector) | |
33 | ||
1f87e9fb | 34 | //______________________________________________________________________________ |
35 | AliAnalysisSelector::AliAnalysisSelector() | |
36 | :TSelector(), | |
37 | fInitialized(kFALSE), | |
38 | fAnalysis(NULL) | |
39 | { | |
40 | // Dummy ctor. | |
41 | fAnalysis = AliAnalysisManager::GetAnalysisManager(); | |
42 | if (fAnalysis) fAnalysis->SetSelector(this); | |
43 | } | |
44 | ||
c52c2132 | 45 | //______________________________________________________________________________ |
46 | AliAnalysisSelector::AliAnalysisSelector(AliAnalysisManager *mgr) | |
47 | :TSelector(), | |
48 | fInitialized(kFALSE), | |
49 | fAnalysis(mgr) | |
50 | { | |
51 | // Constructor. Called by AliAnalysisManager which registers itself on the | |
52 | // selector running on the master. | |
1f87e9fb | 53 | mgr->SetSelector(this); |
c52c2132 | 54 | } |
55 | ||
56 | //______________________________________________________________________________ | |
57 | AliAnalysisSelector::~AliAnalysisSelector() | |
58 | { | |
59 | // Dtor. The analysis manager object is sent in the input list and duplicated | |
60 | // on the workers - it needs to be deleted (?) | |
61 | // if (fAnalysis) delete fAnalysis; | |
62 | } | |
63 | ||
64 | //______________________________________________________________________________ | |
65 | void AliAnalysisSelector::Init(TTree *tree) | |
66 | { | |
67 | // Called after Begin/SlaveBegin, assumes that fAnalysis is already initialized. | |
68 | // Is Init called on workers in case of PROOF. | |
69 | if (!fAnalysis) { | |
70 | Error("Init", "Analysis manager NULL !"); | |
2d626244 | 71 | Abort("Cannot initialize without analysis manager. Aborting."); |
72 | SetStatus(-1); | |
c52c2132 | 73 | return; |
74 | } | |
cd463514 | 75 | if (fAnalysis->GetDebugLevel()>1) { |
981f2614 | 76 | cout << "->AliAnalysisSelector->Init()" << endl; |
77 | } | |
c52c2132 | 78 | if (!tree) { |
79 | Error("Init", "Input tree is NULL !"); | |
2d626244 | 80 | Abort("Cannot initialize without tree. Aborting."); |
81 | SetStatus(-1); | |
c52c2132 | 82 | return; |
83 | } | |
2d626244 | 84 | fInitialized = fAnalysis->Init(tree); |
85 | if (!fInitialized) { | |
86 | Error("Init", "Some error occured during analysis manager initialization. Aborting."); | |
87 | Abort("Error during AliAnalysisManager::Init()"); | |
88 | SetStatus(-1); | |
89 | return; | |
90 | } | |
cd463514 | 91 | if (fAnalysis->GetDebugLevel()>1) { |
981f2614 | 92 | cout << "<-AliAnalysisSelector->Init()" << endl; |
93 | } | |
c52c2132 | 94 | } |
95 | ||
96 | //______________________________________________________________________________ | |
97 | void AliAnalysisSelector::Begin(TTree *) | |
98 | { | |
99 | // Assembly the input list. | |
100 | RestoreAnalysisManager(); | |
cd463514 | 101 | if (fAnalysis && fAnalysis->GetDebugLevel()>1) { |
36e82a52 | 102 | cout << "->AliAnalysisSelector->Begin: Analysis manager restored" << endl; |
103 | } | |
c52c2132 | 104 | } |
105 | ||
106 | //______________________________________________________________________________ | |
107 | void AliAnalysisSelector::SlaveBegin(TTree *tree) | |
108 | { | |
109 | // Called on each worker. We "unpack" analysis manager here and call InitAnalysis. | |
50dfc392 | 110 | TObject::SetObjectStat(kFALSE); |
c52c2132 | 111 | RestoreAnalysisManager(); |
981f2614 | 112 | if (fAnalysis) { |
cd463514 | 113 | if (fAnalysis->GetDebugLevel()>1) { |
981f2614 | 114 | cout << "->AliAnalysisSelector->SlaveBegin() after Restore" << endl; |
115 | } | |
116 | fAnalysis->SlaveBegin(tree); | |
cd463514 | 117 | if (fAnalysis->GetDebugLevel()>1) { |
981f2614 | 118 | cout << "<-AliAnalysisSelector->SlaveBegin()" << endl; |
119 | } | |
120 | } | |
36e82a52 | 121 | } |
c52c2132 | 122 | |
981f2614 | 123 | //______________________________________________________________________________ |
124 | Bool_t AliAnalysisSelector::Notify() | |
125 | { | |
126 | // The Notify() function is called when a new file is opened. This | |
127 | // can be either for a new TTree in a TChain or when when a new TTree | |
128 | // is started when using PROOF. It is normaly not necessary to make changes | |
129 | // to the generated code, but the routine can be extended by the | |
130 | // user if needed. The return value is currently not used. | |
d8a5ee94 | 131 | if (fAnalysis) return fAnalysis->Notify(); |
132 | return kFALSE; | |
36e82a52 | 133 | } |
981f2614 | 134 | |
c52c2132 | 135 | //______________________________________________________________________________ |
136 | Bool_t AliAnalysisSelector::Process(Long64_t entry) | |
137 | { | |
138 | // Event loop. | |
4747b4a7 | 139 | static Int_t count = 0; |
140 | count++; | |
cd463514 | 141 | if (fAnalysis->GetDebugLevel() > 1) { |
981f2614 | 142 | cout << "->AliAnalysisSelector::Process()" << endl; |
36e82a52 | 143 | } |
1ee5d9a3 | 144 | static Bool_t init=kTRUE; |
145 | static Int_t nobjCount = 0; | |
146 | if(init) { | |
147 | nobjCount = TProcessID::GetObjectCount(); | |
148 | init=kFALSE; | |
149 | } | |
150 | TProcessID::SetObjectCount(nobjCount); | |
8e57d9c4 | 151 | Int_t returnCode = fAnalysis->GetEntry(entry); |
f5cbe261 | 152 | if (returnCode <= 0) { |
153 | cout << "Error retrieving event:" << entry << " Skipping ..." << endl; | |
154 | fAnalysis->CountEvent(1,0,1,0); | |
1d0b4d65 | 155 | // Try to skip file |
156 | Abort("Bad stream to file. Trying next image.", kAbortFile); | |
157 | return kFALSE; | |
f5cbe261 | 158 | } else { |
159 | fAnalysis->ExecAnalysis(); | |
c5734a4b | 160 | if (returnCode<100000000) fAnalysis->CountEvent(1,1,0,0); |
f5cbe261 | 161 | } |
cd463514 | 162 | if (fAnalysis->GetDebugLevel() > 1) { |
981f2614 | 163 | cout << "<-AliAnalysisSelector::Process()" << endl; |
164 | } | |
c52c2132 | 165 | return kTRUE; |
166 | } | |
167 | ||
168 | //______________________________________________________________________________ | |
169 | void AliAnalysisSelector::RestoreAnalysisManager() | |
170 | { | |
171 | // Restores analysis manager from the input list. | |
172 | if (!fAnalysis) { | |
173 | TIter next(fInput); | |
174 | TObject *obj; | |
175 | while ((obj=next())) { | |
176 | if (obj->IsA() == AliAnalysisManager::Class()) { | |
177 | fAnalysis = (AliAnalysisManager*)obj; | |
8d7d3b59 | 178 | fAnalysis->SetSelector(this); |
cd463514 | 179 | if (fAnalysis->GetDebugLevel()>1) { |
981f2614 | 180 | cout << "->AliAnalysisSelector->RestoreAnalysisManager: Analysis manager restored" << endl; |
181 | } | |
c52c2132 | 182 | break; |
183 | } | |
184 | } | |
185 | if (!fAnalysis) { | |
186 | Error("SlaveBegin", "Analysis manager not found in the input list"); | |
187 | return; | |
188 | } | |
189 | } | |
190 | } | |
191 | ||
192 | //______________________________________________________________________________ | |
193 | void AliAnalysisSelector::SlaveTerminate() | |
194 | { | |
195 | // The SlaveTerminate() function is called after all entries or objects | |
196 | // have been processed. When running with PROOF SlaveTerminate() is called | |
197 | // on each slave server. | |
2d626244 | 198 | if (fStatus == -1) return; // TSelector won't abort... |
aee5ee44 | 199 | if (fAnalysis->GetAnalysisType() == AliAnalysisManager::kMixingAnalysis) return; |
cd463514 | 200 | if (fAnalysis->GetDebugLevel() > 1) { |
981f2614 | 201 | cout << "->AliAnalysisSelector::SlaveTerminate()" << endl; |
c52c2132 | 202 | } |
203 | fAnalysis->PackOutput(fOutput); | |
cd463514 | 204 | if (fAnalysis->GetDebugLevel() > 1) { |
981f2614 | 205 | cout << "<-AliAnalysisSelector::SlaveTerminate()" << endl; |
206 | } | |
c52c2132 | 207 | } |
208 | ||
209 | //______________________________________________________________________________ | |
210 | void AliAnalysisSelector::Terminate() | |
211 | { | |
212 | // The Terminate() function is the last function to be called during | |
213 | // a query. It always runs on the client, it can be used to present | |
214 | // the results graphically or save the results to file. | |
2d626244 | 215 | if (fStatus == -1) return; // TSelector won't abort... |
c52c2132 | 216 | if (!fAnalysis) { |
36e82a52 | 217 | Error("Terminate","AliAnalysisSelector::Terminate: No analysis manager!!!"); |
c52c2132 | 218 | return; |
219 | } | |
aee5ee44 | 220 | // No Terminate() in case of event mixing |
221 | if (fAnalysis->GetAnalysisType() == AliAnalysisManager::kMixingAnalysis) return; | |
cd463514 | 222 | if (fAnalysis->GetDebugLevel() > 1) { |
981f2614 | 223 | cout << "->AliAnalysisSelector::Terminate()" << endl; |
c52c2132 | 224 | } |
225 | fAnalysis->UnpackOutput(fOutput); | |
226 | fAnalysis->Terminate(); | |
cd463514 | 227 | if (fAnalysis->GetDebugLevel() > 1) { |
981f2614 | 228 | cout << "<-AliAnalysisSelector::Terminate()" << endl; |
229 | } | |
c52c2132 | 230 | } |