ROOT
master
Reference Guide
Loading...
Searching...
No Matches
testMergeCont.C File Reference
Tutorials
»
IO tutorials
Detailed Description
Macro demonstrating the merging of containers.
TFile
*
f
;
TSeqCollection
*GetCollection()
{
TObject
*obj;
#ifndef ClingWorkAroundMissingDynamicScope
# define ClingWorkAroundMissingDynamicScope
#endif
f
=
TFile::Open
(
"hsimple.root"
);
if
( !
f
) {
#ifdef ClingWorkAroundMissingDynamicScope
f
= (
TFile
*)
gROOT
->ProcessLine(
"hsimple(1);"
);
#else
f
=
hsimple
(1);
#endif
}
gROOT
->cd();
TList
*l0 =
new
TList
();
TList
*l01 =
new
TList
();
TH1
*hpx = (
TH1
*)
f
->Get(
"hpx"
);
printf(
"Adding hpx: %d entries\n"
, (
int
)hpx->
GetEntries
());
l01->
Add
(hpx);
TH1
*hpxpy = (
TH1
*)
f
->Get(
"hpxpy"
);
l01->
Add
(hpxpy);
TH1
*hprof = (
TH1
*)
f
->Get(
"hprof"
);
l0->
Add
(hprof);
l0->
Add
(l01);
return
l0;
}
void
testMergeCont()
{
TString
tutdir =
gROOT
->GetTutorialDir();
gROOT
->LoadMacro(tutdir+
"/hsimple.C"
);
TList
*list1 = (
TList
*)GetCollection();
TList
*inputs =
new
TList
();
for
(
Int_t
i=0; i<10; i++) {
inputs->
AddAt
(GetCollection(),0);
list1->
Merge
(inputs);
inputs->
Delete
();
f
->Close();
}
delete
inputs;
TH1F
*hpx = (
TH1F
*)(((
TList
*)list1->
At
(1))->At(0));
printf(
"============================================\n"
);
printf(
"Total hpx: %d entries\n"
, (
int
)hpx->
GetEntries
());
hpx->
Draw
();
list1->
Delete
();
delete
list1;
}
f
#define f(i)
Definition
RSha256.hxx:104
Int_t
int Int_t
Definition
RtypesCore.h:45
gROOT
#define gROOT
Definition
TROOT.h:406
TFile
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition
TFile.h:53
TFile::Open
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition
TFile.cxx:4086
TH1F
1-D histogram with a float per channel (see TH1 documentation)
Definition
TH1.h:623
TH1
TH1 is the base class of all histogram classes in ROOT.
Definition
TH1.h:59
TH1::Draw
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition
TH1.cxx:3068
TH1::GetEntries
virtual Double_t GetEntries() const
Return the current number of entries.
Definition
TH1.cxx:4431
TList
A doubly linked list.
Definition
TList.h:38
TList::AddAt
void AddAt(TObject *obj, Int_t idx) override
Insert object at position idx in the list.
Definition
TList.cxx:304
TList::Add
void Add(TObject *obj) override
Definition
TList.h:81
TList::Delete
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
Definition
TList.cxx:468
TList::At
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition
TList.cxx:355
TObject
Mother of all ROOT objects.
Definition
TObject.h:41
TSeqCollection
Sequenceable collection abstract base class.
Definition
TSeqCollection.h:28
TSeqCollection::Merge
Long64_t Merge(TCollection *list)
Merge this collection with all collections coming in the input list.
Definition
TSeqCollection.cxx:185
TString
Basic string class.
Definition
TString.h:139
hsimple
Definition
hsimple.py:1
Author
The Root Team
Definition in file
testMergeCont.C
.
tutorials
io
testMergeCont.C
ROOT master - Reference Guide Generated on Thu Dec 19 2024 09:47:43 (GVA Time) using Doxygen 1.9.8