]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowCommon/AliStarEventReader.cxx
fixes warnings
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / 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 ClassImp(AliStarEventReader)
38
39 //______________________________________________________________________________
40 AliStarEventReader::AliStarEventReader():
41   TObject(),
42   fFileList(NULL),
43   fEventHeader(NULL),
44   fTracks(NULL),
45   fEvent(NULL)
46 {
47   //ctor
48 }
49 //______________________________________________________________________________
50 AliStarEventReader::AliStarEventReader( const char* inputFileDirectory ):
51   TObject(),
52   fFileList(NULL),
53   fEventHeader( new TNtuple("event","event",
54                             "runId:eventNumber:VtxX:VtxY:VtxZ:BField:refMult:centralityId:numberOfPrimaryTracks:numberOfParticles:h0:h1:h2:h3:h4" )),
55   fTracks( new TNtuple("tracks","tracks",
56                        "ID:Charge:Eta:Phi:Pt:Dca:nHits:nHitsFit:nHitsPoss:nHitsDedx:dEdx:nSigElect:nSigPi:nSigK:nSigProton" )),
57   fEvent(new AliStarEvent(1024))
58 {
59   //ctor
60   MakeFileList ( inputFileDirectory ) ;
61 }
62
63 //______________________________________________________________________________
64 AliStarEventReader::~AliStarEventReader()
65 {
66   //dtor
67   delete fEventHeader;
68   fEventHeader = NULL ;
69   delete fTracks;
70   fTracks = NULL ;
71   delete fFileList;
72   fFileList = NULL ;
73   delete fEvent;
74 }
75
76 //______________________________________________________________________________
77 Bool_t AliStarEventReader::GetNextEvent( )
78 {
79   //gets next event
80   static TFile*    nextFile    = NULL ;
81   static TNtuple*  ntData      = NULL ;
82   static Int_t     DoOnce      = 0 ;
83   static Long64_t  NextEntry   = 0 ;
84   static Long64_t  entries     = 0 ;
85   static Long64_t  FileCounter = 0 ;
86
87   if ( DoOnce == 0 )
88   {
89     DoOnce = 1     ;
90     nextFile = (TFile*)  fFileList->First() ;
91     if ( nextFile == 0 ) return false      ;
92     ntData        = (TNtuple*) ( nextFile->Get("NTtracks") ) ;
93     entries       = ntData->GetEntriesFast() ;
94     Int_t columns = ntData->GetNvar() ;
95     if ( columns != 15 )
96     {
97       cout << "Error in reading Ntuple file: no. columns != 15" << endl ;
98       return false ;
99     }
100     FileCounter++ ;
101     cout << "Start New File " << FileCounter << endl ;
102   }
103
104   while ( nextFile )
105   {
106     while ( NextEntry < entries )
107     {
108       Float_t* header = NULL;
109       Int_t numberOfParticles =  0 ;                   // Number of particle tracks in the next event
110       Long64_t HeaderEntry    =  0 ;                   // Store position of Header and Set Flag in case of EOF or error
111       Long64_t SkipEvent      =  0 ;                   // Flag in case of wrong number of tracks in this event
112
113       ////////////////////////////////////////////
114       fEvent->Reset();           //reset the event
115       ////////////////////////////////////////////
116
117       // Search for the first "Event" record
118
119       delete fEventHeader  ;
120       fEventHeader  = NULL ;               // Note that a first 'event'  ntuple was created in constructor
121       fEventHeader  = new TNtuple("event","event",
122                                   "runId:eventNumber:VtxX:VtxY:VtxZ:BField:refMult:centralityId:numberOfPrimaryTracks:numberOfParticles:h0:h1:h2:h3:h4" ) ;
123
124       for ( Long64_t j = NextEntry ; j < entries ; j++ )
125       {
126         Long64_t BytesRead = ntData->GetEntry(j) ;
127         if ( BytesRead < 60 )
128         {
129           cout << "Warning: error in file or EOF " <<  endl ;
130           HeaderEntry = -1 ;
131           break ;
132         }
133         header = ntData->GetArgs() ;
134         if ( (int) header[10] == -1 && (int) header[11] == -1 && (int) header[12] == -1 &&
135              (int) header[13] == -1 && (int) header[14] == -1 )
136         {
137           fEventHeader->Fill(header) ;
138
139           //////////////////////////////////////////////////
140           fEvent->SetParams(header);  //set the event params
141           //////////////////////////////////////////////////
142
143           numberOfParticles = (int) header[9]  ;   // # of particles passing track cuts, thus in ntuple
144           HeaderEntry = j ;
145           break ;
146         }
147         cout << "Warning: no header entries found in this file" << endl ;
148         HeaderEntry = -1 ;
149       }
150
151       if ( HeaderEntry == -1 ) break ;                 // Break out of main loop if I/O error
152
153       // Get subsequent "track" data
154
155       delete fTracks ;
156       fTracks = NULL ;                 // Note that a first 'tracks' ntuple was created in constructor
157       fTracks = new TNtuple("tracks","tracks",
158                             "ID:Charge:Eta:Phi:Pt:Dca:nHits:nHitsFit:nHitsPoss:nHitsDedx:dEdx:nSigElect:nSigPi:nSigK:nSigProton" ) ;
159
160       for ( Long64_t j = HeaderEntry + 1 ; j < HeaderEntry + 1 + numberOfParticles  ; j++ )
161       {
162         Long64_t BytesRead = ntData->GetEntry(j) ;
163         if ( BytesRead < 60 )
164         {
165           cout << "Warning: error in file sequence or EOF" << endl ;
166           NextEntry = -1 ;
167           break ;
168         }
169         header = ntData->GetArgs() ;
170
171         if ( TMath::IsNaN(header[10]) == 1 )
172         {
173           cout << "IsNan ... dEdx will be zeroed out" << endl ;
174           header[10] = 0 ;
175           header[11] = 999 ;
176           header[12] = 999 ;
177           header[13] = 999 ;
178           header[14] = 999 ;
179           cout << header[0]  << " " << header[1]  << " " << header[2]  << " " << header[3]  << " "
180                << header[4]  << " " << header[5]  << " " << header[6]  << " " << header[7]  << " "
181                << header[8]  << " " << header[9]  << " " << header[10] << " " << header[11] << " "
182                << header[12] << " " << header[13] << " " << header[14] << endl ;  // JT test
183         }
184
185         if ( (int) header[10] == -1 && (int) header[11] == -1 && (int) header[12] == -1 &&
186              (int) header[13] == -1 && (int) header[14] == -1 )
187         {
188           cout << "Warning: Header in the wrong place, skipping event" << endl ;
189           SkipEvent = 1 ;
190           NextEntry = j ;          // Skip event and freeze NextEntry counter
191           break ;
192         }
193
194         fTracks->Fill(header) ;
195
196         ////////////////////////////////////////////////////////////////////
197         fEvent->AddTrack( new AliStarTrack(header) );    //add the new track
198         ////////////////////////////////////////////////////////////////////
199
200         NextEntry = j+1 ;
201       }
202       if ( NextEntry == -1 ) break ;      // Bad record in file, go to next file in fFileList
203       if ( SkipEvent ==  1 ) continue ;   // Bad event, go to next event in this file
204       return true ;                       // Success: Event read OK, note unusual location for a successful return
205     }
206
207     NextEntry = 0 ; // this entry goes before nextFile
208     nextFile = (TFile*) fFileList->After(nextFile) ;
209     if ( nextFile == 0 ) break ;
210     if (ntData) delete ntData;
211     ntData        = (TNtuple*) ( nextFile->Get("NTtracks") ) ;
212     entries       = ntData->GetEntriesFast() ;
213     Int_t columns = ntData->GetNvar() ;
214     if ( columns != 15 )
215     {
216       cout << "Error in reading Ntuple file: no. columns != 15" << endl ;
217       break ;
218     }
219     FileCounter++ ;
220     cout << "Start New File " << FileCounter << endl ;
221   }
222
223   return false ;  // Failure: Error or EOF
224 }
225
226 //______________________________________________________________________________
227 Bool_t AliStarEventReader::AcceptEvent( AliStarEvent* event )
228 {
229   // Cut parameters for each event
230
231   const Float_t VertexXMin  =  -1.0 ;  // cm
232   const Float_t VertexXMax  =   1.0 ;  // cm
233   const Float_t VertexYMin  =  -1.0 ;  // cm
234   const Float_t VertexYMax  =   1.0 ;  // cm
235   const Float_t VertexZMin  = -30.0 ;  // cm
236   const Float_t VertexZMax  =  30.0 ;  // cm
237   const Int_t   MultMin     =    10 ;  // Note: this is a cut on refMult which is not your ordinary multiplicity.
238   const Int_t   MultMax     =  1000 ;  // Refmult corresponds to Zhangbu's ascii table of cuts for centrality bins.
239   const Int_t   BlackEvent  =  3000 ;  // Maximum number of primary tracks allowed in one event (more is an error)
240
241   Int_t refMult               = event->GetRefMult();  // Reference Multiplicity for Centrality determination
242   Int_t NumberOfPrimaryTracks = event->GetNumberOfPrimaryTracks();  // Number of primary tracks in full event (not just those in main TPC)
243
244   if ( refMult < MultMin          || refMult > MultMax )                   return false ;
245   if ( NumberOfPrimaryTracks <= 0 || NumberOfPrimaryTracks > BlackEvent )  return false ;
246
247   // Cut on Vertex location
248   Float_t vertex[3] = { event->GetVtxX(),
249                         event->GetVtxY(),
250                         event->GetVtxZ() };
251
252   if ( vertex[0] < VertexXMin || vertex[0] > VertexXMax )    return false ;  // Skip events that fall outside Vtx cuts
253   if ( vertex[1] < VertexYMin || vertex[1] > VertexYMax )    return false ;
254   if ( vertex[2] < VertexZMin || vertex[2] > VertexZMax )    return false ;
255
256   return true ;
257 }
258
259 //______________________________________________________________________________
260 Bool_t AliStarEventReader::AcceptTrack( AliStarTrack* track )
261 {
262   // Cut Parameters for individual tracks
263
264   const Float_t dcaCut      =   3.0  ;      // cm
265   const Float_t PtMin       =   0.15 ;      // GeV
266   const Float_t PtMax       =   2.0  ;      // GeV
267   const Float_t EtaMin      =  -1.1  ;
268   const Float_t EtaMax      =   1.1  ;
269   const Float_t FitRatio    =  0.52  ;      // Number of hits over number of hits possible
270   const Int_t   nHitMin     =    15  ;      // 15 is typical but sometimes goes as high as 25
271   const Int_t   nHitMax     =   100  ;      // 45 pad rows in the TPC and so anything bigger than 45+Silicon is infinite
272   const Int_t   nHitPossMin =     5  ;      // Don't bother to fit tracks if # possible hits is too low, also protect / 0
273
274   if ( track->GetDCA() >  dcaCut ) return false ;                 // magnitude of 3D DCA for global tracks
275   if ( track->GetNHitsPoss() <  nHitPossMin ||
276        track->GetNHitsPoss() >  nHitMax ) return false ;   // Minimum number of Possible hits, see above.
277   if ( track->GetEta() <  EtaMin  || track->GetEta() >  EtaMax  ) return false ;
278   if ( track->GetPt()  <  PtMin   || track->GetEta() >  PtMax   ) return false ;
279   if ( track->GetNHitsFit() < nHitMin || track->GetNHitsFit() > nHitMax ) return false ;
280   if ( ((1.0*track->GetNHitsFit()) /
281         (1.0*track->GetNHitsPoss())) < FitRatio )  return false ;
282
283   return true ;
284 }
285
286 //______________________________________________________________________________
287 Int_t AliStarEventReader::ParticleID( AliStarTrack* track )
288 {
289   // Note: This is a very simple PID selection scheme.  More elaborate methods (with multiple cuts) may be required.
290   // When you *are* using dEdx information, you must chose a finite number of good Dedx hits ... but the limit should
291   // be about 2/3 of nHitsMin.  This is because some clusters do not form good dEdx hits due to track
292   // merging, etc., and so nHitsDedx is always less than nHitsFit.  A rule of thumb says ~2/3 ratio.
293
294   Int_t ID = 0 ;
295
296   const Int_t   nHitDedxMin =    15  ;       // 10 to 20 is typical.  nHitDedxMin is often chosen to be about 2/3 of nHitMin.
297   const Float_t nSigmaPID   =    2.0 ;       // Number of Sigma cut to apply to PID bands
298
299   // Test on Number of dE/dx hits required, return 0 if not enough hits
300   if ( track->GetNHitsDedx() <  nHitDedxMin ) return ID;
301
302   // Begin PID
303
304   if ( TMath::Abs( track->GetNSigElect() ) >= nSigmaPID )
305   {
306     if ( TMath::Abs( track->GetNSigK()  ) <= nSigmaPID )
307     {
308       ID = 321  ;
309     }
310     if ( TMath::Abs( track->GetNSigProton()  ) <= nSigmaPID )
311     {
312       ID = 2212 ;
313     }
314     if ( TMath::Abs( track->GetNSigPi()  ) <= nSigmaPID )
315     {
316       ID = 211  ;
317     }
318   }
319
320   // Pion is the default in case of ambiguity because it is most abundent. Don't re-arrange order, above.
321
322   return ID ;
323 }
324
325 //______________________________________________________________________________
326 Int_t AliStarEventReader::Centrality( Int_t referenceMultiplicity )
327
328 {
329   // Note Carefully:  Centrality is based on refMult.  This is the 'reference' multiplicity that is measured
330   // independpently from the TPC.  Selecting the centrality bins according to the refMult is something that
331   // is calibrated for each year and each run.  You can get the basic information off the web:
332   // For Example .... http://www.star.bnl.gov/protected/common/common2004/trigger2004/200gev/200gevFaq.html
333   // An index pointing to FAQs, Trigger and Centrality data, for all years, is available at:
334   // http://www.star.bnl.gov/public/all
335   //
336   // Note: Add 0.5 to the (int) that is returned by centrality() when using it as an argument for a histogram
337   // that expects (float) or (double) as input parameters.  This will place the data point in the center of
338   // the bin, avoids ambiguities, and is best for plotting scatter plots and contour plots.
339   // For example histogram2D[1] -> Fill ( (float)CentralityID + 0.5 , SumData )   ;
340   //
341   // The refMult quoted in the Centrality bins array is the lower limit on refMult
342
343
344   Int_t   CentralityBins  [] = { 14 , 31 , 57 , 96 , 150 , 222 , 319 , 441 , 520 , 1000 } ;  // Run4 200 GeV
345   Int_t   MiddleBinID     [] = {  0 ,  1 ,  2 ,  3 ,   4 ,   5 ,   6 ,   7 ,   8 ,    9 } ;  // ID Number
346   //Info  MiddleBinPercent[] = { 85., 75., 65., 55.,  45.,  35.,  25.,  15., 7.5 ,  2.5 } ;  // Percent
347   Int_t   myCentrality  ;
348
349   if      ( referenceMultiplicity < CentralityBins[0] )
350   {
351     myCentrality = MiddleBinID[0] ;
352   }
353   else if ( referenceMultiplicity < CentralityBins[1] )
354   {
355     myCentrality = MiddleBinID[1] ;
356   }
357   else if ( referenceMultiplicity < CentralityBins[2] )
358   {
359     myCentrality = MiddleBinID[2] ;
360   }
361   else if ( referenceMultiplicity < CentralityBins[3] )
362   {
363     myCentrality = MiddleBinID[3] ;
364   }
365   else if ( referenceMultiplicity < CentralityBins[4] )
366   {
367     myCentrality = MiddleBinID[4] ;
368   }
369   else if ( referenceMultiplicity < CentralityBins[5] )
370   {
371     myCentrality = MiddleBinID[5] ;
372   }
373   else if ( referenceMultiplicity < CentralityBins[6] )
374   {
375     myCentrality = MiddleBinID[6] ;
376   }
377   else if ( referenceMultiplicity < CentralityBins[7] )
378   {
379     myCentrality = MiddleBinID[7] ;
380   }
381   else if ( referenceMultiplicity < CentralityBins[8] )
382   {
383     myCentrality = MiddleBinID[8] ;
384   }
385   else
386   {
387     myCentrality = MiddleBinID[9] ;
388   }
389
390   return myCentrality ;
391 }
392
393 //______________________________________________________________________________
394 void AliStarEventReader::PrintEventHeader ( )
395 {
396   // TNtuple* event: names are documented in the next 2 lines
397   // event  = new TNtuple("event","event",
398   //   "runId:eventNumber:VtxX:VtxY:VtxZ:BField:refMult:centralityId:numberOfPrimaryTracks:numberOfParticles:h0:h1:h2:h3:h4" ) ;
399   //
400   Float_t  *eventhd ;
401   eventhd = fEventHeader->GetArgs() ;
402
403   Int_t   runId                  = (int)   eventhd[0]  ;
404   Int_t   eventNumber            = (int)   eventhd[1]  ;
405   Float_t   primaryVertexPosition[3] = { (float) eventhd[2],  // (cm)
406                                          (float) eventhd[3],  // (cm)
407                                          (float) eventhd[4] };  // (cm)
408   Float_t magneticField          = (float) eventhd[5]  ;  // kilogauss
409   Int_t   refMult                = (int)   eventhd[6]  ;  // Raw Mult into 0.5 unit: a relative number, not total Mult.
410   Int_t   centralityId           = (int)   eventhd[7]  ;  // STAR centrality bin # based on refMult
411   Int_t   numberOfPrimaryTracks  = (int)   eventhd[8]  ;  // # of primaries, including FTPC tracks which are not in ntuple
412   Int_t   numberOfParticles      = (int)   eventhd[9]  ;  // # of particles passing track cuts, thus in ntuple
413
414   printf("\n  runId eventNo    VtxX    VtxY    VtxZ  MagFld  refMult  CentBin  #PrimeTrak  #Tracks \n") ;
415   printf("%7d  %6d %7.4f %7.4f %7.3f  %6.3f   %6d   %6d      %6d   %6d \n\n",
416          runId, eventNumber, primaryVertexPosition[0], primaryVertexPosition[1], primaryVertexPosition[2],
417          magneticField, refMult, centralityId, numberOfPrimaryTracks, numberOfParticles ) ;
418
419   //Int_t newCentralityID ;
420   //newCentralityID = Centrality( refMult) ;              // Should be the same as "centralityID" from tape
421   //cout << "Test: should be the same " << newCentralityID << " " << centralityId << endl ;  // JT test
422 }
423
424 //______________________________________________________________________________
425 void AliStarEventReader::PrintTrack ( Int_t counter )
426 {
427   // TNtuple* track: names are documented in the next 2 lines
428   // tracks = new TNtuple("tracks","tracks",
429   //   "ID:Charge:Eta:Phi:Pt:Dca:nHits:nHitsFit:nHitsPoss:nHitsDedx:dEdx:nSigElect:nSigPi:nSigK:nSigProton" ) ;
430   //
431   if ( counter == 0 )
432   {
433     printf(
434       "    id charge     eta     phi      pt     dca  nHits  nFit nPoss ndEdx   dEdx nSigElec nSigPi  nSigK nSigPr\n") ;
435   }
436   Float_t* data = fTracks -> GetArgs()       ;  // Extract data from the track
437   Int_t   id             = (int)   data[0]   ;  // id - a unique integer for each track in this event
438   Int_t   charge         = (int)   data[1]   ;  // +1 or -1
439   Float_t eta            = (float) data[2]   ;  // Pseudo-rapidity at the vertex
440   Float_t phi            = (float) data[3]   ;  // Phi emission angle at the vertexcd
441   Float_t pt             = (float) data[4]   ;  // Pt of the track at the vertex
442   Float_t dca            = (float) data[5]   ;  // Magnitude of 3D DCA to vertex
443   Int_t   nHits          = (int)   data[6]   ;  // Number of clusters available to the Kalman Filter
444   Int_t   nHitsFit       = (int)   data[7]   ;  // Number of clusters used in the fit (after cuts)
445   Int_t   nHitsPoss      = (int)   data[8]   ;  // Number of possible cluster on track (theoretical max)
446   Int_t   nHitsDedx      = (int)   data[9]   ;  // Number of clusters used in the fit (after dEdx cuts)
447   Float_t dEdx           = 1.e6*(float)data[10]  ;  // Measured dEdx (Note: GeV/cm so convert to keV/cm!!)
448   Float_t nSigmaElectron = (float) data[11]  ;  // Number of sigma from electron Bethe-Bloch curve
449   Float_t nSigmaPion     = (float) data[12]  ;  // Number of sigma from pion Bethe-Bloch curve
450   Float_t nSigmaKaon     = (float) data[13]  ;  // Number of sigma from kaon Bethe-Bloch curve
451   Float_t nSigmaProton   = (float) data[14]  ;  // Number of sigma from proton Bethe-Bloch curve
452
453   // Alternative way to access the data
454   nHitsPoss      = (int) ( fTracks->GetLeaf("nHitsPoss")->GetValue() ) ;  // Note alternative method to retrieve data
455   // Using the definition of the original NTuple
456   // TrackTuple      = new TNtuple("NTtracks","NTtracks",
457   // "ID:Charge:Eta:Phi:Pt:Dca:nHits:nHitsFit:nHitsPoss:nHitsDedx:dEdx:nSigElect:nSigPi:nSigK:nSigProton" )
458
459   printf("%6d %4d   %7.3f %7.3f %7.3f %7.4f %6d %5d %5d %5d %6.2f   %6.2f %6.2f %6.2f %6.2f \n",
460          id, charge, eta, phi, pt, dca, nHits, nHitsFit, nHitsPoss, nHitsDedx, dEdx,
461          nSigmaElectron, nSigmaPion, nSigmaKaon, nSigmaProton ) ;
462 }
463
464 //______________________________________________________________________________
465 Bool_t AliStarEventReader::MakeFileList ( const char* input )
466 {
467   //get the files to process
468   TString inputstring(input);
469   inputstring = inputstring.Strip(TString::kBoth);
470   TSystemFile inputfile(inputstring.Data(),"");
471   if (inputfile.IsDirectory())
472     return MakeFileListFromDir(inputstring.Data());
473   else
474     return MakeFileListFromFile(inputstring.Data());
475 }
476
477 //______________________________________________________________________________
478 Bool_t AliStarEventReader::MakeFileListFromDir ( const char* inputFileDirectory )
479 {
480   //get the files to process
481   Int_t  Count        = 0 ;
482   static Int_t DoOnce = 0 ;
483   fFileList =  new TList() ;
484   void*   directory = gSystem->OpenDirectory(inputFileDirectory) ;
485   const char* entry = gSystem->GetDirEntry(directory) ;
486
487   if ( entry == 0 )
488   {
489     cout << endl <<  "Error: \"" << inputFileDirectory << "\" does not exist" << endl << endl ;
490     return false ;
491   }
492   else cout << endl ;
493
494   while(entry != 0)
495   {
496     int len = strlen(entry);
497     if( len >= 5 && strcmp( &entry[len - 5], ".root" ) == 0 )
498     {
499       TString fileName ;
500       fileName = inputFileDirectory ;
501       if( !fileName.EndsWith("/") ) fileName += "/" ;
502       fileName += entry;
503       fFileList->Add ( TFile::Open(fileName) ) ;
504       if ( DoOnce == 0 )
505       {
506         cout << "Add: " << fileName << endl ;
507         DoOnce = 1 ;
508       }
509       Count ++ ;
510     }
511     entry = gSystem->GetDirEntry(directory) ;
512   }
513
514   cout << "Add: " << Count-1 << " more file(s) from this directory for a total of " << Count << " files." << endl ;
515   cout << "Finished creating file list ... preparing to open first file." << endl << endl ;
516   return true ;
517 }
518
519 //______________________________________________________________________________
520 Bool_t AliStarEventReader::MakeFileListFromFile ( const char* inputFile )
521 {
522   //get the files to process, from a text file, one file per line
523   if (!fFileList) fFileList=new TList();
524   ifstream filein;
525   filein.open(inputFile);
526   if (!filein.good()) 
527   {
528     printf("problem reading the file list \"%s\"\n",inputFile);
529     return kFALSE;
530   }
531   TString line;
532   while (filein.good())
533   {
534     printf("opening file: ");
535     line.ReadLine(filein);
536     if (line.Length() == 0) continue;
537     TFile* file = TFile::Open(line.Data());
538     if (!file) 
539     {
540       printf("problem opening file \"%s\"\n",line.Data());
541       continue;
542     }
543     fFileList->Add(file);
544     printf("%s\n",line.Data());
545   }
546   if (fFileList->GetEntries()>0) return kTRUE;
547   return kFALSE;
548 }
549