]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - EVGEN/AliGenReaderTreeK.cxx
Protection in copy constructor
[u/mrichter/AliRoot.git] / EVGEN / AliGenReaderTreeK.cxx
... / ...
CommitLineData
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// Realisation of AliGenReader to be used with AliGenExtFile
19// It reads events from a kinematics TreeK.
20// NextEvent() is used to loop over events
21// and NextParticle() to loop over particles.
22// Author: andreas.morsch@cern.ch
23//
24#include <TFile.h>
25#include <TTree.h>
26#include <TParticle.h>
27#include <TObjString.h>
28#include <TObjArray.h>
29
30#include "AliGenReaderTreeK.h"
31#include "AliHeader.h"
32#include "AliRun.h"
33#include "AliStack.h"
34#include "AliRunLoader.h"
35
36ClassImp(AliGenReaderTreeK)
37
38const TString AliGenReaderTreeK::fgkEventFolderName("GenReaderTreeK");
39
40AliGenReaderTreeK::AliGenReaderTreeK():
41 AliGenReader(),
42 fNcurrent(0),
43 fNparticle(0),
44 fNp(0),
45 fInRunLoader(0),
46 fBaseFile(0),
47 fStack(0),
48 fOnlyPrimaries(kFALSE),
49 fDirs(0x0),
50 fCurrentDir(0)
51{
52// Default constructor
53}
54
55AliGenReaderTreeK::AliGenReaderTreeK(const AliGenReaderTreeK &reader):
56 AliGenReader(reader),
57 fNcurrent(0),
58 fNparticle(0),
59 fNp(0),
60 fInRunLoader(0),
61 fBaseFile(0),
62 fStack(0),
63 fOnlyPrimaries(kFALSE),
64 fDirs(0x0),
65 fCurrentDir(0)
66{
67 reader.Copy(*this);
68}
69
70
71AliGenReaderTreeK::~AliGenReaderTreeK()
72{
73// Destructor
74 delete fInRunLoader;//it cleans all the loaded data
75 delete fDirs;
76}
77
78void AliGenReaderTreeK::Init()
79{
80// Initialization
81// Connect base file and file to read from
82
83 TTree *ali = AliRunLoader::Instance()->TreeE();
84 if (ali) {
85 fBaseFile = ali->GetCurrentFile();
86 } else {
87 printf("\n Warning: Basefile cannot be found !\n");
88 }
89 //if (!fFile) fFile = new TFile(fFileName);
90 if (fInRunLoader == 0x0)
91 {
92 fInRunLoader = AliRunLoader::Open((GetDirName(fCurrentDir++)+"/")+fFileName,fgkEventFolderName);
93 fInRunLoader->LoadHeader();
94 fInRunLoader->LoadKinematics("READ");
95 }
96}
97
98Int_t AliGenReaderTreeK::NextEvent()
99{
100// Read the next event
101// cd to file with old kine tree
102 if (!fBaseFile) Init();
103// Get next event
104
105 if (fNcurrent >= fInRunLoader->GetNumberOfEvents())
106 {
107 if (fCurrentDir >= fDirs->GetEntries())
108 {
109 Warning("NextEvent","No more events");
110 return 0;
111 }
112 delete fInRunLoader;
113 fInRunLoader = AliRunLoader::Open((GetDirName(fCurrentDir++)+"/")+fFileName,fgkEventFolderName);
114 fInRunLoader->LoadHeader();
115 fInRunLoader->LoadKinematics("READ");
116 fNcurrent = 0;
117 }
118 fInRunLoader->GetEvent(fNcurrent);
119 fStack = fInRunLoader->Stack();
120
121// cd back to base file
122 fBaseFile->cd();
123//
124 fNcurrent++;
125 fNparticle = 0;
126 fNp = fStack->GetNtrack();
127 printf("\n Next event contains %d particles", fNp);
128//
129 return fNp;
130}
131
132TParticle* AliGenReaderTreeK::NextParticle()
133{
134//Return next particle
135 TParticle* part = GetParticle(fNparticle++);
136 if (part == 0x0) return 0x0;
137 //if only primaries are to be read, and this particle is not primary enter loop
138 if (fOnlyPrimaries && ( part->GetFirstMother() > -1) )
139 for (;;)
140 { //look for a primary
141 part = GetParticle(fNparticle++);
142 if (part == 0x0) return 0x0;
143 if (part->GetFirstMother() == -1) return part;
144 }
145
146 return part;
147}
148
149void AliGenReaderTreeK::RewindEvent()
150{
151 // Go back to the first particle of the event
152 fNparticle = 0;
153}
154
155
156AliGenReaderTreeK& AliGenReaderTreeK::operator=(const AliGenReaderTreeK& rhs)
157{
158// Assignment operator
159 rhs.Copy(*this);
160 return *this;
161}
162
163void AliGenReaderTreeK::Copy(TObject&) const
164{
165 //
166 // Copy
167 //
168 Fatal("Copy","Not implemented!\n");
169}
170
171
172
173TString& AliGenReaderTreeK::GetDirName(Int_t entry)
174 {
175// Get the current directory name
176
177 TString* retval;//return value
178 if (fDirs == 0x0)
179 {
180 retval = new TString(".");
181 return *retval;
182 }
183
184 if ( (entry>fDirs->GetEntries()) || (entry<0))//if out of bounds return empty string
185 { //note that entry==0 is accepted even if array is empty (size=0)
186 Error("GetDirName","Name out of bounds");
187 retval = new TString();
188 return *retval;
189 }
190
191 if (fDirs->GetEntries() == 0)
192 {
193 retval = new TString(".");
194 return *retval;
195 }
196
197 TObjString *dir = dynamic_cast<TObjString*>(fDirs->At(entry));
198 if(dir == 0x0)
199 {
200 Error("GetDirName","Object in TObjArray is not a TObjString or its descendant");
201 retval = new TString();
202 return *retval;
203 }
204 if (gDebug > 0) Info("GetDirName","Returned ok %s",dir->String().Data());
205 return dir->String();
206 }
207
208void AliGenReaderTreeK::AddDir(const char* dirname)
209{
210 //adds a directory to the list of directories where data are looked for
211 if(fDirs == 0x0)
212 {
213 fDirs = new TObjArray();
214 fDirs->SetOwner(kTRUE);
215 }
216 TObjString *odir= new TObjString(dirname);
217 fDirs->Add(odir);
218}
219
220TParticle* AliGenReaderTreeK::GetParticle(Int_t i)
221 {
222 if (fStack && i<fNp) return fStack->Particle(i);
223 return 0x0;
224 }