StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFileIter.cxx
1 // @(#)root/table:$Name: $:$Id: StFileIter.cxx,v 1.6 2009/03/27 17:40:16 fine Exp $
2 // Author: Valery Fine(fine@bnl.gov) 01/03/2001
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2004, Rene Brun and Fons Rademakers. *
6  * Copyright (C) 2001 [BNL] Brookhaven National Laboratory. *
7  * All rights reserved. *
8  * *
9  * For the licensing terms see $ROOTSYS/LICENSE. *
10  * For the list of contributors see $ROOTSYS/README/CREDITS. *
11  *************************************************************************/
12 
14 // //
15 // Class to iterate (read / write ) the events written to TFile. //
16 // The event is supposed to assign an unique ID in form of //
17 // //
18 // TKey <event Id> ::= eventName "." run_number "." event_number //
19 // //
20 // and stored as the TKey name of the object written //
21 // //
22 // /////// ////////// //////// /////// //////
23 //
24 // void TesStFileIter(){
25 // // This macros tests the various methods of StFileIter class.
26 // gSystem->Load("libTable");
27 //
28 // //First create simple ROOT file
29 // TDataSet *ds = new TDataSet("event");
30 // TObject *nextObject = 0;
31 // TRandom run;
32 // TRandom event;
33 // {
34 // StFileIter outSet("test.root","RECREATE");
35 // UInt_t totalEvent = 10;
36 // UInt_t runNumber = 20010301;
37 // Int_t i=0;
38 // Int_t j=0;
39 // for (;j < 10;j++) {
40 // for (i = 1;i<totalEvent;i++) {
41 // outSet.NextEventPut(ds,UInt_t(i),UInt_t(runNumber+j+10*run.Rndm()-5));
42 // }
43 // }
44 // }
45 // printf(" ----------------------> TFile has been created <--------------------\n");
46 // TFile *f = new TFile("test.root");
47 // StFileIter readObj(f);
48 // // the number of the object available directly from "MyDataSet.root"
49 // Int_t size = readObj.TotalKeys();
50 // printf(" The total number of the objects: %d\n",size);
51 //
52 // //-----------------------------------------------------------------------
53 // // Loop over all objects, read them in to memory one by one
54 //
55 // printf(" -- > Loop over all objects, read them in to memory one by one < -- \n");
56 // for( readObj = 0; int(readObj) < size; ++readObj){
57 // nextObject = *readObj;
58 // printf(" %d bytes of the object \"%s\" of class \"%s\" written with TKey \"%s\" has been read from file\n"
59 // ,readObj.GetObjlen()
60 // ,nextObject->GetName()
61 // ,nextObject->IsA()->GetName()
62 // ,(const char *)readObj
63 // );
64 // delete nextObject;
65 // }
66 // //-----------------------------------------------------------------------
67 // // Now loop over all objects in inverse order
68 // printf(" -- > Now loop over all objects in inverse order < -- \n");
69 // for( readObj = size-1; (int)readObj >= 0; --readObj)
70 // {
71 // nextObject = *readObj;
72 // if (nextObject) {
73 // printf(" Object \"%s\" of class \"%s\" written with TKey \"%s\" has been read from file\n"
74 // ,nextObject->GetName()
75 // , nextObject->IsA()->GetName()
76 // ,(const char *)readObj
77 // );
78 // delete nextObject;
79 // } else {
80 // printf("Error reading file by index\n");
81 // }
82 // }
83 // //-----------------------------------------------------------------------
84 // // Loop over the objects starting from the object with the key name "event.02.01"
85 // printf(" -- > Loop over the objects starting from the object with the key name \"event.02.01\" < -- \n");
86 // for( readObj = "event.02.01"; (const char *)readObj != 0; ++readObj){
87 // nextObject = *readObj;
88 // printf(" Object \"%s\" of class \"%s\" written with Tkey \"%s\" has been read from file\n"
89 // , nextObject->GetName()
90 // , nextObject->IsA()->GetName()
91 // , (const char *)readObj
92 // );
93 // delete nextObject;
94 // }
95 //
96 // printf(" -- > Loop over the objects starting from the 86-th object" < -- \n");
97 // for( readObj = (const char *)(readObj = 86); (const char *)readObj != 0; ++readObj){
98 // nextObject = *readObj;
99 // printf(" Object \"%s\" of class \"%s\" written with Tkey \"%s\" has been read from file\n"
100 // , nextObject->GetName()
101 // , nextObject->IsA()->GetName()
102 // , (const char *)readObj
103 // );
104 // delete nextObject;
105 // }
106 //
107 // }
108 //-----------------------------------------------------------------------
110 
111 
112 #include <assert.h>
113 
114 #include "TEnv.h"
115 #include "TSystem.h"
116 #include "TFile.h"
117 #include "TKey.h"
118 
119 #include "StFileIter.h"
120 #include "TDsKey.h"
121 
122 ClassImp(StFileIter)
123 
124 //__________________________________________________________________________
125 StFileIter::StFileIter(TFile *file) : fFileBackUp(0),fDirectoryBackUp(0), fNestedIterator(0)
126  , fRootFile(file)
127  , fEventName("event"), fRunNumber(UInt_t(-1)),fEventNumber(UInt_t(-1))
128  , fCursorPosition(-1), fOwnTFile(kFALSE)
129 {
130  // Create iterator over all objects from the TFile provided
131  Initialize();
132 }
133 
134 //__________________________________________________________________________
135 StFileIter::StFileIter(TDirectory *directory) : fFileBackUp(0),fDirectoryBackUp(0), fNestedIterator(0)
136  , fRootFile(directory)
137  , fEventName("event"), fRunNumber(UInt_t(-1)),fEventNumber(UInt_t(-1))
138  , fCursorPosition(-1), fOwnTFile(kFALSE)
139 {
140  // Create iterator over all objects from the TDirectory provided
141  Initialize();
142 }
143 //__________________________________________________________________________
144 StFileIter::StFileIter(const char *name, Option_t *option, const char *ftitle
145  , Int_t compress, Int_t /*netopt*/) :fNestedIterator(0),fRootFile (0)
146 {
147  // Open ROOT TFile by the name provided;
148  // This TFile is to be deleted by the StFileIter alone
149  if (name && name[0]) {
150  fOwnTFile = kTRUE;
151  // Map a special file system to rfio
152  // /hpss/in2p3.fr/group/atlas/cppm/data/genz
153  // #setenv HPSSIN bnlhpss:/home/atlasgen/evgen
154  // #example for castor: /castor/cern.ch/user/p/paniccia/evgen
155  fRootFile = TFile::Open(MapName(name),option,ftitle,compress);
156  Initialize();
157  }
158 }
159 
160 //__________________________________________________________________________
161 StFileIter::StFileIter(const StFileIter &dst) : TListIter()
162  ,fFileBackUp(0), fDirectoryBackUp(0), fNestedIterator(0)
163  ,fRootFile(dst.fRootFile),fEventName(dst.fEventName), fRunNumber(dst.fRunNumber)
164  ,fEventNumber(dst.fRunNumber),
165  fCursorPosition(-1), fOwnTFile(dst.fOwnTFile)
166 {
167  // Copy ctor can be used with the "read only" files only.
168  assert(!fRootFile->IsWritable());
169  if (fRootFile && fOwnTFile && !fRootFile->IsWritable()) {
170  // Reopen the file
171  if (fRootFile->InheritsFrom(TFile::Class()))
172  {
173  TFile *thisFile = (TFile *)fRootFile;
174  fRootFile = TFile::Open(MapName(fRootFile->GetName())
175  ,fRootFile->GetOption()
176  ,fRootFile->GetTitle()
177  ,thisFile->GetCompressionLevel());
178  }
179  }
180 
181  Initialize();
182  // Adjust this iterator position
183  SkipObjects(dst.fCursorPosition);
184 }
185 //__________________________________________________________________________
186 StFileIter::~StFileIter()
187 {
188  // StFileIter dtor
189  StFileIter *deleteit = fNestedIterator; fNestedIterator = 0;
190  delete deleteit;
191  if (fRootFile && fOwnTFile ) { // delete own TFile if any
192  if (fRootFile->IsWritable()) fRootFile->Write();
193  fRootFile->Close();
194  delete fRootFile;
195  fRootFile = 0;
196  }
197 }
198 
199 //__________________________________________________________________________
200 void StFileIter::Initialize()
201 {
202  //to be documented
203  if (fRootFile) {
204  fDirection = kIterForward;
205  if (IsOpen()) Reset();
206  else {
207  if (fRootFile && fOwnTFile ) delete fRootFile;
208  fRootFile = 0;
209  }
210  }
211 }
212 //__________________________________________________________________________
213 Bool_t StFileIter::IsOpen() const
214 {
215  Bool_t iOpen = kFALSE;
216  if (fRootFile && !fRootFile->IsZombie() ) {
217  iOpen = kTRUE;
218  if (fRootFile->InheritsFrom(TFile::Class()) && !((TFile*)fRootFile)->IsOpen())
219  iOpen = kFALSE;
220  }
221  return iOpen;
222 }
223 
224 //__________________________________________________________________________
225 TKey *StFileIter::GetCurrentKey() const
226 {
227  // return the pointer to the current TKey
228 
229  return ((StFileIter*)this)->SkipObjects(0);
230 }
231 //__________________________________________________________________________
232 Int_t StFileIter::GetDepth() const
233 {
234  // return the current number of the nested subdirectroies;
235  // = 0 - means there is no subdirectories
236  return fNestedIterator ? fNestedIterator->GetDepth()+1 : 0;
237 }
238 
239 //__________________________________________________________________________
240 const char *StFileIter::GetKeyName() const
241 {
242  // return the name of the current TKey
243  const char *name = 0;
244  TKey *key = GetCurrentKey();
245  if (key) name = key->GetName();
246  return name;
247 }
248 //__________________________________________________________________________
249 TObject *StFileIter::GetObject() const
250 {
251  // read the object from TFile defined by the current TKey
252  //
253  // ATTENTION: memory leak danger !!!
254  // ---------
255  // This method does create a new object and it is the end-user
256  // code responsibility to take care about this object
257  // to avoid memory leak.
258  //
259  return ReadObj(GetCurrentKey());
260 }
261 //__________________________________________________________________________
262 Int_t StFileIter::GetObjlen() const
263 {
264  // Returns the uncompressed length of the current object
265  Int_t lenObj = 0;
266  TKey *key = GetCurrentKey();
267  if (key) lenObj = ((TKey *)key)->GetObjlen();
268  return lenObj;
269 }
270 //__________________________________________________________________________
271 Int_t StFileIter::TotalKeys() const
272 {
273  // The total number of the TKey keys in the current TDirectory only
274  // Usually this means the total number of different objects
275  // those can be read one by one.
276  // It does NOT count the nested sub-TDirectory.
277  // It is too costly and it can be abused.
278 
279  Int_t size = 0;
280  if(fList) size += fList->GetSize();
281  return size;
282 }
283 //__________________________________________________________________________
284 TObject *StFileIter::Next(Int_t nSkip)
285 {
286  // return the pointer to the object defined by next TKey
287  // This method is not recommended. It was done for the sake
288  // of the compatibility with TListIter
289 
290  SkipObjects(nSkip);
291  return GetObject();
292 }
293 
294 //__________________________________________________________________________
295 void StFileIter::PurgeKeys(TList *listOfKeys) {
296  assert(listOfKeys);
297  // Remove the TKey duplication,
298  // leave the keys with highest cycle number only
299  // Sort if first
300  listOfKeys->Sort();
301  TObjLink *lnk = listOfKeys->FirstLink();
302  while(lnk) {
303  TKey *key = (TKey *)lnk->GetObject();
304  Short_t cycle = key->GetCycle();
305  const char *keyName = key->GetName();
306  // Check next object
307  lnk = lnk->Next();
308  if (lnk) {
309  TKey *nextkey = 0;
310  TObjLink *lnkThis = lnk;
311  while ( lnk
312  && (nextkey = (TKey *)lnk->GetObject())
313  && !strcmp(nextkey->GetName(), keyName)
314  ) {
315  // compare the cycles
316  Short_t nextCycle = nextkey->GetCycle() ;
317  //printf(" StFileIter::PurgeKeys found new cycle %s :%d : %d\n",
318  // keyName,cycle ,nextCycle);
319  assert(cycle != nextCycle);
320  TObjLink *lnkNext = lnk->Next();
321  if (cycle > nextCycle ) {
322  delete listOfKeys->Remove(lnk);
323  } else {
324  delete listOfKeys->Remove(lnkThis);
325  cycle = nextCycle;
326  lnkThis = lnk;
327  }
328  lnk = lnkNext;
329  }
330  }
331  }
332 }
333 
334 #include <memory>
335 
336 void AssignPointer(TObjLink* cursor, TObjLink* ptr) { cursor = ptr; }
337 void AssignPointer(std::shared_ptr<TObjLink> cursor, TObjLink* ptr) { cursor.reset(ptr); }
338 
339 //__________________________________________________________________________
340 void StFileIter::Reset()
341 {
342  // Reset the status of the iterator
343  if (fNestedIterator) {
344  StFileIter *it = fNestedIterator;
345  fNestedIterator=0;
346  delete it;
347  }
348  TListIter::Reset();
349  if (!fRootFile->IsWritable()) {
350  TList *listOfKeys = fRootFile->GetListOfKeys();
351  if (listOfKeys) {
352  if (!listOfKeys->IsSorted()) PurgeKeys(listOfKeys);
353  fList = listOfKeys;
354  if (fDirection == kIterForward) {
355  fCursorPosition = 0;
356  AssignPointer(fCurCursor, fList->FirstLink());
357  if (fCurCursor) AssignPointer(fCursor, fCurCursor->Next());
358  } else {
359  fCursorPosition = fList->GetSize()-1;
360  AssignPointer(fCurCursor, fList->LastLink());
361  if (fCurCursor) AssignPointer(fCursor, fCurCursor->Prev());
362  }
363  }
364  }
365 }
366 //__________________________________________________________________________
367 void StFileIter::SetCursorPosition(const char *keyNameToFind)
368 {
369  // Find the key by the name provided
370  Reset();
371  while( (*this != keyNameToFind) && SkipObjects() );
372 }
373 //__________________________________________________________________________
374 TKey *StFileIter::SkipObjects(Int_t nSkip)
375 {
376  //
377  // Returns the TKey pointer to the nSkip TKey object from the current one
378  // nSkip = 0; the state of the iterator is not changed
379  //
380  // nSkip > 0; iterator skips nSkip objects in the container.
381  // the direction of the iteration is
382  // sign(nSkip)*kIterForward
383  //
384  // Returns: TKey that can be used to fetch the object from the TDirectory
385  //
386  TKey *nextObject = fNestedIterator ? fNestedIterator->SkipObjects(nSkip): 0;
387  if (!nextObject) {
388  if (fNestedIterator) {
389  StFileIter *it = fNestedIterator;
390  fNestedIterator = 0;
391  delete it;
392  }
393  Int_t collectionSize = 0;
394  if (fList && (collectionSize = fList->GetSize()) ) {
395  if (fDirection !=kIterForward) nSkip = -nSkip;
396  Int_t newPos = fCursorPosition + nSkip;
397  if (0 <= newPos && newPos < collectionSize) {
398  do {
399  if (fCursorPosition < newPos) {
400  fCursorPosition++;
401  fCurCursor = fCursor;
402  AssignPointer(fCursor, fCursor->Next());
403  } else if (fCursorPosition > newPos) {
404  fCursorPosition--;
405  fCurCursor = fCursor;
406  AssignPointer(fCursor, fCursor->Prev());
407  }
408  } while (fCursorPosition != newPos);
409  if (fCurCursor) nextObject = dynamic_cast<TKey *>(fCurCursor->GetObject());
410  } else {
411  fCurCursor = fCursor = 0;
412  if (newPos < 0) {
413  fCursorPosition = -1;
414  if (fList) AssignPointer(fCursor, fList->FirstLink());
415  } else {
416  fCursorPosition = collectionSize;
417  if (fList) AssignPointer(fCursor, fList->LastLink());
418  }
419  }
420  }
421  }
422  return nextObject;
423 }
424 //__________________________________________________________________________
425 TKey *StFileIter::NextEventKey(UInt_t eventNumber, UInt_t runNumber, const char *name)
426 {
427 
428  // Return the key that name matches the "event" . "run number" . "event number" schema
429  // Attention: Side effect: Been called from the end-user code this method causes
430  // ========= the StFileIter::NextObjectGet() to get the "next" object,
431  // defined by "next" key.
432 
433  Bool_t reset = kFALSE;
434  if (name && name[0] && name[0] != '*') { if (fEventName > name) reset = kTRUE; fEventName = name; }
435  if (runNumber !=UInt_t(-1) ) { if (fRunNumber > runNumber) reset = kTRUE; fRunNumber = runNumber;}
436  if (eventNumber !=UInt_t(-1) ) { if (fEventNumber > eventNumber) reset = kTRUE; fEventNumber = eventNumber;}
437 
438  if (reset) Reset();
439  // TIter &nextKey = *fKeysIterator;
440  TKey *key = 0;
441  TDsKey thisKey;
442  while ( (key = SkipObjects()) ) {
443  if (fDirection==kIterForward) fCursorPosition++;
444  else fCursorPosition--;
445  if ( strcmp(name,"*") ) {
446  thisKey.SetKey(key->GetName());
447  if (thisKey.GetName() < name) continue;
448  if (thisKey.GetName() > name) { key = 0; break; }
449  }
450  // Check "run number"
451  if (runNumber != UInt_t(-1)) {
452  UInt_t thisRunNumber = thisKey.RunNumber();
453  if (thisRunNumber < runNumber) continue;
454  if (thisRunNumber < runNumber) { key = 0; break; }
455  }
456  // Check "event number"
457  if (eventNumber != UInt_t(-1)) {
458  UInt_t thisEventNumber = thisKey.EventNumber();
459  if (thisEventNumber < eventNumber) continue;
460  if (thisEventNumber > eventNumber) {key = 0; break; }
461  }
462  break;
463  }
464  return key;
465 }
466 //__________________________________________________________________________
467 TObject *StFileIter::NextEventGet(UInt_t eventNumber, UInt_t runNumber, const char *name)
468 {
469  // reads, creates and returns the object by TKey name that matches
470  // the "name" ."runNumber" ." eventNumber" schema
471  // Attention: This method does create a new TObject and it is the user
472  // code responsibility to take care (delete) this object to avoid
473  // memory leak.
474 
475  return ReadObj(NextEventKey(eventNumber,runNumber,name));
476 }
477 
478 //__________________________________________________________________________
479 TObject *StFileIter::ReadObj(const TKey *key) const
480 {
481  //Read the next TObject from for the TDirectory by TKey provided
482  TObject *obj = 0;
483  if (fNestedIterator) obj = fNestedIterator->ReadObj(key);
484  else if (key) {
485  obj = ((TKey *)key)->ReadObj();
486  if (obj && obj->InheritsFrom(TDirectory::Class()) )
487  {
488  // create the next iteration level.
489  assert(!fNestedIterator);
490  ((StFileIter*)this)->fNestedIterator = new StFileIter((TDirectory *)obj);
491  // FIXME: needs to set fDirection if needed 02/11/2007 vf
492  }
493  }
494  return obj;
495 }
496 
497 //__________________________________________________________________________
498 Int_t StFileIter::NextEventPut(TObject *obj, UInt_t eventNum, UInt_t runNumber
499  , const char *name)
500 {
501  // Create a special TKey name with obj provided and write it out.
502 
503  Int_t wBytes = 0;
504  if (obj && IsOpen() && fRootFile->IsWritable()) {
505  TDsKey thisKey(runNumber,eventNum);
506  if (name && name[0])
507  thisKey.SetName(name);
508  else
509  thisKey.SetName(obj->GetName());
510 
511  if (fRootFile != gDirectory) {
512  SaveFileScope();
513  fRootFile->cd();
514  }
515  wBytes = obj->Write(thisKey.GetKey());
516  if (fRootFile->InheritsFrom(TFile::Class())) ((TFile*)fRootFile)->Flush();
517  if (fRootFile != gDirectory) RestoreFileScope();
518  }
519  return wBytes;
520 }
521 //__________________________________________________________________________
522 TString StFileIter::MapName(const char *name, const char *localSystemKey,const char *mountedFileSystemKey)
523 {
524  // --------------------------------------------------------------------------------------
525  // MapName(const char *name, const char *localSystemKey,const char *mountedFileSystemKey)
526  // --------------------------------------------------------------------------------------
527  // Substitute the logical name with the real one if any
528  // 1. add a line into system.rootrc or ~/.rootrc or ./.rootrc
529  //
530  // StFileIter.ForeignFileMap mapFile // the name of the file
531  // to map the local name
532  // to the global file service
533  //
534  // If this line is omitted then StFileIter class seeks for
535  // the default mapping file in the current directory "io.config"
536 
537  // 2. If the "io.config" file found then it defines the mapping as follows:
538  //
539  // StFileIter.LocalFileSystem /castor
540  // StFileIter.MountedFileSystem rfio:/castor
541 
542  // If "io.config" doesn't exist then no mapping is to be performed
543  // and all file names are treated "as is"
544 
545  if ( !localSystemKey) localSystemKey = GetLocalFileNameKey();
546  if ( !mountedFileSystemKey) mountedFileSystemKey = GetForeignFileSystemKey();
547  TString newName = name;
548  TString fileMap = gEnv->GetValue(GetResourceName(),GetDefaultMapFileName());
549  const char *localName = 0;
550  const char *foreignName = 0;
551  if ( gSystem->AccessPathName(fileMap) == 0 ){
552  TEnv myMapResource(fileMap);
553  localName = myMapResource.Defined(localSystemKey) ?
554  myMapResource.GetValue(localSystemKey,"") : 0;
555  foreignName = myMapResource.Defined(mountedFileSystemKey) ?
556  myMapResource.GetValue(mountedFileSystemKey,""):0;
557  } else {
558  localName = "/castor"; // This is the default CERN name
559  foreignName = "rfio:/castor"; // and it needs "RFIO"
560  }
561  if (localName && localName[0]
562  && foreignName
563  && foreignName[0]
564  && newName.BeginsWith(localName) )
565  newName.Replace(0,strlen(localName),foreignName);
566  return newName;
567 }
virtual TString GetKey() const
to be documented
Definition: TDsKey.cxx:111
Definition: TDsKey.h:21
virtual void SetKey(const char *key)
to be documented
Definition: TDsKey.cxx:126