Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
hvector.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tree
3 /// \notebook
4 /// Write and read STL vectors in a tree.
5 ///
6 /// \macro_image
7 /// \macro_code
8 ///
9 /// \author The ROOT Team
10 
11 #include <vector>
12 
13 #include "TFile.h"
14 #include "TTree.h"
15 #include "TCanvas.h"
16 #include "TFrame.h"
17 #include "TH1F.h"
18 #include "TBenchmark.h"
19 #include "TRandom.h"
20 #include "TSystem.h"
21 
22 void write()
23 {
24 
25  TFile *f = TFile::Open("hvector.root","RECREATE");
26 
27  if (!f) { return; }
28 
29  // Create one histograms
30  TH1F *hpx = new TH1F("hpx","This is the px distribution",100,-4,4);
31  hpx->SetFillColor(48);
32 
33  std::vector<float> vpx;
34  std::vector<float> vpy;
35  std::vector<float> vpz;
36  std::vector<float> vrand;
37 
38  // Create a TTree
39  TTree *t = new TTree("tvec","Tree with vectors");
40  t->Branch("vpx",&vpx);
41  t->Branch("vpy",&vpy);
42  t->Branch("vpz",&vpz);
43  t->Branch("vrand",&vrand);
44 
45  // Create a new canvas.
46  TCanvas *c1 = new TCanvas("c1","Dynamic Filling Example",200,10,700,500);
47 
48  gRandom->SetSeed();
49  const Int_t kUPDATE = 1000;
50  for (Int_t i = 0; i < 25000; i++) {
51  Int_t npx = (Int_t)(gRandom->Rndm(1)*15);
52 
53  vpx.clear();
54  vpy.clear();
55  vpz.clear();
56  vrand.clear();
57 
58  for (Int_t j = 0; j < npx; ++j) {
59 
60  Float_t px,py,pz;
61  gRandom->Rannor(px,py);
62  pz = px*px + py*py;
63  Float_t random = gRandom->Rndm(1);
64 
65  hpx->Fill(px);
66 
67  vpx.emplace_back(px);
68  vpy.emplace_back(py);
69  vpz.emplace_back(pz);
70  vrand.emplace_back(random);
71 
72  }
73  if (i && (i%kUPDATE) == 0) {
74  if (i == kUPDATE) hpx->Draw();
75  c1->Modified();
76  c1->Update();
77  if (gSystem->ProcessEvents())
78  break;
79  }
80  t->Fill();
81  }
82  f->Write();
83 
84  delete f;
85 }
86 
87 
88 void read()
89 {
90 
91  TFile *f = TFile::Open("hvector.root","READ");
92 
93  if (!f) { return; }
94 
95  TTree *t; f->GetObject("tvec",t);
96 
97  std::vector<float> *vpx = 0;
98 
99  // Create a new canvas.
100  TCanvas *c1 = new TCanvas("c1","Dynamic Filling Example",200,10,700,500);
101 
102  const Int_t kUPDATE = 1000;
103 
104  TBranch *bvpx = 0;
105  t->SetBranchAddress("vpx",&vpx,&bvpx);
106 
107 
108  // Create one histograms
109  TH1F *h = new TH1F("h","This is the px distribution",100,-4,4);
110  h->SetFillColor(48);
111 
112  for (Int_t i = 0; i < 25000; i++) {
113 
114  Long64_t tentry = t->LoadTree(i);
115  bvpx->GetEntry(tentry);
116 
117  for (UInt_t j = 0; j < vpx->size(); ++j) {
118 
119  h->Fill(vpx->at(j));
120 
121  }
122  if (i && (i%kUPDATE) == 0) {
123  if (i == kUPDATE) h->Draw();
124  c1->Modified();
125  c1->Update();
126  if (gSystem->ProcessEvents())
127  break;
128  }
129  }
130 
131  // Since we passed the address of a local variable we need
132  // to remove it.
133  t->ResetBranchAddresses();
134 }
135 
136 
137 void hvector()
138 {
139  gBenchmark->Start("hvector");
140 
141  write();
142  read();
143 
144  gBenchmark->Show("hvector");
145 }