Logo ROOT   6.18/05
Reference Guide
testMergeCont.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_io
3/// \notebook
4/// Macro demonstrating the merging of containers.
5/// \macro_code
6///
7/// \author The Root Team
8
9TFile *f;
10
11
12TSeqCollection *GetCollection()
13{
14 TObject *obj;
15#ifndef ClingWorkAroundMissingDynamicScope
16# define ClingWorkAroundMissingDynamicScope
17#endif
18 f = TFile::Open("hsimple.root");
19 if( !f ) {
20#ifdef ClingWorkAroundMissingDynamicScope
21 f = (TFile*)gROOT->ProcessLine("hsimple(1);");
22#else
23 f = hsimple(1);
24#endif
25 }
26 gROOT->cd();
27 TList *l0 = new TList();
28 TList *l01 = new TList();
29 TH1 *hpx = (TH1*)f->Get("hpx");
30 printf("Adding hpx: %d entries\n", (int)hpx->GetEntries());
31 l01->Add(hpx);
32 TH1 *hpxpy = (TH1*)f->Get("hpxpy");
33 l01->Add(hpxpy);
34 TH1 *hprof = (TH1*)f->Get("hprof");
35 l0->Add(hprof);
36 l0->Add(l01);
37 return l0;
38}
39
40void testMergeCont()
41{
42 TString tutdir = gROOT->GetTutorialDir();
43 gROOT->LoadMacro(tutdir+"/hsimple.C");
44 TList *list1 = (TList *)GetCollection();
45 TList *inputs = new TList();
46 for (Int_t i=0; i<10; i++) {
47 inputs->AddAt(GetCollection(),0);
48 list1->Merge(inputs);
49 inputs->Delete();
50 f->Close();
51 }
52 delete inputs;
53 TH1F *hpx = (TH1F*)(((TList*)list1->At(1))->At(0));
54 printf("============================================\n");
55 printf("Total hpx: %d entries\n", (int)hpx->GetEntries());
56 hpx->Draw();
57 list1->Delete();
58 delete list1;
59}
#define f(i)
Definition: RSha256.hxx:104
int Int_t
Definition: RtypesCore.h:41
#define gROOT
Definition: TROOT.h:414
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseGeneralPurpose, Int_t netopt=0)
Create / open a file.
Definition: TFile.cxx:3980
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
The TH1 histogram class.
Definition: TH1.h:56
virtual Double_t GetEntries() const
Return the current number of entries.
Definition: TH1.cxx:4277
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
A doubly linked list.
Definition: TList.h:44
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual void AddAt(TObject *obj, Int_t idx)
Insert object at position idx in the list.
Definition: TList.cxx:303
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:354
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
Mother of all ROOT objects.
Definition: TObject.h:37
Sequenceable collection abstract base class.
Long64_t Merge(TCollection *list)
Merge this collection with all collections coming in the input list.
Basic string class.
Definition: TString.h:131