Track3Dreco_module.cc
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////
2 //
3 // \file Track3Dreco.cxx
4 //
5 // maddalena.antonello@lngs.infn.it
6 // ornella.palamara@lngs.infn.it
7 // ART port and edits: soderber,echurch@fnal.gov
8 // This algorithm is designed to reconstruct 3D tracks through a simple
9 // 2D-track matching algorithm
10 ////////////////////////////////////////////////////////////////////////
12 
13 // C++ includes
14 #include <algorithm>
15 #include <cmath>
16 #include <iomanip>
17 #include <string>
18 #include <vector>
19 
20 // Framework includes
27 #include "fhiclcpp/ParameterSet.h"
29 
30 // LArSoft includes
31 #include "canvas/Persistency/Common/FindManyP.h"
41 
42 // ROOT includes
43 #include "TF1.h"
44 #include "TGraph.h"
45 #include "TMath.h"
46 #include "TVector2.h"
47 #include "TVector3.h"
48 
49 namespace trkf {
50 
51  class Track3Dreco : public art::EDProducer {
52  public:
53  explicit Track3Dreco(fhicl::ParameterSet const& pset);
54 
55  private:
56  void produce(art::Event& evt) override;
57 
58  int ftmatch; ///< tolerance for time matching (in time samples)
59  double fchi2dof; ///< tolerance for chi2/dof of cluster fit to function
60  std::string fClusterModuleLabel; ///< label for input cluster collection
61  }; // class Track3Dreco
62 
63 }
64 
65 namespace trkf {
66 
67  //-------------------------------------------------
69  {
70  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
71  ftmatch = pset.get<int>("TMatch");
72  fchi2dof = pset.get<double>("Chi2DOFmax");
73 
74  produces<std::vector<recob::Track>>();
75  produces<std::vector<recob::SpacePoint>>();
76  produces<art::Assns<recob::Track, recob::Cluster>>();
77  produces<art::Assns<recob::Track, recob::SpacePoint>>();
78  produces<art::Assns<recob::SpacePoint, recob::Hit>>();
79  produces<art::Assns<recob::Track, recob::Hit>>();
80  }
81 
82  //------------------------------------------------------------------------------------//
83  void
85  {
87  auto const detProp =
89 
90  auto tcol = std::make_unique<std::vector<recob::Track>>();
91  auto spacepoints = std::make_unique<std::vector<recob::SpacePoint>>();
92  auto cassn = std::make_unique<art::Assns<recob::Track, recob::Cluster>>();
93  auto sassn = std::make_unique<art::Assns<recob::Track, recob::SpacePoint>>();
94  auto shassn = std::make_unique<art::Assns<recob::SpacePoint, recob::Hit>>();
95  auto hassn = std::make_unique<art::Assns<recob::Track, recob::Hit>>();
96 
97  // define TPC parameters
98  TString tpcName = geom->GetLArTPCVolumeName();
99 
100  double YC = (geom->DetHalfHeight()) * 2.; // *ArgoNeuT* TPC active-volume height in cm
101  double Angle = geom->Plane(1).Wire(0).ThetaZ(false) -
102  TMath::Pi() / 2.; // wire angle with respect to the vertical direction
103  // Parameters temporary defined here, but possibly to be retrieved somewhere in the code
104  double timetick = 0.198; //time sample in us
105  double presamplings = 60.;
106  const double wireShift =
107  50.; // half the number of wires from the Induction(Collection) plane intersecting with a wire from the Collection(Induction) plane.
108  double plane_pitch = geom->PlanePitch(0, 1); //wire plane pitch in cm
109  double wire_pitch = geom->WirePitch(); //wire pitch in cm
110  double Efield_drift = 0.5; // Electric Field in the drift region in kV/cm
111  double Efield_SI = 0.7; // Electric Field between Shield and Induction planes in kV/cm
112  double Efield_IC = 0.9; // Electric Field between Induction and Collection planes in kV/cm
113  double Temperature = 87.6; // LAr Temperature in K
114 
115  double driftvelocity =
116  detProp.DriftVelocity(Efield_drift, Temperature); //drift velocity in the drift
117  //region (cm/us)
118  double driftvelocity_SI =
119  detProp.DriftVelocity(Efield_SI, Temperature); //drift velocity between shield
120  //and induction (cm/us)
121  double driftvelocity_IC =
122  detProp.DriftVelocity(Efield_IC, Temperature); //drift velocity between induction
123  //and collection (cm/us)
124  double timepitch = driftvelocity * timetick; //time sample (cm)
125  double tSI = plane_pitch / driftvelocity_SI / timetick; //drift time between Shield and
126  //Collection planes (time samples)
127  double tIC = plane_pitch / driftvelocity_IC / timetick; //drift time between Induction and
128  //Collection planes (time samples)
129 
130  // get input Cluster object(s).
131  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
132  evt.getByLabel(fClusterModuleLabel, clusterListHandle);
133 
134  // Declare some vectors..
135  // Induction
136  std::vector<double> Iwirefirsts; // in cm
137  std::vector<double> Iwirelasts; // in cm
138  std::vector<double> Itimefirsts; // in cm
139  std::vector<double> Itimelasts; // in cm
140  std::vector<double> Itimefirsts_line; // in cm
141  std::vector<double> Itimelasts_line; // in cm
142  std::vector<std::vector<art::Ptr<recob::Hit>>> IclusHitlists;
143  std::vector<unsigned int> Icluster_count;
144 
145  // Collection
146  std::vector<double> Cwirefirsts; // in cm
147  std::vector<double> Cwirelasts; // in cm
148  std::vector<double> Ctimefirsts; // in cm
149  std::vector<double> Ctimelasts; // in cm
150  std::vector<double> Ctimefirsts_line; // in cm
151  std::vector<double> Ctimelasts_line; // in cm
152  std::vector<std::vector<art::Ptr<recob::Hit>>> CclusHitlists;
153  std::vector<unsigned int> Ccluster_count;
154 
155  // Some variables for the hit
156  float time; //hit time at maximum
157  unsigned int wire = 0; //hit wire number
158  unsigned int plane = 0; //hit plane number
159 
160  size_t startSPIndex =
161  spacepoints->size(); //index for knowing which spacepoints are with which cluster
162  size_t endSPIndex =
163  spacepoints->size(); //index for knowing which spacepoints are with which cluster
164 
165  art::FindManyP<recob::Hit> fmh(clusterListHandle, evt, fClusterModuleLabel);
166 
167  for (size_t ii = 0; ii < clusterListHandle->size(); ++ii) {
168  art::Ptr<recob::Cluster> cl(clusterListHandle, ii);
169 
170  /////////////////////////
171  //////////// 2D track FIT
172  /////////////////////////
173 
174  // Gaaaaaah! Change me soon!!! But, for now,
175  // let's just chuck one plane's worth of info. EC, 30-Mar-2011.
176  ///\todo: This is very bad practice and should be changed ASAP
177  if (cl->View() == geo::kZ) continue;
178 
179  // Can not be const, cuz we're gonna sort 'em.
180  std::vector<art::Ptr<recob::Hit>> hitlist(fmh.at(ii));
181 
182  if (hitlist.size() == 1)
183  continue; //only one Hit in this Cluster...will cause TGraph fit to fail.
184 
185  // sort the hit list to be sure it is in the correct order
186  // using the Hit < operator
187  std::sort(hitlist.begin(), hitlist.end());
188 
189  TGraph* the2Dtrack = new TGraph(hitlist.size());
190 
191  std::vector<double> wires;
192  std::vector<double> times;
193 
194  int np = 0;
195  //loop over cluster hits
196  for (art::PtrVector<recob::Hit>::const_iterator theHit = hitlist.begin();
197  theHit != hitlist.end();
198  theHit++) {
199  //recover the Hit
200  time = (*theHit)->PeakTime();
201 
202  time -= presamplings;
203 
204  plane = (*theHit)->WireID().Plane;
205  wire = (*theHit)->WireID().Wire;
206 
207  //correct for the distance between wire planes
208  if (plane == 1) time -= tIC; // Collection
209 
210  //transform hit wire and time into cm
211  double wire_cm;
212  if (plane == 0)
213  wire_cm = (double)((wire + 3.95) * wire_pitch);
214  else
215  wire_cm = (double)((wire + 1.84) * wire_pitch);
216 
217  double time_cm;
218  if (time > tSI)
219  time_cm = (double)((time - tSI) * timepitch + tSI * driftvelocity_SI * timetick);
220  else
221  time_cm = time * driftvelocity_SI * timetick;
222 
223  wires.push_back(wire_cm);
224  times.push_back(time_cm);
225 
226  the2Dtrack->SetPoint(np, wire_cm, time_cm);
227  np++;
228  } //end of loop over cluster hits
229 
230  // fit the 2Dtrack and get some info to store
231  try {
232  the2Dtrack->Fit("pol1", "Q");
233  }
234  catch (...) {
235  mf::LogWarning("Track3Dreco") << "The 2D track fit failed";
236  continue;
237  }
238 
239  TF1* pol1 = (TF1*)the2Dtrack->GetFunction("pol1");
240  double par[2];
241  pol1->GetParameters(par);
242  double intercept = par[0];
243  double slope = par[1];
244 
245  double w0 = wires.front(); // first hit wire (cm)
246  double w1 = wires.back(); // last hit wire (cm)
247  double t0 = times.front(); // first hit time (cm)
248  double t1 = times.back(); // last hit time (cm)
249  double t0_line = intercept + (w0)*slope; // time coordinate at wire w0 on the fit line (cm)
250  double t1_line = intercept + (w1)*slope; // time coordinate at wire w1 on the fit line (cm)
251 
252  // actually store the 2Dtrack info
253  switch (plane) {
254  case 0:
255  Iwirefirsts.push_back(w0);
256  Iwirelasts.push_back(w1);
257  Itimefirsts.push_back(t0);
258  Itimelasts.push_back(t1);
259  Itimefirsts_line.push_back(t0_line);
260  Itimelasts_line.push_back(t1_line);
261  IclusHitlists.push_back(hitlist);
262  Icluster_count.push_back(ii);
263  break;
264  case 1:
265  Cwirefirsts.push_back(w0);
266  Cwirelasts.push_back(w1);
267  Ctimefirsts.push_back(t0);
268  Ctimelasts.push_back(t1);
269  Ctimefirsts_line.push_back(t0_line);
270  Ctimelasts_line.push_back(t1_line);
271  CclusHitlists.push_back(hitlist);
272  Ccluster_count.push_back(ii);
273  break;
274  }
275  delete pol1;
276  } // end of loop over all input clusters
277 
278  /////////////////////////////////////////////////////
279  /////// 2D Track Matching and 3D Track Reconstruction
280  /////////////////////////////////////////////////////
281 
282  //loop over Collection view 2D tracks
283  for (size_t collectionIter = 0; collectionIter < CclusHitlists.size(); ++collectionIter) {
284  // Recover previously stored info
285  double Cw0 = Cwirefirsts[collectionIter];
286  double Cw1 = Cwirelasts[collectionIter];
287  //double Ct0 = Ctimefirsts[collectionIter];
288  //double Ct1 = Ctimelasts[collectionIter];
289  double Ct0_line = Ctimefirsts_line[collectionIter];
290  double Ct1_line = Ctimelasts_line[collectionIter];
291  std::vector<art::Ptr<recob::Hit>> hitsCtrk = CclusHitlists[collectionIter];
292 
293  double collLength = std::hypot(Ct1_line - Ct0_line, Cw1 - Cw0);
294 
295  //loop over Induction view 2D tracks
296  for (size_t inductionIter = 0; inductionIter < IclusHitlists.size(); ++inductionIter) {
297  // Recover previously stored info
298  double Iw0 = Iwirefirsts[inductionIter];
299  double Iw1 = Iwirelasts[inductionIter];
300  double It0_line = Itimefirsts_line[inductionIter];
301  double It1_line = Itimelasts_line[inductionIter];
302  std::vector<art::Ptr<recob::Hit>> hitsItrk = IclusHitlists[inductionIter];
303 
304  double indLength = std::hypot(It1_line - It0_line, Iw1 - Iw0);
305 
306  bool forward_match = ((std::abs(Ct0_line - It0_line) < ftmatch * timepitch) &&
307  (std::abs(Ct1_line - It1_line) < ftmatch * timepitch));
308  bool backward_match = ((std::abs(Ct0_line - It1_line) < ftmatch * timepitch) &&
309  (std::abs(Ct1_line - It0_line) < ftmatch * timepitch));
310 
311  // match 2D tracks
312  if (forward_match || backward_match) {
313 
314  // Reconstruct the 3D track
315  TVector3 XYZ0, XYZ1; // track endpoints
316  if (forward_match) {
317  XYZ0.SetXYZ(Ct0_line,
318  (Cw0 - Iw0) / (2. * std::sin(Angle)),
319  (Cw0 + Iw0) / (2. * std::cos(Angle)) - YC / 2. * std::tan(Angle));
320  XYZ1.SetXYZ(Ct1_line,
321  (Cw1 - Iw1) / (2. * std::sin(Angle)),
322  (Cw1 + Iw1) / (2. * std::cos(Angle)) - YC / 2. * std::tan(Angle));
323  }
324  else {
325  XYZ0.SetXYZ(Ct0_line,
326  (Cw0 - Iw1) / (2. * std::sin(Angle)),
327  (Cw0 + Iw1) / (2. * std::cos(Angle)) - YC / 2. * std::tan(Angle));
328  XYZ1.SetXYZ(Ct1_line,
329  (Cw1 - Iw0) / (2. * std::sin(Angle)),
330  (Cw1 + Iw0) / (2. * std::cos(Angle)) - YC / 2. * std::tan(Angle));
331  }
332 
333  //compute track direction in Local co-ordinate system
334  //WARNING: There is an ambiguity introduced here for the case of backwards-going tracks.
335  //If available, vertex info. could sort this out.
336  TVector3 startpointVec, endpointVec;
337  TVector2 collVtx, indVtx;
338  if (XYZ0.Z() <= XYZ1.Z()) {
339  startpointVec.SetXYZ(XYZ0.X(), XYZ0.Y(), XYZ0.Z());
340  endpointVec.SetXYZ(XYZ1.X(), XYZ1.Y(), XYZ1.Z());
341  if (forward_match) {
342  collVtx.Set(Ct0_line, Cw0);
343  indVtx.Set(It0_line, Iw0);
344  }
345  else {
346  collVtx.Set(Ct0_line, Cw0);
347  indVtx.Set(It1_line, Iw1);
348  }
349  }
350  else {
351  startpointVec.SetXYZ(XYZ1.X(), XYZ1.Y(), XYZ1.Z());
352  endpointVec.SetXYZ(XYZ0.X(), XYZ0.Y(), XYZ0.Z());
353  if (forward_match) {
354  collVtx.Set(Ct1_line, Cw1);
355  indVtx.Set(It1_line, Iw1);
356  }
357  else {
358  collVtx.Set(Ct1_line, Cw1);
359  indVtx.Set(It0_line, Iw0);
360  }
361  }
362 
363  //compute track (normalized) cosine directions in the TPC co-ordinate system
364  TVector3 DirCos = endpointVec - startpointVec;
365 
366  //SetMag casues a crash if the magnitude of the vector is zero
367  try {
368  DirCos.SetMag(1.0); //normalize vector
369  }
370  catch (...) {
371  mf::LogWarning("Track3Dreco") << "The Spacepoint is infinitely small";
372  continue;
373  }
374 
375  art::Ptr<recob::Cluster> cl1(clusterListHandle, Icluster_count[inductionIter]);
376  art::Ptr<recob::Cluster> cl2(clusterListHandle, Ccluster_count[collectionIter]);
377  art::PtrVector<recob::Cluster> clustersPerTrack;
378  clustersPerTrack.push_back(cl1);
379  clustersPerTrack.push_back(cl2);
380 
381  /////////////////////////////
382  // Match hits
383  ////////////////////////////
384 
385  std::vector<art::Ptr<recob::Hit>> minhits =
386  hitsCtrk.size() <= hitsItrk.size() ? hitsCtrk : hitsItrk;
387  std::vector<art::Ptr<recob::Hit>> maxhits =
388  hitsItrk.size() <= hitsCtrk.size() ? hitsCtrk : hitsItrk;
389 
390  std::vector<bool> maxhitsMatch(maxhits.size());
391  for (size_t it = 0; it < maxhits.size(); ++it)
392  maxhitsMatch[it] = false;
393 
394  std::vector<recob::Hit*> hits3Dmatched;
395  // For the matching start from the view where the track projection presents less hits
396  unsigned int imaximum = 0;
397  //loop over hits
398 
399  startSPIndex = spacepoints->size();
400 
401  for (size_t imin = 0; imin < minhits.size(); ++imin) {
402  //get wire - time coordinate of the hit
403  unsigned int wire, plane1, plane2;
404  wire = minhits[imin]->WireID().Wire;
405  plane1 = minhits[imin]->WireID().Plane;
406 
407  // get the wire-time co-ordinates of the hit to be matched
408  double w1;
409  if (plane1 == 0)
410  w1 = (double)((wire + 3.95) * wire_pitch);
411  else
412  w1 = (double)((wire + 1.84) * wire_pitch);
413  double temptime1 = minhits[imin]->PeakTime() - presamplings;
414  if (plane1 == 1) temptime1 -= tIC;
415  double t1;
416  if (temptime1 > tSI)
417  t1 = (double)((temptime1 - tSI) * timepitch + tSI * driftvelocity_SI * timetick);
418  else
419  t1 = temptime1 * driftvelocity_SI * timetick;
420 
421  //get the track origin co-ordinates in the two views
422  TVector2 minVtx2D;
423  plane1 == 1 ? minVtx2D.Set(collVtx.X(), collVtx.Y()) :
424  minVtx2D.Set(indVtx.X(), indVtx.Y());
425  TVector2 maxVtx2D;
426  plane1 == 1 ? maxVtx2D.Set(indVtx.X(), indVtx.Y()) :
427  maxVtx2D.Set(collVtx.X(), collVtx.Y());
428 
429  double ratio =
430  (collLength > indLength) ? collLength / indLength : indLength / collLength;
431 
432  //compute the distance of the hit (imin) from the relative track origin
433  double minDistance = ratio * TMath::Sqrt(TMath::Power(t1 - minVtx2D.X(), 2) +
434  TMath::Power(w1 - minVtx2D.Y(), 2));
435 
436  //core matching algorithm
437  double difference = 9999999.;
438 
439  //loop over hits of the other view
440  for (size_t imax = 0; imax < maxhits.size(); ++imax) {
441  if (!maxhitsMatch[imax]) {
442  //get wire - time coordinate of the hit
443  wire = maxhits[imax]->WireID().Wire;
444  plane2 = maxhits[imax]->WireID().Plane;
445 
446  double w2;
447  if (plane2 == 0)
448  w2 = (double)((wire + 3.95) * wire_pitch);
449  else
450  w2 = (double)((wire + 1.84) * wire_pitch);
451  double temptime2 = maxhits[imax]->PeakTime() - presamplings;
452  if (plane2 == 1) temptime2 -= tIC;
453  double t2;
454  if (temptime2 > tSI)
455  t2 = (double)((temptime2 - tSI) * timepitch + tSI * driftvelocity_SI * timetick);
456  else
457  t2 = temptime2 * driftvelocity_SI * timetick;
458 
459  bool timematch = (std::abs(t1 - t2) < ftmatch * timepitch);
460  bool wirematch = (std::abs(w1 - w2) < wireShift * wire_pitch);
461 
462  double maxDistance = TMath::Sqrt(TMath::Power(t2 - maxVtx2D.X(), 2) +
463  TMath::Power(w2 - maxVtx2D.Y(), 2));
464  if (wirematch && timematch && std::abs(maxDistance - minDistance) < difference) {
465  difference = std::abs(maxDistance - minDistance);
466  imaximum = imax;
467  }
468  }
469  } // end loop over max hits
470  maxhitsMatch[imaximum] = true;
471 
473  if (difference != 9999999.) {
474  sp_hits.push_back(minhits[imin]);
475  sp_hits.push_back(maxhits[imaximum]);
476  }
477 
478  // Get the time-wire co-ordinates of the matched hit
479  wire = maxhits[imaximum]->WireID().Wire;
480  plane2 = maxhits[imaximum]->WireID().Plane;
481 
482  double w1_match;
483  if (plane2 == 0)
484  w1_match = (double)((wire + 3.95) * wire_pitch);
485  else
486  w1_match = (double)((wire + 1.84) * wire_pitch);
487  double temptime3 = maxhits[imaximum]->PeakTime() - presamplings;
488  if (plane2 == 1) temptime3 -= tIC;
489  double t1_match;
490  if (temptime3 > tSI)
491  t1_match =
492  (double)((temptime3 - tSI) * timepitch + tSI * driftvelocity_SI * timetick);
493  else
494  t1_match = temptime3 * driftvelocity_SI * timetick;
495 
496  // create the 3D hit, compute its co-ordinates and add it to the 3D hits list
497  double Ct = plane1 == 1 ? t1 : t1_match;
498  double Cw = plane1 == 1 ? w1 : w1_match;
499  double Iw = plane1 == 1 ? w1_match : w1;
500 
501  const TVector3 hit3d(Ct,
502  (Cw - Iw) / (2. * std::sin(Angle)),
503  (Cw + Iw) / (2. * std::cos(Angle)) - YC / 2. * std::tan(Angle));
504  Double_t hitcoord[3];
505  hitcoord[0] = hit3d.X();
506  hitcoord[1] = hit3d.Y();
507  hitcoord[2] = hit3d.Z();
508 
509  Double_t hitcoord_errs[3];
510  for (int i = 0; i < 3; i++)
511  hitcoord_errs[i] = -1.000;
512 
513  //3d point at end of track
514  recob::SpacePoint mysp(hitcoord, hitcoord_errs, -1., spacepoints->size());
515 
516  spacepoints->push_back(mysp);
517 
518  // associate the hits to the space point
519  util::CreateAssn(evt, *spacepoints, sp_hits, *shassn);
520 
521  } //loop over min-hits
522 
523  endSPIndex = spacepoints->size();
524 
525  // Add the 3D track to the vector of the reconstructed tracks
526  if (spacepoints->size() > startSPIndex || clustersPerTrack.size() > 0) {
527 
528  std::vector<TVector3> xyz;
529  xyz.push_back(startpointVec);
530  xyz.push_back(endpointVec);
531  std::vector<TVector3> dir_xyz;
532  dir_xyz.push_back(DirCos);
533  dir_xyz.push_back(DirCos);
534 
535  recob::Track the3DTrack(
538  recob::Track::Flags_t(xyz.size()),
539  false),
540  0,
541  -1.,
542  0,
545  tcol->size());
546  tcol->push_back(the3DTrack);
547 
548  // associate the track with its spacepoints
549  util::CreateAssn(evt, *tcol, *spacepoints, *sassn, startSPIndex, endSPIndex);
550 
551  // associate the track with its clusters
552  util::CreateAssn(evt, *tcol, clustersPerTrack, *cassn);
553 
554  art::FindManyP<recob::Hit> fmhc(clustersPerTrack, evt, fClusterModuleLabel);
555 
556  // get the hits associated with each cluster and associate those with the track
557  for (size_t p = 0; p < clustersPerTrack.size(); ++p) {
558  const std::vector<art::Ptr<recob::Hit>>& hits = fmhc.at(p);
559  util::CreateAssn(evt, *tcol, hits, *hassn);
560  }
561  }
562 
563  } //close match 2D tracks
564 
565  } //close loop over Induction view 2D tracks
566 
567  } //close loop over Collection view 2D tracks
568 
569  mf::LogVerbatim("Summary") << std::setfill('-') << std::setw(175) << "-" << std::setfill(' ');
570  mf::LogVerbatim("Summary") << "Track3Dreco Summary:";
571  for (unsigned int i = 0; i < tcol->size(); ++i) {
572  mf::LogVerbatim("Summary") << tcol->at(i);
573  }
574 
575  evt.put(std::move(tcol));
576  evt.put(std::move(spacepoints));
577  evt.put(std::move(cassn));
578  evt.put(std::move(sassn));
579  evt.put(std::move(hassn));
580  evt.put(std::move(shassn));
581 
582  return;
583  }
584 
586 
587 } // namespace
code to link reconstructed objects back to the MC truth information
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
void produce(art::Event &evt) override
PlaneGeo const & Plane(unsigned int const p, unsigned int const tpc=0, unsigned int const cstat=0) const
Returns the specified wire.
geo::Length_t PlanePitch(geo::TPCID const &tpcid, geo::PlaneID::PlaneID_t p1=0, geo::PlaneID::PlaneID_t p2=1) const
Returns the distance between two planes.
WireGeo const & Wire(unsigned int iwire) const
Definition: PlaneGeo.cxx:506
std::string string
Definition: nybbler.cc:12
ROOT::Math::SMatrix< Double32_t, 5, 5, ROOT::Math::MatRepSym< Double32_t, 5 > > SMatrixSym55
Definition: TrackingTypes.h:85
std::vector< Vector_t > convertCollToVector(std::vector< Vector > const &coll)
Definition: TrackingTypes.h:58
std::string fClusterModuleLabel
label for input cluster collection
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
TrackTrajectory::Flags_t Flags_t
Definition: Track.h:68
Planes which measure Z direction.
Definition: geo_types.h:132
double fchi2dof
tolerance for chi2/dof of cluster fit to function
double ThetaZ() const
Returns angle of wire with respect to z axis in the Y-Z plane in radians.
Definition: WireGeo.h:250
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
QAsciiDict< Entry > cl
art framework interface to geometry description
std::enable_if_t< std::is_arithmetic_v< T >, T > hypot(T x, T y)
Definition: hypot.h:60
T abs(T value)
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
typename data_t::const_iterator const_iterator
Definition: PtrVector.h:55
geo::Length_t DetHalfHeight(geo::TPCID const &tpcid) const
Returns the half height of the active volume of the specified TPC.
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
A trajectory in space reconstructed from hits.
def move(depos, offset)
Definition: depos.py:107
p
Definition: test.py:223
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
std::vector< Point_t > convertCollToPoint(std::vector< Point > const &coll)
Definition: TrackingTypes.h:55
size_type size() const
Definition: PtrVector.h:302
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
int ftmatch
tolerance for time matching (in time samples)
Encapsulate the geometry of a wire.
geo::View_t View() const
Returns the view for this cluster.
Definition: Cluster.h:741
Declaration of signal hit object.
int imax
Definition: tracks.py:195
Encapsulate the construction of a single detector plane.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
Provides recob::Track data product.
Track3Dreco(fhicl::ParameterSet const &pset)
TCEvent evt
Definition: DataStructs.cxx:7
Q_EXPORT QTSManip setfill(int f)
Definition: qtextstream.h:337
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
std::string GetLArTPCVolumeName(geo::TPCID const &tpcid) const
Return the name of specified LAr TPC volume.