]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliFileMerger.cxx
Added loop for extraction of clusters really attached to its track.
[u/mrichter/AliRoot.git] / ANALYSIS / AliFileMerger.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 //////////////////////////////////////////////////////////////////////////
19 //  marian.ivanov@cern.ch
20 //  Utilities for file merging.
21 //  Additional functionality on top of the standard TFileMerger:
22 //
23 //  1. Possibility to Set the reject/accept list.
24 //     1.a)  Only entries selected in accept list are merged. By default all entries are selected
25 //           use AddAccept 0 to specify your desired entry
26 //     1.b)  Entries selected in reject list are not merged. By default the reject list is empty.
27 //
28 //  2. syswatch.log is created diring mergin procedure. 
29 //     Memeory consumption - for reading and for merging can be monitored
30
31 //  RS: Changed merger to respect the structure of files being merged (directories, collections...)
32 //      Additional option: SetNoTrees (default false) to not merge any tree
33 //      The code mostly taken from root's hadd.cxx
34 /*
35   Usage:
36   // Libraries for all classes to be merged should be loaded before using the class
37   gSystem->Load("libANALYSIS");
38   gSystem->Load("libANALYSIScalib");
39   gSystem->Load("libTPCcalib"); 
40   TH1::AddDirectory(0);
41
42   //Example usage starting from the input data list in text file:
43   //
44   AliFileMerger merger;
45   merger.AddReject("esdFriend");
46   merger.IterTXT("calib.list","CalibObjects.root",kFALSE);
47   //
48
49 */
50 //////////////////////////////////////////////////////////////////////////
51  
52
53 #include <fstream>
54 #include <THashList.h>
55 #include <TChain.h>
56 #include <TKey.h>
57 #include <TH1.h>
58 #include <THStack.h>
59 #include "TSystem.h"
60 #include "TFile.h"
61 #include "TGrid.h"
62 #include "TGridResult.h"
63 #include "TObjString.h"
64 #include "TObjArray.h"
65 #include "TMethodCall.h"
66 #include "Riostream.h"
67 #include "AliSysInfo.h"
68 #include "AliFileMerger.h"
69 #include "AliLog.h"
70
71 using std::cerr;
72 using std::endl;
73 using std::cout;
74 using std::ifstream;
75 ClassImp(AliFileMerger)
76
77 ProcInfo_t procInfo;//TMP
78
79 ////////////////////////////////////////////////////////////////////////
80
81 AliFileMerger::AliFileMerger():
82   TNamed(),
83   fRejectMask(0),
84   fAcceptMask(0),
85   fMaxFilesOpen(800),
86   fNoTrees(kFALSE)
87 {
88   //
89   // Default constructor
90   //
91 }
92
93 //______________________________________________________________________
94
95 AliFileMerger::AliFileMerger(const char* name):
96   TNamed(name,name),
97   fRejectMask(0),
98   fAcceptMask(0),
99   fMaxFilesOpen(800),
100   fNoTrees(kFALSE)
101 {
102   //
103   // 
104   //
105 }
106
107
108 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern, Bool_t dontOverwrite){
109
110   //
111   // Merge the files coming out of the calibration job
112   // 
113   TString command;
114   // looking for files to be merged in the output directory
115   command = Form("find %s/ *%s", outputDir, pattern);
116   printf("command: %s\n", command.Data());
117   TGrid::Connect("alien://");
118   TGridResult *res = gGrid->Command(command);
119   if (!res) return;
120   TIter nextmap(res);
121   TMap *map = 0;
122   // loop over the results
123   TList sourcelist;  
124   sourcelist.SetOwner(kTRUE);
125   //
126   while((map=(TMap*)nextmap())) {
127     // getting the turl
128     TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
129     if (!objs || !objs->GetString().Length()) {
130       // Nothing found - skip this output
131       delete res;
132       break;
133     } 
134     printf("looking for file %s\n",(objs->GetString()).Data());
135     AddFile(&sourcelist, (objs->GetString()).Data());;
136   }
137   //
138   IterList(&sourcelist, outputFileName, dontOverwrite);
139   delete res;
140 }
141
142 void AliFileMerger::IterList(const TList* namesList, const char* outputFileName, Bool_t dontOverwrite)
143 {
144   // merge in steps or in one go
145   //
146   gSystem->GetProcInfo(&procInfo);
147   AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
148   //
149   TString outputFile(outputFileName);
150   gSystem->ExpandPathName(outputFile);
151   //
152   int nFiles = namesList->GetEntries();
153   int maxSrcOpen = fMaxFilesOpen - 1;
154   TList filesList;
155   filesList.SetOwner(kTRUE);
156   //
157   TString tmpDest[2] = {outputFile,outputFile}; // names for tmp files
158   int npl = outputFile.Last('.');
159   if (npl<0) npl  = outputFile.Length();
160   for (int i=0;i<2;i++) tmpDest[i].Insert(npl,Form("_TMPMERGE%d_",i));
161   //
162   int nsteps = 0, currTmp = 0, start = 0;
163   for (int ifl=0;ifl<nFiles;ifl++) {
164     int st = ifl%maxSrcOpen;
165     if (st==0 && ifl) { // new chunk should be started, merge what was already accumulated
166       OpenNextChunks(namesList,&filesList,start,ifl-1);
167       start = ifl; // remember where to start next step
168       if (nsteps++) { // if not 1st one, merge the privous chunk with this one
169         filesList.AddFirst(TFile::Open(tmpDest[currTmp].Data()));
170         currTmp = (currTmp==0) ? 1:0;         // swap tmp files
171       }
172       // open temp target
173       TFile* targetTmp = TFile::Open( tmpDest[currTmp].Data(), "RECREATE");
174       if (!targetTmp || targetTmp->IsZombie()) {
175         printf("Error opening temporary file %s\n",tmpDest[currTmp].Data());
176         return;
177       }
178       MergeRootfile(targetTmp, &filesList);
179       targetTmp->Close();
180       delete targetTmp;
181       filesList.Clear(); // close all open files
182     }
183     // nothing to do until needed amount of files is accumulated
184   }
185   // merge last step
186   TFile* target = TFile::Open( outputFile.Data(), (dontOverwrite ? "CREATE":"RECREATE") );
187   if (!target || target->IsZombie()) {
188     cerr << "Error opening target file (does " << outputFileName << " exist?)." << endl;
189     cerr << "Use force = kTRUE to re-creation of output file." << endl;
190     return;
191   }
192   OpenNextChunks(namesList,&filesList,start,nFiles-1);
193   // add result of previous merges
194   if (nsteps) filesList.AddFirst(TFile::Open(tmpDest[currTmp].Data()));
195   MergeRootfile( target, &filesList);
196   target->Close();
197   delete target;
198   filesList.Clear();
199   // 
200   for (int i=0;i<2;i++) gSystem->Exec(Form("if [ -e %s ]; then \nrm %s\nfi",tmpDest[i].Data(),tmpDest[i].Data()));
201   //
202   printf("Merged %d files in %d steps\n",nFiles,++nsteps);
203   //  
204   gSystem->GetProcInfo(&procInfo);
205   AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
206 }
207
208 void AliFileMerger::IterTXT( const char * fileList,  const char* outputFileName, Bool_t dontOverwrite){
209   
210   // Merge the files indicated in the list - fileList
211   // ASCII file option example: 
212   // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
213   
214   // Open the input stream
215   ifstream in;
216   in.open(fileList);
217   // Read the input list of files 
218   TString objfile;
219   Int_t counter=0;
220   TList sourcelist;  
221   sourcelist.SetOwner(kTRUE);
222   while(in.good()) {
223     in >> objfile;
224     if (!objfile.Contains(".root")) continue; // protection
225     gSystem->ExpandPathName(objfile);
226     printf("Add file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
227     AddFile(&sourcelist, objfile.Data());
228   }
229   //
230   IterList(&sourcelist, outputFileName, dontOverwrite);
231   //
232 }
233
234 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
235   //
236   // Storing the results in one single file
237   //
238   TFile *f = new TFile(outputFileName,"recreate");
239   for (Int_t i=0; i<array->GetEntries(); i++){
240     TObject *object0 = array->At(i);
241     if (!object0) continue;
242     object0->Write();
243   }
244   f->Close();
245   delete f;
246 }
247
248
249 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
250   //
251   // Store the results in separate files (one per object)
252   //
253   for (Int_t i=0; i<array->GetEntries(); i++){
254     TObject *object0 = array->At(i);
255     if (!object0) continue;
256     TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
257     object0->Write();
258     f->Close();
259     delete f;
260   }
261 }
262
263 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
264   //
265   // Merging procedure
266   //
267   if (!array) return;
268   static Int_t counter=-1;
269   counter++;
270   TObjArray *carray = new TObjArray;   //array of the objects inside current file
271   carray->SetOwner(kTRUE);
272   
273   // load all objects to  memory
274   
275   TList *farr = fileIn->GetListOfKeys();
276   if (!farr) { 
277     delete carray;
278     return;
279   }
280   for (Int_t ical=0; ical<farr->GetEntries(); ical++){
281     if (!farr->At(ical)) continue;
282     TString name(farr->At(ical)->GetName());
283     if (!IsAccepted(name)) continue;                        // skip not accepted entries
284     TObject *obj = fileIn->Get(name.Data());
285     if (obj) carray->AddLast(obj);
286     AliSysInfo::AddStamp(name.Data(),1,ical,counter);  
287   }
288   
289   if (carray->GetEntries()==0)  { 
290     delete carray;
291     return;
292   }
293   TMethodCall callEnv;
294   Int_t entries =carray->GetEntriesFast();
295   for (Int_t i=0; i<entries; i++){
296     
297     TObjArray *templist = new TObjArray(1);
298     templist->SetOwner(kFALSE);
299     TObject *currentObject = carray->At(i);
300     if (!currentObject) { 
301       delete templist;
302       continue;
303     }
304     printf("%s\n",currentObject->GetName());
305     callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
306     if (!callEnv.IsValid()) {
307       delete templist; 
308       continue;
309     }
310     TString oname=currentObject->GetName();
311     TObject *mergedObject = array->FindObject(currentObject->GetName());
312     if (!mergedObject) {
313       array->AddLast(currentObject);
314       carray->RemoveAt(i);
315       delete templist; 
316       continue;
317     }
318     templist->AddLast(currentObject);
319     callEnv.SetParam((Long_t) templist);
320     callEnv.Execute(mergedObject);
321     AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);  
322     delete templist;
323   }
324   carray->Delete();
325   delete carray;
326 }
327
328 Bool_t AliFileMerger::IsAccepted(TString name){
329   //
330   // Accept/reject logic
331   // name - name of the entry
332   //
333   //  if fAcceptMask specified   - entry has to be in list of selected
334   //  if fRejectMask speciefied  - entry with name speciief in the list are rejected 
335   //
336   Bool_t accept=kTRUE;
337   if (fAcceptMask){
338     //
339     accept=kFALSE;
340     for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
341       if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE;   // entry was selected
342     }
343   }
344   if (!accept) return kFALSE;
345
346   if (fRejectMask){
347     //
348     for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
349       if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE;   // entry was rejected
350     }
351   }
352   return accept;
353 }
354
355
356
357
358 void AliFileMerger::AddReject(const char *reject){
359   //
360   // add reject string to the list of entries to be rejected for merging
361   //
362   if (!fRejectMask) fRejectMask = new TObjArray;
363   fRejectMask->AddLast(new TObjString(reject));
364 }
365 void AliFileMerger::AddAccept(const char *accept){
366   //
367   // add reject string to the list of entries to be rejected for merging
368   //
369   if (!fAcceptMask) fAcceptMask = new TObjArray;
370   fAcceptMask->AddLast(new TObjString(accept));
371
372
373 }
374
375 //___________________________________________________________________________
376 int AliFileMerger::MergeRootfile( TDirectory *target, TList *sourcelist, Bool_t nameFiltering)
377 {
378   // Merge all objects in a directory
379   // modified version of root's hadd.cxx
380   gSystem->GetProcInfo(&procInfo);
381   AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
382   //
383   int status = 0;
384   cout << "Target path: " << target->GetPath() << endl;
385   TString path( (char*)strstr( target->GetPath(), ":" ) );
386   path.Remove( 0, 2 );
387   //
388   // find 1st valid file
389   TDirectory *first_source = (TDirectory*)sourcelist->First();
390   //
391   Int_t nguess = sourcelist->GetSize()+1000;
392   THashList allNames(nguess);
393   ((THashList*)target->GetList())->Rehash(nguess);
394   ((THashList*)target->GetListOfKeys())->Rehash(nguess);
395   TList listH;
396   TString listHargs;
397   listHargs.Form("((TCollection*)0x%lx)", (ULong_t)&listH);
398   //
399   while(first_source) {
400     //
401     TDirectory *current_sourcedir = first_source->GetDirectory(path);
402     if (!current_sourcedir) {
403       first_source = (TDirectory*)sourcelist->After(first_source);
404       continue;
405     }
406     // loop over all keys in this directory
407     TChain *globChain = 0;
408     TIter nextkey( current_sourcedir->GetListOfKeys() );
409     TKey *key, *oldkey=0;
410     //gain time, do not add the objects in the list in memory
411     TH1::AddDirectory(kFALSE);
412     //
413     int counterK = 0;
414     int counterF=0;
415     //
416     while ( (key = (TKey*)nextkey())) {
417       if (current_sourcedir == target) break;
418       //
419       // check if we don't reject this name
420       TString nameK(key->GetName());
421       if (!IsAccepted(nameK) && nameFiltering) {
422         if (!counterF) printf("Object %s is in rejection list, skipping...\n",nameK.Data());
423         continue;
424       }
425       //
426       //keep only the highest cycle number for each key
427       if (oldkey && !strcmp(oldkey->GetName(),key->GetName())) continue;
428       if (!strcmp(key->GetClassName(),"TProcessID")) {key->ReadObj(); continue;}
429       if (allNames.FindObject(key->GetName())) continue;
430       TClass *cl = TClass::GetClass(key->GetClassName());
431       if (!cl || !cl->InheritsFrom(TObject::Class())) {
432         cout << "Cannot merge object type, name: "
433              << key->GetName() << " title: " << key->GetTitle() << endl;
434         continue;
435       }
436       allNames.Add(new TObjString(key->GetName()));
437       AliSysInfo::AddStamp(nameK.Data(),1,++counterK,counterF++); 
438       // read object from first source file
439       //current_sourcedir->cd();
440
441       TObject *obj = key->ReadObj();
442       if (!obj) {
443         AliError(Form("Failed to get the object with key %s from %s",key->GetName(),current_sourcedir->GetFile()->GetName()));
444         continue;
445       }
446
447       if ( obj->IsA()->InheritsFrom( TTree::Class() ) ) {
448         
449         // loop over all source files create a chain of Trees "globChain"
450         if (!fNoTrees) { // 
451           TString obj_name;
452           if (path.Length()) {
453             obj_name = path + "/" + obj->GetName();
454           } else {
455             obj_name = obj->GetName();
456           }
457           globChain = new TChain(obj_name);
458           globChain->Add(first_source->GetName());
459           TFile *nextsource = (TFile*)sourcelist->After( first_source );
460           while ( nextsource ) {
461             //do not add to the list a file that does not contain this Tree
462             TFile *curf = TFile::Open(nextsource->GetName());
463             if (curf) {
464               Bool_t mustAdd = kFALSE;
465               if (curf->FindKey(obj_name)) {
466                 mustAdd = kTRUE;
467               } else {
468                 //we could be more clever here. No need to import the object
469                 //we are missing a function in TDirectory
470                 TObject *aobj = curf->Get(obj_name);
471                 if (aobj) { mustAdd = kTRUE; delete aobj;}
472               }
473               if (mustAdd) {
474                 globChain->Add(nextsource->GetName());
475               }
476             }
477             delete curf;
478             nextsource = (TFile*)sourcelist->After( nextsource );
479           }
480         }
481       } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
482         // it's a subdirectory
483         
484         cout << "Found subdirectory " << obj->GetName() << endl;
485         // create a new subdir of same name and title in the target file
486         target->cd();
487         TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
488         
489         // newdir is now the starting point of another round of merging
490         // newdir still knows its depth within the target file via
491         // GetPath(), so we can still figure out where we are in the recursion
492         status = MergeRootfile( newdir, sourcelist, kFALSE);
493         if (status) return status;
494         
495       } else if ( obj->InheritsFrom(TObject::Class())
496                   && obj->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
497         // object implements Merge(TCollection*)
498         
499         // loop over all source files and merge same-name object
500         TFile *nextsource = (TFile*)sourcelist->After( first_source );
501         while ( nextsource ) {
502           // make sure we are at the correct directory level by cd'ing to path
503           TDirectory *ndir = nextsource->GetDirectory(path);
504           if (ndir) {
505             ndir->cd();
506             TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
507             if (key2) {
508               TObject *hobj = key2->ReadObj();
509               if (!hobj) {
510                 cout << "Failed to get the object with key " << key2->GetName() << " from " << 
511                   ndir->GetFile()->GetName() << "/" << ndir->GetName() << endl;
512                 nextsource = (TFile*)sourcelist->After( nextsource );
513                 continue;
514               }
515               //
516               hobj->ResetBit(kMustCleanup);
517               listH.Add(hobj);
518               Int_t error = 0;
519               obj->Execute("Merge", listHargs.Data(), &error); // RS Probleme here
520               if (error) {
521                 cerr << "Error calling Merge() on " << obj->GetName()
522                      << " with the corresponding object in " << nextsource->GetName() << endl;
523               }
524               listH.Delete();
525         // get the number of processed entries to be put in the syswatch.log
526         Double_t numberOfEntries = -1;
527         if (obj->IsA()->GetMethodAllAny("GetEntries"))
528         {
529           TMethodCall getEntries(obj->IsA(), "GetEntries", "");
530           getEntries.Execute(obj, numberOfEntries);
531         }
532               AliSysInfo::AddStamp(nameK.Data(),1,counterK,counterF++,numberOfEntries); 
533             }
534           }
535           nextsource = (TFile*)sourcelist->After( nextsource );
536         }
537       } else if ( obj->IsA()->InheritsFrom( THStack::Class() ) ) {
538         THStack *hstack1 = (THStack*) obj;
539         TList* l = new TList();
540         
541         // loop over all source files and merge the histos of the
542         // corresponding THStacks with the one pointed to by "hstack1"
543         TFile *nextsource = (TFile*)sourcelist->After( first_source );
544         while ( nextsource ) {
545           // make sure we are at the correct directory level by cd'ing to path
546           TDirectory *ndir = nextsource->GetDirectory(path);
547           if (ndir) {
548             ndir->cd();
549             TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(hstack1->GetName());
550             if (key2) {
551               THStack *hstack2 = (THStack*) key2->ReadObj();
552               l->Add(hstack2->GetHists()->Clone());
553               delete hstack2;
554               AliSysInfo::AddStamp(nameK.Data(),1,counterK,counterF++); 
555             }
556           }
557           
558           nextsource = (TFile*)sourcelist->After( nextsource );
559         }
560         hstack1->GetHists()->Merge(l);
561         l->Delete();
562       } else {
563         // object is of no type that we can merge
564         cout << "Cannot merge object type, name: "
565              << obj->GetName() << " title: " << obj->GetTitle() << endl;
566         
567         // loop over all source files and write similar objects directly to the output file
568         TFile *nextsource = (TFile*)sourcelist->After( first_source );
569         while ( nextsource ) {
570           // make sure we are at the correct directory level by cd'ing to path
571           TDirectory *ndir = nextsource->GetDirectory(path);
572           if (ndir) {
573             ndir->cd();
574             TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
575             if (key2) {
576               TObject *nobj = key2->ReadObj();
577               nobj->ResetBit(kMustCleanup);
578               int nbytes1 = target->WriteTObject(nobj, key2->GetName(), "SingleKey" );
579               if (nbytes1 <= 0) status = -1;
580               delete nobj;
581             }
582           }
583           nextsource = (TFile*)sourcelist->After( nextsource );
584         }
585       }
586       
587       // now write the merged histogram (which is "in" obj) to the target file
588       // note that this will just store obj in the current directory level,
589       // which is not persistent until the complete directory itself is stored
590       // by "target->Write()" below
591       target->cd();
592       
593       //!!if the object is a tree, it is stored in globChain...
594       if(obj->IsA()->InheritsFrom( TDirectory::Class() )) {
595         //printf("cas d'une directory\n");
596       } else if(obj->IsA()->InheritsFrom( TTree::Class() )) {
597         if (!fNoTrees) {
598           globChain->ls();
599           globChain->Merge(target->GetFile(),0,"keep fast");
600           delete globChain;
601         }
602       } else {
603         int nbytes2 = obj->Write( key->GetName(), TObject::kSingleKey );
604         if (nbytes2 <= 0) status = -1;
605       }
606       oldkey = key;
607       delete obj;
608     } // while ( ( TKey *key = (TKey*)nextkey() ) )
609     first_source = (TDirectory*)sourcelist->After(first_source);
610   }
611   // save modifications to target file
612   target->SaveSelf(kTRUE);
613   //
614   gSystem->GetProcInfo(&procInfo);
615   AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
616
617   return status;
618 }
619
620 //___________________________________________________________________________
621 int AliFileMerger::OpenNextChunks(const TList* namesList, TList* filesList, Int_t from, Int_t to)
622 {
623   gSystem->GetProcInfo(&procInfo);
624   AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
625
626   filesList->Clear();
627   int nEnt = namesList->GetEntries();
628   from = from<nEnt ? from : nEnt;
629   to   = to<nEnt ? to : nEnt;
630   int count = 0;
631   for (int i=from;i<=to;i++) {
632     TNamed* fnam = (TNamed*)namesList->At(i);
633     if (!fnam) continue;
634     TString fnamS(fnam->GetName());
635     gSystem->ExpandPathName(fnamS);
636     if (fnamS.BeginsWith("alien://") && !gGrid) TGrid::Connect("alien");
637     TFile* source = TFile::Open(fnam->GetName());
638     if( source==0 ) { printf("Failed to open file %s, will skip\n",fnam->GetName()); continue; }
639     filesList->Add(source);
640     printf("Opened file %s\n",fnam->GetName());
641     count++;
642   }
643   gSystem->GetProcInfo(&procInfo);
644   AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
645
646   return count;
647 }
648
649
650 //___________________________________________________________________________
651 int AliFileMerger::AddFile(TList* namesList, std::string entry)
652 {
653   // add a new file to the list of files
654   //  static int count(0);
655   if( entry.empty() ) return 0;
656   size_t j =entry.find_first_not_of(' ');
657   if( j==std::string::npos ) return 0;
658   entry = entry.substr(j);
659   if( entry.substr(0,1)=="@") {
660     std::ifstream indirect_file(entry.substr(1).c_str() );
661     if( ! indirect_file.is_open() ) {
662       std::cerr<< "Could not open indirect file " << entry.substr(1) << std::endl;
663       return 1;
664     }
665     while( indirect_file ){
666       std::string line;
667       std::getline(indirect_file, line);
668       if( AddFile(namesList, line)!=0 ) return 1;;
669     }
670     return 0;
671   }
672   //  cout << "Source file " << (++count) << ": " << entry << endl;
673   namesList->Add(new TNamed(entry,""));
674   return 0;
675 }