Merge Selective


Merge only part of the content of a set of files. This macro demonstrates how to merge only a part of the content of a set of input files, specified via the interface.

TFileMerger::AddObjectNames(const char *names)

The method can be called several times to add object names, or using a single string with names separated by a blank. Directory names contained in the files to be merged are accepted.

Two modes are supported:

  1. kOnlyListed: via TFileMerger::PartialMerge(kOnlyListed) This will merge only the objects in the files having the names in the specified list. If a folder is specified, its whole content will be merged

  2. kSkipListed: via TFileMerger::PartialMerge(kSkipListed) This will skip merging of specified objects. If a folder is specified, its whole content will be skipped.

Important note: the kOnlyListed and kSkipListed flags have to be bitwise OR-ed on top of the merging defaults: kAll | kIncremental (as in the example)

The files to be merged have the following structure:

  • hpx (TH1F)
  • hpxpy (TH2F)
  • hprof (TProfile)
  • ntuple (TNtuple)
  • folder (TDirectory)
    • hpx1 (TH1F)

The example first merges exclusively hprof and the content of "folder", producing the file exclusive.root, then merges all content but skipping hprof and the content of "folder". The result can be inspected in the browser.

Author: The Root Team
This notebook tutorial was automatically generated with ROOTBOOK-izer (Beta) from the macro found in the ROOT repository on Thursday, January 19, 2017 at 04:31 PM.


In [1]:
void CreateFile(const char *);



A helper function is created:


In [2]:
%%cpp -d
void CreateFile(const char *fname)
{
   TFile *example = (TFile*)gROOT->ProcessLineFast("hsimple(1)");
   if (!example) return;
   TH1F *hpx = (TH1F*)example->Get("hpx");
   hpx->SetName("hpx1");
   TFile::Cp(example->GetName(), fname);
   TFile *file = TFile::Open(fname, "UPDATE");
   file->mkdir("folder")->cd();
   hpx->Write();
   file->Close();
   example->Close();
   TString sname(fname);
   if (sname.Contains("000")) {
      TFile::Cp(fname, "original.root");
      TFile::Open("original.root");
   }
}

Arguments are defined.


In [3]:
Int_t nfiles=5;

Create the files to be merged


In [4]:
TStopwatch timer;
timer.Start();
TString tutdir = gROOT->GetTutorialsDir();
if (gROOT->LoadMacro(tutdir + "/hsimple.C")) return;
Int_t i;
for (i=0; i<nfiles; i++) CreateFile(Form("tomerge%03d.root",i));


hsimple   : Real Time =   0.05 seconds Cpu Time =   0.03 seconds
[TFile::Cp] Total 0.40 MB	|====================| 100.00 % [1431.5 MB/s]
[TFile::Cp] Total 0.40 MB	|====================| 100.00 % [1952.0 MB/s]
[TFile::Cp] Total 0.40 MB	|====================| 100.00 % [1788.2 MB/s]
[TFile::Cp] Total 0.40 MB	|====================| 100.00 % [1612.5 MB/s]
[TFile::Cp] Total 0.40 MB	|====================| 100.00 % [1945.7 MB/s]
[TFile::Cp] Total 0.40 MB	|====================| 100.00 % [2006.9 MB/s]

Merge only the listed objects


In [5]:
TFileMerger *fm;
fm = new TFileMerger(kFALSE);
fm->OutputFile("exclusive.root");
fm->AddObjectNames("hprof folder");
for (i=0; i<nfiles; i++) fm->AddFile(Form("tomerge%03d.root",i));

Must add new merging flag on top of the the default ones


In [6]:
Int_t default_mode = TFileMerger::kAll | TFileMerger::kIncremental;
Int_t mode = default_mode | TFileMerger::kOnlyListed;
fm->PartialMerge(mode);
fm->Reset();

Skip merging of the listed objects


In [7]:
fm->OutputFile("skipped.root");
fm->AddObjectNames("hprof folder");
for (i=0; i<nfiles; i++) fm->AddFile(Form("tomerge%03d.root",i));

Must add new merging flag on top of the the default ones


In [8]:
mode = default_mode | TFileMerger::kSkipListed;
fm->PartialMerge(mode);
delete fm;

Cleanup initial files


In [9]:
for (i=0; i<nfiles; i++) gSystem->Unlink(Form("tomerge%03d.root",i));

Open files to inspect in the browser


In [10]:
TFile::Open("exclusive.root");
TFile::Open("skipped.root");
new TBrowser();
timer.Stop();
timer.Print();


Real time 0:00:00, CP time 0.390