Public Member Functions | Private Member Functions | Private Attributes | List of all members
trkf::CosmicTracker Class Reference
Inheritance diagram for trkf::CosmicTracker:
art::EDProducer art::detail::Producer art::detail::LegacyModule art::Modifier art::ModuleBase art::ProductRegistryHelper

Public Member Functions

 CosmicTracker (fhicl::ParameterSet const &pset)
 
- Public Member Functions inherited from art::EDProducer
 EDProducer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDProducer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Producer
virtual ~Producer () noexcept
 
 Producer (fhicl::ParameterSet const &)
 
 Producer (Producer const &)=delete
 
 Producer (Producer &&)=delete
 
Produceroperator= (Producer const &)=delete
 
Produceroperator= (Producer &&)=delete
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Modifier
 ~Modifier () noexcept
 
 Modifier ()
 
 Modifier (Modifier const &)=delete
 
 Modifier (Modifier &&)=delete
 
Modifieroperator= (Modifier const &)=delete
 
Modifieroperator= (Modifier &&)=delete
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

Private Member Functions

void produce (art::Event &evt)
 

Private Attributes

cluster::ClusterMatchTQ fClusterMatch
 
trkf::CosmicTrackerAlg fCTAlg
 
std::string fClusterModuleLabel
 label for input cluster collection More...
 
std::string fSortDir
 sort space points More...
 
bool fStitchTracks
 Stitch tracks from different TPCs. More...
 
double fDisCut
 Distance cut for track merging. More...
 
double fAngCut
 Angle cut for track merging. More...
 
bool fTrajOnly
 Only use trajectory points from TrackTrajectoryAlg for debugging. More...
 

Additional Inherited Members

- Public Types inherited from art::EDProducer
using ModuleType = EDProducer
 
using WorkerType = WorkerT< EDProducer >
 
- Public Types inherited from art::detail::Producer
template<typename UserConfig , typename KeysToIgnore = void>
using Table = Modifier::Table< UserConfig, KeysToIgnore >
 
- Public Types inherited from art::Modifier
template<typename UserConfig , typename UserKeysToIgnore = void>
using Table = ProducerTable< UserConfig, detail::ModuleConfig, UserKeysToIgnore >
 
- Static Public Member Functions inherited from art::EDProducer
static void commitEvent (EventPrincipal &ep, Event &e)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 216 of file CosmicTracker_module.cc.

Constructor & Destructor Documentation

trkf::CosmicTracker::CosmicTracker ( fhicl::ParameterSet const &  pset)
explicit

Definition at line 243 of file CosmicTracker_module.cc.

244  : EDProducer{pset}
245  , fClusterMatch(pset.get<fhicl::ParameterSet>("ClusterMatch"))
246  , fCTAlg(pset.get<fhicl::ParameterSet>("CTAlg"))
247  {
248  fClusterModuleLabel = pset.get<std::string>("ClusterModuleLabel");
249  fSortDir = pset.get<std::string>("SortDirection", "+z");
250  fStitchTracks = pset.get<bool>("StitchTracks");
251  fDisCut = pset.get<double>("DisCut");
252  fAngCut = pset.get<double>("AngCut");
253  fTrajOnly = pset.get<bool>("TrajOnly");
254 
255  produces<std::vector<recob::Track>>();
256  produces<std::vector<recob::SpacePoint>>();
257  produces<art::Assns<recob::Track, recob::Cluster>>();
258  produces<art::Assns<recob::Track, recob::SpacePoint>>();
259  produces<art::Assns<recob::SpacePoint, recob::Hit>>();
260  produces<art::Assns<recob::Track, recob::Hit>>();
261  }
std::string string
Definition: nybbler.cc:12
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
bool fTrajOnly
Only use trajectory points from TrackTrajectoryAlg for debugging.
trkf::CosmicTrackerAlg fCTAlg
bool fStitchTracks
Stitch tracks from different TPCs.
double fAngCut
Angle cut for track merging.
cluster::ClusterMatchTQ fClusterMatch
double fDisCut
Distance cut for track merging.
std::string fClusterModuleLabel
label for input cluster collection
std::string fSortDir
sort space points

Member Function Documentation

void trkf::CosmicTracker::produce ( art::Event evt)
privatevirtual

Implements art::EDProducer.

Definition at line 265 of file CosmicTracker_module.cc.

