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"
75 ClassImp(AliFileMerger)
77 ProcInfo_t procInfo;//TMP
79 ////////////////////////////////////////////////////////////////////////
81 AliFileMerger::AliFileMerger():
89 // Default constructor
93 //______________________________________________________________________
95 AliFileMerger::AliFileMerger(const char* name):
108 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern, Bool_t dontOverwrite){
111 // Merge the files coming out of the calibration job
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);
122 // loop over the results
124 sourcelist.SetOwner(kTRUE);
126 while((map=(TMap*)nextmap())) {
128 TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
129 if (!objs || !objs->GetString().Length()) {
130 // Nothing found - skip this output
134 printf("looking for file %s\n",(objs->GetString()).Data());
135 AddFile(&sourcelist, (objs->GetString()).Data());;
138 IterList(&sourcelist, outputFileName, dontOverwrite);
142 void AliFileMerger::IterList(const TList* namesList, const char* outputFileName, Bool_t dontOverwrite)
144 // merge in steps or in one go
146 gSystem->GetProcInfo(&procInfo);
147 AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
149 TString outputFile(outputFileName);
150 gSystem->ExpandPathName(outputFile);
152 int nFiles = namesList->GetEntries();
153 int maxSrcOpen = fMaxFilesOpen - 1;
155 filesList.SetOwner(kTRUE);
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));
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
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());
178 MergeRootfile(targetTmp, &filesList);
181 filesList.Clear(); // close all open files
183 // nothing to do until needed amount of files is accumulated
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;
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);
200 for (int i=0;i<2;i++) gSystem->Exec(Form("if [ -e %s ]; then \nrm %s\nfi",tmpDest[i].Data(),tmpDest[i].Data()));
202 printf("Merged %d files in %d steps\n",nFiles,++nsteps);
204 gSystem->GetProcInfo(&procInfo);
205 AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
208 void AliFileMerger::IterTXT( const char * fileList, const char* outputFileName, Bool_t dontOverwrite){
210 // Merge the files indicated in the list - fileList
211 // ASCII file option example:
212 // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
214 // Open the input stream
217 // Read the input list of files
221 sourcelist.SetOwner(kTRUE);
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());
230 IterList(&sourcelist, outputFileName, dontOverwrite);
234 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
236 // Storing the results in one single file
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;
249 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
251 // Store the results in separate files (one per object)
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");
263 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
268 static Int_t counter=-1;
270 TObjArray *carray = new TObjArray; //array of the objects inside current file
271 carray->SetOwner(kTRUE);
273 // load all objects to memory
275 TList *farr = fileIn->GetListOfKeys();
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);
289 if (carray->GetEntries()==0) {
294 Int_t entries =carray->GetEntriesFast();
295 for (Int_t i=0; i<entries; i++){
297 TObjArray *templist = new TObjArray(1);
298 templist->SetOwner(kFALSE);
299 TObject *currentObject = carray->At(i);
300 if (!currentObject) {
304 printf("%s\n",currentObject->GetName());
305 callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
306 if (!callEnv.IsValid()) {
310 TString oname=currentObject->GetName();
311 TObject *mergedObject = array->FindObject(currentObject->GetName());
313 array->AddLast(currentObject);
318 templist->AddLast(currentObject);
319 callEnv.SetParam((Long_t) templist);
320 callEnv.Execute(mergedObject);
321 AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);
328 Bool_t AliFileMerger::IsAccepted(TString name){
330 // Accept/reject logic
331 // name - name of the entry
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
340 for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
341 if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE; // entry was selected
344 if (!accept) return kFALSE;
348 for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
349 if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE; // entry was rejected
358 void AliFileMerger::AddReject(const char *reject){
360 // add reject string to the list of entries to be rejected for merging
362 if (!fRejectMask) fRejectMask = new TObjArray;
363 fRejectMask->AddLast(new TObjString(reject));
365 void AliFileMerger::AddAccept(const char *accept){
367 // add reject string to the list of entries to be rejected for merging
369 if (!fAcceptMask) fAcceptMask = new TObjArray;
370 fAcceptMask->AddLast(new TObjString(accept));
375 //___________________________________________________________________________
376 int AliFileMerger::MergeRootfile( TDirectory *target, TList *sourcelist, Bool_t nameFiltering)
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));
384 cout << "Target path: " << target->GetPath() << endl;
385 TString path( (char*)strstr( target->GetPath(), ":" ) );
388 // find 1st valid file
389 TDirectory *first_source = (TDirectory*)sourcelist->First();
391 Int_t nguess = sourcelist->GetSize()+1000;
392 THashList allNames(nguess);
393 ((THashList*)target->GetList())->Rehash(nguess);
394 ((THashList*)target->GetListOfKeys())->Rehash(nguess);
397 listHargs.Form("((TCollection*)0x%lx)", (ULong_t)&listH);
399 while(first_source) {
401 TDirectory *current_sourcedir = first_source->GetDirectory(path);
402 if (!current_sourcedir) {
403 first_source = (TDirectory*)sourcelist->After(first_source);
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);
416 while ( (key = (TKey*)nextkey())) {
417 if (current_sourcedir == target) break;
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());
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;
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();
441 TObject *obj = key->ReadObj();
443 AliError(Form("Failed to get the object with key %s from %s",key->GetName(),current_sourcedir->GetFile()->GetName()));
447 if ( obj->IsA()->InheritsFrom( TTree::Class() ) ) {
449 // loop over all source files create a chain of Trees "globChain"
453 obj_name = path + "/" + obj->GetName();
455 obj_name = obj->GetName();
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());
464 Bool_t mustAdd = kFALSE;
465 if (curf->FindKey(obj_name)) {
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;}
474 globChain->Add(nextsource->GetName());
478 nextsource = (TFile*)sourcelist->After( nextsource );
481 } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
482 // it's a subdirectory
484 cout << "Found subdirectory " << obj->GetName() << endl;
485 // create a new subdir of same name and title in the target file
487 TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
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;
495 } else if ( obj->InheritsFrom(TObject::Class())
496 && obj->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
497 // object implements Merge(TCollection*)
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);
506 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
508 TObject *hobj = key2->ReadObj();
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 );
516 hobj->ResetBit(kMustCleanup);
519 obj->Execute("Merge", listHargs.Data(), &error); // RS Probleme here
521 cerr << "Error calling Merge() on " << obj->GetName()
522 << " with the corresponding object in " << nextsource->GetName() << endl;
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"))
529 TMethodCall getEntries(obj->IsA(), "GetEntries", "");
530 getEntries.Execute(obj, numberOfEntries);
532 AliSysInfo::AddStamp(nameK.Data(),1,counterK,counterF++,numberOfEntries);
535 nextsource = (TFile*)sourcelist->After( nextsource );
537 } else if ( obj->IsA()->InheritsFrom( THStack::Class() ) ) {
538 THStack *hstack1 = (THStack*) obj;
539 TList* l = new TList();
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);
549 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(hstack1->GetName());
551 THStack *hstack2 = (THStack*) key2->ReadObj();
552 l->Add(hstack2->GetHists()->Clone());
554 AliSysInfo::AddStamp(nameK.Data(),1,counterK,counterF++);
558 nextsource = (TFile*)sourcelist->After( nextsource );
560 hstack1->GetHists()->Merge(l);
563 // object is of no type that we can merge
564 cout << "Cannot merge object type, name: "
565 << obj->GetName() << " title: " << obj->GetTitle() << endl;
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);
574 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
576 TObject *nobj = key2->ReadObj();
577 nobj->ResetBit(kMustCleanup);
578 int nbytes1 = target->WriteTObject(nobj, key2->GetName(), "SingleKey" );
579 if (nbytes1 <= 0) status = -1;
583 nextsource = (TFile*)sourcelist->After( nextsource );
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
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() )) {
599 globChain->Merge(target->GetFile(),0,"keep fast");
603 int nbytes2 = obj->Write( key->GetName(), TObject::kSingleKey );
604 if (nbytes2 <= 0) status = -1;
608 } // while ( ( TKey *key = (TKey*)nextkey() ) )
609 first_source = (TDirectory*)sourcelist->After(first_source);
611 // save modifications to target file
612 target->SaveSelf(kTRUE);
614 gSystem->GetProcInfo(&procInfo);
615 AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
620 //___________________________________________________________________________
621 int AliFileMerger::OpenNextChunks(const TList* namesList, TList* filesList, Int_t from, Int_t to)
623 gSystem->GetProcInfo(&procInfo);
624 AliInfo(Form(">> memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
627 int nEnt = namesList->GetEntries();
628 from = from<nEnt ? from : nEnt;
629 to = to<nEnt ? to : nEnt;
631 for (int i=from;i<=to;i++) {
632 TNamed* fnam = (TNamed*)namesList->At(i);
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());
643 gSystem->GetProcInfo(&procInfo);
644 AliInfo(Form("<< memory usage %ld %ld", procInfo.fMemResident, procInfo.fMemVirtual));
650 //___________________________________________________________________________
651 int AliFileMerger::AddFile(TList* namesList, std::string entry)
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;
665 while( indirect_file ){
667 std::getline(indirect_file, line);
668 if( AddFile(namesList, line)!=0 ) return 1;;
672 // cout << "Source file " << (++count) << ": " << entry << endl;
673 namesList->Add(new TNamed(entry,""));