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