Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
mergeSelective.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_io
3 /// \notebook -nodraw
4 /// Merge only part of the content of a set of files.
5 /// This macro demonstrates how to merge only a part of the content of a set
6 /// of input files, specified via the interface.
7 /// ~~~{.cpp}
8 /// TFileMerger::AddObjectNames(const char *names)
9 /// ~~~
10 /// The method can be called several times to add object names, or using a single
11 /// string with names separated by a blank. Directory names contained in the files
12 /// to be merged are accepted.
13 ///
14 /// Two modes are supported:
15 /// 1. kOnlyListed: via <tt>TFileMerger::PartialMerge(kOnlyListed)</tt>
16 /// This will merge only the objects in the files having the names in the
17 /// specified list. If a folder is specified, its whole content will be merged
18 ///
19 /// 2. kSkipListed: via <tt>TFileMerger::PartialMerge(kSkipListed)</tt>
20 /// This will skip merging of specified objects. If a folder is specified, its
21 /// whole content will be skipped.
22 ///
23 /// Important note: the kOnlyListed and kSkipListed flags have to be bitwise OR-ed
24 /// on top of the merging defaults: kAll | kIncremental (as in the example)
25 ///
26 /// The files to be merged have the following structure:
27 /// - hpx (TH1F)
28 /// - hpxpy (TH2F)
29 /// - hprof (TProfile)
30 /// - ntuple (TNtuple)
31 /// - folder (TDirectory)
32 /// - hpx1 (TH1F)
33 ///
34 /// The example first merges exclusively hprof and the content of "folder",
35 /// producing the file exclusive.root, then merges all content but skipping
36 /// hprof and the content of "folder". The result can be inspected in the
37 /// browser.
38 ///
39 /// \macro_code
40 ///
41 /// \author The Root Team
42 
43 
44 void CreateFile(const char *);
45 
46 void mergeSelective(Int_t nfiles=5)
47 {
48 
49  // Create the files to be merged
50  TStopwatch timer;
51  timer.Start();
52  TString tutdir = gROOT->GetTutorialDir();
53  if (gROOT->LoadMacro(tutdir + "/hsimple.C")) return;
54  Int_t i;
55  for (i=0; i<nfiles; i++) CreateFile(Form("tomerge%03d.root",i));
56 
57  //------------------------------------
58  // Merge only the listed objects
59  //------------------------------------
60  TFileMerger *fm;
61  fm = new TFileMerger(kFALSE);
62  fm->OutputFile("exclusive.root");
63  fm->AddObjectNames("hprof folder");
64  for (i=0; i<nfiles; i++) fm->AddFile(Form("tomerge%03d.root",i));
65  // Must add new merging flag on top of the the default ones
66  Int_t default_mode = TFileMerger::kAll | TFileMerger::kIncremental;
67  Int_t mode = default_mode | TFileMerger::kOnlyListed;
68  fm->PartialMerge(mode);
69  fm->Reset();
70 
71  //------------------------------------
72  // Skip merging of the listed objects
73  //------------------------------------
74  fm->OutputFile("skipped.root");
75  fm->AddObjectNames("hprof folder");
76  for (i=0; i<nfiles; i++) fm->AddFile(Form("tomerge%03d.root",i));
77  // Must add new merging flag on top of the the default ones
78  mode = default_mode | TFileMerger::kSkipListed;
79  fm->PartialMerge(mode);
80  delete fm;
81 
82 
83  // Cleanup initial files
84  for (i=0; i<nfiles; i++) gSystem->Unlink(Form("tomerge%03d.root",i));
85  // Open files to inspect in the browser
86  TFile::Open("exclusive.root");
87  TFile::Open("skipped.root");
88  new TBrowser();
89  timer.Stop();
90  timer.Print();
91 }
92 
93 void CreateFile(const char *fname)
94 {
95  TFile *example = (TFile*)gROOT->ProcessLineFast("hsimple(1)");
96  if (!example) return;
97  TH1F *hpx = (TH1F*)example->Get("hpx");
98  hpx->SetName("hpx1");
99  TFile::Cp(example->GetName(), fname);
100  TFile *file = TFile::Open(fname, "UPDATE");
101  file->mkdir("folder")->cd();
102  hpx->Write();
103  file->Close();
104  example->Close();
105  TString sname(fname);
106  if (sname.Contains("000")) {
107  TFile::Cp(fname, "original.root");
108  TFile::Open("original.root");
109  }
110 }