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