]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/AliSelector.cxx
fixed the way the event bias selection is done
[u/mrichter/AliRoot.git] / PWG0 / AliSelector.cxx
1 /* $Id$ */
2
3 // The class definition in esdV0.h has been generated automatically
4 // by the ROOT utility TTree::MakeSelector(). This class is derived
5 // from the ROOT class TSelector. For more information on the TSelector
6 // framework see $ROOTSYS/README/README.SELECTOR or the ROOT User Manual.
7
8 // The following methods are defined in this file:
9 //    Begin():        called everytime a loop on the tree starts,
10 //                    a convenient place to create your histograms.
11 //    SlaveBegin():   called after Begin(), when on PROOF called only on the
12 //                    slave servers.
13 //    Process():      called for each event, in this function you decide what
14 //                    to read and fill your histograms.
15 //    SlaveTerminate: called at the end of the loop on the tree, when on PROOF
16 //                    called only on the slave servers.
17 //    Terminate():    called at the end of the loop on the tree,
18 //                    a convenient place to draw/fit your histograms.
19 //
20 // To use this file, try the following session on your Tree T:
21 //
22 // Root > T->Process("AliSelector.C")
23 // Root > T->Process("AliSelector.C","some options")
24 // Root > T->Process("AliSelector.C+")
25 //
26
27 #include "AliSelector.h"
28
29 #include <TStyle.h>
30 #include <TSystem.h>
31 #include <TCanvas.h>
32 #include <TRegexp.h>
33 #include <TTime.h>
34 #include <TParticle.h>
35 #include <TParticlePDG.h>
36 #include <TFriendElement.h>
37 #include <TTree.h>
38 #include <TChain.h>
39 #include <TFile.h>
40
41 #include <AliLog.h>
42 #include <AliESD.h>
43 #include <AliHeader.h>
44
45 ClassImp(AliSelector)
46
47 AliSelector::AliSelector() :
48   TSelector(),
49   fChain(0),
50   fESD(0),
51   fCountFiles(0),
52   fKineFile(0),
53   fHeaderFile(0),
54   fHeaderTree(0),
55   fHeader(0)
56 {
57   //
58   // Constructor. Initialization of pointers
59   //
60 }
61
62 AliSelector::~AliSelector()
63 {
64   //
65   // Destructor
66   //
67
68   // histograms are in the output list and deleted when the output
69   // list is deleted by the TSelector dtor
70 }
71
72 void AliSelector::Begin(TTree *)
73 {
74   // The Begin() function is called at the start of the query.
75   // When running with PROOF Begin() is only called on the client.
76   // The tree argument is deprecated (on PROOF 0 is passed).
77 }
78
79 void AliSelector::SlaveBegin(TTree * tree)
80 {
81   // The SlaveBegin() function is called after the Begin() function.
82   // When running with PROOF SlaveBegin() is called on each slave server.
83   // The tree argument is deprecated (on PROOF 0 is passed).
84
85   Init(tree);
86
87   AliDebug(AliLog::kDebug, "=======SLAVEBEGIN========");
88   AliDebug(AliLog::kDebug, Form("Hostname: %s", gSystem->HostName()));
89   AliDebug(AliLog::kDebug, Form("Time: %s", gSystem->Now().AsString()));
90
91   TString option = GetOption();
92 }
93
94 void AliSelector::Init(TTree *tree)
95 {
96   // The Init() function is called when the selector needs to initialize
97   // a new tree or chain. Typically here the branch addresses of the tree
98   // will be set. It is normaly not necessary to make changes to the
99   // generated code, but the routine can be extended by the user if needed.
100   // Init() will be called many times when running with PROOF.
101
102   AliDebug(AliLog::kDebug, "=========Init==========");
103
104   // Set branch addresses
105   if (tree == 0)
106   {
107     AliDebug(AliLog::kError, "ERROR: tree argument is 0.");
108     return;
109   }
110
111   fChain = dynamic_cast<TChain*> (tree);
112   if (fChain == 0)
113   {
114     AliDebug(AliLog::kDebug, "ERROR: tree argument could not be casted to TChain.");
115     return;
116   }
117
118   fChain->SetBranchAddress("ESD", &fESD);
119   if (fESD != 0)
120     AliDebug(AliLog::kInfo, "INFO: Found ESD branch in chain.");
121
122   /*fChain->SetBranchAddress("Header", &fHeader);
123   if (fHeader != 0)
124     AliDebug(AliLog::kInfo, "INFO: Found event header branch in chain.");*/
125 }
126
127 Bool_t AliSelector::Notify()
128 {
129   // The Notify() function is called when a new file is opened. This
130   // can be either for a new TTree in a TChain or when when a new TTree
131   // is started when using PROOF. Typically here the branch pointers
132   // will be retrieved. 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.
135
136   AliDebug(AliLog::kDebug, "=========NOTIFY==========");
137   AliDebug(AliLog::kDebug, Form("Hostname: %s", gSystem->HostName()));
138   AliDebug(AliLog::kDebug, Form("Time: %s", gSystem->Now().AsString()));
139
140   ++fCountFiles;
141   TFile *f = fChain->GetCurrentFile();
142   AliDebug(AliLog::kInfo, Form("Processing %d. file %s", fCountFiles, f->GetName()));
143
144   DeleteKinematicsFile();
145   DeleteHeaderFile();
146
147   return kTRUE;
148 }
149
150 Bool_t AliSelector::Process(Long64_t entry)
151 {
152   // The Process() function is called for each entry in the tree (or possibly
153   // keyed object in the case of PROOF) to be processed. The entry argument
154   // specifies which entry in the currently loaded tree is to be processed.
155   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
156   // to read either all or the required parts of the data. When processing
157   // keyed objects with PROOF, the object is already loaded and is available
158   // via the fObject pointer.
159   //
160   // This function should contain the "body" of the analysis. It can contain
161   // simple or elaborate selection criteria, run algorithms on the data
162   // of the event and typically fill histograms.
163
164   // WARNING when a selector is used with a TChain, you must use
165   //  the pointer to the current TTree to call GetEntry(entry).
166   //  The entry is always the local entry number in the current tree.
167   //  Assuming that fChain is the pointer to the TChain being processed,
168   //  use fChain->GetTree()->GetEntry(entry).
169
170   AliDebug(AliLog::kDebug, Form("=========PROCESS========== Entry %lld", entry));
171
172   if (!fChain)
173   {
174     AliDebug(AliLog::kError, "ERROR: fChain is 0.");
175     return kFALSE;
176   }
177
178   fChain->GetTree()->GetEntry(entry);
179
180   /*
181   // debugging
182   if (fESD)
183     AliDebug(AliLog::kDebug, Form("ESD: We have %d tracks.", fESD->GetNumberOfTracks()));
184
185   if (fHeader)
186     AliDebug(AliLog::kDebug, Form("Header: We have %d primaries.", fHeader->GetNprimary()));
187
188   TTree* kinematics = GetKinematics();
189   if (kinematics)
190     AliDebug(AliLog::kDebug, Form("Kinematics: We have %lld particles.", kinematics->GetEntries()));
191   */
192
193   return kTRUE;
194 }
195
196 void AliSelector::SlaveTerminate()
197 {
198   // The SlaveTerminate() function is called after all entries or objects
199   // have been processed. When running with PROOF SlaveTerminate() is called
200   // on each slave server.
201
202   DeleteKinematicsFile();
203   DeleteHeaderFile();
204 }
205
206 void AliSelector::Terminate()
207 {
208   // The Terminate() function is the last function to be called during
209   // a query. It always runs on the client, it can be used to present
210   // the results graphically or save the results to file.
211
212   AliDebug(AliLog::kDebug, "=========TERMINATE==========");
213 }
214
215 TTree* AliSelector::GetKinematics()
216 {
217   // Returns kinematics tree corresponding to current ESD active in fChain
218   // Loads the kinematics from the kinematics file, the file is identified by replacing "AliESDs" to
219   // "Kinematics" in the file path of the ESD file. This is a hack, to be changed!
220
221   if (!fKineFile)
222   {
223     if (!fChain->GetCurrentFile())
224       return 0;
225
226     TString fileName(fChain->GetCurrentFile()->GetName());
227     fileName.ReplaceAll("AliESDs", "Kinematics");
228
229     AliDebug(AliLog::kInfo, Form("Opening %s", fileName.Data()));
230
231     fKineFile = TFile::Open(fileName);
232     if (!fKineFile)
233       return 0;
234   }
235
236   return dynamic_cast<TTree*> (fKineFile->Get(Form("Event%d/TreeK", fChain->GetTree()->GetReadEntry())));
237 }
238
239 void AliSelector::DeleteKinematicsFile()
240 {
241   //
242   // Closes the kinematics file and deletes the pointer.
243   //
244
245   if (fKineFile)
246   {
247     fKineFile->Close();
248     delete fKineFile;
249     fKineFile = 0;
250   }
251 }
252
253 AliHeader* AliSelector::GetHeader()
254 {
255   // Returns header corresponding to current ESD active in fChain
256   // Loads the header from galice.root, the file is identified by replacing "AliESDs" to
257   // "galice" in the file path of the ESD file. This is a hack, to be changed!
258
259   if (!fHeaderFile || !fHeaderTree)
260   {
261     if (!fChain->GetCurrentFile())
262       return 0;
263
264     TString fileName(fChain->GetCurrentFile()->GetName());
265     fileName.ReplaceAll("AliESDs", "galice");
266
267     AliDebug(AliLog::kInfo, Form("Opening %s", fileName.Data()));
268
269     fHeaderFile = TFile::Open(fileName);
270     if (!fHeaderFile)
271       return 0;
272
273     fHeaderTree = dynamic_cast<TTree*> (fHeaderFile->Get("TE"));
274     if (!fHeaderTree)
275       return 0;
276
277     fHeaderTree->SetBranchAddress("Header", &fHeader);
278   }
279
280   fHeaderTree->GetEntry(fChain->GetTree()->GetReadEntry());
281
282   return fHeader;
283 }
284
285 void AliSelector::DeleteHeaderFile()
286 {
287   //
288   // Closes the kinematics file and deletes the pointer.
289   //
290
291   if (fHeaderFile)
292   {
293     fHeaderFile->Close();
294     delete fHeaderFile;
295     fHeaderTree = 0;
296     fHeader = 0;
297   }
298 }
299
300 Bool_t AliSelector::IsPrimaryCharged(TParticle* aParticle, Int_t aTotalPrimaries) const
301 {
302   //
303   // Returns if the given particle is a primary particle
304   // This function or a equivalent should be available in some common place of AliRoot
305   //
306
307   // if the particle has a daughter primary, we do not want to count it
308   if (aParticle->GetFirstDaughter() != -1 && aParticle->GetFirstDaughter() < aTotalPrimaries)
309   {
310     AliDebug(AliLog::kDebug+1, "Dropping particle because it has a daughter among the primaries.");
311     return kFALSE;
312   }
313
314   Int_t pdgCode = TMath::Abs(aParticle->GetPdgCode());
315
316   // skip quarks and gluon
317   if (pdgCode <= 10 || pdgCode == 21)
318   {
319     AliDebug(AliLog::kDebug+1, "Dropping particle because it is a quark or gluon.");
320     return kFALSE;
321   }
322
323   if (strcmp(aParticle->GetName(),"XXX") == 0)
324   {
325     AliDebug(AliLog::kDebug, Form("WARNING: There is a particle named XXX."));
326     return kFALSE;
327   }
328
329   TParticlePDG* pdgPart = aParticle->GetPDG();
330
331   if (strcmp(pdgPart->ParticleClass(),"Unknown") == 0)
332   {
333     AliDebug(AliLog::kDebug, Form("WARNING: There is a particle with an unknown particle class (pdg code %d).", pdgCode));
334     return kFALSE;
335   }
336
337   if (pdgPart->Charge() == 0)
338   {
339     return kFALSE;
340     AliDebug(AliLog::kDebug+1, "Dropping particle because it is not charged.");
341   }
342
343   return kTRUE;
344 }