EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4TrackFastSimEval.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4TrackFastSimEval.cc
1 
8 #include "PHG4TrackFastSimEval.h"
9 
13 #include <trackbase_historic/SvtxVertex.h> // for SvtxVertex
15 
16 #include <g4main/PHG4Hit.h>
18 #include <g4main/PHG4Particle.h>
20 #include <g4main/PHG4VtxPoint.h>
21 
22 #include <pdbcalbase/PdbParameterMap.h>
23 #include <phparameter/PHParameters.h>
24 
26 #include <fun4all/PHTFileServer.h>
27 #include <fun4all/SubsysReco.h> // for SubsysReco
28 
29 #include <phool/getClass.h>
30 #include <phool/phool.h>
31 
32 #include <TH2.h>
33 #include <TSystem.h>
34 #include <TTree.h>
35 #include <TVector3.h>
36 
37 #include <cassert>
38 #include <cmath>
39 #include <iostream>
40 #include <map> // for _Rb_tree_const_ite...
41 #include <utility> // for pair
42 
43 #define LogDebug(exp) std::cout << "DEBUG: " << __FILE__ << ": " << __LINE__ << ": " << exp << "\n"
44 #define LogError(exp) std::cout << "ERROR: " << __FILE__ << ": " << __LINE__ << ": " << exp << "\n"
45 #define LogWarning(exp) std::cout << "WARNING: " << __FILE__ << ": " << __LINE__ << ": " << exp << "\n"
46 
47 using namespace std;
48 
49 const string xyzt[] = {"x", "y", "z", "t"};
50 
51 //----------------------------------------------------------------------------//
52 //-- Constructor:
53 //-- simple initialization
54 //----------------------------------------------------------------------------//
55 PHG4TrackFastSimEval::PHG4TrackFastSimEval(const string &name, const string &filename, const string &trackmapname)
56  : SubsysReco(name)
57  , m_TruthInfoContainer(nullptr)
58  , m_TrackMap(nullptr)
59  , m_VertexMap(nullptr)
60  , m_TracksEvalTree(nullptr)
61  , m_VertexEvalTree(nullptr)
62  , m_H2D_DeltaMomVsTruthMom(nullptr)
63  , m_H2D_DeltaMomVsTruthEta(nullptr)
64  , m_EventCounter(0)
65  , m_OutFileName(filename)
66  , m_TrackMapName(trackmapname)
67 {
69 }
70 
71 //----------------------------------------------------------------------------//
72 //-- Init():
73 //-- Intialize all histograms, trees, and ntuples
74 //----------------------------------------------------------------------------//
76 {
78 }
79 
80 //----------------------------------------------------------------------------//
81 //-- InitRun():
82 //-- add related hit object
83 //----------------------------------------------------------------------------//
85 {
86  if (Verbosity())
87  cout << PHWHERE << " Openning file " << m_OutFileName << endl;
88  PHTFileServer::get().open(m_OutFileName, "RECREATE");
89 
90  // create TTree
91  m_TracksEvalTree = new TTree("tracks", "FastSim Eval => tracks");
92  m_TracksEvalTree->Branch("event", &m_TTree_Event, "event/I");
93  m_TracksEvalTree->Branch("gtrackID", &m_TTree_gTrackID, "gtrackID/I");
94  m_TracksEvalTree->Branch("gflavor", &m_TTree_gFlavor, "gflavor/I");
95  m_TracksEvalTree->Branch("gpx", &m_TTree_gpx, "gpx/F");
96  m_TracksEvalTree->Branch("gpy", &m_TTree_gpy, "gpy/F");
97  m_TracksEvalTree->Branch("gpz", &m_TTree_gpz, "gpz/F");
98  m_TracksEvalTree->Branch("gvx", &m_TTree_gvx, "gvx/F");
99  m_TracksEvalTree->Branch("gvy", &m_TTree_gvy, "gvy/F");
100  m_TracksEvalTree->Branch("gvz", &m_TTree_gvz, "gvz/F");
101  m_TracksEvalTree->Branch("gvt", &m_TTree_gvt, "gvt/F");
102  m_TracksEvalTree->Branch("trackID", &m_TTree_TrackID, "trackID/I");
103  m_TracksEvalTree->Branch("charge", &m_TTree_Charge, "charge/I");
104  m_TracksEvalTree->Branch("nhits", &m_TTree_nHits, "nhits/I");
105  m_TracksEvalTree->Branch("px", &m_TTree_px, "px/F");
106  m_TracksEvalTree->Branch("py", &m_TTree_py, "py/F");
107  m_TracksEvalTree->Branch("pz", &m_TTree_pz, "pz/F");
108  m_TracksEvalTree->Branch("pcax", &m_TTree_pcax, "pcax/F");
109  m_TracksEvalTree->Branch("pcay", &m_TTree_pcay, "pcay/F");
110  m_TracksEvalTree->Branch("pcaz", &m_TTree_pcaz, "pcaz/F");
111  m_TracksEvalTree->Branch("dca2d", &m_TTree_dca2d, "dca2d/F");
112 
113  // next a stat. on hits
114  PHParameters PHG4TrackFastSim_Parameter("PHG4TrackFastSim");
115 
116  PdbParameterMap *nodeparams = findNode::getClass<PdbParameterMap>(topNode,
117  "PHG4TrackFastSim_Parameter");
118  if (not nodeparams)
119  {
120  cout << __PRETTY_FUNCTION__ << " : Warning, missing PHG4TrackFastSim_Parameter node and skip saving hits"
121  << endl;
122  }
123  else
124  {
125  PHG4TrackFastSim_Parameter.FillFrom(nodeparams);
126  if (Verbosity())
127  {
128  cout << __PRETTY_FUNCTION__ << " PHG4TrackFastSim_Parameter : ";
129  PHG4TrackFastSim_Parameter.Print();
130  }
131 
132  auto range = PHG4TrackFastSim_Parameter.get_all_int_params();
133  for (auto iter = range.first; iter != range.second; ++iter)
134  {
135  const string &phg4hit_node_name = iter->first;
136  const int &phg4hit_node_id = iter->second;
137 
138  cout << __PRETTY_FUNCTION__ << " Prepare PHG4Hit node name " << phg4hit_node_name
139  << " with ID = " << phg4hit_node_id << endl;
140 
141  string branch_name = string("nHit_") + phg4hit_node_name;
142  m_TracksEvalTree->Branch(branch_name.c_str(),
143  &m_TTree_HitContainerID_nHits_map[phg4hit_node_id],
144  (branch_name + "/I").c_str());
145  }
146  }
147 
148  m_H2D_DeltaMomVsTruthEta = new TH2D("DeltaMomVsTruthEta",
149  "#frac{#Delta p}{truth p} vs. truth #eta", 54, -4.5, +4.5, 1000, -1,
150  1);
151 
152  m_H2D_DeltaMomVsTruthMom = new TH2D("DeltaMomVsTruthMom",
153  "#frac{#Delta p}{truth p} vs. truth p", 41, -0.5, 40.5, 1000, -1,
154  1);
155 
156  // create TTree - vertex
157  m_VertexEvalTree = new TTree("vertex", "FastSim Eval => vertces");
158  m_VertexEvalTree->Branch("event", &m_TTree_Event, "event/I");
159  m_VertexEvalTree->Branch("gvx", &m_TTree_gvx, "gvx/F");
160  m_VertexEvalTree->Branch("gvy", &m_TTree_gvy, "gvy/F");
161  m_VertexEvalTree->Branch("gvz", &m_TTree_gvz, "gvz/F");
162  m_VertexEvalTree->Branch("gvt", &m_TTree_gvt, "gvt/F");
163  m_VertexEvalTree->Branch("vx", &m_TTree_vx, "vx/F");
164  m_VertexEvalTree->Branch("vy", &m_TTree_vy, "vy/F");
165  m_VertexEvalTree->Branch("vz", &m_TTree_vz, "vz/F");
166  m_VertexEvalTree->Branch("deltavx", &m_TTree_DeltaVx, "deltavx/F");
167  m_VertexEvalTree->Branch("deltavy", &m_TTree_DeltaVy, "deltavy/F");
168  m_VertexEvalTree->Branch("deltavz", &m_TTree_DeltaVz, "deltavz/F");
169  m_VertexEvalTree->Branch("gID", &m_TTree_gTrackID, "gID/I");
170  m_VertexEvalTree->Branch("ID", &m_TTree_TrackID, "ID/I");
171  m_VertexEvalTree->Branch("ntracks", &m_TTree_nTracks, "ntracks/I");
172  m_VertexEvalTree->Branch("n_from_truth", &m_TTree_nFromTruth, "n_from_truth/I");
173 
174  for (map<string, unsigned int>::const_iterator iter = m_ProjectionNameMap.begin(); iter != m_ProjectionNameMap.end(); ++iter)
175  {
176  for (int i = 0; i < 4; i++)
177  {
178  string bname = iter->first + "_proj_" + xyzt[i];
179  string bdef = bname + "/F";
180 
181  // fourth element is the path length
182  if (i == 3)
183  {
184  bdef = iter->first + "_proj_path_length" + "/F";
185  }
186 
187  m_TracksEvalTree->Branch(bname.c_str(), &m_TTree_proj_vec[iter->second][i], bdef.c_str());
188  }
189 
190  for (int i = 0; i < 3; i++)
191  {
192  string bname = iter->first + "_proj_p" + xyzt[i];
193  string bdef = bname + "/F";
194  m_TracksEvalTree->Branch(bname.c_str(), &m_TTree_proj_p_vec[iter->second][i], bdef.c_str());
195  }
196  string nodename = "G4HIT_" + iter->first;
197  PHG4HitContainer *hits = findNode::getClass<PHG4HitContainer>(topNode, nodename);
198  if (hits)
199  {
200  for (int i = 0; i < 4; i++)
201  {
202  string bname = iter->first + "_" + xyzt[i];
203  string bdef = bname + "/F";
204  m_TracksEvalTree->Branch(bname.c_str(), &m_TTree_ref_vec[iter->second][i], bdef.c_str());
205  }
206  for (int i = 0; i < 3; i++)
207  {
208  string bname = iter->first + "_p" + xyzt[i];
209  string bdef = bname + "/F";
210 
211  m_TracksEvalTree->Branch(bname.c_str(), &m_TTree_ref_p_vec[iter->second][i], bdef.c_str());
212  }
213  }
214  if (!hits && Verbosity() > 0)
215  {
216  cout << "InitRun: could not find " << nodename << endl;
217  }
218  }
219 
221 }
222 
223 //----------------------------------------------------------------------------//
224 //-- process_event():
225 //-- Call user instructions for every event.
226 //-- This function contains the analysis structure.
227 //----------------------------------------------------------------------------//
229 {
230  m_EventCounter++;
231  if (Verbosity() >= 2 and m_EventCounter % 1000 == 0)
232  cout << PHWHERE << "Events processed: " << m_EventCounter << endl;
233 
234  //std::cout << "Opening nodes" << std::endl;
235  GetNodes(topNode);
236 
237  //std::cout << "Filling trees" << std::endl;
238  fill_track_tree(topNode);
239  fill_vertex_tree(topNode);
240  //std::cout << "DONE" << std::endl;
241 
243 }
244 
245 //----------------------------------------------------------------------------//
246 //-- End():
247 //-- End method, wrap everything up
248 //----------------------------------------------------------------------------//
250 {
252 
253  m_TracksEvalTree->Write();
254  m_VertexEvalTree->Write();
255 
256  m_H2D_DeltaMomVsTruthEta->Write();
257  m_H2D_DeltaMomVsTruthMom->Write();
258 
259  //PHTFileServer::get().close();
260 
262 }
263 
264 //----------------------------------------------------------------------------//
265 //-- fill_tree():
266 //-- Fill the trees with truth, track fit, and cluster information
267 //----------------------------------------------------------------------------//
269 {
270  // Make sure to reset all the TTree variables before trying to set them.
271 
273  {
274  LogError("m_TruthInfoContainer not found!");
275  return;
276  }
277 
278  if (!m_TrackMap)
279  {
280  LogError("m_TrackMap not found!");
281  return;
282  }
283 
286  for (PHG4TruthInfoContainer::ConstIterator truth_itr = range.first;
287  truth_itr != range.second; ++truth_itr)
288  {
289  reset_variables();
291 
292  PHG4Particle *g4particle = truth_itr->second;
293  if (!g4particle)
294  {
295  LogDebug("");
296  continue;
297  }
298 
299  SvtxTrack_FastSim *track = nullptr;
300 
301  if (Verbosity()) cout << __PRETTY_FUNCTION__ << "TRACKmap size " << m_TrackMap->size() << std::endl;
302  for (SvtxTrackMap::ConstIter track_itr = m_TrackMap->begin();
303  track_itr != m_TrackMap->end();
304  track_itr++)
305  {
306  //std::cout << "TRACK * " << track_itr->first << std::endl;
307  SvtxTrack_FastSim *temp = dynamic_cast<SvtxTrack_FastSim *>(track_itr->second);
308  if (!temp)
309  {
310  if (Verbosity())
311  {
312  std::cout << "PHG4TrackFastSimEval::fill_track_tree - ignore track that is not a SvtxTrack_FastSim:";
313  track_itr->second->identify();
314  }
315  continue;
316  }
317  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " PARTICLE!" << std::endl;
318 
319  if ((temp->get_truth_track_id() - g4particle->get_track_id()) == 0)
320  {
321  track = temp;
322  }
323  }
324 
325  //std::cout << "B2" << std::endl;
326  m_TTree_gTrackID = g4particle->get_track_id();
327  m_TTree_gFlavor = g4particle->get_pid();
328 
329  m_TTree_gpx = g4particle->get_px();
330  m_TTree_gpy = g4particle->get_py();
331  m_TTree_gpz = g4particle->get_pz();
332 
333  m_TTree_gvx = NAN;
334  m_TTree_gvy = NAN;
335  m_TTree_gvz = NAN;
336  m_TTree_gvt = NAN;
338  if (vtx)
339  {
340  m_TTree_gvx = vtx->get_x();
341  m_TTree_gvy = vtx->get_y();
342  m_TTree_gvz = vtx->get_z();
343  m_TTree_gvt = vtx->get_t();
344  }
345 
346  if (track)
347  {
348  //std::cout << "C1" << std::endl;
349  m_TTree_TrackID = track->get_id();
350  m_TTree_Charge = track->get_charge();
351  m_TTree_nHits = track->size_clusters();
352 
353  m_TTree_px = track->get_px();
354  m_TTree_py = track->get_py();
355  m_TTree_pz = track->get_pz();
356  m_TTree_pcax = track->get_x();
357  m_TTree_pcay = track->get_y();
358  m_TTree_pcaz = track->get_z();
359  m_TTree_dca2d = track->get_dca2d();
360 
361  TVector3 truth_mom(m_TTree_gpx, m_TTree_gpy, m_TTree_gpz);
362  TVector3 reco_mom(m_TTree_px, m_TTree_py, m_TTree_pz);
363  //std::cout << "C2" << std::endl;
364 
365  m_H2D_DeltaMomVsTruthMom->Fill(truth_mom.Mag(), (reco_mom.Mag() - truth_mom.Mag()) / truth_mom.Mag());
366  m_H2D_DeltaMomVsTruthEta->Fill(truth_mom.Eta(), (reco_mom.Mag() - truth_mom.Mag()) / truth_mom.Mag());
367  // find projections
368  for (SvtxTrack::ConstStateIter trkstates = track->begin_states();
369  trkstates != track->end_states();
370  ++trkstates)
371  {
372  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " checking " << trkstates->second->get_name() << endl;
373  map<string, unsigned int>::const_iterator iter = m_ProjectionNameMap.find(trkstates->second->get_name());
374  if (iter != m_ProjectionNameMap.end())
375  {
376  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " found " << trkstates->second->get_name() << endl;
377  // setting the projection (xyz and pxpypz)
378  for (int i = 0; i < 3; i++)
379  {
380  m_TTree_proj_vec[iter->second][i] = trkstates->second->get_pos(i);
381  m_TTree_proj_p_vec[iter->second][i] = trkstates->second->get_mom(i);
382  }
383  // fourth element is the path length
384  m_TTree_proj_vec[iter->second][3] = trkstates->first;
385 
386  string nodename = "G4HIT_" + trkstates->second->get_name();
387  PHG4HitContainer *hits = findNode::getClass<PHG4HitContainer>(topNode, nodename);
388  if (!hits)
389  {
390  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " could not find " << nodename << endl;
391  continue;
392  }
393  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " number of hits: " << hits->size() << endl;
394  PHG4HitContainer::ConstRange hit_range = hits->getHits();
395  for (PHG4HitContainer::ConstIterator hit_iter = hit_range.first; hit_iter != hit_range.second; hit_iter++)
396  {
397  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " checking hit id " << hit_iter->second->get_trkid() << " against " << track->get_truth_track_id() << endl;
398  if (hit_iter->second->get_trkid() - track->get_truth_track_id() == 0)
399  {
400  if (Verbosity()) cout << __PRETTY_FUNCTION__ << " found hit with id " << hit_iter->second->get_trkid() << endl;
401  if (iter->second > m_ProjectionNameMap.size())
402  {
403  cout << "bad index: " << iter->second << endl;
404  gSystem->Exit(1);
405  }
406  m_TTree_ref_vec[iter->second][0] = hit_iter->second->get_x(0);
407  m_TTree_ref_vec[iter->second][1] = hit_iter->second->get_y(0);
408  m_TTree_ref_vec[iter->second][2] = hit_iter->second->get_z(0);
409  m_TTree_ref_vec[iter->second][3] = hit_iter->second->get_t(0);
410 
411  m_TTree_ref_p_vec[iter->second][0] = hit_iter->second->get_px(0);
412  m_TTree_ref_p_vec[iter->second][1] = hit_iter->second->get_py(0);
413  m_TTree_ref_p_vec[iter->second][2] = hit_iter->second->get_pz(0);
414  }
415  }
416  }
417  } // find projections
418 
419  // find number of hits
420  for (const auto &g4hit_id_hitset : track->g4hit_ids())
421  {
422  const int &g4hit_id = g4hit_id_hitset.first;
423  const set<PHG4HitDefs::keytype> &g4hit_set = g4hit_id_hitset.second;
424  //
425  auto nhit_iter = m_TTree_HitContainerID_nHits_map.find(g4hit_id);
426  assert(nhit_iter != m_TTree_HitContainerID_nHits_map.end());
427  //
428  nhit_iter->second = g4hit_set.size();
429  }
430 
431  } // if (track)
432 
433  m_TracksEvalTree->Fill();
434  } // PHG4TruthInfoContainer::ConstRange range = m_TruthInfoContainer->GetPrimaryParticleRange();
435 
436  return;
437 }
438 
439 //----------------------------------------------------------------------------//
440 //-- fill_tree():
441 //-- Fill the trees with truth, track fit, and cluster information
442 //----------------------------------------------------------------------------//
444 {
446  {
447  LogError("m_TruthInfoContainer not found!");
448  return;
449  }
450 
451  if (!m_TrackMap)
452  {
453  LogError("m_TrackMap not found!");
454  return;
455  }
456 
457  if (!m_VertexMap)
458  {
459  return;
460  }
461 
462  for (SvtxVertexMap::Iter iter = m_VertexMap->begin();
463  iter != m_VertexMap->end();
464  ++iter)
465  {
466  SvtxVertex *vertex = iter->second;
467 
468  // Make sure to reset all the TTree variables before trying to set them.
469  reset_variables();
470  //std::cout << "A1" << std::endl;
472 
473  if (!vertex)
474  {
475  LogDebug("");
476  continue;
477  }
478 
479  //std::cout << "C1" << std::endl;
480  m_TTree_TrackID = vertex->get_id();
481  m_TTree_nTracks = vertex->size_tracks();
482 
483  m_TTree_vx = vertex->get_x();
484  m_TTree_vy = vertex->get_y();
485  m_TTree_vz = vertex->get_z();
486  m_TTree_DeltaVx = sqrt(vertex->get_error(1, 1));
487  m_TTree_DeltaVy = sqrt(vertex->get_error(2, 2));
488  m_TTree_DeltaVz = sqrt(vertex->get_error(3, 3));
489 
490  // best matched vertex
491  PHG4VtxPoint *best_vtx = nullptr;
492  int best_n_match = -1;
493  map<PHG4VtxPoint *, int> vertex_match_map;
494  for (auto iter = vertex->begin_tracks(); iter != vertex->end_tracks(); ++iter)
495  {
496  const auto &trackid = *iter;
497  const auto trackIter = m_TrackMap->find(trackid);
498 
499  if (trackIter == m_TrackMap->end()) continue;
500 
501  SvtxTrack_FastSim *temp = dynamic_cast<SvtxTrack_FastSim *>(trackIter->second);
502 
503  if (!temp) continue;
504 
505  const auto g4trackID = temp->get_truth_track_id();
506  const PHG4Particle *g4particle = m_TruthInfoContainer->GetParticle(g4trackID);
507  assert(g4particle);
509 
510  int n_match = ++vertex_match_map[vtx];
511 
512  if (n_match > best_n_match)
513  {
514  best_n_match = n_match;
515  best_vtx = vtx;
516  }
517  }
518  if (best_vtx)
519  {
520  m_TTree_gvx = best_vtx->get_x();
521  m_TTree_gvy = best_vtx->get_y();
522  m_TTree_gvz = best_vtx->get_z();
523  m_TTree_gvt = best_vtx->get_t();
524 
525  m_TTree_nFromTruth = best_n_match;
526  m_TTree_gTrackID = best_vtx->get_id();
527  }
528  m_VertexEvalTree->Fill();
529  }
530  //std::cout << "B3" << std::endl;
531 
532  return;
533 }
534 
535 //----------------------------------------------------------------------------//
536 //-- reset_variables():
537 //-- Reset all the tree variables to their default values.
538 //-- Needs to be called at the start of every event
539 //----------------------------------------------------------------------------//
541 {
542  m_TTree_Event = -9999;
543 
544  //-- truth
545  m_TTree_gTrackID = -9999;
546  m_TTree_gFlavor = -9999;
547  m_TTree_gpx = NAN;
548  m_TTree_gpy = NAN;
549  m_TTree_gpz = NAN;
550 
551  m_TTree_gvx = NAN;
552  m_TTree_gvy = NAN;
553  m_TTree_gvz = NAN;
554  m_TTree_gvt = NAN;
555 
556  //-- reco
557  m_TTree_TrackID = -9999;
558  m_TTree_Charge = -9999;
559  m_TTree_nHits = -9999;
560  m_TTree_px = NAN;
561  m_TTree_py = NAN;
562  m_TTree_pz = NAN;
563  m_TTree_pcax = NAN;
564  m_TTree_pcay = NAN;
565  m_TTree_pcaz = NAN;
566  m_TTree_dca2d = NAN;
567 
568  m_TTree_vx = NAN;
569  m_TTree_vy = NAN;
570  m_TTree_vz = NAN;
571  m_TTree_DeltaVx = NAN;
572  m_TTree_DeltaVy = NAN;
573  m_TTree_DeltaVz = NAN;
574  m_TTree_nTracks = -9999;
575  m_TTree_nFromTruth = -9999;
576  for (auto &elem : m_TTree_proj_vec) std::fill(elem.begin(), elem.end(), -9999);
577  for (auto &elem : m_TTree_proj_p_vec) std::fill(elem.begin(), elem.end(), -9999);
578  for (auto &elem : m_TTree_ref_vec) std::fill(elem.begin(), elem.end(), -9999);
579  for (auto &elem : m_TTree_ref_p_vec) std::fill(elem.begin(), elem.end(), -9999);
580  for (auto &pair : m_TTree_HitContainerID_nHits_map) pair.second=0;
581 }
582 
583 //----------------------------------------------------------------------------//
584 //-- GetNodes():
585 //-- Get all the all the required nodes off the node tree
586 //----------------------------------------------------------------------------//
588 {
589  //DST objects
590  //Truth container
591  m_TruthInfoContainer = findNode::getClass<PHG4TruthInfoContainer>(topNode,
592  "G4TruthInfo");
594  {
595  cout << PHWHERE << " PHG4TruthInfoContainer node not found on node tree"
596  << endl;
598  }
599 
600  m_TrackMap = findNode::getClass<SvtxTrackMap>(topNode,
602  //std::cout << m_TrackMapName << std::endl;
603  if (!m_TrackMap)
604  {
605  cout << PHWHERE << "SvtxTrackMap node with name "
606  << m_TrackMapName
607  << " not found on node tree"
608  << endl;
610  }
611 
612  m_VertexMap = findNode::getClass<SvtxVertexMap>(topNode, "SvtxVertexMap");
613  if (!m_VertexMap && Verbosity())
614  {
615  cout << PHWHERE << "SvtxTrackMap node with name SvtxVertexMap not found on node tree. Will not build the vertex eval tree"
616  << endl;
617  }
618 
620 }
621 
623 {
624  vector<float> floatvec{-9999, -9999, -9999, -9999};
625  m_TTree_proj_vec.push_back(floatvec);
626  m_TTree_proj_p_vec.push_back(floatvec);
627  m_TTree_ref_vec.push_back(floatvec);
628  m_TTree_ref_p_vec.push_back(floatvec);
629  // using m_ProjectionNameMap.size() makes sure it starts with 0 and then increments by 1
630  // for each additional projection
631  m_ProjectionNameMap.insert(make_pair(name, m_ProjectionNameMap.size()));
632  return;
633 }