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