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"
70 ClassImp(AliFileMerger)
72 ////////////////////////////////////////////////////////////////////////
74 AliFileMerger::AliFileMerger():
81 // Default constructor
85 //______________________________________________________________________
87 AliFileMerger::AliFileMerger(const char* name):
99 void AliFileMerger::IterAlien(const char* outputDir, const char* outputFileName, const char* pattern, Bool_t dontOverwrite){
102 // Merge the files coming out of the calibration job
104 TString outputFile(outputFileName);
105 gSystem->ExpandPathName(outputFile);
107 // looking for files to be merged in the output directory
108 command = Form("find %s/ *%s", outputDir, pattern);
109 printf("command: %s\n", command.Data());
110 TGrid::Connect("alien://");
111 TGridResult *res = gGrid->Command(command);
115 // loop over the results
117 sourcelist.SetOwner(kTRUE);
119 while((map=(TMap*)nextmap())) {
121 TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
122 if (!objs || !objs->GetString().Length()) {
123 // Nothing found - skip this output
127 printf("looking for file %s\n",(objs->GetString()).Data());
128 AddFile(&sourcelist, (objs->GetString()).Data());;
130 printf("List of files to be merged:\n");
133 TFile* target = TFile::Open( outputFile.Data(), (dontOverwrite ? "CREATE":"RECREATE") );
134 if (!target || target->IsZombie()) {
135 cerr << "Error opening target file (does " << outputFileName << " exist?)." << endl;
136 cerr << "Use force = kTRUE to re-creation of output file." << endl;
139 MergeRootfile( target, &sourcelist);
144 void AliFileMerger::IterTXT( const char * fileList, const char* outputFileName, Bool_t dontOverwrite){
146 // Merge the files indicated in the list - fileList
147 // ASCII file option example:
148 // find `pwd`/ | grep AliESDfriends_v1.root > calib.list
150 // Open the input stream
153 // Read the input list of files
157 sourcelist.SetOwner(kTRUE);
160 if (!objfile.Contains(".root")) continue; // protection
161 gSystem->ExpandPathName(objfile);
162 printf("Add file:Counter\t%d\tMerging file %s\n",counter++,objfile.Data());
163 AddFile(&sourcelist, objfile.Data());
166 printf("List of files to be merged:\n");
169 TString outputFile(outputFileName);
170 gSystem->ExpandPathName(outputFile);
171 TFile* target = TFile::Open( outputFile.Data(), (dontOverwrite ? "CREATE":"RECREATE") );
172 if (!target || target->IsZombie()) {
173 cerr << "Error opening target file (does " << outputFileName << " exist?)." << endl;
174 cerr << "Use force = kTRUE to re-creation of output file." << endl;
178 MergeRootfile( target, &sourcelist);
182 void AliFileMerger::StoreResults(TObjArray * array, const char* outputFileName){
184 // Storing the results in one single file
186 TFile *f = new TFile(outputFileName,"recreate");
187 for (Int_t i=0; i<array->GetEntries(); i++){
188 TObject *object0 = array->At(i);
189 if (!object0) continue;
197 void AliFileMerger::StoreSeparateResults(TObjArray * array, const char* outputFileName){
199 // Store the results in separate files (one per object)
201 for (Int_t i=0; i<array->GetEntries(); i++){
202 TObject *object0 = array->At(i);
203 if (!object0) continue;
204 TFile *f = new TFile(Form("%s_%s.root",outputFileName,object0->GetName()),"recreate");
211 void AliFileMerger::Merge(TFile* fileIn, TObjArray * array){
216 static Int_t counter=-1;
218 TObjArray *carray = new TObjArray; //array of the objects inside current file
219 carray->SetOwner(kTRUE);
221 // load all objects to memory
223 TList *farr = fileIn->GetListOfKeys();
228 for (Int_t ical=0; ical<farr->GetEntries(); ical++){
229 if (!farr->At(ical)) continue;
230 TString name(farr->At(ical)->GetName());
231 if (!IsAccepted(name)) continue; // skip not accepted entries
232 TObject *obj = fileIn->Get(name.Data());
233 if (obj) carray->AddLast(obj);
234 AliSysInfo::AddStamp(name.Data(),1,ical,counter);
237 if (carray->GetEntries()==0) {
242 Int_t entries =carray->GetEntriesFast();
243 for (Int_t i=0; i<entries; i++){
245 TObjArray *templist = new TObjArray(1);
246 templist->SetOwner(kFALSE);
247 TObject *currentObject = carray->At(i);
248 if (!currentObject) {
252 printf("%s\n",currentObject->GetName());
253 callEnv.InitWithPrototype(currentObject->IsA(), "Merge", "TCollection*");
254 if (!callEnv.IsValid()) {
258 TString oname=currentObject->GetName();
259 TObject *mergedObject = array->FindObject(currentObject->GetName());
261 array->AddLast(currentObject);
266 templist->AddLast(currentObject);
267 callEnv.SetParam((Long_t) templist);
268 callEnv.Execute(mergedObject);
269 AliSysInfo::AddStamp(currentObject->GetName(),2,i,counter);
276 Bool_t AliFileMerger::IsAccepted(TString name){
278 // Accept/reject logic
279 // name - name of the entry
281 // if fAcceptMask specified - entry has to be in list of selected
282 // if fRejectMask speciefied - entry with name speciief in the list are rejected
288 for (Int_t iaccept=0; iaccept<fAcceptMask->GetEntries(); iaccept++){
289 if (name.Contains(fAcceptMask->At(iaccept)->GetName())) accept=kTRUE; // entry was selected
292 if (!accept) return kFALSE;
296 for (Int_t ireject=0; ireject<fRejectMask->GetEntries(); ireject++){
297 if (name.Contains(fRejectMask->At(ireject)->GetName())) accept=kFALSE; // entry was rejected
306 void AliFileMerger::AddReject(const char *reject){
308 // add reject string to the list of entries to be rejected for merging
310 if (!fRejectMask) fRejectMask = new TObjArray;
311 fRejectMask->AddLast(new TObjString(reject));
313 void AliFileMerger::AddAccept(const char *accept){
315 // add reject string to the list of entries to be rejected for merging
317 if (!fAcceptMask) fAcceptMask = new TObjArray;
318 fAcceptMask->AddLast(new TObjString(accept));
323 //___________________________________________________________________________
324 int AliFileMerger::MergeRootfile( TDirectory *target, TList *sourcelist)
326 // Merge all objects in a directory
327 // modified version of root's hadd.cxx
330 cout << "Target path: " << target->GetPath() << endl;
331 TString path( (char*)strstr( target->GetPath(), ":" ) );
334 // find 1st valid file
335 TDirectory *first_source = (TDirectory*)sourcelist->First();
337 Int_t nguess = sourcelist->GetSize()+1000;
338 THashList allNames(nguess);
339 ((THashList*)target->GetList())->Rehash(nguess);
340 ((THashList*)target->GetListOfKeys())->Rehash(nguess);
343 listHargs.Form("((TCollection*)0x%lx)", (ULong_t)&listH);
345 while(first_source) {
347 TDirectory *current_sourcedir = first_source->GetDirectory(path);
348 if (!current_sourcedir) {
349 first_source = (TDirectory*)sourcelist->After(first_source);
352 // loop over all keys in this directory
353 TChain *globChain = 0;
354 TIter nextkey( current_sourcedir->GetListOfKeys() );
355 TKey *key, *oldkey=0;
356 //gain time, do not add the objects in the list in memory
357 TH1::AddDirectory(kFALSE);
361 while ( (key = (TKey*)nextkey())) {
362 if (current_sourcedir == target) break;
364 // check if we don't reject this name
365 TString nameK(key->GetName());
366 if (!IsAccepted(nameK)) {
367 if (!counterF) printf("Object %s is in rejection list, skipping...\n",nameK.Data());
371 //keep only the highest cycle number for each key
372 if (oldkey && !strcmp(oldkey->GetName(),key->GetName())) continue;
373 if (!strcmp(key->GetClassName(),"TProcessID")) {key->ReadObj(); continue;}
374 if (allNames.FindObject(key->GetName())) continue;
375 TClass *cl = TClass::GetClass(key->GetClassName());
376 if (!cl || !cl->InheritsFrom(TObject::Class())) {
377 cout << "Cannot merge object type, name: "
378 << key->GetName() << " title: " << key->GetTitle() << endl;
381 allNames.Add(new TObjString(key->GetName()));
382 AliSysInfo::AddStamp(nameK.Data(),1,counterK++,counterF-1);
383 // read object from first source file
384 //current_sourcedir->cd();
385 TObject *obj = key->ReadObj();
386 //printf("keyname=%s, obj=%x\n",key->GetName(),obj);
388 if ( obj->IsA()->InheritsFrom( TTree::Class() ) ) {
390 // loop over all source files create a chain of Trees "globChain"
394 obj_name = path + "/" + obj->GetName();
396 obj_name = obj->GetName();
398 globChain = new TChain(obj_name);
399 globChain->Add(first_source->GetName());
400 TFile *nextsource = (TFile*)sourcelist->After( first_source );
401 while ( nextsource ) {
402 //do not add to the list a file that does not contain this Tree
403 TFile *curf = TFile::Open(nextsource->GetName());
405 Bool_t mustAdd = kFALSE;
406 if (curf->FindKey(obj_name)) {
409 //we could be more clever here. No need to import the object
410 //we are missing a function in TDirectory
411 TObject *aobj = curf->Get(obj_name);
412 if (aobj) { mustAdd = kTRUE; delete aobj;}
415 globChain->Add(nextsource->GetName());
419 nextsource = (TFile*)sourcelist->After( nextsource );
422 } else if ( obj->IsA()->InheritsFrom( TDirectory::Class() ) ) {
423 // it's a subdirectory
425 cout << "Found subdirectory " << obj->GetName() << endl;
426 // create a new subdir of same name and title in the target file
428 TDirectory *newdir = target->mkdir( obj->GetName(), obj->GetTitle() );
430 // newdir is now the starting point of another round of merging
431 // newdir still knows its depth within the target file via
432 // GetPath(), so we can still figure out where we are in the recursion
433 status = MergeRootfile( newdir, sourcelist);
434 if (status) return status;
436 } else if ( obj->InheritsFrom(TObject::Class())
437 && obj->IsA()->GetMethodWithPrototype("Merge", "TCollection*") ) {
438 // object implements Merge(TCollection*)
440 // loop over all source files and merge same-name object
441 TFile *nextsource = (TFile*)sourcelist->After( first_source );
442 while ( nextsource ) {
443 // make sure we are at the correct directory level by cd'ing to path
444 TDirectory *ndir = nextsource->GetDirectory(path);
447 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
449 TObject *hobj = key2->ReadObj();
450 hobj->ResetBit(kMustCleanup);
453 obj->Execute("Merge", listHargs.Data(), &error);
455 cerr << "Error calling Merge() on " << obj->GetName()
456 << " with the corresponding object in " << nextsource->GetName() << endl;
461 nextsource = (TFile*)sourcelist->After( nextsource );
463 } else if ( obj->IsA()->InheritsFrom( THStack::Class() ) ) {
464 THStack *hstack1 = (THStack*) obj;
465 TList* l = new TList();
467 // loop over all source files and merge the histos of the
468 // corresponding THStacks with the one pointed to by "hstack1"
469 TFile *nextsource = (TFile*)sourcelist->After( first_source );
470 while ( nextsource ) {
471 // make sure we are at the correct directory level by cd'ing to path
472 TDirectory *ndir = nextsource->GetDirectory(path);
475 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(hstack1->GetName());
477 THStack *hstack2 = (THStack*) key2->ReadObj();
478 l->Add(hstack2->GetHists()->Clone());
483 nextsource = (TFile*)sourcelist->After( nextsource );
485 hstack1->GetHists()->Merge(l);
488 // object is of no type that we can merge
489 cout << "Cannot merge object type, name: "
490 << obj->GetName() << " title: " << obj->GetTitle() << endl;
492 // loop over all source files and write similar objects directly to the output file
493 TFile *nextsource = (TFile*)sourcelist->After( first_source );
494 while ( nextsource ) {
495 // make sure we are at the correct directory level by cd'ing to path
496 TDirectory *ndir = nextsource->GetDirectory(path);
499 TKey *key2 = (TKey*)gDirectory->GetListOfKeys()->FindObject(key->GetName());
501 TObject *nobj = key2->ReadObj();
502 nobj->ResetBit(kMustCleanup);
503 int nbytes1 = target->WriteTObject(nobj, key2->GetName(), "SingleKey" );
504 if (nbytes1 <= 0) status = -1;
508 nextsource = (TFile*)sourcelist->After( nextsource );
512 // now write the merged histogram (which is "in" obj) to the target file
513 // note that this will just store obj in the current directory level,
514 // which is not persistent until the complete directory itself is stored
515 // by "target->Write()" below
518 //!!if the object is a tree, it is stored in globChain...
519 if(obj->IsA()->InheritsFrom( TDirectory::Class() )) {
520 //printf("cas d'une directory\n");
521 } else if(obj->IsA()->InheritsFrom( TTree::Class() )) {
524 globChain->Merge(target->GetFile(),0,"keep fast");
528 int nbytes2 = obj->Write( key->GetName(), TObject::kSingleKey );
529 if (nbytes2 <= 0) status = -1;
533 } // while ( ( TKey *key = (TKey*)nextkey() ) )
534 first_source = (TDirectory*)sourcelist->After(first_source);
536 // save modifications to target file
537 target->SaveSelf(kTRUE);
542 //___________________________________________________________________________
543 int AliFileMerger::AddFile(TList* sourcelist, std::string entry)
545 // add a new file to the list of files
546 // static int count(0);
547 if( entry.empty() ) return 0;
548 size_t j =entry.find_first_not_of(' ');
549 if( j==std::string::npos ) return 0;
550 entry = entry.substr(j);
551 if( entry.substr(0,1)=="@") {
552 std::ifstream indirect_file(entry.substr(1).c_str() );
553 if( ! indirect_file.is_open() ) {
554 std::cerr<< "Could not open indirect file " << entry.substr(1) << std::endl;
557 while( indirect_file ){
559 std::getline(indirect_file, line);
560 if( AddFile(sourcelist, line)!=0 )return 1;;
564 // cout << "Source file " << (++count) << ": " << entry << endl;
566 TFile* source = TFile::Open( entry.c_str());
568 cout << "Failed to open " << entry << " will skip" << endl;
571 sourcelist->Add(source);