hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
TrackHistos.cxx
Go to the documentation of this file.
1#include "TrackHistos.h"
2#include "TLorentzVector.h"
3#include "TVector3.h"
4#include <iostream>
5
7
9
10 std::vector<std::string> trkTypes;
11 trkTypes.push_back("topEle");
12 trkTypes.push_back("botEle");
13 trkTypes.push_back("topPos");
14 trkTypes.push_back("botPos");
15 std::string h_name = "";
16 for (auto trkType : trkTypes)
17 {
18 for (auto hist : _h_configs.items()) {
19
20 h_name = m_name+"_"+trkType+"_"+hist.key();
21
22 //Get the extension of the name to decide the histogram to create
23 //i.e. _h = TH1D, _hh = TH2D, _ge = TGraphErrors, _p = TProfile ...
24
25 std::size_t found = (hist.key()).find_last_of("_");
26 std::string extension = hist.key().substr(found+1);
27
28 if (extension == "h") {
29 histos1d[h_name] = plot1D(h_name, trkType+" "+std::string(hist.value().at("xtitle")),
30 hist.value().at("bins"),
31 hist.value().at("minX"),
32 hist.value().at("maxX"));
33
34 std::string ytitle = hist.value().at("ytitle");
35
36 histos1d[h_name]->GetYaxis()->SetTitle(ytitle.c_str());
37
38 if (hist.value().contains("labels")) {
39 std::vector<std::string> labels = hist.value().at("labels").get<std::vector<std::string> >();
40
41 if (labels.size() < hist.value().at("bins")) {
42 std::cout<<"Cannot apply labels to histogram:"<<h_name<<std::endl;
43 }
44 else {
45 for (int i = 1; i<=hist.value().at("bins");++i)
46 histos1d[h_name]->GetXaxis()->SetBinLabel(i,labels[i-1].c_str());
47 }//bins
48 }//labels
49 }//1D histo
50
51 else if (extension == "hh") {
52 histos2d[h_name] = plot2D(h_name,
53 trkType+" "+std::string(hist.value().at("xtitle")),hist.value().at("binsX"),hist.value().at("minX"),hist.value().at("maxX"),
54 hist.value().at("ytitle"),hist.value().at("binsY"),hist.value().at("minY"),hist.value().at("maxY"));
55 }
56 else
57 std::cout<<"Error in histo definition "<<h_name<<std::endl;
58 }//loop on config
59 }//loop on trkTypes
60}
61
63
64 //TODO improve naming
65 std::string h_name = "";
66
67 //TODO improve binning
68 if (doTrkCompPlots) {
69 /*
70 for (unsigned int itp = 0; itp<tPs.size(); ++itp){
71
72 if (debug_)
73 std::cout<<"Bulding:: TH2::" + m_name+"_"+tPs[itp]+"_vs_"+tPs[itp] << std::endl;
74
75 histos2d[m_name+"_"+tPs[itp]+"_vs_"+tPs[itp] ] = plot2D(m_name+tPs[itp]+"_vs_"+tPs[itp],
76 tPs[itp],axes[tPs[itp]][0],axes[tPs[itp]][1],axes[tPs[itp]][2],
77 tPs[itp],axes[tPs[itp]][0],axes[tPs[itp]][1],axes[tPs[itp]][2]);
78
79 }//loop on vars
80 */
81 }//do track comparison
82 /*
83 histos2d[m_name+"_vtxY_vs_vtxX"] = plot2D(m_name+"_vtxY_vs_vtxX",
84 "vtxX",axes["vtx_X"][0],axes["vtx_X"][1],axes["vtx_X"][2],
85 "vtxY",axes["vtx_Y"][0],axes["vtx_Y"][1],axes["vtx_Y"][2]);
86 */
87
88
89}//define 2dhistos
90
91
93 Particle* ele,
94 Particle* pos,
95 Track* ele_trk,
96 Track* pos_trk,
97 float weight) {
98
99 Fill1DVertex(vtx,weight);
100
101 CalCluster eleClus = ele->getCluster();
102 CalCluster posClus = pos->getCluster();
103
104 //TODO remove hardcode!
105 if (ele_trk)
106 Fill1DTrack(ele_trk,weight,"ele_");
107 if (pos_trk)
108 Fill1DTrack(pos_trk,weight,"pos_");
109
110
111 TLorentzVector p_ele;
112 //p_ele.SetPxPyPzE(ele->getMomentum()[0], ele->getMomentum()[1],ele->getMomentum()[2],ele->getEnergy());
113 p_ele.SetPxPyPzE(ele_trk->getMomentum()[0],ele_trk->getMomentum()[1],ele_trk->getMomentum()[2],ele->getEnergy());
114
115 TLorentzVector p_pos;
116 //p_pos.SetPxPyPzE(pos->getMomentum()[0], pos->getMomentum()[1],pos->getMomentum()[2],pos->getEnergy());
117 p_pos.SetPxPyPzE(pos_trk->getMomentum()[0],pos_trk->getMomentum()[1],pos_trk->getMomentum()[2],pos->getEnergy());
118
119 //Fill ele and pos information
120 //Fill1DHisto("ele_p_h",p_ele.P(),weight);
121 //Fill1DHisto("pos_p_h",p_pos.P(),weight);
122 Fill1DHisto("ele_clusE_h",eleClus.getEnergy(),weight);
123 Fill1DHisto("pos_clusE_h",posClus.getEnergy(),weight);
124 Fill1DHisto("ele_EoP_h",eleClus.getEnergy()/p_ele.P(),weight);
125 Fill1DHisto("ele_trkClusX_h",ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0],weight);
126 Fill1DHisto("ele_trkClusY_h",ele_trk->getPositionAtEcal()[1] - eleClus.getPosition()[1],weight);
127 Fill1DHisto("pos_trkClusX_h",pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0],weight);
128 Fill1DHisto("pos_trkClusY_h",pos_trk->getPositionAtEcal()[1] - posClus.getPosition()[1],weight);
129 Fill1DHisto("pos_EoP_h",posClus.getEnergy()/p_pos.P(),weight);
130 Fill2DHisto("EoP_hh", eleClus.getEnergy()/p_ele.P(), posClus.getEnergy()/p_pos.P(),weight);
131 Fill2DHisto("ele_trkClusX_p_hh", ele_trk->getP(),ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
132 Fill2DHisto("pos_trkClusX_p_hh", pos_trk->getP(),pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
133 Fill2DHisto("ele_trkClusX_phi_hh", ele_trk->getPhi(),ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
134 Fill2DHisto("pos_trkClusX_phi_hh", pos_trk->getPhi(),pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
135 Fill2DHisto("ele_trkClusX_z0_hh", ele_trk->getZ0(),ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
136 Fill2DHisto("pos_trkClusX_z0_hh", pos_trk->getZ0(),pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
137 Fill2DHisto("ele_trkClusX_d0_hh", ele_trk->getD0(),ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
138 Fill2DHisto("pos_trkClusX_d0_hh", pos_trk->getD0(),pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
139 Fill2DHisto("ele_trkClusX_tanL_hh", ele_trk->getTanLambda(),ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
140 Fill2DHisto("pos_trkClusX_tanL_hh", pos_trk->getTanLambda(), pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
141 Fill2DHisto("ele_trkClusX_clusX_hh", eleClus.getPosition()[0],ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
142 Fill2DHisto("pos_trkClusX_clusX_hh", posClus.getPosition()[0], pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
143 for(int i = 0; i < ele_trk->getHitLayers().size(); i++)
144 {
145 Fill2DHisto("ele_trkClusX_hitLay_hh", (float) ele_trk->getHitLayers().at(i), ele_trk->getPositionAtEcal()[0] - eleClus.getPosition()[0], weight);
146 }
147 for(int i = 0; i < pos_trk->getHitLayers().size(); i++)
148 {
149 Fill2DHisto("pos_trkClusX_hitLay_hh", (float) pos_trk->getHitLayers().at(i), pos_trk->getPositionAtEcal()[0] - posClus.getPosition()[0], weight);
150 }
151
152
153 //Compute some extra variables
154
155 //TODO::Rotate them
156 p_ele.RotateY(-0.0305);
157 p_pos.RotateY(-0.0305);
158
159 //Massless electrons. TODO fix initialization
160 TLorentzVector p_beam(0.,0.,2.3,2.3);
161 TLorentzVector p_v0 = p_ele+p_pos;
162 TLorentzVector p_miss = p_beam - p_v0;
163
164 double thetax_v0_val = TMath::ATan2(p_v0.X(),p_v0.Z());
165 double thetax_pos_val = TMath::ATan2(p_pos.X(),p_pos.Z());
166
167 double thetay_miss_val = TMath::ATan2(p_miss.Y(),p_miss.Z());
168 double thetay_pos_val = TMath::ATan2(p_pos.Y(),p_pos.Z());
169
170 double pt_ele = p_ele.Pt();
171 double pt_pos = p_pos.Pt();
172
173 double pt_asym_val = (pt_ele-pt_pos) / (pt_ele+pt_pos);
174
175 double thetay_diff_val;
176
177 if (thetay_pos_val>0) {
178 thetay_diff_val = thetay_miss_val - thetay_pos_val;
179 }
180 else {
181 thetay_diff_val = thetay_pos_val - thetay_miss_val;
182 }
183
184 //Fill event information
185
186 //Esum
187 Fill1DHisto("Pmiss_h", p_miss.P(),weight);
188 Fill1DHisto("Esum_h",ele->getEnergy() + pos->getEnergy(),weight);
189 Fill1DHisto("EsumClus_h",eleClus.getEnergy() + posClus.getEnergy(),weight);
190 Fill2DHisto("EClus_hh", eleClus.getEnergy() , posClus.getEnergy(),weight);
191 Fill2DHisto("InvM_eleP_hh", p_ele.P(), vtx->getInvMass(),weight);
192 Fill2DHisto("InvM_posP_hh", p_pos.P(), vtx->getInvMass(), weight);
193 Fill1DHisto("Psum_h",p_ele.P() + p_pos.P(),weight);
194 Fill1DHisto("PtAsym_h",pt_asym_val,weight);
195 Fill1DHisto("thetax_v0_h",thetax_v0_val,weight);
196 Fill1DHisto("thetax_pos_h",thetax_pos_val,weight);
197 Fill1DHisto("thetay_pos_h",thetay_pos_val,weight);
198 Fill1DHisto("thetay_miss_h",thetay_miss_val,weight);
199 Fill1DHisto("thetay_diff_h",thetay_diff_val,weight);
200}
201
202
203void TrackHistos::Fill2DTrack(Track* track, float weight, const std::string& trkname) {
204
205
206 if (track) {
207
208 double d0 = track->getD0();
209 double z0 = track->getZ0();
210 //Fill2DHisto(trkname+"tanlambda_vs_phi0_hh",track->getPhi(),track->getTanLambda(), weight);
211 Fill2DHisto(trkname+"d0_vs_p_hh",track->getP(),d0,weight);
212 Fill2DHisto(trkname+"d0_vs_phi0_hh",track->getPhi(),d0,weight);
213 Fill2DHisto(trkname+"d0_vs_tanlambda_hh",track->getTanLambda(),d0,weight);
214
215 Fill2DHisto(trkname+"z0_vs_p_hh",track->getP(),z0,weight);
216 Fill2DHisto(trkname+"phi0_vs_p_hh",track->getP(),track->getPhi(),weight);
217 Fill2DHisto(trkname+"z0_vs_phi0_hh",track->getPhi(),z0,weight);
218 Fill2DHisto(trkname+"z0_vs_tanlambda_hh",track->getTanLambda(),z0,weight);
219
220 Fill2DHisto(trkname+"TanLambda_vs_Phi_hh" ,track->getPhi() , track->getTanLambda() ,weight);
221 Fill2DHisto(trkname+"p_vs_Phi_hh" ,track->getPhi() , track->getP() ,weight);
222 Fill2DHisto(trkname+"p_vs_TanLambda_hh" ,track->getTanLambda() , track->getP() ,weight);
223
224
225 }
226}
227
228void TrackHistos::Fill1DTrack(Track* track, float weight, const std::string& trkname) {
229
230 double charge = (double) track->getCharge();
231
232 //2D hits
233 int n_hits_2d = track->getTrackerHitCount();
234 if (!track->isKalmanTrack())
235 n_hits_2d*=2;
236
237 Fill1DHisto(trkname+"d0_h" ,track->getD0() ,weight);
238 Fill1DHisto(trkname+"Phi_h" ,track->getPhi() ,weight);
239 Fill1DHisto(trkname+"Omega_h" ,track->getOmega() ,weight);
240 Fill1DHisto(trkname+"pT_h" ,-1*charge*track->getPt(),weight);
241 Fill1DHisto(trkname+"p_h" ,track->getP() ,weight);
242 Fill1DHisto(trkname+"invpT_h" ,-1*charge/track->getPt(),weight);
243 Fill1DHisto(trkname+"TanLambda_h",track->getTanLambda() ,weight);
244 Fill1DHisto(trkname+"Z0_h" ,track->getZ0() ,weight);
245 Fill1DHisto(trkname+"Z0oTanLambda_h",track->getZ0()/track->getTanLambda() ,weight);
246 Fill1DHisto(trkname+"time_h" ,track->getTrackTime() ,weight);
247 Fill1DHisto(trkname+"chi2_h" ,track->getChi2() ,weight);
248 Fill1DHisto(trkname+"chi2ndf_h" ,track->getChi2Ndf() ,weight);
249 Fill1DHisto(trkname+"nShared_h" ,track->getNShared() ,weight);
250 Fill1DHisto(trkname+"nHits_2d_h" ,n_hits_2d ,weight);
251 Fill1DHisto(trkname+"track_xpos_h",track->getPosition().at(0) ,weight);
252 Fill1DHisto(trkname+"track_ypos_h",track->getPosition().at(1) ,weight);
253 Fill1DHisto(trkname+"track_zpos_h",track->getPosition().at(2) ,weight);
254 Fill1DHisto(trkname+"xpos_at_ecal_h",track->getPositionAtEcal().at(0) ,weight);
255 Fill1DHisto(trkname+"ypos_at_ecal_h",track->getPositionAtEcal().at(1) ,weight);
256 Fill1DHisto(trkname+"zpos_at_ecal_h",track->getPositionAtEcal().at(2) ,weight);
257
258 //Top vs Bot
259 if(track->getTanLambda() > 0.0)
260 Fill1DHisto(trkname+"top_track_z0_h", track->getZ0(), weight);
261 else
262 Fill1DHisto(trkname+"bot_track_z0_h", track->getZ0(), weight);
263
264 //Track param errors
265 Fill1DHisto(trkname+"d0_err_h" ,track->getD0Err() ,weight);
266 Fill1DHisto(trkname+"Phi_err_h" ,track->getPhiErr() ,weight);
267 Fill1DHisto(trkname+"Omega_err_h" ,track->getOmegaErr() ,weight);
268 Fill1DHisto(trkname+"TanLambda_err_h",track->getTanLambdaErr() ,weight);
269 Fill1DHisto(trkname+"Z0_err_h" ,track->getZ0Err() ,weight);
270
271 for (int ihit=0; ihit<track->getHitLayers().size();++ihit)
272 {
273 int hit2d = track->getHitLayers().at(ihit);
274 Fill1DHisto(trkname+"hit_lay_h",(float) hit2d ,weight);
275 }
276
277 //All Tracks
278 Fill1DHisto(trkname+"sharingHits_h",0,weight);
279 if (track->getNShared() == 0)
280 Fill1DHisto(trkname+"sharingHits_h",1.,weight);
281 else {
282 //track has shared hits
283 if (track->getSharedLy0())
284 Fill1DHisto(trkname+"sharingHits_h",2.,weight);
285 if (track->getSharedLy1())
286 Fill1DHisto(trkname+"sharingHits_h",3.,weight);
287 if (track->getSharedLy0() && track->getSharedLy1())
288 Fill1DHisto(trkname+"sharingHits_h",4.,weight);
289 if (!track->getSharedLy0() && !track->getSharedLy1())
290 Fill1DHisto(trkname+"sharingHits_h",5.,weight);
291 }
292
293 if (track -> is345Seed())
294 Fill1DHisto(trkname+"strategy_h",0,weight);
295 if (track-> is456Seed())
296 Fill1DHisto(trkname+"strategy_h",1,weight);
297 if (track-> is123SeedC4())
298 Fill1DHisto(trkname+"strategy_h",2,weight);
299 if (track->is123SeedC5())
300 Fill1DHisto(trkname+"strategy_h",3,weight);
301 if (track->isMatchedTrack())
302 Fill1DHisto(trkname+"strategy_h",4,weight);
303 if (track->isGBLTrack())
304 Fill1DHisto(trkname+"strategy_h",5,weight);
305
306
307 Fill1DHisto(trkname+"type_h",track->getType(),weight);
308}
309
310void TrackHistos::Fill1DVertex(Vertex* vtx, float weight) {
311
312 Fill1DHisto("vtx_chi2_h" ,vtx->getChi2(),weight);
313 Fill2DHisto("vtx_XY_hh",vtx->getX(),vtx->getY(),weight);
314 Fill1DHisto("vtx_X_h" ,vtx->getX(),weight);
315 Fill1DHisto("vtx_Y_h" ,vtx->getY(),weight);
316 Fill1DHisto("vtx_Z_h" ,vtx->getZ(),weight);
317
318 TVector3 vtxPosSvt;
319 vtxPosSvt.SetX(vtx->getX());
320 vtxPosSvt.SetY(vtx->getY());
321 vtxPosSvt.SetZ(vtx->getZ());
322
323 vtxPosSvt.RotateY(-0.0305);
324
325 Fill2DHisto("vtx_XY_svt_hh",vtxPosSvt.X(),vtxPosSvt.Y(),weight);
326 Fill1DHisto("vtx_X_svt_h",vtxPosSvt.X(),weight);
327 Fill1DHisto("vtx_Y_svt_h",vtxPosSvt.Y(),weight);
328 Fill1DHisto("vtx_Z_svt_h",vtxPosSvt.Z(),weight);
329
330
331 // 0 xx 1 xy 2 xz 3 yy 4 yz 5 zz
332 Fill1DHisto("vtx_sigma_X_h",sqrt(vtx->getCovariance()[0]),weight);
333 Fill1DHisto("vtx_sigma_Y_h",sqrt(vtx->getCovariance()[3]),weight);
334 Fill1DHisto("vtx_sigma_Z_h",sqrt(vtx->getCovariance()[5]),weight);
335 Fill1DHisto("vtx_InvM_h" ,vtx->getInvMass(),weight);
336 Fill1DHisto("vtx_InvMErr_Z_h",vtx->getInvMassErr(),weight);
337 Fill1DHisto("vtx_px_h",vtx->getP().X());
338 Fill1DHisto("vtx_py_h",vtx->getP().Y());
339 Fill1DHisto("vtx_pz_h",vtx->getP().Z());
340 Fill1DHisto("vtx_p_h" ,vtx->getP().Mag());
341}
342
343void TrackHistos::Fill1DHistograms(Track *track, Vertex* vtx, float weight ) {
344
345 if (track) {
346 Fill1DTrack(track);
347 }
348
349 //Vertices
350
351 if (vtx) {
352 Fill1DVertex(vtx);
353 }
354}
355
356
357void TrackHistos::Fill1DTrackTruth(Track *track, Track* truth_track, float weight, const std::string& trkname) {
358
359 if (!track || !truth_track)
360 return;
361
362 //Momentum
363 std::vector<double> trk_mom = track->getMomentum();
364 std::vector<double> trk_truth_mom = truth_track->getMomentum();
365
366 double d0 = track->getD0();
367 double d0err = track->getD0Err();
368 double d0_truth = truth_track->getD0();
369 double phi = track->getPhi();
370 double phi_truth = truth_track->getPhi();
371 double phierr = track->getPhiErr();
372 double omega = track->getOmega();
373 double omega_truth = truth_track->getOmega();
374 double omegaerr = track->getOmegaErr();
375 double tanLambda = track->getTanLambda();
376 double tanLambda_truth = truth_track->getTanLambda();
377 double tanLambdaerr = track->getTanLambdaErr();
378 double z0 = track->getZ0();
379 double z0_truth = truth_track->getZ0();
380 double z0err = track->getZ0Err();
381 double p = track->getP();
382 //Charge different wrt Robert's plots.
383 double invPt = -1.*(double) track->getCharge()/track->getPt();
384 double p_truth = truth_track->getP();
385 double invPt_truth = -1*(double) track->getCharge()/truth_track->getPt();
386
387 double diff_percent_invpT = ((invPt - invPt_truth) / invPt_truth) * 100.;
388
389 // truth residuals
390 Fill1DHisto(trkname+"d0_truth_res_h", d0 - d0_truth , weight);
391 Fill1DHisto(trkname+"Phi_truth_res_h", phi - phi_truth , weight);
392 Fill1DHisto(trkname+"Omega_truth_res_h", omega - omega_truth , weight);
393 Fill1DHisto(trkname+"TanLambda_truth_res_h",tanLambda - tanLambda_truth , weight);
394 Fill1DHisto(trkname+"Z0_truth_res_h", z0 - z0_truth , weight);
395 Fill1DHisto(trkname+"p_truth_res_h", p - p_truth , weight);
396 Fill1DHisto(trkname+"invpT_truth_res_h", invPt - invPt_truth , weight);
397 Fill1DHisto(trkname+"invpT_truth_res_percent_h", diff_percent_invpT , weight);
398 Fill1DHisto(trkname+"px_truth_res_h", trk_mom[0] - trk_truth_mom[0] , weight);
399 Fill1DHisto(trkname+"py_truth_res_h", trk_mom[1] - trk_truth_mom[1] , weight);
400 Fill1DHisto(trkname+"pz_truth_res_h", trk_mom[2] - trk_truth_mom[2] , weight);
401
402 // truth pulls
403 Fill1DHisto(trkname+"d0_truth_pull_h", (d0 - d0_truth) / d0err, weight);
404 Fill1DHisto(trkname+"Phi_truth_pull_h", (phi - phi_truth) / phierr , weight);
405 Fill1DHisto(trkname+"Omega_truth_pull_h", (omega - omega_truth) / omegaerr, weight);
406 Fill1DHisto(trkname+"TanLambda_truth_pull_h",(tanLambda - tanLambda_truth) / tanLambdaerr, weight);
407 Fill1DHisto(trkname+"Z0_truth_pull_h", (z0 - z0_truth) / z0err, weight);
408
409}
410
411
412
413void TrackHistos::Fill2DHistograms(Vertex* vtx, float weight) {
414
415 if (vtx) {
416
417 //TODO Improve this.
418 TVector3 vtxPosSvt;
419 vtxPosSvt.SetX(vtx->getX());
420 vtxPosSvt.SetY(vtx->getY());
421 vtxPosSvt.SetZ(vtx->getZ());
422
423 vtxPosSvt.RotateY(-0.0305);
424
425
426 double vtxP = vtx->getP().Mag();
427
428 Fill2DHisto("vtx_InvM_vtx_z_hh",vtx->getInvMass(),vtx->getZ(),weight);
429 Fill2DHisto("vtx_InvM_vtx_svt_z_hh",vtx->getInvMass(),vtxPosSvt.Z(),weight);
430 Fill2DHisto("vtx_p_svt_z_hh",vtxP,vtxPosSvt.Z(),weight);
431 Fill2DHisto("vtx_p_svt_x_hh",vtxP,vtxPosSvt.X(),weight);
432 Fill2DHisto("vtx_p_svt_y_hh",vtxP,vtxPosSvt.Y(),weight);
433
434 Fill2DHisto("vtx_svt_y_svt_z_hh",vtxPosSvt.Y(),vtxPosSvt.Z(),weight);
435
436 Fill2DHisto("vtx_p_sigmaZ_hh",vtxP,vtx->getCovariance()[5],weight);
437 Fill2DHisto("vtx_p_sigmaX_hh",vtxP,vtx->getCovariance()[3],weight);
438 Fill2DHisto("vtx_p_sigmaY_hh",vtxP,vtx->getCovariance()[0],weight);
439 }
440}
441
442void TrackHistos::FillTrackComparisonHistograms(Track* track_x, Track* track_y, float weight) {
443
444 if (doTrkCompPlots) {
445 /*
446 histos2d[m_name+"_d0_vs_d0" ]->Fill(track_x->getD0(),track_y->getD0(),weight);
447 histos2d[m_name+"_Phi_vs_Phi" ]->Fill(track_x->getPhi(),track_y->getPhi(),weight);
448 histos2d[m_name+"_Omega_vs_Omega" ]->Fill(track_x->getOmega(),track_y->getOmega(),weight);
449 histos2d[m_name+"_TanLambda_vs_TanLambda"]->Fill(track_x->getTanLambda(),track_y->getTanLambda(),weight);
450 histos2d[m_name+"_Z0_vs_Z0" ]->Fill(track_x->getZ0(),track_y->getZ0(),weight);
451 histos2d[m_name+"_time_vs_time" ]->Fill(track_x->getTrackTime(),track_y->getTrackTime(),weight);
452 histos2d[m_name+"_chi2_vs_chi2" ]->Fill(track_x->getChi2Ndf(),
453 track_y->getChi2Ndf(),
454 weight);
455 */
456 }
457
458}
459
460
461//Residual Plots ============ They should probably go somewhere else ====================
462
463
464void TrackHistos::FillResidualHistograms(Track* track, int ly, double res, double sigma) {
465
466 double trk_mom = track->getP();
467 std::string lyr = std::to_string(ly);
468
469 TrackerHit* hit = nullptr;
470 //Get the hits on track
471 for (int ihit = 0; ihit<track->getSvtHits().GetEntries();++ihit) {
472 TrackerHit* tmphit = (TrackerHit*) track->getSvtHits().At(ihit);
473 if (tmphit->getLayer() == ly) {
474 hit = tmphit;
475 break;
476 }
477 }
478
479 if (!hit) {
480 std::cout<<"Hit-on-track residual infos not found on hit on track list for ly="<<ly<<std::endl;
481 }
482
483 double hit_y = -9999.;
484 if (hit) {
485 hit_y = hit->getPosition()[1];
486 }
487
488 //General Plots
489 Fill1DHisto("u_res_ly_"+lyr+"_h",res);
490 Fill2DHisto("u_res_ly_"+lyr+"_vsp_hh",trk_mom,res);
491 Fill2DHisto("u_res_ly_"+lyr+"_vsy_hh",hit_y,res);
492
493 //Top = 0 bottom=1 - Per Volume
494 std::string vol = track->getTanLambda()>0 ? "top" : "bot";
495 Fill1DHisto("u_res_ly_"+lyr+"_"+vol+"_h",res);
496 Fill2DHisto("u_res_ly_"+lyr+"_"+vol+"_vsp_hh",trk_mom,res);
497
498}
double getEnergy() const
Definition CalCluster.h:73
std::vector< double > getPosition() const
Definition CalCluster.h:63
void Fill2DHisto(const std::string &histoName, float valuex, float valuey, float weight=1.)
description
std::map< std::string, TH2F * > histos2d
description
TH2F * plot2D(std::string name, std::string xtitle, int nbinsX, float xmin, float xmax, std::string ytitle, int nbinsY, float ymin, float ymax)
description
std::map< std::string, TH1F * > histos1d
description
void Fill1DHisto(const std::string &histoName, float value, float weight=1.)
description
TH1F * plot1D(const std::string &name, const std::string &xtitle, int nbinsX, float xmin, float xmax)
description
json _h_configs
description
std::string m_name
description
CalCluster getCluster() const
Definition Particle.h:62
double getEnergy() const
Definition Particle.h:143
void FillTrackComparisonHistograms(Track *track_x, Track *track_y, float weight=1.)
description
void Fill2DTrack(Track *track, float weight=1., const std::string &trkname="")
Fill 2D track.
void FillResidualHistograms(Track *track, int ly, double res, double sigma)
Fill residual histograms.
void Fill1DHistograms(Track *track=nullptr, Vertex *vtx=nullptr, float weight=1.)
Fill 1D histograms.
void BuildAxes()
description
bool doTrkCompPlots
void Fill1DVertex(Vertex *vtx, float weight=1.)
description
virtual void Define2DHistos()
description
void Fill1DTrackTruth(Track *track, Track *truth_track, float weight=1., const std::string &="")
Truth comparison.
void Fill2DHistograms(Vertex *vtx=nullptr, float weight=1.)
Fill 2D histograms.
void Fill1DTrack(Track *track, float weight=1., const std::string &trkname="")
Fill 1D track.
void DefineTrkHitHistos()
description
Definition Track.h:32
std::vector< int > getHitLayers()
Definition Track.h:99
double getPhi() const
Definition Track.h:88
double getOmega() const
Definition Track.h:90
std::vector< double > getPositionAtEcal()
Definition Track.cxx:53
double getD0() const
Definition Track.h:87
double getTanLambda() const
Definition Track.h:91
std::vector< double > getMomentum()
Definition Track.h:285
double getZ0() const
Definition Track.h:92
double getPt() const
Definition Track.h:293
double getP() const
Definition Track.h:291
int getLayer() const
Definition TrackerHit.h:109
std::vector< double > getPosition() const
Definition TrackerHit.h:47