]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCEventHandler.cxx
Additional protection added.
[u/mrichter/AliRoot.git] / STEER / AliMCEventHandler.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 //---------------------------------------------------------------------------------
18 //                          Class AliMCEventHandler
19 // This class gives access to MC truth during the analysis.
20 // Monte Carlo truth is containe in the kinematics tree (produced particles) and 
21 // the tree of reference hits.
22 //      
23 // Origin: Andreas Morsch, CERN, andreas.morsch@cern.ch 
24 //---------------------------------------------------------------------------------
25
26
27
28 #include "AliMCEventHandler.h"
29 #include "AliMCEvent.h"
30 #include "AliPDG.h"
31 #include "AliTrackReference.h"
32 #include "AliHeader.h"
33 #include "AliStack.h"
34 #include "AliLog.h"
35
36 #include <TTree.h>
37 #include <TFile.h>
38 #include <TParticle.h>
39 #include <TString.h>
40 #include <TClonesArray.h>
41 #include <TDirectoryFile.h>
42
43 ClassImp(AliMCEventHandler)
44
45 AliMCEventHandler::AliMCEventHandler() :
46     AliVEventHandler(),
47     fMCEvent(new AliMCEvent()),
48     fFileE(0),
49     fFileK(0),
50     fFileTR(0),
51     fTreeE(0),
52     fTreeK(0),
53     fTreeTR(0),
54     fDirK(0),
55     fDirTR(0),
56     fNEvent(-1),
57     fEvent(-1),
58     fPathName(new TString("./")),
59     fExtension(""),
60     fFileNumber(0),
61     fEventsPerFile(0),
62     fReadTR(kTRUE),
63     fInitOk(kFALSE)
64 {
65   //
66   // Default constructor
67   //
68   // Be sure to add all particles to the PDG database
69   AliPDG::AddParticlesToPdgDataBase();
70 }
71
72 AliMCEventHandler::AliMCEventHandler(const char* name, const char* title) :
73     AliVEventHandler(name, title),
74     fMCEvent(new AliMCEvent()),
75     fFileE(0),
76     fFileK(0),
77     fFileTR(0),
78     fTreeE(0),
79     fTreeK(0),
80     fTreeTR(0),
81     fDirK(0),
82     fDirTR(0),
83     fNEvent(-1),
84     fEvent(-1),
85     fPathName(new TString("./")),
86     fExtension(""),
87     fFileNumber(0),
88     fEventsPerFile(0),
89     fReadTR(kTRUE),
90     fInitOk(kFALSE)
91 {
92   //
93   // Constructor
94   //
95   // Be sure to add all particles to the PDG database
96   AliPDG::AddParticlesToPdgDataBase();
97 }
98 AliMCEventHandler::~AliMCEventHandler()
99
100     // Destructor
101     delete fMCEvent;
102     delete fFileE;
103     delete fFileK;
104     delete fFileTR;
105 }
106
107 Bool_t AliMCEventHandler::Init(Option_t* opt)
108
109     // Initialize input
110     //
111     if (!(strcmp(opt, "proof")) || !(strcmp(opt, "local"))) return kTRUE;
112     //
113     fFileE = TFile::Open(Form("%sgalice.root", fPathName->Data()));
114     if (!fFileE) {
115         AliError(Form("AliMCEventHandler:galice.root not found in directory %s ! \n", fPathName->Data()));
116         fInitOk = kFALSE;
117         return kFALSE;
118     }
119     
120     //
121     // Tree E
122     fFileE->GetObject("TE", fTreeE);
123     // Connect Tree E to the MCEvent
124     fMCEvent->ConnectTreeE(fTreeE);
125     fNEvent = fTreeE->GetEntries();
126     //
127     // Tree K
128     fFileK = TFile::Open(Form("%sKinematics%s.root", fPathName->Data(), fExtension));
129     if (!fFileK) {
130         AliError(Form("AliMCEventHandler:Kinematics.root not found in directory %s ! \n", fPathName));
131         fInitOk = kFALSE;
132         return kFALSE;
133     }
134     
135     fEventsPerFile = fFileK->GetNkeys() - fFileK->GetNProcessIDs();
136     //
137     // Tree TR
138     if (fReadTR) {
139         fFileTR = TFile::Open(Form("%sTrackRefs%s.root", fPathName->Data(), fExtension));
140         if (!fFileTR) {
141             AliError(Form("AliMCEventHandler:TrackRefs.root not found in directory %s ! \n", fPathName->Data()));
142             fInitOk = kFALSE;
143             return kFALSE;
144         }
145     }
146     //
147     // Reset the event number
148     fEvent      = -1;
149     fFileNumber =  0;
150     AliInfo(Form("Number of events in this directory %5d \n", fNEvent));
151     fInitOk = kTRUE;
152     return kTRUE;
153 }
154
155 Bool_t AliMCEventHandler::GetEvent(Int_t iev)
156 {
157     // Load the event number iev
158     //
159     // Calculate the file number
160     if (!fInitOk) return kFALSE;
161     
162     Int_t inew  = iev / fEventsPerFile;
163     if (inew != fFileNumber) {
164         fFileNumber = inew;
165         if (!OpenFile(fFileNumber)){
166             return kFALSE;
167         }
168     }
169     // Folder name
170     char folder[20];
171     sprintf(folder, "Event%d", iev);
172     // TreeE
173     fTreeE->GetEntry(iev);
174     // Tree K
175     fFileK->GetObject(folder, fDirK);
176     if (!fDirK) {
177         AliWarning(Form("AliMCEventHandler: Event #%5d not found\n", iev));
178         return kFALSE;
179     }
180     fDirK ->GetObject("TreeK", fTreeK);
181     // Connect TreeK to MCEvent
182     fMCEvent->ConnectTreeK(fTreeK);
183     //Tree TR 
184     if (fFileTR) {
185         // Check which format has been read
186         fFileTR->GetObject(folder, fDirTR);
187         fDirTR->GetObject("TreeTR", fTreeTR);
188         //
189         // Connect TR to MCEvent
190         fMCEvent->ConnectTreeTR(fTreeTR);
191     }
192     //
193     return kTRUE;
194 }
195
196 Bool_t AliMCEventHandler::OpenFile(Int_t i)
197 {
198     // Open file i
199     if (i > 0) {
200         fExtension = Form("%d", i);
201     } else {
202         fExtension = "";
203     }
204     
205     
206     delete fFileK;
207     fFileK = TFile::Open(Form("%sKinematics%s.root", fPathName->Data(), fExtension));
208     if (!fFileK) {
209         AliError(Form("AliMCEventHandler:Kinematics%s.root not found in directory %s ! \n", fExtension, fPathName->Data()));
210         fInitOk = kFALSE;
211         return kFALSE;
212     }
213     
214     if (fReadTR) {
215         delete fFileTR;
216         fFileTR = TFile::Open(Form("%sTrackRefs%s.root", fPathName->Data(), fExtension));
217         if (!fFileTR) {
218             AliWarning(Form("AliMCEventHandler:TrackRefs%s.root not found in directory %s ! \n", fExtension, fPathName->Data()));
219             fInitOk = kFALSE;
220             return kFALSE;
221         }
222     }
223     
224     fInitOk = kTRUE;
225     return kTRUE;
226 }
227
228 Bool_t AliMCEventHandler::BeginEvent(Long64_t entry)
229
230     // Read the next event
231     if (entry == -1) {
232         fEvent++;
233         entry = fEvent;
234     } else {
235         fEvent = entry;
236     }
237
238     if (entry >= fNEvent) {
239         AliWarning(Form("AliMCEventHandler: Event number out of range %5d %5d\n", entry,fNEvent));
240         return kFALSE;
241     }
242     return GetEvent(entry);
243 }
244
245 Int_t AliMCEventHandler::GetParticleAndTR(Int_t i, TParticle*& particle, TClonesArray*& trefs)
246 {
247     // Retrieve entry i
248     if (!fInitOk) {
249         return 0;
250     } else {
251         return (fMCEvent->GetParticleAndTR(i, particle, trefs));
252     }
253 }
254
255 void AliMCEventHandler::DrawCheck(Int_t i, Int_t search)
256 {
257     // Retrieve entry i and draw momentum vector and hits
258     fMCEvent->DrawCheck(i, search);
259 }
260
261 Bool_t AliMCEventHandler::Notify(const char *path)
262 {
263   // Notify about directory change
264   // The directory is taken from the 'path' argument
265   // Reconnect trees
266     TString fileName(path);
267     if(fileName.Contains("AliESDs.root")){
268         fileName.ReplaceAll("AliESDs.root", "");
269     }
270     else if(fileName.Contains("AliAOD.root")){
271         fileName.ReplaceAll("AliAOD.root", "");
272     }
273     else if(fileName.Contains("galice.root")){
274         // for running with galice and kinematics alone...
275         fileName.ReplaceAll("galice.root", "");
276     }
277     
278     *fPathName = fileName;
279     AliInfo(Form("Notify() Path: %s\n", fPathName->Data()));
280     
281     ResetIO();
282     InitIO("");
283
284     return kTRUE;
285 }
286
287 void AliMCEventHandler::ResetIO()
288 {
289 //  Clear header and stack
290     if (fInitOk) fMCEvent->Clean();
291     
292 // Delete Tree E
293     delete fTreeE; fTreeE = 0;
294
295 // Reset files
296     if (fFileE)  {delete fFileE;  fFileE  = 0;}
297     if (fFileK)  {delete fFileK;  fFileK  = 0;}
298     if (fFileTR) {delete fFileTR; fFileTR = 0;}
299     fExtension="";
300     fInitOk = kFALSE;
301 }
302
303                             
304 Bool_t AliMCEventHandler::FinishEvent()
305 {
306     // Clean-up after each event
307     delete fDirTR;  fDirTR = 0;
308     delete fDirK;   fDirK  = 0;    
309     if (fInitOk) fMCEvent->FinishEvent();
310     return kTRUE;
311 }
312
313 Bool_t AliMCEventHandler::Terminate()
314
315     // Dummy 
316     return kTRUE;
317 }
318
319 Bool_t AliMCEventHandler::TerminateIO()
320
321     // Dummy
322     return kTRUE;
323 }
324     
325
326 void AliMCEventHandler::SetInputPath(const char* fname)
327 {
328     // Set the input path name
329     delete fPathName;
330     fPathName = new TString(fname);
331 }