LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
CosmicTracker_module.cc
Go to the documentation of this file.
1 //
3 // CosmicTracker
4 //
5 // Tracker to reconstruct cosmic ray muons and neutrino interactions
6 //
7 // tjyang@fnal.gov
8 //
9 //
10 // ** Modified by Muhammad Elnimr to check multiple TPCs for the 35t prototype.
11 //
12 // mmelnimr@as.ua.edu
13 // April 2014
14 //
15 // Major modifications to separate algorithms from the module and
16 // use Bruce's TrackTrajectoryAlg to get trajectory points. T. Yang
17 // June 2015
19 
20 // C++ includes
21 #include <algorithm>
22 #include <iomanip>
23 #include <iostream>
24 #include <string>
25 
26 // Framework includes
35 #include "fhiclcpp/ParameterSet.h"
37 
38 // LArSoft includes
49 
50 // ROOT includes
51 #include "TMath.h"
52 #include "TVector3.h"
53 
54 struct trkPoint {
55 public:
56  TVector3 pos;
57  TVector3 dir; //direction cosines
59 };
60 
61 bool AnglesConsistent(const TVector3& p1,
62  const TVector3& p2,
63  const TVector3& a1,
64  const TVector3& a2,
65  double angcut)
66 {
67  double angle1 = a1.Angle(a2);
68  if (angle1 > TMath::PiOver2()) angle1 = TMath::Pi() - angle1;
69  double angle2 = a1.Angle(p1 - p2);
70  if (angle2 > TMath::PiOver2()) angle2 = TMath::Pi() - angle2;
71  if (angle1 < angcut && angle2 < angcut)
72  return true;
73  else
74  return false;
75 }
76 
77 // compare end points and directions
78 bool MatchTrack(const std::vector<trkPoint>& trkpts1,
79  const std::vector<trkPoint>& trkpts2,
80  double angcut)
81 {
82  bool match = false;
83  if (!trkpts1.size()) return match;
84  if (!trkpts2.size()) return match;
85  if ((trkpts1[0].hit)->WireID().Cryostat == (trkpts2[0].hit)->WireID().Cryostat &&
86  (trkpts1[0].hit)->WireID().TPC == (trkpts2[0].hit)->WireID().TPC)
87  return match;
88  // art::ServiceHandle<geo::Geometry const> geom;
89  // const geo::TPCGeo &thetpc1 = geom->TPC((trkpts1[0].hit)->WireID().TPC, (trkpts1[0].hit)->WireID().Cryostat);
90  // const geo::TPCGeo &thetpc2 = geom->TPC((trkpts2[0].hit)->WireID().TPC, (trkpts2[0].hit)->WireID().Cryostat);
91 
92  //std::cout<<trkpts1[0].pos.Y()<<" "<<trkpts1.back().pos.Y()<<" "<<trkpts1[0].pos.Z()<<" "<<trkpts1.back().pos.Z()<<std::endl;
93  //std::cout<<trkpts2[0].pos.Y()<<" "<<trkpts2.back().pos.Y()<<" "<<trkpts2[0].pos.Z()<<" "<<trkpts2.back().pos.Z()<<std::endl;
94 
95  // if (thetpc1.InFiducialY(trkpts1[0].pos.Y(),5)&&
96  // thetpc1.InFiducialY(trkpts1.back().pos.Y(),5)&&
97  // thetpc1.InFiducialZ(trkpts1[0].pos.Z(),5)&&
98  // thetpc1.InFiducialZ(trkpts1.back().pos.Z(),5)) return match;
99  // //std::cout<<"pass 1"<<std::endl;
100  // if (thetpc2.InFiducialY(trkpts2[0].pos.Y(),5)&&
101  // thetpc2.InFiducialY(trkpts2.back().pos.Y(),5)&&
102  // thetpc2.InFiducialZ(trkpts2[0].pos.Z(),5)&&
103  // thetpc2.InFiducialZ(trkpts2.back().pos.Z(),5)) return match;
104  // //std::cout<<"pass 2"<<std::endl;
105 
106  if (AnglesConsistent(trkpts1[0].pos, trkpts2[0].pos, trkpts1[0].dir, trkpts2[0].dir, angcut))
107  match = true;
108  if (AnglesConsistent(
109  trkpts1.back().pos, trkpts2[0].pos, trkpts1.back().dir, trkpts2[0].dir, angcut))
110  match = true;
111  if (AnglesConsistent(
112  trkpts1[0].pos, trkpts2.back().pos, trkpts1[0].dir, trkpts2.back().dir, angcut))
113  match = true;
114  if (AnglesConsistent(
115  trkpts1.back().pos, trkpts2.back().pos, trkpts1.back().dir, trkpts2.back().dir, angcut))
116  match = true;
117 
118  return match;
119 }
120 
122 {
123  return h1->WireID().Wire < h2->WireID().Wire;
124 }
125 
126 bool sp_sort_x0(const trkPoint& tp1, const trkPoint& tp2)
127 {
128  return tp1.pos.X() < tp2.pos.X();
129 }
130 
131 bool sp_sort_x1(const trkPoint& tp1, const trkPoint& tp2)
132 {
133  return tp1.pos.X() > tp2.pos.X();
134 }
135 
136 bool sp_sort_y0(const trkPoint& tp1, const trkPoint& tp2)
137 {
138  return tp1.pos.Y() < tp2.pos.Y();
139 }
140 
141 bool sp_sort_y1(const trkPoint& tp1, const trkPoint& tp2)
142 {
143  return tp1.pos.Y() > tp2.pos.Y();
144 }
145 
146 bool sp_sort_z0(const trkPoint& tp1, const trkPoint& tp2)
147 {
148  return tp1.pos.Z() < tp2.pos.Z();
149 }
150 
151 bool sp_sort_z1(const trkPoint& tp1, const trkPoint& tp2)
152 {
153  return tp1.pos.Z() > tp2.pos.Z();
154 }
155 
157 {
158  const double* xyz1 = h1.XYZ();
159  const double* xyz2 = h2.XYZ();
160  return xyz1[0] < xyz2[0];
161 }
162 
164 {
165  const double* xyz1 = h1.XYZ();
166  const double* xyz2 = h2.XYZ();
167  return xyz1[0] > xyz2[0];
168 }
169 
171 {
172  const double* xyz1 = h1.XYZ();
173  const double* xyz2 = h2.XYZ();
174  return xyz1[1] < xyz2[1];
175 }
176 
178 {
179  const double* xyz1 = h1.XYZ();
180  const double* xyz2 = h2.XYZ();
181  return xyz1[1] > xyz2[1];
182 }
183 
185 {
186  const double* xyz1 = h1.XYZ();
187  const double* xyz2 = h2.XYZ();
188  return xyz1[2] < xyz2[2];
189 }
190 
192 {
193  const double* xyz1 = h1.XYZ();
194  const double* xyz2 = h2.XYZ();
195  return xyz1[2] > xyz2[2];
196 }
197 
198 namespace trkf {
199 
201  public:
202  explicit CosmicTracker(fhicl::ParameterSet const& pset);
203 
204  private:
205  void produce(art::Event& evt);
206 
209 
210  std::string fClusterModuleLabel;
211 
212  std::string fSortDir;
213 
215  double fAngCut;
216 
217  bool fTrajOnly;
218 
219  }; // class CosmicTracker
220 
221 }
222 
223 namespace trkf {
224 
225  //-------------------------------------------------
226  CosmicTracker::CosmicTracker(fhicl::ParameterSet const& pset)
227  : EDProducer{pset}
228  , fClusterMatch(pset.get<fhicl::ParameterSet>("ClusterMatch"))
229  , fCTAlg(pset.get<fhicl::ParameterSet>("CTAlg"))
230  {
231  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
232  fSortDir = pset.get<std::string>("SortDirection", "+z");
233  fStitchTracks = pset.get<bool>("StitchTracks");
234  fAngCut = pset.get<double>("AngCut");
235  fTrajOnly = pset.get<bool>("TrajOnly");
236 
237  produces<std::vector<recob::Track>>();
238  produces<std::vector<recob::SpacePoint>>();
239  produces<art::Assns<recob::Track, recob::Cluster>>();
240  produces<art::Assns<recob::Track, recob::SpacePoint>>();
241  produces<art::Assns<recob::SpacePoint, recob::Hit>>();
242  produces<art::Assns<recob::Track, recob::Hit>>();
243  }
244 
245  //------------------------------------------------------------------------------------//
247  {
248 
249  // get services
251 
252  std::unique_ptr<std::vector<recob::Track>> tcol(new std::vector<recob::Track>);
253  std::unique_ptr<std::vector<recob::SpacePoint>> spcol(new std::vector<recob::SpacePoint>);
254  std::unique_ptr<art::Assns<recob::Track, recob::SpacePoint>> tspassn(
256  std::unique_ptr<art::Assns<recob::Track, recob::Cluster>> tcassn(
258  std::unique_ptr<art::Assns<recob::Track, recob::Hit>> thassn(
260  std::unique_ptr<art::Assns<recob::SpacePoint, recob::Hit>> shassn(
262 
263  //double timetick = detprop->SamplingRate()*1e-3; //time sample in us
264  //double presamplings = detprop->TriggerOffset(); // presamplings in ticks
265  //double plane_pitch = geom->PlanePitch(0,1); //wire plane pitch in cm
266  //double wire_pitch = geom->WirePitch(); //wire pitch in cm
267  //double Efield_drift = larprop->Efield(0); // Electric Field in the drift region in kV/cm
268  //double Temperature = detprop->Temperature(); // LAr Temperature in K
269 
270  //double driftvelocity = larprop->DriftVelocity(Efield_drift,Temperature); //drift velocity in the drift region (cm/us)
271  //double timepitch = driftvelocity*timetick; //time sample (cm)
272 
273  // get input Cluster object(s).
274  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
275  std::vector<art::Ptr<recob::Cluster>> clusterlist;
276  if (evt.getByLabel(fClusterModuleLabel, clusterListHandle))
277  art::fill_ptr_vector(clusterlist, clusterListHandle);
278 
279  art::FindManyP<recob::Hit> fm(clusterListHandle, evt, fClusterModuleLabel);
280 
281  // find matched clusters
282  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
283  auto const detProp =
285  auto const matchedclusters = fClusterMatch.MatchedClusters(detProp, clusterlist, fm);
286 
287  // get track space points
288  std::vector<std::vector<trkPoint>> trkpts(matchedclusters.size());
289  for (size_t itrk = 0; itrk < matchedclusters.size(); ++itrk) { //loop over tracks
290 
291  std::vector<art::Ptr<recob::Hit>> hitlist;
292  for (size_t iclu = 0; iclu < matchedclusters[itrk].size(); ++iclu) { //loop over clusters
293 
294  std::vector<art::Ptr<recob::Hit>> hits = fm.at(matchedclusters[itrk][iclu]);
295  for (size_t ihit = 0; ihit < hits.size(); ++ihit) {
296  hitlist.push_back(hits[ihit]);
297  }
298  }
299  // reconstruct space points and directions
300  fCTAlg.SPTReco(clockData, detProp, hitlist);
301  if (!fTrajOnly) {
302  if (!fCTAlg.trkPos.size()) continue;
303  for (size_t i = 0; i < hitlist.size(); ++i) {
304  trkPoint trkpt;
305  trkpt.pos = fCTAlg.trkPos[i];
306  trkpt.dir = fCTAlg.trkDir[i];
307  trkpt.hit = hitlist[i];
308  trkpts[itrk].push_back(trkpt);
309  }
310  if (fSortDir == "+x") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_x0);
311  if (fSortDir == "-x") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_x1);
312  if (fSortDir == "+y") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_y0);
313  if (fSortDir == "-y") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_y1);
314  if (fSortDir == "+z") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_z0);
315  if (fSortDir == "-z") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_z1);
316  }
317 
318  if (fTrajOnly) { //debug only
319  if (!fCTAlg.trajPos.size()) continue;
320  size_t spStart = spcol->size();
321  std::vector<recob::SpacePoint> spacepoints;
322  for (size_t ipt = 0; ipt < fCTAlg.trajPos.size(); ++ipt) {
324  double hitcoord[3];
325  hitcoord[0] = fCTAlg.trajPos[ipt].X();
326  hitcoord[1] = fCTAlg.trajPos[ipt].Y();
327  hitcoord[2] = fCTAlg.trajPos[ipt].Z();
328  double err[6] = {util::kBogusD};
329  recob::SpacePoint mysp(hitcoord,
330  err,
332  spStart + spacepoints.size()); //3d point at end of track
333  spacepoints.push_back(mysp);
334  spcol->push_back(mysp);
335  util::CreateAssn(evt, *spcol, fCTAlg.trajHit[ipt], *shassn);
336  } // ihit
337  size_t spEnd = spcol->size();
338  if (fSortDir == "+x") {
339  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_x0);
340  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_x0);
341  }
342  if (fSortDir == "-x") {
343  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_x1);
344  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_x1);
345  }
346  if (fSortDir == "+y") {
347  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_y0);
348  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_y0);
349  }
350  if (fSortDir == "-y") {
351  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_y1);
352  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_y1);
353  }
354  if (fSortDir == "+z") {
355  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_z0);
356  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_z0);
357  }
358  if (fSortDir == "-z") {
359  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_z1);
360  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_z1);
361  }
362 
363  if (spacepoints.size() > 0) {
364 
365  // make a vector of the trajectory points along the track
366  std::vector<TVector3> xyz(spacepoints.size());
367  for (size_t s = 0; s < spacepoints.size(); ++s) {
368  xyz[s] = TVector3(spacepoints[s].XYZ());
369  }
370  //Calculate track direction cosines
371  TVector3 startpointVec, endpointVec, DirCos;
372  startpointVec = xyz[0];
373  endpointVec = xyz.back();
374  DirCos = endpointVec - startpointVec;
375  //SetMag casues a crash if the magnitude of the vector is zero
376  try {
377  DirCos.SetMag(1.0); //normalize vector
378  }
379  catch (...) {
380  std::cout << "The Spacepoint is infinitely small" << std::endl;
381  continue;
382  }
383  std::vector<TVector3> dircos(spacepoints.size(), DirCos);
384 
385  tcol->push_back(
388  recob::Track::Flags_t(xyz.size()),
389  false),
390  0,
391  -1.,
392  0,
395  tcol->size()));
396 
397  // make associations between the track and space points
398  util::CreateAssn(evt, *tcol, *spcol, *tspassn, spStart, spEnd);
399 
400  // and the hits and track
401  std::vector<art::Ptr<recob::Hit>> trkhits;
402  for (size_t ihit = 0; ihit < hitlist.size(); ++ihit) {
403  trkhits.push_back(hitlist[ihit]);
404  }
405  util::CreateAssn(evt, *tcol, trkhits, *thassn);
406  }
407  }
408  } //itrk
409 
410  // by default creates a spacepoint and direction for each hit
411  if (!fTrajOnly) {
412  std::vector<std::vector<unsigned int>> trkidx;
413  //stitch tracks from different TPCs
414  if (fStitchTracks) {
415  for (size_t itrk1 = 0; itrk1 < trkpts.size(); ++itrk1) {
416  for (size_t itrk2 = itrk1 + 1; itrk2 < trkpts.size(); ++itrk2) {
417  if (MatchTrack(trkpts[itrk1], trkpts[itrk2], fAngCut)) {
418  int found1 = -1;
419  int found2 = -1;
420  for (size_t i = 0; i < trkidx.size(); ++i) {
421  for (size_t j = 0; j < trkidx[i].size(); ++j) {
422  if (trkidx[i][j] == itrk1) found1 = i;
423  if (trkidx[i][j] == itrk2) found2 = i;
424  }
425  }
426  if (found1 == -1 && found2 == -1) {
427  std::vector<unsigned int> tmp;
428  tmp.push_back(itrk1);
429  tmp.push_back(itrk2);
430  trkidx.push_back(tmp);
431  }
432  else if (found1 == -1 && found2 != -1) {
433  trkidx[found2].push_back(itrk1);
434  }
435  else if (found1 != -1 && found2 == -1) {
436  trkidx[found1].push_back(itrk2);
437  }
438  else if (found1 != found2) { //merge two vectors
439  trkidx[found1].insert(
440  trkidx[found1].end(), trkidx[found2].begin(), trkidx[found2].end());
441  trkidx.erase(trkidx.begin() + found2);
442  }
443  } //found match
444  } //itrk2
445  } //itrk1
446  for (size_t itrk = 0; itrk < trkpts.size(); ++itrk) {
447  bool found = false;
448  for (size_t i = 0; i < trkidx.size(); ++i) {
449  for (size_t j = 0; j < trkidx[i].size(); ++j) {
450  if (trkidx[i][j] == itrk) found = true;
451  }
452  }
453  if (!found) {
454  std::vector<unsigned int> tmp;
455  tmp.push_back(itrk);
456  trkidx.push_back(tmp);
457  }
458  }
459  } //stitch
460  else {
461  trkidx.resize(trkpts.size());
462  for (size_t i = 0; i < trkpts.size(); ++i) {
463  trkidx[i].push_back(i);
464  }
465  }
466 
467  //make recob::track and associations
468  for (size_t i = 0; i < trkidx.size(); ++i) {
469  //all track points
470  std::vector<trkPoint> finaltrkpts;
471  //all the clusters associated with the current track
472  std::vector<art::Ptr<recob::Cluster>> clustersPerTrack;
473  //all hits
474  std::vector<art::Ptr<recob::Hit>> hitlist;
475  for (size_t j = 0; j < trkidx[i].size(); ++j) {
476  for (size_t k = 0; k < trkpts[trkidx[i][j]].size(); ++k) {
477  finaltrkpts.push_back(trkpts[trkidx[i][j]][k]);
478  hitlist.push_back(trkpts[trkidx[i][j]][k].hit);
479  } //k
480  for (size_t iclu = 0; iclu < matchedclusters[trkidx[i][j]].size(); ++iclu) {
481  art::Ptr<recob::Cluster> cluster(clusterListHandle,
482  matchedclusters[trkidx[i][j]][iclu]);
483  clustersPerTrack.push_back(cluster);
484  }
485 
486  } //j
487  if (fStitchTracks) {
488  if (fSortDir == "+x") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_x0);
489  if (fSortDir == "-x") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_x1);
490  if (fSortDir == "+y") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_y0);
491  if (fSortDir == "-y") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_y1);
492  if (fSortDir == "+z") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_z0);
493  if (fSortDir == "-z") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_z1);
494  }
495  size_t spStart = spcol->size();
496  std::vector<recob::SpacePoint> spacepoints;
497  for (size_t ipt = 0; ipt < finaltrkpts.size(); ++ipt) {
499  sp_hits.push_back(finaltrkpts[ipt].hit);
500  double hitcoord[3];
501  hitcoord[0] = finaltrkpts[ipt].pos.X();
502  hitcoord[1] = finaltrkpts[ipt].pos.Y();
503  hitcoord[2] = finaltrkpts[ipt].pos.Z();
504  double err[6] = {util::kBogusD};
505  recob::SpacePoint mysp(hitcoord,
506  err,
508  spStart + spacepoints.size()); //3d point at end of track
509  spacepoints.push_back(mysp);
510  spcol->push_back(mysp);
511  util::CreateAssn(evt, *spcol, sp_hits, *shassn);
512  } // ipt
513  size_t spEnd = spcol->size();
514  if (spacepoints.size() > 0) {
515  // make a vector of the trajectory points along the track
516  std::vector<TVector3> xyz(spacepoints.size());
517  std::vector<TVector3> dircos(spacepoints.size());
518  for (size_t s = 0; s < spacepoints.size(); ++s) {
519  xyz[s] = TVector3(spacepoints[s].XYZ());
520  dircos[s] = finaltrkpts[s].dir;
521  //flip direction if needed.
522  if (spacepoints.size() > 1) {
523  if (s == 0) {
524  TVector3 xyz1 = TVector3(spacepoints[s + 1].XYZ());
525  TVector3 dir = xyz1 - xyz[s];
526  if (dir.Angle(dircos[s]) > 0.8 * TMath::Pi()) { dircos[s] = -dircos[s]; }
527  }
528  else {
529  TVector3 dir = xyz[s] - xyz[s - 1];
530  if (dir.Angle(dircos[s]) > 0.8 * TMath::Pi()) { dircos[s] = -dircos[s]; }
531  }
532  }
533  }
534  tcol->push_back(
537  recob::Track::Flags_t(xyz.size()),
538  false),
539  0,
540  -1.,
541  0,
544  tcol->size()));
545 
546  // make associations between the track and space points
547  util::CreateAssn(evt, *tcol, *spcol, *tspassn, spStart, spEnd);
548 
549  // now the track and clusters
550  util::CreateAssn(evt, *tcol, clustersPerTrack, *tcassn);
551 
552  // and the hits and track
553  std::vector<art::Ptr<recob::Hit>> trkhits;
554  for (size_t ihit = 0; ihit < hitlist.size(); ++ihit) {
555  trkhits.push_back(hitlist[ihit]);
556  }
557  util::CreateAssn(evt, *tcol, trkhits, *thassn);
558  }
559 
560  } //i
561  }
562  mf::LogVerbatim("Summary") << std::setfill('-') << std::setw(175) << "-" << std::setfill(' ');
563  mf::LogVerbatim("Summary") << "CosmicTracker Summary:";
564  for (unsigned int i = 0; i < tcol->size(); ++i)
565  mf::LogVerbatim("Summary") << tcol->at(i);
566  mf::LogVerbatim("Summary") << "CosmicTracker Summary End:";
567 
568  evt.put(std::move(tcol));
569  evt.put(std::move(spcol));
570  evt.put(std::move(tspassn));
571  evt.put(std::move(tcassn));
572  evt.put(std::move(thassn));
573  evt.put(std::move(shassn));
574 
575  return;
576  }
577 
579 
580 } // namespace
bool spt_sort_x1(const recob::SpacePoint h1, const recob::SpacePoint h2)
bool sp_sort_x1(const trkPoint &tp1, const trkPoint &tp2)
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
bool sp_sort_y1(const trkPoint &tp1, const trkPoint &tp2)
std::vector< Vector_t > convertCollToVector(std::vector< Vector > const &coll)
Definition: TrackingTypes.h:78
Declaration of signal hit object.
TrackTrajectory::Flags_t Flags_t
Definition: Track.h:67
std::vector< TVector3 > trajPos
bool sp_sort_z1(const trkPoint &tp1, const trkPoint &tp2)
ROOT::Math::SMatrix< Double32_t, 5, 5, ROOT::Math::MatRepSym< Double32_t, 5 >> SMatrixSym55
art::Ptr< recob::Hit > hit
std::vector< TVector3 > trkDir
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
Float_t tmp
Definition: plot.C:35
bool fTrajOnly
Only use trajectory points from TrackTrajectoryAlg for debugging.
Cluster finding and building.
trkf::CosmicTrackerAlg fCTAlg
bool fStitchTracks
Stitch tracks from different TPCs.
bool MatchTrack(const std::vector< trkPoint > &trkpts1, const std::vector< trkPoint > &trkpts2, double angcut)
double fAngCut
Angle cut for track merging.
bool sp_sort_z0(const trkPoint &tp1, const trkPoint &tp2)
bool sp_sort_y0(const trkPoint &tp1, const trkPoint &tp2)
std::vector< std::vector< unsigned int > > MatchedClusters(const detinfo::DetectorPropertiesData &detProp, const std::vector< art::Ptr< recob::Cluster >> &clusterlist, const art::FindManyP< recob::Hit > &fm) const
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
bool sp_sort_x0(const trkPoint &tp1, const trkPoint &tp2)
#define a2
bool AnglesConsistent(const TVector3 &p1, const TVector3 &p2, const TVector3 &a1, const TVector3 &a2, double angcut)
void hits()
Definition: readHits.C:15
bool SortByWire(art::Ptr< recob::Hit > const &h1, art::Ptr< recob::Hit > const &h2)
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
A trajectory in space reconstructed from hits.
bool spt_sort_y0(const recob::SpacePoint h1, const recob::SpacePoint h2)
bool spt_sort_z1(const recob::SpacePoint h1, const recob::SpacePoint h2)
Provides recob::Track data product.
void produce(art::Event &evt)
const Double32_t * XYZ() const
Definition: SpacePoint.h:78
std::vector< Point_t > convertCollToPoint(std::vector< Point > const &coll)
Definition: TrackingTypes.h:68
cluster::ClusterMatchTQ fClusterMatch
Declaration of cluster object.
Detector simulation of raw signals on wires.
TH1F * h2
Definition: plot.C:44
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
std::vector< TVector3 > trkPos
std::string fClusterModuleLabel
label for input cluster collection
bool spt_sort_z0(const recob::SpacePoint h1, const recob::SpacePoint h2)
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Utility object to perform functions of association.
TDirectory * dir
Definition: macro.C:5
void SPTReco(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, std::vector< art::Ptr< recob::Hit >> &fHits)
bool spt_sort_x0(const recob::SpacePoint h1, const recob::SpacePoint h2)
bool spt_sort_y1(const recob::SpacePoint h1, const recob::SpacePoint h2)
TH1F * h1
Definition: plot.C:41
decltype(auto) constexpr begin(T &&obj)
ADL-aware version of std::begin.
Definition: StdUtils.h:69
constexpr double kBogusD
obviously bogus double value
TCEvent evt
Definition: DataStructs.cxx:8
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
std::string fSortDir
sort space points
art framework interface to geometry description
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:49
#define a1
std::vector< std::vector< art::Ptr< recob::Hit > > > trajHit