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