]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Base/AliStarEventReader.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWG / FLOW / Base / AliStarEventReader.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 // Title:   Class for ccessing data from STAR NTuples
17 //          produces data encapsulated in AliStarEvent and AliStarTrack classes
18 //
19 // Origin:  Jim Thomas,        jhthomas@lbl.gov
20 //          Mikolaj Krzewicki, mikolaj.krzewicki@cern.ch
21
22 #include <Riostream.h>
23
24 #include <TSystem.h>
25 #include <TSystemFile.h>
26 #include <TFile.h>
27 #include <TList.h>
28 #include <TLeaf.h>
29 #include <TMath.h>
30 #include <TNtuple.h>
31 #include <TString.h>
32
33 #include "AliStarEventReader.h"
34 #include "AliStarEvent.h"
35 #include "AliStarTrack.h"
36
37 using std::cout;
38 using std::endl;
39 using std::ifstream;
40 ClassImp(AliStarEventReader)
41
42 //______________________________________________________________________________
43 AliStarEventReader::AliStarEventReader():
44   TObject(),
45   fFileList(NULL),
46   fEvent(NULL)
47 {
48   //ctor
49 }
50 //______________________________________________________________________________
51 AliStarEventReader::AliStarEventReader( const char* inputFileDirectory ):
52   TObject(),
53   fFileList(NULL),
54   fEvent(new AliStarEvent(1024))
55 {
56   //ctor
57   MakeFileList ( inputFileDirectory ) ;
58 }
59
60 //______________________________________________________________________________
61 AliStarEventReader::~AliStarEventReader()
62 {
63   //dtor
64   delete fFileList;
65   delete fEvent;
66 }
67
68 //______________________________________________________________________________
69 Bool_t AliStarEventReader::GetNextEvent( )
70 {
71   //gets next event
72   static TFile*    nextFile    = NULL ;
73   static TNtuple*  ntData      = NULL ;
74   static Int_t     doOnce      = 0 ;
75   static Long64_t  nextEntry   = 0 ;
76   static Long64_t  entries     = 0 ;
77   static Long64_t  fileCounter = 0 ;
78
79   if ( doOnce == 0 )
80   {
81     doOnce = 1     ;
82     nextFile = (TFile*)  fFileList->First() ;
83     if ( nextFile == 0 ) return false      ;
84     ntData        = (TNtuple*) ( nextFile->Get("NTtracks") ) ;
85     entries       = ntData->GetEntriesFast() ;
86     Int_t columns = ntData->GetNvar() ;
87     if ( columns != 15 )
88     {
89       cout << "Error in reading Ntuple file: no. columns != 15" << endl ;
90       return false ;
91     }
92     fileCounter++ ;
93     cout << "Start New File " << fileCounter << endl ;
94   }
95
96   while ( nextFile )
97   {
98     while ( nextEntry < entries )
99     {
100       Float_t* header = NULL;
101       Int_t numberOfParticles =  0 ;                   // Number of particle tracks in the next event
102       Long64_t headerEntry    =  0 ;                   // Store position of Header and Set Flag in case of EOF or error
103       Long64_t skipEvent      =  0 ;                   // Flag in case of wrong number of tracks in this event
104
105       fEvent->Reset();           //reset the event
106
107       // Search for the first "Event" record
108
109       for ( Long64_t j = nextEntry ; j < entries ; j++ )
110       {
111         Long64_t BytesRead = ntData->GetEntry(j) ;
112         if ( BytesRead < 60 )
113         {
114           cout << "Warning: error in file or EOF " <<  endl ;
115           headerEntry = -1 ;
116           break ;
117         }
118         header = ntData->GetArgs() ;
119         if ( (int) header[10] == -1 && (int) header[11] == -1 && (int) header[12] == -1 &&
120              (int) header[13] == -1 && (int) header[14] == -1 )
121         {
122           fEvent->SetParams(header);  //set the event params
123
124           numberOfParticles = (int) header[9]  ;   // # of particles passing track cuts, thus in ntuple
125           headerEntry = j ;
126           break ;
127         }
128         cout << "Warning: no header entries found in this file" << endl ;
129         headerEntry = -1 ;
130       }
131
132       if ( headerEntry == -1 ) break ;                 // Break out of main loop if I/O error
133
134       // Get subsequent "track" data
135       for ( Long64_t j = headerEntry + 1 ; j < headerEntry + 1 + numberOfParticles  ; j++ )
136       {
137         Long64_t BytesRead = ntData->GetEntry(j) ;
138         if ( BytesRead < 60 )
139         {
140           cout << "Warning: error in file sequence or EOF" << endl ;
141           nextEntry = -1 ;
142           break ;
143         }
144         header = ntData->GetArgs() ;
145
146         if ( TMath::IsNaN(header[10]) == 1 )
147         {
148           cout << "IsNan ... dEdx will be zeroed out" << endl ;
149           header[10] = 0 ;
150           header[11] = 999 ;
151           header[12] = 999 ;
152           header[13] = 999 ;
153           header[14] = 999 ;
154           cout << header[0]  << " " << header[1]  << " " << header[2]  << " " << header[3]  << " "
155                << header[4]  << " " << header[5]  << " " << header[6]  << " " << header[7]  << " "
156                << header[8]  << " " << header[9]  << " " << header[10] << " " << header[11] << " "
157                << header[12] << " " << header[13] << " " << header[14] << endl ;  // JT test
158         }
159
160         if ( (int) header[10] == -1 && (int) header[11] == -1 && (int) header[12] == -1 &&
161              (int) header[13] == -1 && (int) header[14] == -1 )
162         {
163           cout << "Warning: Header in the wrong place, skipping event" << endl ;
164           skipEvent = 1 ;
165           nextEntry = j ;          // Skip event and freeze nextEntry counter
166           break ;
167         }
168
169         fEvent->AddTrack( new AliStarTrack(header) );    //add the new track
170
171         nextEntry = j+1 ;
172       }
173       if ( nextEntry == -1 ) break ;      // Bad record in file, go to next file in fFileList
174       if ( skipEvent ==  1 ) continue ;   // Bad event, go to next event in this file
175       return true ;                       // Success: Event read OK, note unusual location for a successful return
176     }
177
178     nextEntry = 0 ; // this entry goes before nextFile
179     nextFile = (TFile*) fFileList->After(nextFile) ;
180     if ( nextFile == 0 ) break ;
181     if (ntData) delete ntData;
182     ntData        = (TNtuple*) ( nextFile->Get("NTtracks") ) ;
183     entries       = ntData->GetEntriesFast() ;
184     Int_t columns = ntData->GetNvar() ;
185     if ( columns != 15 )
186     {
187       cout << "Error in reading Ntuple file: no. columns != 15" << endl ;
188       break ;
189     }
190     fileCounter++ ;
191     cout << "Start New File " << fileCounter << endl ;
192   }
193
194   return false ;  // Failure: Error or EOF
195 }
196
197 //______________________________________________________________________________
198 Bool_t AliStarEventReader::MakeFileList ( const char* input )
199 {
200   //get the files to process
201   TString inputstring(input);
202   inputstring = inputstring.Strip(TString::kBoth);
203   TSystemFile inputfile(inputstring.Data(),"");
204   if (inputfile.IsDirectory())
205     return MakeFileListFromDir(inputstring.Data());
206   else
207     return MakeFileListFromFile(inputstring.Data());
208 }
209
210 //______________________________________________________________________________
211 Bool_t AliStarEventReader::MakeFileListFromDir ( const char* inputFileDirectory )
212 {
213   //get the files to process
214   Int_t  Count        = 0 ;
215   static Int_t doOnce = 0 ;
216   fFileList =  new TList() ;
217   void*   directory = gSystem->OpenDirectory(inputFileDirectory) ;
218   const char* entry = gSystem->GetDirEntry(directory) ;
219
220   if ( entry == 0 )
221   {
222     cout << endl <<  "Error: \"" << inputFileDirectory << "\" does not exist" << endl << endl ;
223     return false ;
224   }
225   else cout << endl ;
226
227   while(entry != 0)
228   {
229     int len = strlen(entry);
230     if( len >= 5 && strcmp( &entry[len - 5], ".root" ) == 0 )
231     {
232       TString fileName ;
233       fileName = inputFileDirectory ;
234       if( !fileName.EndsWith("/") ) fileName += "/" ;
235       fileName += entry;
236       fFileList->Add ( TFile::Open(fileName) ) ;
237       if ( doOnce == 0 )
238       {
239         cout << "Add: " << fileName << endl ;
240         doOnce = 1 ;
241       }
242       Count ++ ;
243     }
244     entry = gSystem->GetDirEntry(directory) ;
245   }
246
247   cout << "Add: " << Count-1 << " more file(s) from this directory for a total of " << Count << " files." << endl ;
248   cout << "Finished creating file list ... preparing to open first file." << endl << endl ;
249   return true ;
250 }
251
252 //______________________________________________________________________________
253 Bool_t AliStarEventReader::MakeFileListFromFile ( const char* inputFile )
254 {
255   //get the files to process, from a text file, one file per line
256   if (!fFileList) fFileList=new TList();
257   ifstream filein;
258   filein.open(inputFile);
259   if (!filein.good()) 
260   {
261     printf("problem reading the file list \"%s\"\n",inputFile);
262     return kFALSE;
263   }
264   TString line;
265   while (filein.good())
266   {
267     printf("opening file: ");
268     line.ReadLine(filein);
269     if (line.Length() == 0) continue;
270     TFile* file = TFile::Open(line.Data());
271     if (!file) 
272     {
273       printf("problem opening file \"%s\"\n",line.Data());
274       continue;
275     }
276     fFileList->Add(file);
277     printf("%s\n",line.Data());
278   }
279   if (fFileList->GetEntries()>0) return kTRUE;
280   return kFALSE;
281 }
282