-
Notifications
You must be signed in to change notification settings - Fork 1.3k
/
Copy pathmathcoreVectorIO.C
168 lines (131 loc) · 3.78 KB
/
mathcoreVectorIO.C
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
/// \file
/// \ingroup tutorial_math
/// \notebook -nodraw
/// Example of I/O of a GenVector Lorentz Vectors in a Tree and comparison with legacy TLorentzVector.
///
/// A ROOT tree is written and read in both using either a XYZTVector or a TLorentzVector.
///
/// To execute the macro type in:
///
/// ~~~{.cpp}
/// root[0] .x mathcoreVectorIO.C
/// ~~~
///
/// \macro_output
/// \macro_code
///
/// \author Lorenzo Moneta
#include "TRandom2.h"
#include "TStopwatch.h"
#include "TSystem.h"
#include "TFile.h"
#include "TTree.h"
#include "TH1D.h"
#include "TCanvas.h"
#include <iostream>
#include "TLorentzVector.h"
#include "Math/Vector4D.h"
using namespace ROOT::Math;
void write(int n) {
TRandom2 R;
TStopwatch timer;
R.SetSeed(1);
timer.Start();
double s = 0;
for (int i = 0; i < n; ++i) {
s += R.Gaus(0,10);
s += R.Gaus(0,10);
s += R.Gaus(0,10);
s += R.Gaus(100,10);
}
timer.Stop();
std::cout << s/double(n) << std::endl;
std::cout << " Time for Random gen " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
TFile f1("mathcoreVectorIO_1.root","RECREATE");
// create tree
TTree t1("t1","Tree with new LorentzVector");
XYZTVector *v1 = new XYZTVector();
t1.Branch("LV branch","ROOT::Math::XYZTVector",&v1);
R.SetSeed(1);
timer.Start();
for (int i = 0; i < n; ++i) {
double Px = R.Gaus(0,10);
double Py = R.Gaus(0,10);
double Pz = R.Gaus(0,10);
double E = R.Gaus(100,10);
v1->SetCoordinates(Px,Py,Pz,E);
t1.Fill();
}
f1.Write();
timer.Stop();
std::cout << " Time for new Vector " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
t1.Print();
// create tree with old LV
TFile f2("mathcoreVectorIO_2.root","RECREATE");
TTree t2("t2","Tree with TLorentzVector");
TLorentzVector * v2 = new TLorentzVector();
TLorentzVector::Class()->IgnoreTObjectStreamer();
TVector3::Class()->IgnoreTObjectStreamer();
t2.Branch("TLV branch","TLorentzVector",&v2,16000,2);
R.SetSeed(1);
timer.Start();
for (int i = 0; i < n; ++i) {
double Px = R.Gaus(0,10);
double Py = R.Gaus(0,10);
double Pz = R.Gaus(0,10);
double E = R.Gaus(100,10);
v2->SetPxPyPzE(Px,Py,Pz,E);
t2.Fill();
}
f2.Write();
timer.Stop();
std::cout << " Time for old Vector " << timer.RealTime() << " " << timer.CpuTime() << endl;
t2.Print();
}
void read() {
TRandom R;
TStopwatch timer;
TFile f1("mathcoreVectorIO_1.root");
// create tree
TTree *t1 = (TTree*)f1.Get("t1");
XYZTVector *v1 = 0;
t1->SetBranchAddress("LV branch",&v1);
timer.Start();
int n = (int) t1->GetEntries();
std::cout << " Tree Entries " << n << std::endl;
double etot=0;
for (int i = 0; i < n; ++i) {
t1->GetEntry(i);
etot += v1->Px();
etot += v1->Py();
etot += v1->Pz();
etot += v1->E();
}
timer.Stop();
std::cout << " Time for new Vector " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
std::cout << " TOT average : n = " << n << "\t " << etot/double(n) << endl;
// create tree with old LV
TFile f2("mathcoreVectorIO_2.root");
TTree *t2 = (TTree*)f2.Get("t2");
TLorentzVector * v2 = 0;
t2->SetBranchAddress("TLV branch",&v2);
timer.Start();
n = (int) t2->GetEntries();
std::cout << " Tree Entries " << n << std::endl;
etot = 0;
for (int i = 0; i < n; ++i) {
t2->GetEntry(i);
etot += v2->Px();
etot += v2->Py();
etot += v2->Pz();
etot += v2->E();
}
timer.Stop();
std::cout << " Time for old Vector " << timer.RealTime() << " " << timer.CpuTime() << endl;
std::cout << " TOT average:\t" << etot/double(n) << endl;
}
void mathcoreVectorIO() {
int nEvents = 100000;
write(nEvents);
read();
}