-
Notifications
You must be signed in to change notification settings - Fork 1.3k
/
Copy pathmathcoreVectorCollection.C
170 lines (143 loc) · 4.34 KB
/
mathcoreVectorCollection.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
169
/// \file
/// \ingroup tutorial_math
/// \notebook -js
/// Example showing how to write and read a std vector of ROOT::Math LorentzVector in a ROOT tree.
///
/// In the write() function a variable number of track Vectors is generated
/// according to a Poisson distribution with random momentum uniformly distributed
/// in phi and eta.
/// In the read() the vectors are read back and the content analysed and
/// some information such as number of tracks per event or the track pt
/// distributions are displayed in a canvas.
///
/// To execute the macro type in:
///
/// ~~~{.cpp}
/// root[0]: .x mathcoreVectorCollection.C
/// ~~~
///
/// \macro_image
/// \macro_output
/// \macro_code
///
/// \author Andras Zsenei
#include "TRandom.h"
#include "TStopwatch.h"
#include "TSystem.h"
#include "TFile.h"
#include "TTree.h"
#include "TH1D.h"
#include "TCanvas.h"
#include "TMath.h"
#include <iostream>
// CLING does not understand some files included by LorentzVector
#include "Math/Vector3D.h"
#include "Math/Vector4D.h"
using namespace ROOT::Math;
double write(int n) {
TRandom R;
TStopwatch timer;
TFile f1("mathcoreLV.root","RECREATE");
// create tree
TTree t1("t1","Tree with new LorentzVector");
std::vector<ROOT::Math::XYZTVector> tracks;
std::vector<ROOT::Math::XYZTVector> * pTracks = &tracks;
t1.Branch("tracks","std::vector<ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > >",&pTracks);
double M = 0.13957; // set pi+ mass
timer.Start();
double sum = 0;
for (int i = 0; i < n; ++i) {
int nPart = R.Poisson(5);
pTracks->clear();
pTracks->reserve(nPart);
for (int j = 0; j < nPart; ++j) {
double px = R.Gaus(0,10);
double py = R.Gaus(0,10);
double pt = sqrt(px*px +py*py);
double eta = R.Uniform(-3,3);
double phi = R.Uniform(0.0 , 2*TMath::Pi() );
RhoEtaPhiVector vcyl( pt, eta, phi);
// set energy
double E = sqrt( vcyl.R()*vcyl.R() + M*M);
XYZTVector q( vcyl.X(), vcyl.Y(), vcyl.Z(), E);
// fill track vector
pTracks->push_back(q);
// evaluate sum of components to check
sum += q.x()+q.y()+q.z()+q.t();
}
t1.Fill();
}
f1.Write();
timer.Stop();
std::cout << " Time for new Vector " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
t1.Print();
return sum;
}
double read() {
TRandom R;
TStopwatch timer;
TH1D * h1 = new TH1D("h1","total event energy ",100,0,1000.);
TH1D * h2 = new TH1D("h2","Number of track per event",21,-0.5,20.5);
TH1D * h3 = new TH1D("h3","Track Energy",100,0,200);
TH1D * h4 = new TH1D("h4","Track Pt",100,0,100);
TH1D * h5 = new TH1D("h5","Track Eta",100,-5,5);
TH1D * h6 = new TH1D("h6","Track Cos(theta)",100,-1,1);
TFile f1("mathcoreLV.root");
// create tree
TTree *t1 = (TTree*)f1.Get("t1");
std::vector<ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > > * pTracks = nullptr;
t1->SetBranchAddress("tracks",&pTracks);
timer.Start();
int n = (int) t1->GetEntries();
std::cout << " Tree Entries " << n << std::endl;
double sum=0;
for (int i = 0; i < n; ++i) {
t1->GetEntry(i);
int ntrk = pTracks->size();
h3->Fill(ntrk);
XYZTVector q;
for (int j = 0; j < ntrk; ++j) {
XYZTVector v = (*pTracks)[j];
q += v;
h3->Fill(v.E());
h4->Fill(v.Pt());
h5->Fill(v.Eta());
h6->Fill(cos(v.Theta()));
sum += v.x() + v.y() + v.z() + v.t();
}
h1->Fill(q.E() );
h2->Fill(ntrk);
}
timer.Stop();
std::cout << " Time for new Vector " << timer.RealTime() << " " << timer.CpuTime() << std::endl;
TCanvas *c1 = new TCanvas("c1","demo of Trees",10,10,600,800);
c1->Divide(2,3);
c1->cd(1);
h1->Draw();
c1->cd(2);
h2->Draw();
c1->cd(3);
h3->Draw();
c1->cd(3);
h3->Draw();
c1->cd(4);
h4->Draw();
c1->cd(5);
h5->Draw();
c1->cd(6);
h6->Draw();
return sum;
}
int mathcoreVectorCollection() {
int nEvents = 10000;
double s1 = write(nEvents);
double s2 = read();
if (fabs(s1-s2) > s1*1.E-15 ) {
std::cout << "ERROR: Found difference in Vector when reading ( " << s1 << " != " << s2 << " diff = " << fabs(s1-s2) << " ) " << std::endl;
return -1;
}
return 0;
}
int main() {
return mathcoreVectorCollection();
}