1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //////////////////////////////////////////////////////////////////////////
19 // marian.ivanov@cern.ch
20 // Utilities for file merging.
21 // Additional functionality on top of the standard TFileMerger:
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.
28 // 2. syswatch.log is created diring mergin procedure.
29 // Memeory consumption - for reading and for merging can be monitored
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
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");
42 //Example usage starting from the input data list in text file:
45 merger.AddReject("esdFriend");
46 merger.IterTXT("calib.list","CalibObjects.root",kFALSE);
50 //////////////////////////////////////////////////////////////////////////
54 #include <THashList.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"
71 ClassImp(AliFileMerger)
73 ProcInfo_t procInfo;//TMP
75 ////////////////////////////////////////////////////////////////////////
77 AliFileMerger::AliFileMerger():
85 // Default constructor
89 //______________________________________________________________________
91 AliFileMerger::AliFileMerger(const char* name):
104 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern, Bool_t dontOverwrite){
107 // Merge the files coming out of the calibration job
110 // looking for files to be merged in the output directory
111 command = Form("find %s/ *%s", outputDir, pattern);
112 printf("command: %s\n", command.Data());
113 TGrid::Connect("alien://");
114 TGridResult *res = gGrid->Command(command);
118 // loop over the results
120 sourcelist.SetOwner(kTRUE);
122 while((map=(TMap*)nextmap())) {
124 TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
125 if (!objs || !objs->GetString().Length()) {
126 // Nothing found - skip this output
130 printf("looking for file %s\n",(objs->GetString()).Data());
131 AddFile(&sourcelist, (objs->GetString()).Data());;
134 IterList(&sourcelist, outputFileName, dontOverwrite);
138 void AliFileMerger::IterList(const TList* namesList, const char* outputFileName, Bool_t dontOverwrite)
140 // merge in steps or in one go
142 gSystem->GetProcInfo(&procInfo);
143 AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
145 TString outputFile(outputFileName);
146 gSystem->ExpandPathName(outputFile);
148 int nFiles = namesList->GetEntries();
149 int maxSrcOpen = fMaxFilesOpen - 1;
151 filesList.SetOwner(kTRUE);
153 TString tmpDest[2] = {outputFile,outputFile}; // names for tmp files
154 int npl = outputFile.Last('.');
155 if (npl<0) npl = outputFile.Length();
156 for (int i=0;i<2;i++) tmpDest[i].Insert(npl,Form("_TMPMERGE%d_",i));
158 int nsteps = 0, currTmp = 0, start = 0;
159 for (int ifl=0;ifl<nFiles;ifl++) {
160 int st = ifl%maxSrcOpen;
161 if (st==0 && ifl) { // new chunk should be started, merge what was already accumulated
162 OpenNextChunks(namesList,&filesList,start,ifl-1);
163 start = ifl; // remember where to start next step
164 if (nsteps++) { // if not 1st one, merge the privous chunk with this one
165 filesList.AddFirst(TFile::Open(tmpDest[currTmp].Data()));
166 currTmp = (currTmp==0) ? 1:0; // swap tmp files
169 TFile* targetTmp = TFile::Open( tmpDest[currTmp].Data(), "RECREATE");
170 if (!targetTmp || targetTmp->IsZombie()) {
171 printf("Error opening temporary file %s\n",tmpDest[currTmp].Data());
174 MergeRootfile(targetTmp, &filesList);
177 filesList.Clear(); // close all open files
179 // nothing to do until needed amount of files is accumulated
182 TFile* target = TFile::Open( outputFile.Data(), (dontOverwrite ? "CREATE":"RECREATE") );
183 if (!target || target->IsZombie()) {
184 cerr << "Error opening target file (does " << outputFileName << " exist?)." << endl;
185 cerr << "Use force = kTRUE to re-creation of output file." << endl;
188 OpenNextChunks(namesList,&filesList,start,nFiles-1);
189 // add result of previous merges
190 if (nsteps) filesList.AddFirst(TFile::Open(tmpDest[currTmp].Data()));
191 MergeRootfile( target, &filesList);
196 for (int i=0;i<2;i++) gSystem->Exec(Form("if [ -e %s ]; then \nrm %s\nfi",tmpDest[i].Data(),tmpDest[i].Data()));
198 printf("Merged %d files in %d steps\n",nFiles,++nsteps);
200 gSystem->GetProcInfo(&procInfo);
201 AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
204 void AliFileMerger::IterTXT( const char * fileList, const char* outputFileName, Bool_t dontOverwrite){
206 // Merge the files indicated in the list - fileList
207 // ASCII file option example:
208 // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
210 // Open the input stream
213 // Read the input list of files
217 sourcelist.SetOwner(kTRUE);
220 if (!objfile.Contains(".root")) continue; // protection
221 gSystem->ExpandPathName(objfile);
222 printf("Add file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
223 AddFile(&sourcelist, objfile.Data());
226 IterList(&sourcelist, outputFileName, dontOverwrite);
230 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
232 // Storing the results in one single file
234 TFile *f = new TFile(outputFileName,"recreate");
235 for (Int_t i=0; i<array->GetEntries(); i++){
236 TObject *object0 = array->At(i);
237 if (!object0) continue;
245 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
247 // Store the results in separate files (one per object)
249 for (Int_t i=0; i<array->GetEntries(); i++){
250 TObject *object0 = array->At(i);
251 if (!object0) continue;
252 TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
259 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
264 static Int_t counter=-1;
266 TObjArray *carray = new TObjArray; //array of the objects inside current file
267 carray->SetOwner(kTRUE);
269 // load all objects to memory
271 TList *farr = fileIn->GetListOfKeys();
276 for (Int_t ical=0; ical<farr->GetEntries(); ical++){
277 if (!farr->At(ical)) continue;
278 TString name(farr->At(ical)->GetName());
279 if (!IsAccepted(name)) continue; // skip not accepted entries
280 TObject *obj = fileIn->Get(name.Data());
281 if (obj) carray->AddLast(obj);
282 AliSysInfo::AddStamp(name.Data(),1,ical,counter);
285 if (carray->GetEntries()==0) {
290 Int_t entries =carray->GetEntriesFast();
291 for (Int_t i=0; i<entries; i++){
293 TObjArray *templist = new TObjArray(1);
294 templist->SetOwner(kFALSE);
295 TObject *currentObject = carray->At(i);
296 if (!currentObject) {
300 printf("%s\n",currentObject->GetName());
301 callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
302 if (!callEnv.IsValid()) {
306 TString oname=currentObject->GetName();
307 TObject *mergedObject = array->FindObject(currentObject->GetName());
309 array->AddLast(currentObject);
314 templist->AddLast(currentObject);
315 callEnv.SetParam((Long_t) templist);
316 callEnv.Execute(mergedObject);
317 AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);
324 Bool_t AliFileMerger::IsAccepted(TString name){
326 // Accept/reject logic
327 // name - name of the entry
329 // if fAcceptMask specified - entry has to be in list of selected
330 // if fRejectMask speciefied - entry with name speciief in the list are rejected
336 for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
337 if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE; // entry was selected
340 if (!accept) return kFALSE;
344 for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
345 if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE; // entry was rejected
354 void AliFileMerger::AddReject(const char *reject){
356 // add reject string to the list of entries to be rejected for merging
358 if (!fRejectMask) fRejectMask = new TObjArray;
359 fRejectMask->AddLast(new TObjString(reject));
361 void AliFileMerger::AddAccept(const char *accept){
363 // add reject string to the list of entries to be rejected for merging
365 if (!fAcceptMask) fAcceptMask = new TObjArray;
366 fAcceptMask->AddLast(new TObjString(accept));
371 //___________________________________________________________________________
372 int AliFileMerger::MergeRootfile( TDirectory *target, TList *sourcelist)
374 // Merge all objects in a directory
375 // modified version of root's hadd.cxx
376 gSystem->GetProcInfo(&procInfo);
377 AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
382 cout << "Target path: " << target->GetPath() << endl;
383 TString path( (char*)strstr( target->GetPath(), ":" ) );
386 // find 1st valid file
387 TDirectory *first_source = (TDirectory*)sourcelist->First();
389 Int_t nguess = sourcelist->GetSize()+1000;
390 THashList allNames(nguess);
391 ((THashList*)target->GetList())->Rehash(nguess);
392 ((THashList*)target->GetListOfKeys())->Rehash(nguess);
395 listHargs.Form("((TCollection*)0x%lx)", (ULong_t)&listH);
397 while(first_source) {
399 TDirectory *current_sourcedir = first_source->GetDirectory(path);
400 if (!current_sourcedir) {
401 first_source = (TDirectory*)sourcelist->After(first_source);
404 // loop over all keys in this directory
405 TChain *globChain = 0;
406 TIter nextkey( current_sourcedir->GetListOfKeys() );
407 TKey *key, *oldkey=0;
408 //gain time, do not add the objects in the list in memory
409 TH1::AddDirectory(kFALSE);
413 while ( (key = (TKey*)nextkey())) {
414 if (current_sourcedir == target) break;
416 // check if we don't reject this name
417 TString nameK(key->GetName());
418 if (!IsAccepted(nameK)) {
419 if (!counterF) printf("Object %s is in rejection list, skipping...\n",nameK.Data());
423 //keep only the highest cycle number for each key
424 if (oldkey && !strcmp(oldkey->GetName(),key->GetName())) continue;
425 if (!strcmp(key->GetClassName(),"TProcessID")) {key->ReadObj(); continue;}
426 if (allNames.FindObject(key->GetName())) continue;
427 TClass *cl = TClass::GetClass(key->GetClassName());
428 if (!cl || !cl->InheritsFrom(TObject::Class())) {
429 cout << "Cannot merge object type, name: "
430 << key->GetName() << " title: " << key->GetTitle() << endl;
433 allNames.Add(new TObjString(key->GetName()));
434 AliSysInfo::AddStamp(nameK.Data(),1,counterK++,counterF-1);
435 // read object from first source file
436 //current_sourcedir->cd();
438 TObject *obj = key->ReadObj();
440 AliError(Form("Failed to get the object with key %s from %s",key->GetName(),current_sourcedir->GetFile()->GetName()));
444 if ( obj->IsA()->InheritsFrom( TTree::Class() ) ) {
446 // loop over all source files create a chain of Trees "globChain"
450 obj_name = path + "/" + obj->GetName();
452 obj_name = obj->GetName();
454 globChain = new TChain(obj_name);
455 globChain->Add(first_source->GetName());
456 TFile *nextsource = (TFile*)sourcelist->After( first_source );
457 while ( nextsource ) {
458 //do not add to the list a file that does not contain this Tree
459 TFile *curf = TFile::Open(nextsource->GetName());
461 Bool_t mustAdd = kFALSE;
462 if (curf->FindKey(obj_name)) {
465 //we could be more clever here. No need to import the object
466 //we are missing a function in TDirectory
467 TObject *aobj = curf->Get(obj_name);
468 if (aobj) { mustAdd = kTRUE; delete aobj;}
471 globChain->Add(nextsource->GetName());
475 nextsource = (TFile*)sourcelist->After( nextsource );
478 } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
479 // it's a subdirectory
481 cout << "Found subdirectory " << obj->GetName() << endl;
482 // create a new subdir of same name and title in the target file
484 TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
486 // newdir is now the starting point of another round of merging
487 // newdir still knows its depth within the target file via
488 // GetPath(), so we can still figure out where we are in the recursion
489 status = MergeRootfile( newdir, sourcelist);
490 if (status) return status;
492 } else if ( obj->InheritsFrom(TObject::Class())
493 && obj->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
494 // object implements Merge(TCollection*)
496 // loop over all source files and merge same-name object
497 TFile *nextsource = (TFile*)sourcelist->After( first_source );
498 while ( nextsource ) {
499 // make sure we are at the correct directory level by cd'ing to path
500 TDirectory *ndir = nextsource->GetDirectory(path);
503 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
505 TObject *hobj = key2->ReadObj();
507 cout << "Failed to get the object with key " << key2->GetName() << " from " <<
508 ndir->GetFile()->GetName() << "/" << ndir->GetName() << endl;
509 nextsource = (TFile*)sourcelist->After( nextsource );
513 hobj->ResetBit(kMustCleanup);
516 obj->Execute("Merge", listHargs.Data(), &error); // RS Probleme here
518 cerr << "Error calling Merge() on " << obj->GetName()
519 << " with the corresponding object in " << nextsource->GetName() << endl;
524 nextsource = (TFile*)sourcelist->After( nextsource );
526 } else if ( obj->IsA()->InheritsFrom( THStack::Class() ) ) {
527 THStack *hstack1 = (THStack*) obj;
528 TList* l = new TList();
530 // loop over all source files and merge the histos of the
531 // corresponding THStacks with the one pointed to by "hstack1"
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);
538 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(hstack1->GetName());
540 THStack *hstack2 = (THStack*) key2->ReadObj();
541 l->Add(hstack2->GetHists()->Clone());
546 nextsource = (TFile*)sourcelist->After( nextsource );
548 hstack1->GetHists()->Merge(l);
551 // object is of no type that we can merge
552 cout << "Cannot merge object type, name: "
553 << obj->GetName() << " title: " << obj->GetTitle() << endl;
555 // loop over all source files and write similar objects directly to the output file
556 TFile *nextsource = (TFile*)sourcelist->After( first_source );
557 while ( nextsource ) {
558 // make sure we are at the correct directory level by cd'ing to path
559 TDirectory *ndir = nextsource->GetDirectory(path);
562 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
564 TObject *nobj = key2->ReadObj();
565 nobj->ResetBit(kMustCleanup);
566 int nbytes1 = target->WriteTObject(nobj, key2->GetName(), "SingleKey" );
567 if (nbytes1 <= 0) status = -1;
571 nextsource = (TFile*)sourcelist->After( nextsource );
575 // now write the merged histogram (which is "in" obj) to the target file
576 // note that this will just store obj in the current directory level,
577 // which is not persistent until the complete directory itself is stored
578 // by "target->Write()" below
581 //!!if the object is a tree, it is stored in globChain...
582 if(obj->IsA()->InheritsFrom( TDirectory::Class() )) {
583 //printf("cas d'une directory\n");
584 } else if(obj->IsA()->InheritsFrom( TTree::Class() )) {
587 globChain->Merge(target->GetFile(),0,"keep fast");
591 int nbytes2 = obj->Write( key->GetName(), TObject::kSingleKey );
592 if (nbytes2 <= 0) status = -1;
596 } // while ( ( TKey *key = (TKey*)nextkey() ) )
597 first_source = (TDirectory*)sourcelist->After(first_source);
599 // save modifications to target file
600 target->SaveSelf(kTRUE);
602 gSystem->GetProcInfo(&procInfo);
603 AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
608 //___________________________________________________________________________
609 int AliFileMerger::OpenNextChunks(const TList* namesList, TList* filesList, Int_t from, Int_t to)
611 gSystem->GetProcInfo(&procInfo);
612 AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
615 int nEnt = namesList->GetEntries();
616 from = from<nEnt ? from : nEnt;
617 to = to<nEnt ? to : nEnt;
619 for (int i=from;i<=to;i++) {
620 TNamed* fnam = (TNamed*)namesList->At(i);
622 TString fnamS(fnam->GetName());
623 gSystem->ExpandPathName(fnamS);
624 if (fnamS.BeginsWith("alien://") && !gGrid) TGrid::Connect("alien");
625 TFile* source = TFile::Open(fnam->GetName());
626 if( source==0 ) { printf("Failed to open file %s, will skip\n",fnam->GetName()); continue; }
627 filesList->Add(source);
628 printf("Opened file %s\n",fnam->GetName());
631 gSystem->GetProcInfo(&procInfo);
632 AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
638 //___________________________________________________________________________
639 int AliFileMerger::AddFile(TList* namesList, std::string entry)
641 // add a new file to the list of files
642 // static int count(0);
643 if( entry.empty() ) return 0;
644 size_t j =entry.find_first_not_of(' ');
645 if( j==std::string::npos ) return 0;
646 entry = entry.substr(j);
647 if( entry.substr(0,1)=="@") {
648 std::ifstream indirect_file(entry.substr(1).c_str() );
649 if( ! indirect_file.is_open() ) {
650 std::cerr<< "Could not open indirect file " << entry.substr(1) << std::endl;
653 while( indirect_file ){
655 std::getline(indirect_file, line);
656 if( AddFile(namesList, line)!=0 ) return 1;;
660 // cout << "Source file " << (++count) << ": " << entry << endl;
661 namesList->Add(new TNamed(entry,""));