266  {
267 
268  // get services
270 
271  std::unique_ptr<std::vector<recob::Track>> tcol(new std::vector<recob::Track>);
272  std::unique_ptr<std::vector<recob::SpacePoint>> spcol(new std::vector<recob::SpacePoint>);
273  std::unique_ptr<art::Assns<recob::Track, recob::SpacePoint>> tspassn(
275  std::unique_ptr<art::Assns<recob::Track, recob::Cluster>> tcassn(
277  std::unique_ptr<art::Assns<recob::Track, recob::Hit>> thassn(
279  std::unique_ptr<art::Assns<recob::SpacePoint, recob::Hit>> shassn(
281 
282  //double timetick = detprop->SamplingRate()*1e-3; //time sample in us
283  //double presamplings = detprop->TriggerOffset(); // presamplings in ticks
284  //double plane_pitch = geom->PlanePitch(0,1); //wire plane pitch in cm
285  //double wire_pitch = geom->WirePitch(); //wire pitch in cm
286  //double Efield_drift = larprop->Efield(0); // Electric Field in the drift region in kV/cm
287  //double Temperature = detprop->Temperature(); // LAr Temperature in K
288 
289  //double driftvelocity = larprop->DriftVelocity(Efield_drift,Temperature); //drift velocity in the drift region (cm/us)
290  //double timepitch = driftvelocity*timetick; //time sample (cm)
291 
292  // get input Cluster object(s).
293  art::Handle<std::vector<recob::Cluster>> clusterListHandle;
294  std::vector<art::Ptr<recob::Cluster>> clusterlist;
295  if (evt.getByLabel(fClusterModuleLabel, clusterListHandle))
296  art::fill_ptr_vector(clusterlist, clusterListHandle);
297 
298  art::FindManyP<recob::Hit> fm(clusterListHandle, evt, fClusterModuleLabel);
299 
300  // find matched clusters
301  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
302  auto const detProp =
304  auto const matchedclusters = fClusterMatch.MatchedClusters(clockData, detProp, clusterlist, fm);
305 
306  // get track space points
307  std::vector<std::vector<trkPoint>> trkpts(matchedclusters.size());
308  for (size_t itrk = 0; itrk < matchedclusters.size(); ++itrk) { //loop over tracks
309 
310  std::vector<art::Ptr<recob::Hit>> hitlist;
311  for (size_t iclu = 0; iclu < matchedclusters[itrk].size(); ++iclu) { //loop over clusters
312 
313  std::vector<art::Ptr<recob::Hit>> hits = fm.at(matchedclusters[itrk][iclu]);
314  for (size_t ihit = 0; ihit < hits.size(); ++ihit) {
315  hitlist.push_back(hits[ihit]);
316  }
317  }
318  // reconstruct space points and directions
319  fCTAlg.SPTReco(clockData, detProp, hitlist);
320  if (!fTrajOnly) {
321  if (!fCTAlg.trkPos.size()) continue;
322  for (size_t i = 0; i < hitlist.size(); ++i) {
323  trkPoint trkpt;
324  trkpt.pos = fCTAlg.trkPos[i];
325  trkpt.dir = fCTAlg.trkDir[i];
326  trkpt.hit = hitlist[i];
327  trkpts[itrk].push_back(trkpt);
328  }
329  if (fSortDir == "+x") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_x0);
330  if (fSortDir == "-x") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_x1);
331  if (fSortDir == "+y") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_y0);
332  if (fSortDir == "-y") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_y1);
333  if (fSortDir == "+z") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_z0);
334  if (fSortDir == "-z") std::sort(trkpts[itrk].begin(), trkpts[itrk].end(), sp_sort_z1);
335  }
336 
337  if (fTrajOnly) { //debug only
338  if (!fCTAlg.trajPos.size()) continue;
339  size_t spStart = spcol->size();
340  std::vector<recob::SpacePoint> spacepoints;
341  for (size_t ipt = 0; ipt < fCTAlg.trajPos.size(); ++ipt) {
343  double hitcoord[3];
344  hitcoord[0] = fCTAlg.trajPos[ipt].X();
345  hitcoord[1] = fCTAlg.trajPos[ipt].Y();
346  hitcoord[2] = fCTAlg.trajPos[ipt].Z();
347  double err[6] = {util::kBogusD};
348  recob::SpacePoint mysp(hitcoord,
349  err,
351  spStart + spacepoints.size()); //3d point at end of track
352  spacepoints.push_back(mysp);
353  spcol->push_back(mysp);
354  util::CreateAssn(evt, *spcol, fCTAlg.trajHit[ipt], *shassn);
355  } // ihit
356  size_t spEnd = spcol->size();
357  if (fSortDir == "+x") {
358  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_x0);
359  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_x0);
360  }
361  if (fSortDir == "-x") {
362  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_x1);
363  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_x1);
364  }
365  if (fSortDir == "+y") {
366  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_y0);
367  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_y0);
368  }
369  if (fSortDir == "-y") {
370  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_y1);
371  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_y1);
372  }
373  if (fSortDir == "+z") {
374  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_z0);
375  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_z0);
376  }
377  if (fSortDir == "-z") {
378  std::sort(spacepoints.begin(), spacepoints.end(), spt_sort_z1);
379  std::sort(spcol->begin() + spStart, spcol->begin() + spEnd, spt_sort_z1);
380  }
381 
382  if (spacepoints.size() > 0) {
383 
384  // make a vector of the trajectory points along the track
385  std::vector<TVector3> xyz(spacepoints.size());
386  for (size_t s = 0; s < spacepoints.size(); ++s) {
387  xyz[s] = TVector3(spacepoints[s].XYZ());
388  }
389  //Calculate track direction cosines
390  TVector3 startpointVec, endpointVec, DirCos;
391  startpointVec = xyz[0];
392  endpointVec = xyz.back();
393  DirCos = endpointVec - startpointVec;
394  //SetMag casues a crash if the magnitude of the vector is zero
395  try {
396  DirCos.SetMag(1.0); //normalize vector
397  }
398  catch (...) {
399  std::cout << "The Spacepoint is infinitely small" << std::endl;
400  continue;
401  }
402  std::vector<TVector3> dircos(spacepoints.size(), DirCos);
403 
404  tcol->push_back(
407  recob::Track::Flags_t(xyz.size()),
408  false),
409  0,
410  -1.,
411  0,
414  tcol->size()));
415 
416  // make associations between the track and space points
417  util::CreateAssn(evt, *tcol, *spcol, *tspassn, spStart, spEnd);
418 
419  // and the hits and track
420  std::vector<art::Ptr<recob::Hit>> trkhits;
421  for (size_t ihit = 0; ihit < hitlist.size(); ++ihit) {
422  trkhits.push_back(hitlist[ihit]);
423  }
424  util::CreateAssn(evt, *tcol, trkhits, *thassn);
425  }
426  }
427  } //itrk
428 
429  // by default creates a spacepoint and direction for each hit
430  if (!fTrajOnly) {
431  std::vector<std::vector<unsigned int>> trkidx;
432  //stitch tracks from different TPCs
433  if (fStitchTracks) {
434  for (size_t itrk1 = 0; itrk1 < trkpts.size(); ++itrk1) {
435  for (size_t itrk2 = itrk1 + 1; itrk2 < trkpts.size(); ++itrk2) {
436  if (MatchTrack(trkpts[itrk1], trkpts[itrk2], fDisCut, fAngCut)) {
437  int found1 = -1;
438  int found2 = -1;
439  for (size_t i = 0; i < trkidx.size(); ++i) {
440  for (size_t j = 0; j < trkidx[i].size(); ++j) {
441  if (trkidx[i][j] == itrk1) found1 = i;
442  if (trkidx[i][j] == itrk2) found2 = i;
443  }
444  }
445  if (found1 == -1 && found2 == -1) {
446  std::vector<unsigned int> tmp;
447  tmp.push_back(itrk1);
448  tmp.push_back(itrk2);
449  trkidx.push_back(tmp);
450  }
451  else if (found1 == -1 && found2 != -1) {
452  trkidx[found2].push_back(itrk1);
453  }
454  else if (found1 != -1 && found2 == -1) {
455  trkidx[found1].push_back(itrk2);
456  }
457  else if (found1 != found2) { //merge two vectors
458  trkidx[found1].insert(
459  trkidx[found1].end(), trkidx[found2].begin(), trkidx[found2].end());
460  trkidx.erase(trkidx.begin() + found2);
461  }
462  } //found match
463  } //itrk2
464  } //itrk1
465  for (size_t itrk = 0; itrk < trkpts.size(); ++itrk) {
466  bool found = false;
467  for (size_t i = 0; i < trkidx.size(); ++i) {
468  for (size_t j = 0; j < trkidx[i].size(); ++j) {
469  if (trkidx[i][j] == itrk) found = true;
470  }
471  }
472  if (!found) {
473  std::vector<unsigned int> tmp;
474  tmp.push_back(itrk);
475  trkidx.push_back(tmp);
476  }
477  }
478  } //stitch
479  else {
480  trkidx.resize(trkpts.size());
481  for (size_t i = 0; i < trkpts.size(); ++i) {
482  trkidx[i].push_back(i);
483  }
484  }
485 
486  //make recob::track and associations
487  for (size_t i = 0; i < trkidx.size(); ++i) {
488  //all track points
489  std::vector<trkPoint> finaltrkpts;
490  //all the clusters associated with the current track
491  std::vector<art::Ptr<recob::Cluster>> clustersPerTrack;
492  //all hits
493  std::vector<art::Ptr<recob::Hit>> hitlist;
494  for (size_t j = 0; j < trkidx[i].size(); ++j) {
495  for (size_t k = 0; k < trkpts[trkidx[i][j]].size(); ++k) {
496  finaltrkpts.push_back(trkpts[trkidx[i][j]][k]);
497  hitlist.push_back(trkpts[trkidx[i][j]][k].hit);
498  } //k
499  for (size_t iclu = 0; iclu < matchedclusters[trkidx[i][j]].size(); ++iclu) {
500  art::Ptr<recob::Cluster> cluster(clusterListHandle,
501  matchedclusters[trkidx[i][j]][iclu]);
502  clustersPerTrack.push_back(cluster);
503  }
504 
505  } //j
506  if (fStitchTracks) {
507  if (fSortDir == "+x") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_x0);
508  if (fSortDir == "-x") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_x1);
509  if (fSortDir == "+y") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_y0);
510  if (fSortDir == "-y") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_y1);
511  if (fSortDir == "+z") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_z0);
512  if (fSortDir == "-z") std::sort(finaltrkpts.begin(), finaltrkpts.end(), sp_sort_z1);
513  }
514  size_t spStart = spcol->size();
515  std::vector<recob::SpacePoint> spacepoints;
516  for (size_t ipt = 0; ipt < finaltrkpts.size(); ++ipt) {
518  sp_hits.push_back(finaltrkpts[ipt].hit);
519  double hitcoord[3];
520  hitcoord[0] = finaltrkpts[ipt].pos.X();
521  hitcoord[1] = finaltrkpts[ipt].pos.Y();
522  hitcoord[2] = finaltrkpts[ipt].pos.Z();
523  double err[6] = {util::kBogusD};
524  recob::SpacePoint mysp(hitcoord,
525  err,
527  spStart + spacepoints.size()); //3d point at end of track
528  spacepoints.push_back(mysp);
529  spcol->push_back(mysp);
530  util::CreateAssn(evt, *spcol, sp_hits, *shassn);
531  } // ipt
532  size_t spEnd = spcol->size();
533  if (spacepoints.size() > 0) {
534  // make a vector of the trajectory points along the track
535  std::vector<TVector3> xyz(spacepoints.size());
536  std::vector<TVector3> dircos(spacepoints.size());
537  for (size_t s = 0; s < spacepoints.size(); ++s) {
538  xyz[s] = TVector3(spacepoints[s].XYZ());
539  dircos[s] = finaltrkpts[s].dir;
540  //flip direction if needed.
541  if (spacepoints.size() > 1) {
542  if (s == 0) {
543  TVector3 xyz1 = TVector3(spacepoints[s + 1].XYZ());
544  TVector3 dir = xyz1 - xyz[s];
545  if (dir.Angle(dircos[s]) > 0.8 * TMath::Pi()) { dircos[s] = -dircos[s]; }
546  }
547  else {
548  TVector3 dir = xyz[s] - xyz[s - 1];
549  if (dir.Angle(dircos[s]) > 0.8 * TMath::Pi()) { dircos[s] = -dircos[s]; }
550  }
551  }
552  }
553  tcol->push_back(
556  recob::Track::Flags_t(xyz.size()),
557  false),
558  0,
559  -1.,
560  0,
563  tcol->size()));
564 
565  // make associations between the track and space points
566  util::CreateAssn(evt, *tcol, *spcol, *tspassn, spStart, spEnd);
567 
568  // now the track and clusters
569  util::CreateAssn(evt, *tcol, clustersPerTrack, *tcassn);
570 
571  // and the hits and track
572  std::vector<art::Ptr<recob::Hit>> trkhits;
573  for (size_t ihit = 0; ihit < hitlist.size(); ++ihit) {
574  trkhits.push_back(hitlist[ihit]);
575  }
576  util::CreateAssn(evt, *tcol, trkhits, *thassn);
577  }
578 
579  } //i
580  }
581  mf::LogVerbatim("Summary") << std::setfill('-') << std::setw(175) << "-" << std::setfill(' ');
582  mf::LogVerbatim("Summary") << "CosmicTracker Summary:";
583  for (unsigned int i = 0; i < tcol->size(); ++i)
584  mf::LogVerbatim("Summary") << tcol->at(i);
585  mf::LogVerbatim("Summary") << "CosmicTracker Summary End:";
586 
587  evt.put(std::move(tcol));
588  evt.put(std::move(spcol));
589  evt.put(std::move(tspassn));
590  evt.put(std::move(tcassn));
591  evt.put(std::move(thassn));
592  evt.put(std::move(shassn));
593 
594  return;
595  }
bool spt_sort_x1(const recob::SpacePoint h1, const recob::SpacePoint h2)
end
while True: pbar.update(maxval-len(onlies[E][S])) #print iS, "/", len(onlies[E][S]) found = False for...
bool sp_sort_x1(const trkPoint &tp1, const trkPoint &tp2)
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
void cluster(In first, In last, Out result, Pred *pred)
Definition: NNClusters.h:41
bool sp_sort_y1(const trkPoint &tp1, const trkPoint &tp2)
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
TrackTrajectory::Flags_t Flags_t
Definition: Track.h:68
std::vector< TVector3 > trajPos
bool sp_sort_z1(const trkPoint &tp1, const trkPoint &tp2)
bool MatchTrack(const std::vector< trkPoint > &trkpts1, const std::vector< trkPoint > &trkpts2, double discut, double angcut)
art::Ptr< recob::Hit > hit
std::vector< TVector3 > trkDir
bool fTrajOnly
Only use trajectory points from TrackTrajectoryAlg for debugging.
string dir
Cluster finding and building.
trkf::CosmicTrackerAlg fCTAlg
bool fStitchTracks
Stitch tracks from different TPCs.
std::vector< std::vector< unsigned int > > MatchedClusters(const detinfo::DetectorClocksData &clockdata, const detinfo::DetectorPropertiesData &detProp, const std::vector< art::Ptr< recob::Cluster >> &clusterlist, const art::FindManyP< recob::Hit > &fm) const
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)
bool sp_sort_x0(const trkPoint &tp1, const trkPoint &tp2)
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
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
bool spt_sort_y0(const recob::SpacePoint h1, const recob::SpacePoint h2)
bool spt_sort_z1(const recob::SpacePoint h1, const recob::SpacePoint h2)
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
string tmp
Definition: languages.py:63
cluster::ClusterMatchTQ fClusterMatch
void err(const char *fmt,...)
Definition: message.cpp:226
double fDisCut
Distance cut for track merging.
Detector simulation of raw signals on wires.
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
std::vector< TVector3 > trkPos
std::string fClusterModuleLabel
label for input cluster collection
bool spt_sort_z0(const recob::SpacePoint h1, const recob::SpacePoint h2)
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)
static constexpr double fm
Definition: Units.h:75
decltype(auto) constexpr begin(T &&obj)
ADL-aware version of std::begin.
Definition: StdUtils.h:72
constexpr double kBogusD
obviously bogus double value
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:297
std::string fSortDir
sort space points
Q_EXPORT QTSManip setfill(int f)
Definition: qtextstream.h:337
static QCString * s
Definition: config.cpp:1042
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
QTextStream & endl(QTextStream &s)
std::vector< std::vector< art::Ptr< recob::Hit > > > trajHit

Member Data Documentation

double trkf::CosmicTracker::fAngCut
private

Angle cut for track merging.

Definition at line 232 of file CosmicTracker_module.cc.

cluster::ClusterMatchTQ trkf::CosmicTracker::fClusterMatch
private

Definition at line 223 of file CosmicTracker_module.cc.

std::string trkf::CosmicTracker::fClusterModuleLabel
private

label for input cluster collection

Definition at line 226 of file CosmicTracker_module.cc.

trkf::CosmicTrackerAlg trkf::CosmicTracker::fCTAlg
private

Definition at line 224 of file CosmicTracker_module.cc.

double trkf::CosmicTracker::fDisCut
private

Distance cut for track merging.

Definition at line 231 of file CosmicTracker_module.cc.

std::string trkf::CosmicTracker::fSortDir
private

sort space points

Definition at line 228 of file CosmicTracker_module.cc.

bool trkf::CosmicTracker::fStitchTracks
private

Stitch tracks from different TPCs.

Definition at line 230 of file CosmicTracker_module.cc.

bool trkf::CosmicTracker::fTrajOnly
private

Only use trajectory points from TrackTrajectoryAlg for debugging.

Definition at line 234 of file CosmicTracker_module.cc.


The documentation for this class was generated from the following file: