All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Cluster3D_module.cc
Go to the documentation of this file.
1 /**
2  * @file Cluster3D_module.cc
3  *
4  * Class: Cluster3D
5  * Module Type: Producer
6  *
7  * @brief Producer module to create 3D clusters from input recob::Hit objects
8  *
9  * This producer module will drive the 3D association of recob::Hit objects
10  * to form 3D clusters. This information will be output as:
11  * 1) a PFParticle to anchor all the other objects (as associations)
12  * 2) three recob::Cluster objects representing 2D hit clusterswith associations
13  * to the 2D hits comprising each of them
14  * 3) One or more recob::PCAxis objects representing the Principal Components
15  * Analysis output of the space points associated to the 3D objects
16  * 4) recob::SpacePoints representing the accepted 3D points for each PFParticle
17  * 5) recob::Seed objects and associated seed hits representing candidate straight
18  * line segments in the space point collection.
19  *
20  * The module has two main sections
21  * 1) Find the 3D clusters of 3D hits
22  * 2) Get the output objects for each of the 3D clusters
23  * See the code below for more detail on these steps
24  *
25  * Note that the general 3D cluster suite of algorithms make extensive use of a set of data objects
26  * which contain volatile data members. At the end of the routine these are used to make the output
27  * LArSoft data products described above. See LarData/RecoObjects/Cluster3D.h
28  *
29  * Configuration parameters:
30  * HitFinderModuleLabel: the producer module responsible for making the recob:Hits to use
31  * EnableMonitoring: if true then basic monitoring of the module performed
32  * ClusterAlg: Parameter block required by the 3D clustering algorithm
33  * PrincipalComponentsAlg: Parameter block required by the Principal Components Analysis Algorithm
34  * SkeletonAlg: Parameter block required by the 3D skeletonization algorithm
35  * SeedFinderAlg: Parameter block required by the Hough Seed Finder algorithm
36  * PCASeedFinderAlg: Parameter block required by the PCA Seed Finder algorithm
37  * ParrallelHitsAlg: Parameter block required by the parallel hits algorithm
38  *
39  * The current producer module does not try to analyze or break apart PFParticles
40  * so, for example, all tracks emanating from a common vertex will be associated
41  * to a single PFParticle
42  *
43  * @author usher@slac.stanford.edu
44  */
45 
46 // Framework Includes
47 #include "art/Framework/Core/EDProducer.h"
48 #include "art/Framework/Core/ModuleMacros.h"
49 #include "art/Framework/Services/Registry/ServiceHandle.h"
50 #include "art/Framework/Principal/Event.h"
51 #include "art/Persistency/Common/PtrMaker.h"
52 #include "art/Utilities/make_tool.h"
53 #include "art_root_io/TFileService.h"
54 #include "cetlib/cpu_timer.h"
55 
56 // LArSoft includes
70 
85 
86 // ROOT includes
87 #include "TTree.h"
88 #include "TVector3.h"
89 
90 // std includes
91 #include <iostream>
92 #include <memory>
93 #include <string>
94 
95 //------------------------------------------------------------------------------------------------------------------------------------------
96 
97 namespace lar_cluster3d {
98  using Hit3DToSPPtrMap =
99  std::unordered_map<const reco::ClusterHit3D*, art::Ptr<recob::SpacePoint>>;
100  using RecobHitVector = art::PtrVector<recob::Hit>;
101  using RecobSpacePointVector = art::PtrVector<recob::SpacePoint>;
102 
103  //------------------------------------------------------------------------------------------------------------------------------------------
104  // Definition of the producer module here
105 
106  /**
107  * @brief Definition of the Cluster3D class
108  */
109  class Cluster3D : public art::EDProducer {
110  public:
111  explicit Cluster3D(fhicl::ParameterSet const& pset);
112 
113  private:
114  void beginJob() override;
115  void produce(art::Event& evt) override;
116 
118  public:
119  ArtOutputHandler(art::Event& evt,
120  std::string& pathName,
121  std::string& vertexName,
122  std::string& extremeName)
123  : artPCAxisVector(new std::vector<recob::PCAxis>)
124  , artPFParticleVector(new std::vector<recob::PFParticle>)
125  , artClusterVector(new std::vector<recob::Cluster>)
126  , artSpacePointVector(new std::vector<recob::SpacePoint>)
127  , artPathPointVector(new std::vector<recob::SpacePoint>)
128  , artVertexPointVector(new std::vector<recob::SpacePoint>)
129  , artExtremePointVector(new std::vector<recob::SpacePoint>)
130  , artSeedVector(new std::vector<recob::Seed>)
131  , artEdgeVector(new std::vector<recob::Edge>)
132  , artPathEdgeVector(new std::vector<recob::Edge>)
133  , artVertexEdgeVector(new std::vector<recob::Edge>)
134  , artClusterAssociations(new art::Assns<recob::Cluster, recob::Hit>)
135  , artPFPartAxisAssociations(new art::Assns<recob::PFParticle, recob::PCAxis>)
136  , artPFPartClusAssociations(new art::Assns<recob::PFParticle, recob::Cluster>)
137  , artPFPartSPAssociations(new art::Assns<recob::SpacePoint, recob::PFParticle>)
138  , artPFPartPPAssociations(new art::Assns<recob::SpacePoint, recob::PFParticle>)
139  , artPFPartSeedAssociations(new art::Assns<recob::PFParticle, recob::Seed>)
140  , artPFPartEdgeAssociations(new art::Assns<recob::Edge, recob::PFParticle>)
141  , artPFPartPathEdgeAssociations(new art::Assns<recob::Edge, recob::PFParticle>)
142  , artSeedHitAssociations(new art::Assns<recob::Seed, recob::Hit>)
143  , artSPHitAssociations(new art::Assns<recob::Hit, recob::SpacePoint>)
144  , artPPHitAssociations(new art::Assns<recob::Hit, recob::SpacePoint>)
145  , artEdgeSPAssociations(new art::Assns<recob::SpacePoint, recob::Edge>)
146  , artEdgePPAssociations(new art::Assns<recob::SpacePoint, recob::Edge>)
147  , fEvt(evt)
148  , fSPPtrMaker(evt)
149  , fSPPtrMakerPath(evt, pathName)
150  , fEdgePtrMaker(evt)
151  , fEdgePtrMakerPath(evt, pathName)
152  , fPathName(pathName)
153  , fVertexName(vertexName)
154  , fExtremeName(extremeName)
155  {}
156 
157  void
159  {
161  }
162 
163  void
164  makeSpacePointHitAssns(std::vector<recob::SpacePoint>& spacePointVector,
165  RecobHitVector& recobHits,
166  art::Assns<recob::Hit, recob::SpacePoint>& spHitAssns,
167  const std::string& path = "")
168  {
169  for (auto& hit : recobHits)
170  util::CreateAssn(fEvt, spacePointVector, hit, spHitAssns, path);
171  }
172 
173  void
175  {
180  artPCAxisVector->size() - 2,
181  artPCAxisVector->size());
182  }
183 
184  void
185  makePFPartSeedAssns(size_t numSeedsStart)
186  {
189  *artSeedVector,
191  numSeedsStart,
192  artSeedVector->size());
193  }
194 
195  void
196  makePFPartClusterAssns(size_t clusterStart)
197  {
202  clusterStart,
203  artClusterVector->size());
204  }
205 
206  void
208  std::vector<recob::SpacePoint>& spacePointVector,
209  art::Assns<recob::SpacePoint, recob::PFParticle>& pfPartSPAssociations,
210  size_t spacePointStart,
211  const std::string& instance = "")
212  {
213  for (size_t idx = spacePointStart; idx < spacePointVector.size(); idx++) {
214  art::Ptr<recob::SpacePoint> spacePoint = makeSpacePointPtr(idx, instance);
215  util::CreateAssn(fEvt, *artPFParticleVector, spacePoint, pfPartSPAssociations);
216  }
217  }
218 
219  void
220  makePFPartEdgeAssns(std::vector<recob::Edge>& edgeVector,
221  art::Assns<recob::Edge, recob::PFParticle>& pfPartEdgeAssociations,
222  size_t edgeStart,
223  const std::string& instance = "")
224  {
225  for (size_t idx = edgeStart; idx < edgeVector.size(); idx++) {
226  art::Ptr<recob::Edge> edge = makeEdgePtr(idx, instance);
227 
228  util::CreateAssn(fEvt, *artPFParticleVector, edge, pfPartEdgeAssociations);
229  }
230  }
231 
232  void
233  makeEdgeSpacePointAssns(std::vector<recob::Edge>& edgeVector,
234  RecobSpacePointVector& spacePointVector,
235  art::Assns<recob::SpacePoint, recob::Edge>& edgeSPAssociations,
236  const std::string& path = "")
237  {
238  for (auto& spacePoint : spacePointVector)
239  util::CreateAssn(fEvt, edgeVector, spacePoint, edgeSPAssociations, path);
240  }
241 
242  void
244  {
245  fEvt.put(std::move(artPCAxisVector));
246  fEvt.put(std::move(artPFParticleVector));
247  fEvt.put(std::move(artClusterVector));
248  fEvt.put(std::move(artSpacePointVector));
249  fEvt.put(std::move(artSeedVector));
250  fEvt.put(std::move(artEdgeVector));
251  fEvt.put(std::move(artPFPartAxisAssociations));
252  fEvt.put(std::move(artPFPartClusAssociations));
253  fEvt.put(std::move(artClusterAssociations));
254  fEvt.put(std::move(artPFPartSPAssociations));
255  fEvt.put(std::move(artPFPartSeedAssociations));
256  fEvt.put(std::move(artPFPartEdgeAssociations));
257  fEvt.put(std::move(artSeedHitAssociations));
258  fEvt.put(std::move(artSPHitAssociations));
259  fEvt.put(std::move(artEdgeSPAssociations));
260  fEvt.put(std::move(artPathEdgeVector), fPathName);
261  fEvt.put(std::move(artPathPointVector), fPathName);
262  fEvt.put(std::move(artEdgePPAssociations), fPathName);
263  fEvt.put(std::move(artPFPartPPAssociations), fPathName);
265  fEvt.put(std::move(artPPHitAssociations), fPathName);
266  fEvt.put(std::move(artVertexEdgeVector), fVertexName);
267  fEvt.put(std::move(artVertexPointVector), fVertexName);
268  fEvt.put(std::move(artExtremePointVector), fExtremeName);
269  }
270 
271  art::Ptr<recob::SpacePoint>
272  makeSpacePointPtr(size_t index, const std::string& instance = "")
273  {
274  if (empty(instance)) { return fSPPtrMaker(index); }
275  return fSPPtrMakerPath(index);
276  };
277 
278  art::Ptr<recob::Edge>
279  makeEdgePtr(size_t index, const std::string& instance = "")
280  {
281  if (empty(instance)) { return fEdgePtrMaker(index); }
282  return fEdgePtrMakerPath(index);
283  };
284 
285  std::unique_ptr<std::vector<recob::PCAxis>> artPCAxisVector;
286  std::unique_ptr<std::vector<recob::PFParticle>> artPFParticleVector;
287  std::unique_ptr<std::vector<recob::Cluster>> artClusterVector;
288  std::unique_ptr<std::vector<recob::SpacePoint>> artSpacePointVector;
289  std::unique_ptr<std::vector<recob::SpacePoint>> artPathPointVector;
290  std::unique_ptr<std::vector<recob::SpacePoint>> artVertexPointVector;
291  std::unique_ptr<std::vector<recob::SpacePoint>> artExtremePointVector;
292  std::unique_ptr<std::vector<recob::Seed>> artSeedVector;
293  std::unique_ptr<std::vector<recob::Edge>> artEdgeVector;
294  std::unique_ptr<std::vector<recob::Edge>> artPathEdgeVector;
295  std::unique_ptr<std::vector<recob::Edge>> artVertexEdgeVector;
296 
297  std::unique_ptr<art::Assns<recob::Cluster, recob::Hit>> artClusterAssociations;
298  std::unique_ptr<art::Assns<recob::PFParticle, recob::PCAxis>> artPFPartAxisAssociations;
299  std::unique_ptr<art::Assns<recob::PFParticle, recob::Cluster>> artPFPartClusAssociations;
300  std::unique_ptr<art::Assns<recob::SpacePoint, recob::PFParticle>> artPFPartSPAssociations;
301  std::unique_ptr<art::Assns<recob::SpacePoint, recob::PFParticle>> artPFPartPPAssociations;
302  std::unique_ptr<art::Assns<recob::PFParticle, recob::Seed>> artPFPartSeedAssociations;
303  std::unique_ptr<art::Assns<recob::Edge, recob::PFParticle>> artPFPartEdgeAssociations;
304  std::unique_ptr<art::Assns<recob::Edge, recob::PFParticle>> artPFPartPathEdgeAssociations;
305  std::unique_ptr<art::Assns<recob::Seed, recob::Hit>> artSeedHitAssociations;
306  std::unique_ptr<art::Assns<recob::Hit, recob::SpacePoint>> artSPHitAssociations;
307  std::unique_ptr<art::Assns<recob::Hit, recob::SpacePoint>> artPPHitAssociations;
308  std::unique_ptr<art::Assns<recob::SpacePoint, recob::Edge>> artEdgeSPAssociations;
309  std::unique_ptr<art::Assns<recob::SpacePoint, recob::Edge>> artEdgePPAssociations;
310 
311  private:
312  art::Event& fEvt;
313  art::PtrMaker<recob::SpacePoint> fSPPtrMaker;
314  art::PtrMaker<recob::SpacePoint> fSPPtrMakerPath;
315  art::PtrMaker<recob::Edge> fEdgePtrMaker;
316  art::PtrMaker<recob::Edge> fEdgePtrMakerPath;
317  std::string& fPathName;
318  std::string& fVertexName;
319  std::string& fExtremeName;
320  };
321 
322  /**
323  * @brief Event Preparation
324  *
325  * @param evt the ART event
326  */
327  void PrepareEvent(const art::Event& evt);
328 
329  /**
330  * @brief Initialize the internal monitoring
331  */
332  void InitializeMonitoring();
333 
334  /**
335  * @brief An interface to the seed finding algorithm
336  *
337  * @param evt the ART event
338  * @param cluster structure of information representing a single cluster
339  * @param hitToPtrMap This maps our Cluster2D hits back to art Ptr's to reco Hits
340  * @param seedVec the output vector of candidate seeds
341  * @param seedHitAssns the associations between the seeds and the 2D hits making them
342  */
343  void findTrackSeeds(art::Event& evt,
345  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
346  std::vector<recob::Seed>& seedVec,
347  art::Assns<recob::Seed, recob::Hit>& seedHitAssns) const;
348 
349  /**
350  * @brief Attempt to split clusters by using a minimum spanning tree
351  *
352  * @param clusterParameters The given cluster parameters object to try to split
353  * @param clusterParametersList The list of clusters
354  */
355  void splitClustersWithMST(reco::ClusterParameters& clusterParameters,
356  reco::ClusterParametersList& clusterParametersList) const;
357 
358  /**
359  * @brief Attempt to split clusters using the output of the Hough Filter
360  *
361  * @param clusterParameters The given cluster parameters object to try to split
362  * @param clusterParametersList The list of clusters
363  */
364  void splitClustersWithHough(reco::ClusterParameters& clusterParameters,
365  reco::ClusterParametersList& clusterParametersList) const;
366 
367  /**
368  * @brief Special routine to handle creating and saving space points
369  *
370  * @param output the object containting the art output
371  * @param clusterParameters Cluster info to output (in internal format)
372  * @param pfParticleParent The parent ID reference for the output PFParticle
373  * @param hitToPtrMap This maps our Cluster2D hits back to art Ptr's to reco Hits
374  */
376  std::vector<recob::SpacePoint>& spacePointVec,
377  art::Assns<recob::Hit, recob::SpacePoint>& spHitAssns,
378  reco::HitPairListPtr& clusHitPairVector,
379  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
380  Hit3DToSPPtrMap& hit3DToSPPtrMap,
381  const std::string& path = "") const;
382 
383  /**
384  * @brief Special routine to handle creating and saving space points
385  *
386  * @param output the object containting the art output
387  * @param clusHitPairVector List of 3D hits to output as "extreme" space points
388  */
390  reco::ConvexHullKinkTupleList& clusHitPairVector) const;
391 
392  /**
393  * @brief Special routine to handle creating and saving space points & edges associated to voronoi diagrams
394  *
395  * @param output the object containting the art output
396  * @param vertexList list of vertices in the diagram
397  * @param HalfEdgeList list of half edges in the diagram
398  */
401  dcel2d::HalfEdgeList&) const;
402 
403  /**
404  * @brief Special routine to handle creating and saving space points & edges PCA points
405  *
406  * @param output the object containting the art output
407  * @param clusterParamsList List of clusters to get PCA's from
408  */
409  using IdxToPCAMap = std::map<size_t, const reco::PrincipalComponents*>;
410 
413  IdxToPCAMap&) const;
414 
415  /**
416  * @brief This will produce art output for daughters noting that it needs to be done recursively
417  *
418  * @param output the object containting the art output
419  * @param clusterParameters Cluster info to output (in internal format)
420  * @param pfParticleParent The parent ID reference for the output PFParticle
421  * @param daughterList List of PFParticle indices for stored daughters
422  * @param hitToPtrMap This maps our Cluster2D hits back to art Ptr's to reco Hits
423  */
426  reco::ClusterParameters& clusterParameters,
427  size_t pfParticleParent,
428  IdxToPCAMap& idxToPCAMap,
429  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
430  Hit3DToSPPtrMap& hit3DToSPPtrMap,
431  Hit3DToSPPtrMap& best3DToSPPtrMap) const;
432 
433  /**
434  * @brief Top level output routine, allows checking cluster status
435  *
436  * @param hitPairList List of all 3D Hits in internal Cluster3D format
437  * @param clusterParametersList Data structure containing the cluster information to output
438  * @param hitToPtrMap This maps our Cluster2D hits back to art Ptr's to reco Hits
439  */
442  reco::HitPairList& hitPairList,
443  reco::ClusterParametersList& clusterParametersList,
444  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap) const;
445 
446  /**
447  * @brief Produces the art output from all the work done in this producer module
448  *
449  * @param output the object containting the art output
450  * @param clusterParameters Cluster info to output (in internal format)
451  * @param pfParticleParent The parent ID reference for the output PFParticle
452  * @param hitToPtrMap This maps our Cluster2D hits back to art Ptr's to reco Hits
453  */
454  size_t ConvertToArtOutput(util::GeometryUtilities const& gser,
456  reco::ClusterParameters& clusterParameters,
457  size_t pfParticleParent,
458  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
459  Hit3DToSPPtrMap& hit3DToSPPtrMap,
460  Hit3DToSPPtrMap& best3DToSPPtrMap) const;
461 
462  /**
463  * @brief There are several places we will want to know if a candidate cluster is a
464  * "parallel hits" type of cluster so encapsulate that here.
465  *
466  * @param pca The Principal Components Analysis parameters for the cluster
467  */
468  bool
470  {
471  return fabs(pca.getEigenVectors().row(2)(0)) > m_parallelHitsCosAng &&
472  3. * sqrt(pca.getEigenValues()(1)) > m_parallelHitsTransWid;
473  }
474 
475  /**
476  * @brief Count number of end of line daughters
477  *
478  * @param clusterParams input cluster parameters to look at
479  */
480  size_t countUltimateDaughters(reco::ClusterParameters& clusterParameters) const;
481 
482  /**
483  * Algorithm parameters
484  */
485  bool m_onlyMakSpacePoints; ///< If true we don't do the full cluster 3D processing
486  bool m_enableMonitoring; ///< Turn on monitoring of this algorithm
487  float m_parallelHitsCosAng; ///< Cut for PCA 3rd axis angle to X axis
488  float m_parallelHitsTransWid; ///< Cut on transverse width of cluster (PCA 2nd eigenvalue)
489 
490  /**
491  * Tree variables for output
492  */
493  TTree* m_pRecoTree; ///<
494  int m_run; ///<
495  int m_event; ///<
496  int m_hits; ///< Keeps track of the number of hits seen
497  int m_hits3D; ///< Keeps track of the number of 3D hits made
498  float m_totalTime; ///< Keeps track of total execution time
499  float m_artHitsTime; ///< Keeps track of time to recover hits
500  float m_makeHitsTime; ///< Keeps track of time to build 3D hits
501  float m_buildNeighborhoodTime; ///< Keeps track of time to build epsilon neighborhood
502  float m_dbscanTime; ///< Keeps track of time to run DBScan
503  float m_clusterMergeTime; ///< Keeps track of the time to merge clusters
504  float m_pathFindingTime; ///< Keeps track of the path finding time
505  float m_finishTime; ///< Keeps track of time to run output module
506  std::string m_pathInstance; ///< Special instance for path points
507  std::string m_vertexInstance; ///< Special instance name for vertex points
508  std::string m_extremeInstance; ///< Instance name for the extreme points
509 
510  // Algorithms
511  std::unique_ptr<lar_cluster3d::IHit3DBuilder>
512  m_hit3DBuilderAlg; ///< Builds the 3D hits to operate on
513  std::unique_ptr<lar_cluster3d::IClusterAlg>
514  m_clusterAlg; ///< Algorithm to do 3D space point clustering
515  std::unique_ptr<lar_cluster3d::IClusterModAlg>
516  m_clusterMergeAlg; ///< Algorithm to do cluster merging
517  std::unique_ptr<lar_cluster3d::IClusterModAlg>
518  m_clusterPathAlg; ///< Algorithm to do cluster path finding
519  std::unique_ptr<lar_cluster3d::IClusterParametersBuilder>
520  m_clusterBuilder; ///< Common cluster builder tool
521  PrincipalComponentsAlg m_pcaAlg; ///< Principal Components algorithm
522  SkeletonAlg m_skeletonAlg; ///< Skeleton point finder
524  PCASeedFinderAlg m_pcaSeedFinderAlg; ///< Use PCA axis to find seeds
525  ParallelHitsSeedFinderAlg m_parallelHitsAlg; ///< Deal with parallel hits clusters
526  };
527 
528  DEFINE_ART_MODULE(Cluster3D)
529 
530 } // namespace lar_cluster3d
531 
532 //------------------------------------------------------------------------------------------------------------------------------------------
533 // implementation follows
534 
535 namespace lar_cluster3d {
536 
537  Cluster3D::Cluster3D(fhicl::ParameterSet const& pset)
538  : EDProducer{pset}
539  , m_pcaAlg(pset.get<fhicl::ParameterSet>("PrincipalComponentsAlg"))
540  , m_skeletonAlg(pset.get<fhicl::ParameterSet>("SkeletonAlg"))
541  , m_seedFinderAlg(pset.get<fhicl::ParameterSet>("SeedFinderAlg"))
542  , m_pcaSeedFinderAlg(pset.get<fhicl::ParameterSet>("PCASeedFinderAlg"))
543  , m_parallelHitsAlg(pset.get<fhicl::ParameterSet>("ParallelHitsAlg"))
544  {
545  m_onlyMakSpacePoints = pset.get<bool>("MakeSpacePointsOnly", false);
546  m_enableMonitoring = pset.get<bool>("EnableMonitoring", false);
547  m_parallelHitsCosAng = pset.get<float>("ParallelHitsCosAng", 0.999);
548  m_parallelHitsTransWid = pset.get<float>("ParallelHitsTransWid", 25.0);
549  m_pathInstance = pset.get<std::string>("PathPointsName", "Path");
550  m_vertexInstance = pset.get<std::string>("VertexPointsName", "Vertex");
551  m_extremeInstance = pset.get<std::string>("ExtremePointsName", "Extreme");
552 
553  m_hit3DBuilderAlg = art::make_tool<lar_cluster3d::IHit3DBuilder>(
554  pset.get<fhicl::ParameterSet>("Hit3DBuilderAlg"));
555  m_clusterAlg =
556  art::make_tool<lar_cluster3d::IClusterAlg>(pset.get<fhicl::ParameterSet>("ClusterAlg"));
557  m_clusterMergeAlg = art::make_tool<lar_cluster3d::IClusterModAlg>(
558  pset.get<fhicl::ParameterSet>("ClusterMergeAlg"));
559  m_clusterPathAlg = art::make_tool<lar_cluster3d::IClusterModAlg>(
560  pset.get<fhicl::ParameterSet>("ClusterPathAlg"));
561  m_clusterBuilder = art::make_tool<lar_cluster3d::IClusterParametersBuilder>(
562  pset.get<fhicl::ParameterSet>("ClusterParamsBuilder"));
563 
564  // Handle special case of Space Point building outputting a new hit collection
565  m_hit3DBuilderAlg->produces(producesCollector());
566 
567  produces<std::vector<recob::PCAxis>>();
568  produces<std::vector<recob::PFParticle>>();
569  produces<std::vector<recob::Cluster>>();
570  produces<std::vector<recob::Seed>>();
571  produces<std::vector<recob::Edge>>();
572 
573  produces<art::Assns<recob::PFParticle, recob::PCAxis>>();
574  produces<art::Assns<recob::PFParticle, recob::Cluster>>();
575  produces<art::Assns<recob::PFParticle, recob::Seed>>();
576  produces<art::Assns<recob::Edge, recob::PFParticle>>();
577  produces<art::Assns<recob::Seed, recob::Hit>>();
578  produces<art::Assns<recob::Cluster, recob::Hit>>();
579 
580  produces<std::vector<recob::SpacePoint>>();
581  produces<art::Assns<recob::SpacePoint, recob::PFParticle>>();
582  produces<art::Assns<recob::Hit, recob::SpacePoint>>();
583  produces<art::Assns<recob::SpacePoint, recob::Edge>>();
584 
585  produces<std::vector<recob::SpacePoint>>(m_pathInstance);
586  produces<std::vector<recob::Edge>>(m_pathInstance);
587  produces<art::Assns<recob::SpacePoint, recob::PFParticle>>(m_pathInstance);
588  produces<art::Assns<recob::Edge, recob::PFParticle>>(m_pathInstance);
589  produces<art::Assns<recob::Hit, recob::SpacePoint>>(m_pathInstance);
590  produces<art::Assns<recob::SpacePoint, recob::Edge>>(m_pathInstance);
591 
592  produces<std::vector<recob::Edge>>(m_vertexInstance);
593  produces<std::vector<recob::SpacePoint>>(m_vertexInstance);
594 
595  produces<std::vector<recob::SpacePoint>>(m_extremeInstance);
596  }
597 
598  //------------------------------------------------------------------------------------------------------------------------------------------
599 
600  void
602  {
603  /**
604  * @brief beginJob will be tasked with initializing monitoring, in necessary, but also to init the
605  * geometry and detector services (and this probably needs to go in a "beginEvent" method?)
606  */
608  }
609 
610  //------------------------------------------------------------------------------------------------------------------------------------------
611 
612  void
613  Cluster3D::produce(art::Event& evt)
614  {
615  mf::LogInfo("Cluster3D") << " *** Cluster3D::produce(...) [Run=" << evt.run()
616  << ", Event=" << evt.id().event() << "] Starting Now! *** "
617  << std::endl;
618 
619  // Set up for monitoring the timing... at some point this should be removed in favor of
620  // external profilers
621  cet::cpu_timer theClockTotal;
622  cet::cpu_timer theClockFinish;
623 
624  if (m_enableMonitoring) theClockTotal.start();
625 
626  // This really only does anything if we are monitoring since it clears our tree variables
627  this->PrepareEvent(evt);
628 
629  // Get instances of the primary data structures needed
630  reco::ClusterParametersList clusterParametersList;
631  IHit3DBuilder::RecobHitToPtrMap clusterHitToArtPtrMap;
632  std::unique_ptr<reco::HitPairList> hitPairList(
633  new reco::HitPairList); // Potentially lots of hits, use heap instead of stack
634 
635  // Call the algorithm that builds 3D hits and stores the hit collection
636  m_hit3DBuilderAlg->Hit3DBuilder(evt, *hitPairList, clusterHitToArtPtrMap);
637 
638  // Only do the rest if we are not in the mode of only building space points (requested by ML folks)
639  if (!m_onlyMakSpacePoints) {
640  // Call the main workhorse algorithm for building the local version of candidate 3D clusters
641  m_clusterAlg->Cluster3DHits(*hitPairList, clusterParametersList);
642 
643  // Try merging clusters
644  m_clusterMergeAlg->ModifyClusters(clusterParametersList);
645 
646  // Run the path finding
647  m_clusterPathAlg->ModifyClusters(clusterParametersList);
648  }
649 
650  if (m_enableMonitoring) theClockFinish.start();
651 
652  // Get the art ouput object
654 
655  // Call the module that does the end processing (of which there is quite a bit of work!)
656  // This goes here to insure that something is always written to the data store
657  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
658  auto const detProp =
659  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt, clockData);
660  util::GeometryUtilities const gser{*lar::providerFrom<geo::Geometry>(), clockData, detProp};
661  ProduceArtClusters(gser, output, *hitPairList, clusterParametersList, clusterHitToArtPtrMap);
662 
663  // Output to art
664  output.outputObjects();
665 
666  // If monitoring then deal with the fallout
667  if (m_enableMonitoring) {
668  theClockFinish.stop();
669  theClockTotal.stop();
670 
671  m_run = evt.run();
672  m_event = evt.id().event();
673  m_totalTime = theClockTotal.accumulated_real_time();
677  m_dbscanTime = m_clusterAlg->getTimeToExecute(IClusterAlg::RUNDBSCAN) +
678  m_clusterAlg->getTimeToExecute(IClusterAlg::BUILDCLUSTERINFO);
679  m_clusterMergeTime = m_clusterMergeAlg->getTimeToExecute();
680  m_pathFindingTime = m_clusterPathAlg->getTimeToExecute();
681  m_finishTime = theClockFinish.accumulated_real_time();
682  m_hits = static_cast<int>(clusterHitToArtPtrMap.size());
683  m_hits3D = static_cast<int>(hitPairList->size());
684  m_pRecoTree->Fill();
685 
686  mf::LogDebug("Cluster3D") << "*** Cluster3D total time: " << m_totalTime
687  << ", art: " << m_artHitsTime << ", make: " << m_makeHitsTime
688  << ", build: " << m_buildNeighborhoodTime
689  << ", clustering: " << m_dbscanTime
690  << ", merge: " << m_clusterMergeTime
691  << ", path: " << m_pathFindingTime << ", finish: " << m_finishTime
692  << std::endl;
693  }
694 
695  // Will we ever get here? ;-)
696  return;
697  }
698 
699  //------------------------------------------------------------------------------------------------------------------------------------------
700 
701  void
703  {
704  art::ServiceHandle<art::TFileService> tfs;
705  m_pRecoTree = tfs->make<TTree>("monitoring", "LAr Reco");
706  m_pRecoTree->Branch("run", &m_run, "run/I");
707  m_pRecoTree->Branch("event", &m_event, "event/I");
708  m_pRecoTree->Branch("hits", &m_hits, "hits/I");
709  m_pRecoTree->Branch("hits3D", &m_hits3D, "hits3D/I");
710  m_pRecoTree->Branch("totalTime", &m_totalTime, "time/F");
711  m_pRecoTree->Branch("artHitsTime", &m_artHitsTime, "time/F");
712  m_pRecoTree->Branch("makeHitsTime", &m_makeHitsTime, "time/F");
713  m_pRecoTree->Branch("buildneigborhoodTime", &m_buildNeighborhoodTime, "time/F");
714  m_pRecoTree->Branch("dbscanTime", &m_dbscanTime, "time/F");
715  m_pRecoTree->Branch("clusterMergeTime", &m_clusterMergeTime, "time/F");
716  m_pRecoTree->Branch("pathfindingtime", &m_pathFindingTime, "time/F");
717  m_pRecoTree->Branch("finishTime", &m_finishTime, "time/F");
718 
719  m_clusterPathAlg->initializeHistograms(*tfs.get());
720  }
721 
722  //------------------------------------------------------------------------------------------------------------------------------------------
723 
724  void
725  Cluster3D::PrepareEvent(const art::Event& evt)
726  {
727  m_run = evt.run();
728  m_event = evt.id().event();
729  m_hits = 0;
730  m_hits3D = 0;
731  m_totalTime = 0.f;
732  m_artHitsTime = 0.f;
733  m_makeHitsTime = 0.f;
735  m_dbscanTime = 0.f;
736  m_pathFindingTime = 0.f;
737  m_finishTime = 0.f;
738  }
739 
740  //------------------------------------------------------------------------------------------------------------------------------------------
741 
742  void
743  Cluster3D::findTrackSeeds(art::Event& evt,
745  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
746  std::vector<recob::Seed>& seedVec,
747  art::Assns<recob::Seed, recob::Hit>& seedHitAssns) const
748  {
749  /**
750  * @brief This method provides an interface to various algorithms for finding candiate
751  * recob::Seed objects and, as well, their candidate related seed hits
752  */
753 
754  // Make sure we are using the right pca
755  reco::PrincipalComponents& fullPCA = cluster.getFullPCA();
756  reco::PrincipalComponents& skeletonPCA = cluster.getSkeletonPCA();
757  reco::HitPairListPtr& hitPairListPtr = cluster.getHitPairListPtr();
758  reco::HitPairListPtr skeletonListPtr;
759 
760  // We want to work with the "skeleton" hits so first step is to call the algorithm to
761  // recover only these hits from the entire input collection
762  m_skeletonAlg.GetSkeletonHits(hitPairListPtr, skeletonListPtr);
763 
764  // Skeleton hits are nice but we can do better if we then make a pass through to "average"
765  // the skeleton hits position in the Y-Z plane
766  m_skeletonAlg.AverageSkeletonPositions(skeletonListPtr);
767 
768  SeedHitPairListPairVec seedHitPairVec;
769 
770  // Some combination of the elements below will be used to determine which seed finding algorithm
771  // to pursue below
772  float eigenVal0 = 3. * sqrt(skeletonPCA.getEigenValues()[0]);
773  float eigenVal1 = 3. * sqrt(skeletonPCA.getEigenValues()[1]);
774  float eigenVal2 = 3. * sqrt(skeletonPCA.getEigenValues()[2]);
775  float transRMS = std::hypot(eigenVal0, eigenVal1);
776 
777  bool foundGoodSeed(false);
778 
779  // Choose a method for finding the seeds based on the PCA that was run...
780  // Currently we have an ad hoc if-else block which I hope will be improved soon!
781  if (aParallelHitsCluster(fullPCA)) {
782  // In this case we have a track moving relatively parallel to the wire plane with lots of
783  // ambiguous 3D hits. Your best bet here is to use the "parallel hits" algorithm to get the
784  // best axis and seeds
785  // This algorithm does not fail (foundGoodSeed will always return true)
786  foundGoodSeed = m_parallelHitsAlg.findTrackSeeds(hitPairListPtr, skeletonPCA, seedHitPairVec);
787  }
788  else if (eigenVal2 > 40. && transRMS < 5.) {
789  // If the input cluster is relatively "straight" then chances are it is a single straight track,
790  // probably a CR muon, and we can simply use the PCA to determine the seed
791  // This algorithm will check both "ends" of the input hits and if the angles become inconsistent
792  // then it will "fail"
793  foundGoodSeed =
794  m_pcaSeedFinderAlg.findTrackSeeds(skeletonListPtr, skeletonPCA, seedHitPairVec);
795  }
796 
797  // In the event the above two methods failed then we hit it with the real seed finder
798  if (!foundGoodSeed) {
799  // If here then we have a complicated 3D cluster and we'll use the hough transform algorithm to
800  // return a list of candidate seeds and seed hits
801  m_seedFinderAlg.findTrackSeeds(skeletonListPtr, skeletonPCA, seedHitPairVec);
802  }
803 
804  // Go through the returned lists and build out the art friendly seeds and hits
805  for (const auto& seedHitPair : seedHitPairVec) {
806  seedVec.push_back(seedHitPair.first);
807 
808  // We use a set here because our 3D hits can share 2D hits
809  // The set will make sure we get unique combinations of 2D hits
810  std::set<art::Ptr<recob::Hit>> seedHitSet;
811  for (const auto& hit3D : seedHitPair.second) {
812  for (const auto& hit2D : hit3D->getHits()) {
813  if (!hit2D) continue;
814 
815  const recob::Hit* recobHit = hit2D->getHit();
816  seedHitSet.insert(hitToPtrMap[recobHit]);
817  }
818  }
819 
820  RecobHitVector seedHitVec;
821  for (const auto& hit2D : seedHitSet)
822  seedHitVec.push_back(hit2D);
823 
824  util::CreateAssn(evt, seedVec, seedHitVec, seedHitAssns);
825  }
826  }
827 
828  struct Hit3DDistanceOrder {
829  bool
830  operator()(const std::pair<float, const reco::ClusterHit3D*>& left,
831  const std::pair<float, const reco::ClusterHit3D*>& right)
832  {
833  return left.first < right.first;
834  }
835  };
836 
837  void
839  reco::ClusterParametersList& clusterParametersList) const
840  {
841  // This is being left in place for future development. Essentially, it was an attempt to implement
842  // a Minimum Spanning Tree as a way to split a particular cluster topology, one where two straight
843  // tracks cross closely enought to appear as one cluster. As of Feb 2, 2015 I think the idea is still
844  // worth merit so am leaving this module in place for now.
845  //
846  // If this routine is called then we believe we have a cluster which needs splitting.
847  // The way we will do this is to use a Minimum Spanning Tree algorithm to associate all
848  // hits together by their distance apart. In theory, we should be able to split the cluster
849  // by finding the largest distance and splitting at that point.
850  //
851  // Typedef some data structures that we will use.
852  // Start with the adjacency map
853  typedef std::pair<float, const reco::ClusterHit3D*> DistanceHit3DPair;
854  typedef std::list<DistanceHit3DPair> DistanceHit3DPairList;
855  typedef std::map<const reco::ClusterHit3D*, DistanceHit3DPairList> Hit3DToDistanceMap;
856 
857  // Now typedef the lists we'll keep
858  typedef std::list<const reco::ClusterHit3D*> Hit3DList;
859  typedef std::pair<Hit3DList::iterator, Hit3DList::iterator> Hit3DEdgePair;
860  typedef std::pair<float, Hit3DEdgePair> DistanceEdgePair;
861  typedef std::list<DistanceEdgePair> DistanceEdgePairList;
862 
863  struct DistanceEdgePairOrder {
864  bool
865  operator()(const DistanceEdgePair& left, const DistanceEdgePair& right) const
866  {
867  return left.first > right.first;
868  }
869  };
870 
871  // Recover the hits we'll work on.
872  // Note that we use on the skeleton hits so will need to recover them
873  reco::HitPairListPtr& hitPairListPtr = clusterParameters.getHitPairListPtr();
874  reco::HitPairListPtr skeletonListPtr;
875 
876  // We want to work with the "skeleton" hits so first step is to call the algorithm to
877  // recover only these hits from the entire input collection
878  m_skeletonAlg.GetSkeletonHits(hitPairListPtr, skeletonListPtr);
879 
880  // Skeleton hits are nice but we can do better if we then make a pass through to "average"
881  // the skeleton hits position in the Y-Z plane
882  m_skeletonAlg.AverageSkeletonPositions(skeletonListPtr);
883 
884  // First task is to define and build the adjacency map
885  Hit3DToDistanceMap hit3DToDistanceMap;
886 
887  for (reco::HitPairListPtr::const_iterator hit3DOuterItr = skeletonListPtr.begin();
888  hit3DOuterItr != skeletonListPtr.end();) {
889  const reco::ClusterHit3D* hit3DOuter = *hit3DOuterItr++;
890  DistanceHit3DPairList& outerHitList = hit3DToDistanceMap[hit3DOuter];
891  TVector3 outerPos(
892  hit3DOuter->getPosition()[0], hit3DOuter->getPosition()[1], hit3DOuter->getPosition()[2]);
893 
894  for (reco::HitPairListPtr::const_iterator hit3DInnerItr = hit3DOuterItr;
895  hit3DInnerItr != skeletonListPtr.end();
896  hit3DInnerItr++) {
897  const reco::ClusterHit3D* hit3DInner = *hit3DInnerItr;
898  TVector3 innerPos(
899  hit3DInner->getPosition()[0], hit3DInner->getPosition()[1], hit3DInner->getPosition()[2]);
900  TVector3 deltaPos = innerPos - outerPos;
901  float hitDistance(float(deltaPos.Mag()));
902 
903  if (hitDistance > 20.) continue;
904 
905  hit3DToDistanceMap[hit3DInner].emplace_back(hitDistance, hit3DOuter);
906  outerHitList.emplace_back(hitDistance, hit3DInner);
907  }
908 
909  // Make sure our membership bit is clear
911  }
912 
913  // Make pass through again to order each of the lists
914  for (auto& mapPair : hit3DToDistanceMap) {
915  mapPair.second.sort(Hit3DDistanceOrder());
916  }
917 
918  // Get the containers for the MST to operate on/with
919  Hit3DList hit3DList;
920  DistanceEdgePairList distanceEdgePairList;
921 
922  // Initialize with first element
923  hit3DList.emplace_back(skeletonListPtr.front());
924  distanceEdgePairList.emplace_back(
925  DistanceEdgePair(0., Hit3DEdgePair(hit3DList.begin(), hit3DList.begin())));
926 
927  skeletonListPtr.front()->setStatusBit(reco::ClusterHit3D::SELECTEDBYMST);
928 
929  float largestDistance(0.);
930  float averageDistance(0.);
931 
932  // Now run the MST
933  // Basically, we loop until the MST list is the same size as the input list
934  while (hit3DList.size() < skeletonListPtr.size()) {
935  Hit3DList::iterator bestHit3DIter = hit3DList.begin();
936  float bestDist = 10000000.;
937 
938  // Loop through all hits currently in the list and look for closest hit not in the list
939  for (Hit3DList::iterator hit3DIter = hit3DList.begin(); hit3DIter != hit3DList.end();
940  hit3DIter++) {
941  const reco::ClusterHit3D* hit3D = *hit3DIter;
942 
943  // For the given 3D hit, find the closest to it that is not already in the list
944  DistanceHit3DPairList& nearestList = hit3DToDistanceMap[hit3D];
945 
946  while (!nearestList.empty()) {
947  const reco::ClusterHit3D* hit3DToCheck = nearestList.front().second;
948 
949  if (!hit3DToCheck->bitsAreSet(reco::ClusterHit3D::SELECTEDBYMST)) {
950  if (nearestList.front().first < bestDist) {
951  bestHit3DIter = hit3DIter;
952  bestDist = nearestList.front().first;
953  }
954  break;
955  }
956 
957  nearestList.pop_front();
958  }
959  }
960 
961  if (bestDist > largestDistance) largestDistance = bestDist;
962 
963  averageDistance += bestDist;
964 
965  // Now we add the best hit not in the list to our list, keep track of the distance
966  // to the object it was closest to
967  const reco::ClusterHit3D* bestHit3D = *bestHit3DIter; // "best" hit already in the list
968  const reco::ClusterHit3D* nextHit3D =
969  hit3DToDistanceMap[bestHit3D].front().second; // "next" hit we are adding to the list
970 
971  Hit3DList::iterator nextHit3DIter = hit3DList.insert(hit3DList.end(), nextHit3D);
972 
973  distanceEdgePairList.emplace_back(bestDist, Hit3DEdgePair(bestHit3DIter, nextHit3DIter));
974 
976  }
977 
978  averageDistance /= float(hit3DList.size());
979 
980  float thirdDist = 2. * sqrt(clusterParameters.getSkeletonPCA().getEigenValues()[2]);
981 
982  // Ok, find the largest distance in the iterator map
983  distanceEdgePairList.sort(DistanceEdgePairOrder());
984 
985  DistanceEdgePairList::iterator largestDistIter = distanceEdgePairList.begin();
986 
987  for (DistanceEdgePairList::iterator edgeIter = distanceEdgePairList.begin();
988  edgeIter != distanceEdgePairList.end();
989  edgeIter++) {
990  if (edgeIter->first < thirdDist) break;
991 
992  largestDistIter = edgeIter;
993  }
994 
995  reco::HitPairListPtr::iterator breakIter = largestDistIter->second.second;
996  reco::HitPairListPtr bestList;
997 
998  bestList.resize(std::distance(hit3DList.begin(), breakIter));
999 
1000  std::copy(hit3DList.begin(), breakIter, bestList.begin());
1001 
1002  // Remove from the grand hit list and see what happens...
1003  // The pieces below are incomplete and were really for testing only.
1004  hitPairListPtr.sort();
1005  bestList.sort();
1006 
1007  reco::HitPairListPtr::iterator newListEnd = std::set_difference(hitPairListPtr.begin(),
1008  hitPairListPtr.end(),
1009  bestList.begin(),
1010  bestList.end(),
1011  hitPairListPtr.begin());
1012 
1013  hitPairListPtr.erase(newListEnd, hitPairListPtr.end());
1014  }
1015 
1016  class CopyIfInRange {
1017  public:
1018  CopyIfInRange(float maxRange) : m_maxRange(maxRange) {}
1019 
1020  bool
1021  operator()(const reco::ClusterHit3D* hit3D) const
1022  {
1023  return hit3D->getDocaToAxis() < m_maxRange;
1024  }
1025 
1026  private:
1027  float m_maxRange;
1028  };
1029 
1030  void
1032  reco::ClusterParametersList& clusterParametersList) const
1033  {
1034  // @brief A method for splitted "crossed tracks" clusters into separate clusters
1035  //
1036  // If this routine is called then we believe we have a cluster which needs splitting.
1037  // The specific topology we are looking for is two long straight tracks which cross at some
1038  // point in close proximity so their hits were joined into a single 3D cluster. The method
1039  // to split this topology is to let the hough transform algorithm find the two leading candidates
1040  // and then to see if we use those to build two clusters instead of one.
1041 
1042  // Recover the hits we'll work on.
1043  // Note that we use on the skeleton hits so will need to recover them
1044  reco::HitPairListPtr& hitPairListPtr = clusterParameters.getHitPairListPtr();
1045  reco::HitPairListPtr skeletonListPtr;
1046 
1047  // We want to work with the "skeleton" hits so first step is to call the algorithm to
1048  // recover only these hits from the entire input collection
1049  m_skeletonAlg.GetSkeletonHits(hitPairListPtr, skeletonListPtr);
1050 
1051  // Skeleton hits are nice but we can do better if we then make a pass through to "average"
1052  // the skeleton hits position in the Y-Z plane
1053  m_skeletonAlg.AverageSkeletonPositions(skeletonListPtr);
1054 
1055  // Define the container for our lists of hits
1056  reco::HitPairListPtrList hitPairListPtrList;
1057 
1058  // Now feed this to the Hough Transform to find candidate straight lines
1060  skeletonListPtr, clusterParameters.getSkeletonPCA(), hitPairListPtrList);
1061 
1062  // We need at least two lists or else there is nothing to do
1063  if (hitPairListPtrList.size() < 2) return;
1064 
1065  // The game plan will be the following:
1066  // 1) Take the first list of hits and run the PCA on this to get an axis
1067  // - Then calculate the 3d doca for ALL hits in the cluster to this axis
1068  // - Move all hits within "3 sigam" of the axis to a new list
1069  // 2) run the PCA on the second list of hits to get that axis
1070  // - Then calculate the 3d doca for all hits in our first list
1071  // - Copy hits in the first list which are within 3 sigma of the new axis
1072  // back into our original cluster - these are shared hits
1073  reco::HitPairListPtrList::iterator hitPairListIter = hitPairListPtrList.begin();
1074  reco::HitPairListPtr& firstHitList = *hitPairListIter++;
1075  reco::PrincipalComponents firstHitListPCA;
1076 
1077  m_pcaAlg.PCAAnalysis_3D(firstHitList, firstHitListPCA);
1078 
1079  // Make sure we have a successful calculation.
1080  if (firstHitListPCA.getSvdOK()) {
1081  // The fill routines below will expect to see unused 2D hits so we need to clear the
1082  // status bits... and I am not sure of a better way...
1083  for (const auto& hit3D : hitPairListPtr) {
1084  for (const auto& hit2D : hit3D->getHits())
1085  if (hit2D) hit2D->clearStatusBits(0x1);
1086  }
1087 
1088  // Calculate the 3D doca's for the hits which were used to make this PCA
1089  m_pcaAlg.PCAAnalysis_calc3DDocas(firstHitList, firstHitListPCA);
1090 
1091  // Divine from the ether some maximum allowed range for transfering hits
1092  float allowedHitRange = 6. * firstHitListPCA.getAveHitDoca();
1093 
1094  // Now go through and calculate the 3D doca's for ALL the hits in the original cluster
1095  m_pcaAlg.PCAAnalysis_calc3DDocas(hitPairListPtr, firstHitListPCA);
1096 
1097  // Let's make a new cluster to hold the hits
1098  clusterParametersList.push_back(reco::ClusterParameters());
1099 
1100  // Can we get a reference to what we just created?
1101  reco::ClusterParameters& newClusterParams = clusterParametersList.back();
1102 
1103  reco::HitPairListPtr& newClusterHitList = newClusterParams.getHitPairListPtr();
1104 
1105  newClusterHitList.resize(hitPairListPtr.size());
1106 
1107  // Do the actual copy of the hits we want
1108  reco::HitPairListPtr::iterator newListEnd = std::copy_if(hitPairListPtr.begin(),
1109  hitPairListPtr.end(),
1110  newClusterHitList.begin(),
1111  CopyIfInRange(allowedHitRange));
1112 
1113  // Shrink to fit
1114  newClusterHitList.resize(std::distance(newClusterHitList.begin(), newListEnd));
1115 
1116  // And now remove these hits from the original cluster
1117  hitPairListPtr.remove_if(CopyIfInRange(allowedHitRange));
1118 
1119  // Get an empty hit to cluster map...
1120  reco::Hit2DToClusterMap hit2DToClusterMap;
1121 
1122  // Now "fill" the cluster parameters but turn off the hit rejection
1123  m_clusterBuilder->FillClusterParams(newClusterParams, hit2DToClusterMap, 0., 1.);
1124 
1125  // Set the skeleton pca to the value calculated above on input
1126  clusterParameters.getSkeletonPCA() = firstHitListPCA;
1127 
1128  // We are done with splitting out one track. Because the two tracks cross in
1129  // close proximity, this is the one case where we might consider sharing 3D hits
1130  // So let's make a little detour here to try to copy some of those hits back into
1131  // the main hit list
1132  reco::HitPairListPtr& secondHitList = *hitPairListIter;
1133  reco::PrincipalComponents secondHitListPCA;
1134 
1135  m_pcaAlg.PCAAnalysis_3D(secondHitList, secondHitListPCA);
1136 
1137  // Make sure we have a successful calculation.
1138  if (secondHitListPCA.getSvdOK()) {
1139  // Calculate the 3D doca's for the hits which were used to make this PCA
1140  m_pcaAlg.PCAAnalysis_calc3DDocas(secondHitList, secondHitListPCA);
1141 
1142  // Since this is the "other" cluster, we'll be a bit more generous in adding back hits
1143  float newAllowedHitRange = 6. * secondHitListPCA.getAveHitDoca();
1144 
1145  // Go through and calculate the 3D doca's for the hits in our new candidate cluster
1146  m_pcaAlg.PCAAnalysis_calc3DDocas(newClusterHitList, secondHitListPCA);
1147 
1148  // Create a temporary list to fill with the hits we might want to save
1149  reco::HitPairListPtr tempHitList(newClusterHitList.size());
1150 
1151  // Do the actual copy of the hits we want...
1152  reco::HitPairListPtr::iterator tempListEnd =
1153  std::copy_if(newClusterHitList.begin(),
1154  newClusterHitList.end(),
1155  tempHitList.begin(),
1156  CopyIfInRange(newAllowedHitRange));
1157 
1158  hitPairListPtr.insert(hitPairListPtr.end(), tempHitList.begin(), tempListEnd);
1159  }
1160 
1161  // Of course, now we need to modify the original cluster parameters
1162  reco::ClusterParameters originalParams(hitPairListPtr);
1163 
1164  // Now "fill" the cluster parameters but turn off the hit rejection
1165  m_clusterBuilder->FillClusterParams(originalParams, hit2DToClusterMap, 0., 1.);
1166 
1167  // Overwrite original cluster parameters with our new values
1168  clusterParameters.getClusterParams() = originalParams.getClusterParams();
1169  clusterParameters.getFullPCA() = originalParams.getFullPCA();
1170  clusterParameters.getSkeletonPCA() = secondHitListPCA;
1171  }
1172  }
1173 
1174  void
1176  ArtOutputHandler& output,
1177  reco::HitPairList& hitPairVector,
1178  reco::ClusterParametersList& clusterParametersList,
1179  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap) const
1180  {
1181  /**
1182  * @brief The workhorse to take the candidate 3D clusters and produce all of the necessary art output
1183  */
1184 
1185  mf::LogDebug("Cluster3D") << " *** Cluster3D::ProduceArtClusters() *** " << std::endl;
1186 
1187  // Make sure there is something to do here!
1188  if (!clusterParametersList.empty()) {
1189  // This is the loop over candidate 3D clusters
1190  // Note that it might be that the list of candidate clusters is modified by splitting
1191  // So we use the following construct to make sure we get all of them
1192  for (auto& clusterParameters : clusterParametersList) {
1193  // It should be straightforward at this point to transfer information from our vector of clusters
1194  // to the larsoft objects... of course we still have some work to do first, in particular to
1195  // find the candidate seeds and their seed hits
1196 
1197  // The chances of getting here and this condition not being true are probably zero... but check anyway
1198  if (!clusterParameters.getFullPCA().getSvdOK()) {
1199  mf::LogDebug("Cluster3D")
1200  << "--> no feature extraction done on this cluster!!" << std::endl;
1201  continue;
1202  }
1203 
1204  // Keep track of hit 3D to SP for when we do edges
1205  Hit3DToSPPtrMap hit3DToSPPtrMap;
1206  Hit3DToSPPtrMap best3DToSPPtrMap;
1207 
1208  // Do a special output of voronoi vertices here...
1209  dcel2d::VertexList& vertexList = clusterParameters.getVertexList();
1210  dcel2d::HalfEdgeList& halfEdgeList = clusterParameters.getHalfEdgeList();
1211 
1212  MakeAndSaveVertexPoints(output, vertexList, halfEdgeList);
1213 
1214  // Special case handling... if no daughters then call standard conversion routine to make sure space points
1215  // created, etc.
1216  if (clusterParameters.daughterList().empty()) {
1217  ConvertToArtOutput(gser,
1218  output,
1219  clusterParameters,
1221  hitToPtrMap,
1222  hit3DToSPPtrMap,
1223  best3DToSPPtrMap);
1224 
1225  // Get the extreme points
1226  MakeAndSaveKinkPoints(output,
1227  clusterParameters.getConvexHull().getConvexHullKinkPoints());
1228  }
1229  // Otherwise, the cluster has daughters so we handle specially
1230  else {
1231  // Set up to keep track of parent/daughters
1232  IdxToPCAMap idxToPCAMap;
1233  size_t numTotalDaughters = countUltimateDaughters(clusterParameters);
1234  size_t pfParticleIdx(output.artPFParticleVector->size() + numTotalDaughters);
1235 
1236  FindAndStoreDaughters(gser,
1237  output,
1238  clusterParameters,
1239  pfParticleIdx,
1240  idxToPCAMap,
1241  hitToPtrMap,
1242  hit3DToSPPtrMap,
1243  best3DToSPPtrMap);
1244 
1245  // Need to make a daughter vec from our map
1246  std::vector<size_t> daughterVec;
1247 
1248  for (auto& idxToPCA : idxToPCAMap)
1249  daughterVec.emplace_back(idxToPCA.first);
1250 
1251  // Now create/handle the parent PFParticle
1252  recob::PFParticle pfParticle(
1253  13, pfParticleIdx, recob::PFParticle::kPFParticlePrimary, daughterVec);
1254  output.artPFParticleVector->push_back(pfParticle);
1255 
1256  recob::PCAxis::EigenVectors eigenVecs;
1257  double eigenVals[] = {0., 0., 0.};
1258  double avePosition[] = {0., 0., 0.};
1259 
1260  eigenVecs.resize(3);
1261 
1262  reco::PrincipalComponents& skeletonPCA = clusterParameters.getSkeletonPCA();
1263 
1264  for (size_t outerIdx = 0; outerIdx < 3; outerIdx++) {
1265  avePosition[outerIdx] = skeletonPCA.getAvePosition()[outerIdx];
1266  eigenVals[outerIdx] = skeletonPCA.getEigenValues()[outerIdx];
1267 
1268  eigenVecs[outerIdx].resize(3);
1269 
1270  // Be careful here... eigen stores in column major order buy default
1271  for (size_t innerIdx = 0; innerIdx < 3; innerIdx++)
1272  eigenVecs[outerIdx][innerIdx] = skeletonPCA.getEigenVectors().row(outerIdx)[innerIdx];
1273  }
1274 
1275  recob::PCAxis skelPcAxis(skeletonPCA.getSvdOK(),
1276  skeletonPCA.getNumHitsUsed(),
1277  eigenVals,
1278  eigenVecs,
1279  avePosition,
1280  skeletonPCA.getAveHitDoca(),
1281  output.artPCAxisVector->size());
1282 
1283  output.artPCAxisVector->push_back(skelPcAxis);
1284 
1285  reco::PrincipalComponents& fullPCA = clusterParameters.getFullPCA();
1286 
1287  for (size_t outerIdx = 0; outerIdx < 3; outerIdx++) {
1288  avePosition[outerIdx] = fullPCA.getAvePosition()[outerIdx];
1289  eigenVals[outerIdx] = fullPCA.getEigenValues()[outerIdx];
1290 
1291  for (size_t innerIdx = 0; innerIdx < 3; innerIdx++)
1292  eigenVecs[outerIdx][innerIdx] = fullPCA.getEigenVectors().row(outerIdx)(innerIdx);
1293  }
1294 
1295  recob::PCAxis fullPcAxis(fullPCA.getSvdOK(),
1296  fullPCA.getNumHitsUsed(),
1297  eigenVals,
1298  eigenVecs,
1299  avePosition,
1300  fullPCA.getAveHitDoca(),
1301  output.artPCAxisVector->size());
1302 
1303  output.artPCAxisVector->push_back(fullPcAxis);
1304 
1305  // Create associations to the PFParticle
1306  output.makePFPartPCAAssns();
1307 
1308  // Make associations to all space points for this cluster
1309  MakeAndSaveSpacePoints(output,
1310  *output.artSpacePointVector.get(),
1311  *output.artSPHitAssociations.get(),
1312  clusterParameters.getHitPairListPtr(),
1313  hitToPtrMap,
1314  hit3DToSPPtrMap);
1315 
1316  // Get the extreme points
1317  MakeAndSaveKinkPoints(output,
1318  clusterParameters.getConvexHull().getConvexHullKinkPoints());
1319 
1320  // Build the edges now
1321  size_t edgeStart(output.artEdgeVector->size());
1322 
1323  for (const auto& edge : clusterParameters.getConvexHull().getConvexHullEdgeList()) {
1324  RecobSpacePointVector spacePointVec;
1325 
1326  try {
1327  spacePointVec.push_back(hit3DToSPPtrMap.at(std::get<0>(edge)));
1328  spacePointVec.push_back(hit3DToSPPtrMap.at(std::get<1>(edge)));
1329  }
1330  catch (...) {
1331  std::cout << "Caught exception in looking up space point ptr... " << std::get<0>(edge)
1332  << ", " << std::get<1>(edge) << std::endl;
1333  continue;
1334  }
1335 
1336  output.artEdgeVector->emplace_back(std::get<2>(edge),
1337  spacePointVec[0].key(),
1338  spacePointVec[1].key(),
1339  output.artEdgeVector->size());
1340 
1341  output.makeEdgeSpacePointAssns(
1342  *output.artEdgeVector.get(), spacePointVec, *output.artEdgeSPAssociations.get());
1343  }
1344 
1345  output.makePFPartEdgeAssns(
1346  *output.artEdgeVector.get(), *output.artPFPartEdgeAssociations.get(), edgeStart);
1347  }
1348  }
1349  }
1350 
1351  // Right now error matrix is uniform...
1352  int nFreePoints(0);
1353 
1354  // Run through the HitPairVector and add any unused hit pairs to the list
1355  for (auto& hitPair : hitPairVector) {
1356  if (hitPair.bitsAreSet(reco::ClusterHit3D::MADESPACEPOINT)) continue;
1357 
1358  double spacePointPos[] = {
1359  hitPair.getPosition()[0], hitPair.getPosition()[1], hitPair.getPosition()[2]};
1360  double spacePointErr[] = {1., 0., 0., 1., 0., 1.};
1361  double chisq = hitPair.getHitChiSquare();
1362 
1363  RecobHitVector recobHits;
1364 
1365  for (const auto hit : hitPair.getHits()) {
1366  if (!hit) {
1367  chisq = -1000.;
1368  continue;
1369  }
1370 
1371  art::Ptr<recob::Hit> hitPtr = hitToPtrMap[hit->getHit()];
1372  recobHits.push_back(hitPtr);
1373  }
1374 
1375  nFreePoints++;
1376 
1377  spacePointErr[1] = hitPair.getTotalCharge();
1378  spacePointErr[3] = hitPair.getChargeAsymmetry();
1379 
1380  output.artSpacePointVector->push_back(
1381  recob::SpacePoint(spacePointPos, spacePointErr, chisq, output.artSpacePointVector->size()));
1382 
1383  if (!recobHits.empty())
1384  output.makeSpacePointHitAssns(
1385  *output.artSpacePointVector.get(), recobHits, *output.artSPHitAssociations.get());
1386  }
1387 
1388  std::cout << "++++>>>> total num hits: " << hitPairVector.size()
1389  << ", num free: " << nFreePoints << std::endl;
1390  }
1391 
1392  size_t
1394  {
1395  size_t localCount(0);
1396 
1397  if (!clusterParameters.daughterList().empty()) {
1398  for (auto& clusterParams : clusterParameters.daughterList())
1399  localCount += countUltimateDaughters(clusterParams);
1400  }
1401  else
1402  localCount++;
1403 
1404  return localCount;
1405  }
1406 
1407  size_t
1409  ArtOutputHandler& output,
1410  reco::ClusterParameters& clusterParameters,
1411  size_t pfParticleParent,
1412  IdxToPCAMap& idxToPCAMap,
1413  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
1414  Hit3DToSPPtrMap& hit3DToSPPtrMap,
1415  Hit3DToSPPtrMap& best3DToSPPtrMap) const
1416  {
1417  // This is a recursive routine, we keep calling ourself as long as the daughter list is non empty
1418  if (!clusterParameters.daughterList().empty()) {
1419  for (auto& clusterParams : clusterParameters.daughterList())
1420  FindAndStoreDaughters(gser,
1421  output,
1422  clusterParams,
1423  pfParticleParent,
1424  idxToPCAMap,
1425  hitToPtrMap,
1426  hit3DToSPPtrMap,
1427  best3DToSPPtrMap);
1428  }
1429  // Otherwise we want to store the information
1430  else {
1431  size_t daughterIdx = ConvertToArtOutput(gser,
1432  output,
1433  clusterParameters,
1434  pfParticleParent,
1435  hitToPtrMap,
1436  hit3DToSPPtrMap,
1437  best3DToSPPtrMap);
1438 
1439  idxToPCAMap[daughterIdx] = &clusterParameters.getFullPCA();
1440  }
1441 
1442  return idxToPCAMap.size();
1443  }
1444 
1445  size_t
1447  ArtOutputHandler& output,
1448  reco::ClusterParameters& clusterParameters,
1449  size_t pfParticleParent,
1450  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
1451  Hit3DToSPPtrMap& hit3DToSPPtrMap,
1452  Hit3DToSPPtrMap& best3DToSPPtrMap) const
1453  {
1454  // prepare the algorithm to compute the cluster characteristics;
1455  // we use the "standard" one here, except that we override selected items
1456  // (so, thanks to metaprogramming, we finally have wrappers of wrappers);
1457  // configuration would happen here, but we are using the default
1458  // configuration for that algorithm
1459  using OverriddenClusterParamsAlg_t =
1461 
1463 
1464  // It should be straightforward at this point to transfer information from our vector of clusters
1465  // to the larsoft objects... of course we still have some work to do first, in particular to
1466  // find the candidate seeds and their seed hits
1467 
1468  // We keep track of 2 PCA axes, the first is the "full" PCA run over all the 3D hits in the
1469  // candidate cluster. The second will be that derived from just using the "skeleton" hits.
1470  // Make a copy of the full PCA to keep that, then get a reference for the skeleton PCA
1471  reco::PrincipalComponents& fullPCA = clusterParameters.getFullPCA();
1472  reco::PrincipalComponents& skeletonPCA = clusterParameters.getSkeletonPCA();
1473 
1474  size_t clusterStart = output.artClusterVector->size();
1475 
1476  // Start loop over views to build out the hit lists and the 2D cluster objects
1477  for (const auto& clusParametersPair : clusterParameters.getClusterParams()) {
1478  const reco::RecobClusterParameters& clusParams = clusParametersPair.second;
1479 
1480  // Protect against a missing view
1481  if (clusParams.m_view == geo::kUnknown) continue;
1482 
1483  // Get a vector of hit pointers to seed the cluster algorithm
1484  std::vector<const recob::Hit*> recobHitVec(clusParams.m_hitVector.size());
1485 
1486  std::transform(clusParams.m_hitVector.begin(),
1487  clusParams.m_hitVector.end(),
1488  recobHitVec.begin(),
1489  [](const auto& hit2D) { return hit2D->getHit(); });
1490 
1491  // Get the tdc/wire slope... from the unit vector...
1492  double startWire(clusParams.m_startWire);
1493  double endWire(clusParams.m_endWire);
1494  double startTime(clusParams.m_startTime);
1495  double endTime(clusParams.m_endTime);
1496 
1497  // feed the algorithm with all the cluster hits
1498  ClusterParamAlgo.ImportHits(gser, recobHitVec);
1499 
1500  // create the recob::Cluster directly in the vector
1501  cluster::ClusterCreator artCluster(gser,
1502  ClusterParamAlgo, // algo
1503  startWire, // start_wire
1504  0., // sigma_start_wire
1505  startTime, // start_tick
1506  clusParams.m_sigmaStartTime, // sigma_start_tick
1507  endWire, // end_wire
1508  0., // sigma_end_wire,
1509  endTime, // end_tick
1510  clusParams.m_sigmaEndTime, // sigma_end_tick
1511  output.artClusterVector->size(), // ID
1512  clusParams.m_view, // view
1513  clusParams.m_plane, // plane
1514  recob::Cluster::Sentry // sentry
1515  );
1516 
1517  output.artClusterVector->emplace_back(artCluster.move());
1518 
1519  // We love looping. In this case, our list of hits is comprised of "ClusterHits" and we need to get a RecobHitVector instead...
1520  RecobHitVector recobHits;
1521 
1522  for (const auto& hit2D : clusParams.m_hitVector) {
1523  if (hit2D == nullptr || hit2D->getHit() == nullptr) continue;
1524 
1525  art::Ptr<recob::Hit> hitPtr = hitToPtrMap[hit2D->getHit()];
1526  recobHits.push_back(hitPtr);
1527  }
1528 
1529  output.makeClusterHitAssns(recobHits);
1530  }
1531 
1532  // Last, let's try to get seeds for tracking..
1533  // Keep track of how many we have so far
1534  size_t numSeedsStart = output.artSeedVector->size();
1535 
1536  // Empty daughter vector for now
1537  std::vector<size_t> nullVector;
1538  size_t pfParticleIdx(output.artPFParticleVector->size());
1539 
1540  recob::PFParticle pfParticle(13, pfParticleIdx, pfParticleParent, nullVector);
1541  output.artPFParticleVector->push_back(pfParticle);
1542 
1543  // Assume that if we are a daughter particle then there is a set of "best" 3D points and the convex hull will have been calculated from them
1544  // If there is no best list then assume the convex hull is from all of the points...
1545  std::string instance("");
1546  reco::HitPairListPtr* hit3DListPtr = &clusterParameters.getHitPairListPtr();
1547  Hit3DToSPPtrMap* hit3DToPtrMap = &hit3DToSPPtrMap;
1548 
1549  auto spVector = output.artSpacePointVector.get();
1550  auto edgeVector = output.artEdgeVector.get();
1551  auto pfPartEdgeAssns = output.artPFPartEdgeAssociations.get();
1552  auto spEdgeAssns = output.artEdgeSPAssociations.get();
1553  auto spHitAssns = output.artSPHitAssociations.get();
1554  auto pfPartSPAssns = output.artPFPartSPAssociations.get();
1555 
1556  // Make associations to all space points for this cluster
1557  int spacePointStart = spVector->size();
1558 
1560  output, *spVector, *spHitAssns, *hit3DListPtr, hitToPtrMap, *hit3DToPtrMap);
1561 
1562  // And make sure to have associations to the PFParticle
1563  output.makePFPartSpacePointAssns(*spVector, *pfPartSPAssns, spacePointStart);
1564 
1565  // If they exist, make space points for the Path points
1566  if (!clusterParameters.getBestHitPairListPtr().empty()) {
1567  instance = m_pathInstance;
1568  hit3DListPtr = &clusterParameters.getBestHitPairListPtr();
1569  hit3DToPtrMap = &best3DToSPPtrMap;
1570 
1571  spVector = output.artPathPointVector.get();
1572  edgeVector = output.artPathEdgeVector.get();
1573  pfPartEdgeAssns = output.artPFPartPathEdgeAssociations.get();
1574  spEdgeAssns = output.artEdgePPAssociations.get();
1575  spHitAssns = output.artPPHitAssociations.get();
1576  pfPartSPAssns = output.artPFPartPPAssociations.get();
1577 
1578  spacePointStart = spVector->size();
1579 
1581  output, *spVector, *spHitAssns, *hit3DListPtr, hitToPtrMap, *hit3DToPtrMap, instance);
1582 
1583  // And make sure to have associations to the PFParticle
1584  output.makePFPartSpacePointAssns(*spVector, *pfPartSPAssns, spacePointStart, instance);
1585  }
1586 
1587  // Now handle the edges according to whether associated with regular or "best" space points
1588  if (!clusterParameters.getBestEdgeList().empty()) {
1589  size_t edgeStart = edgeVector->size();
1590 
1591  for (const auto& edge : clusterParameters.getBestEdgeList()) {
1592  RecobSpacePointVector spacePointVec;
1593 
1594  try {
1595  spacePointVec.push_back(hit3DToPtrMap->at(std::get<0>(edge)));
1596  spacePointVec.push_back(hit3DToPtrMap->at(std::get<1>(edge)));
1597  }
1598  catch (...) {
1599  std::cout << "Caught exception in looking up space point ptr... " << std::get<0>(edge)
1600  << ", " << std::get<1>(edge) << std::endl;
1601  continue;
1602  }
1603 
1604  edgeVector->emplace_back(
1605  std::get<2>(edge), spacePointVec[0].key(), spacePointVec[1].key(), edgeVector->size());
1606 
1607  output.makeEdgeSpacePointAssns(*edgeVector, spacePointVec, *spEdgeAssns, instance);
1608  }
1609 
1610  output.makePFPartEdgeAssns(*edgeVector, *pfPartEdgeAssns, edgeStart, instance);
1611  }
1612 
1613  // Look at making the PCAxis and associations - for both the skeleton (the first) and the full
1614  // First need some float to double conversion containers
1615  recob::PCAxis::EigenVectors eigenVecs;
1616  double eigenVals[] = {0., 0., 0.};
1617  double avePosition[] = {0., 0., 0.};
1618 
1619  eigenVecs.resize(3);
1620 
1621  for (size_t outerIdx = 0; outerIdx < 3; outerIdx++) {
1622  avePosition[outerIdx] = skeletonPCA.getAvePosition()[outerIdx];
1623  eigenVals[outerIdx] = skeletonPCA.getEigenValues()[outerIdx];
1624 
1625  eigenVecs[outerIdx].resize(3);
1626 
1627  for (size_t innerIdx = 0; innerIdx < 3; innerIdx++)
1628  eigenVecs[outerIdx][innerIdx] = skeletonPCA.getEigenVectors().row(outerIdx)(innerIdx);
1629  }
1630 
1631  recob::PCAxis skelPcAxis(skeletonPCA.getSvdOK(),
1632  skeletonPCA.getNumHitsUsed(),
1633  eigenVals,
1634  eigenVecs,
1635  avePosition,
1636  skeletonPCA.getAveHitDoca(),
1637  output.artPCAxisVector->size());
1638 
1639  output.artPCAxisVector->push_back(skelPcAxis);
1640 
1641  for (size_t outerIdx = 0; outerIdx < 3; outerIdx++) {
1642  avePosition[outerIdx] = fullPCA.getAvePosition()[outerIdx];
1643  eigenVals[outerIdx] = fullPCA.getEigenValues()[outerIdx];
1644 
1645  for (size_t innerIdx = 0; innerIdx < 3; innerIdx++)
1646  eigenVecs[outerIdx][innerIdx] = fullPCA.getEigenVectors().row(outerIdx)(innerIdx);
1647  }
1648 
1649  recob::PCAxis fullPcAxis(fullPCA.getSvdOK(),
1650  fullPCA.getNumHitsUsed(),
1651  eigenVals, //fullPCA.getEigenValues(),
1652  eigenVecs, //fullPCA.getEigenVectors(),
1653  avePosition, //fullPCA.getAvePosition(),
1654  fullPCA.getAveHitDoca(),
1655  output.artPCAxisVector->size());
1656 
1657  output.artPCAxisVector->push_back(fullPcAxis);
1658 
1659  // Create associations to the PFParticle
1660  output.makePFPartPCAAssns();
1661  output.makePFPartSeedAssns(numSeedsStart);
1662  output.makePFPartClusterAssns(clusterStart);
1663 
1664  return pfParticleIdx;
1665  }
1666 
1667  void
1669  std::vector<recob::SpacePoint>& spacePointVec,
1670  art::Assns<recob::Hit, recob::SpacePoint>& spHitAssns,
1671  reco::HitPairListPtr& clusHitPairVector,
1672  IHit3DBuilder::RecobHitToPtrMap& hitToPtrMap,
1673  Hit3DToSPPtrMap& hit3DToSPPtrMap,
1674  const std::string& instance) const
1675  {
1676  // Reserve space...
1677  spacePointVec.reserve(spacePointVec.size() + clusHitPairVector.size());
1678 
1679  // Right now error matrix is uniform...
1680  double spError[] = {1., 0., 1., 0., 0., 1.};
1681 
1682  // Copy these hits to the vector to be stored with the event
1683  for (auto& hitPair : clusHitPairVector) {
1684  // Skip those space points that have already been created
1685  if (hit3DToSPPtrMap.find(hitPair) != hit3DToSPPtrMap.end()) continue;
1686 
1687  // Don't make space point if this hit was "rejected"
1688  if (hitPair->bitsAreSet(reco::ClusterHit3D::REJECTEDHIT)) continue;
1689 
1690  double chisq = hitPair->getHitChiSquare(); // secret handshake...
1691 
1692  // Mark this hit pair as in use
1693  hitPair->setStatusBit(reco::ClusterHit3D::MADESPACEPOINT);
1694 
1695  // Create and store the space point
1696  size_t spacePointID = spacePointVec.size();
1697  double spacePointPos[] = {
1698  hitPair->getPosition()[0], hitPair->getPosition()[1], hitPair->getPosition()[2]};
1699 
1700  spError[1] = hitPair->getTotalCharge();
1701  spError[3] = hitPair->getChargeAsymmetry();
1702 
1703  spacePointVec.emplace_back(spacePointPos, spError, chisq, spacePointID);
1704 
1705  // Update mapping
1706  hit3DToSPPtrMap[hitPair] = output.makeSpacePointPtr(spacePointID, instance);
1707 
1708  // space point hits associations
1709  RecobHitVector recobHits;
1710 
1711  for (const auto& hit : hitPair->getHits()) {
1712  if (!hit) continue;
1713 
1714  art::Ptr<recob::Hit> hitPtr = hitToPtrMap[hit->getHit()];
1715  recobHits.push_back(hitPtr);
1716  }
1717 
1718  if (!recobHits.empty())
1719  output.makeSpacePointHitAssns(spacePointVec, recobHits, spHitAssns, instance);
1720  }
1721 
1722  return;
1723  }
1724 
1725  void
1727  reco::ConvexHullKinkTupleList& kinkTupleVec) const
1728  {
1729  // Right now error matrix is uniform...
1730  double spError[] = {1., 0., 1., 0., 0., 1.};
1731 
1732  // Copy these hits to the vector to be stored with the event
1733  for (auto& kinkTuple : kinkTupleVec) {
1734  const reco::ClusterHit3D* hit = std::get<2>(std::get<0>(kinkTuple));
1735 
1736  double chisq = hit->getHitChiSquare(); // secret handshake...
1737 
1738  // Create and store the space point
1739  double spacePointPos[] = {
1740  hit->getPosition()[0], hit->getPosition()[1], hit->getPosition()[2]};
1741 
1742  spError[1] = hit->getTotalCharge();
1743  spError[3] = hit->getChargeAsymmetry();
1744 
1745  output.artExtremePointVector->push_back(
1746  recob::SpacePoint(spacePointPos, spError, chisq, output.artExtremePointVector->size()));
1747  }
1748 
1749  return;
1750  }
1751 
1752  void
1754  dcel2d::VertexList& vertexList,
1755  dcel2d::HalfEdgeList& halfEdgeList) const
1756  {
1757  // We actually do two things here:
1758  // 1) Create space points to represent the vertex locations of the voronoi diagram
1759  // 2) Create the edges that link the space points together
1760 
1761  // Set up the space point creation
1762  // Right now error matrix is uniform...
1763  double spError[] = {1., 0., 1., 0., 0., 1.};
1764  double chisq = 1.;
1765 
1766  // Keep track of the vertex to space point association
1767  std::map<const dcel2d::Vertex*, size_t> vertexToSpacePointMap;
1768 
1769  // Copy these hits to the vector to be stored with the event
1770  for (auto& vertex : vertexList) {
1771  const dcel2d::Coords& coords = vertex.getCoords();
1772 
1773  // Create and store the space point
1774  double spacePointPos[] = {coords[0], coords[1], coords[2]};
1775 
1776  vertexToSpacePointMap[&vertex] = output.artVertexPointVector->size();
1777 
1778  output.artVertexPointVector->emplace_back(
1779  spacePointPos, spError, chisq, output.artVertexPointVector->size());
1780  }
1781 
1782  // Try to avoid double counting
1783  std::set<const dcel2d::HalfEdge*> halfEdgeSet;
1784 
1785  // Build the edges now
1786  for (const auto& halfEdge : halfEdgeList) {
1787  // Recover twin
1788  const dcel2d::HalfEdge* twin = halfEdge.getTwinHalfEdge();
1789 
1790  // It can happen that we have no twin... and also check that we've not been here before
1791  if (twin && halfEdgeSet.find(twin) == halfEdgeSet.end()) {
1792  // Recover the vertices
1793  const dcel2d::Vertex* fromVertex = twin->getTargetVertex();
1794  const dcel2d::Vertex* toVertex = halfEdge.getTargetVertex();
1795 
1796  // It can happen for the open edges that there is no target vertex
1797  if (!toVertex || !fromVertex) continue;
1798 
1799  if (vertexToSpacePointMap.find(fromVertex) == vertexToSpacePointMap.end() ||
1800  vertexToSpacePointMap.find(toVertex) == vertexToSpacePointMap.end())
1801  continue;
1802 
1803  // Need the distance between vertices
1804  Eigen::Vector3f distVec = toVertex->getCoords() - fromVertex->getCoords();
1805 
1806  output.artVertexEdgeVector->emplace_back(distVec.norm(),
1807  vertexToSpacePointMap[fromVertex],
1808  vertexToSpacePointMap[toVertex],
1809  output.artEdgeVector->size());
1810 
1811  halfEdgeSet.insert(&halfEdge);
1812  }
1813  }
1814 
1815  return;
1816  }
1817 
1818  void
1820  const reco::PrincipalComponents& fullPCA,
1821  IdxToPCAMap& idxToPCAMap) const
1822  {
1823  // We actually do two things here:
1824  // 1) Create space points from the centroids of the PCA for each cluster
1825  // 2) Create the edges that link the space points together
1826 
1827  // The first task is to put the list of PCA's into some semblance of order... they may be
1828  // preordered by likely they are piecewise ordered so fix that here
1829 
1830  // We'll need the current PCA axis to determine doca and arclen
1831  Eigen::Vector3f avePosition(
1832  fullPCA.getAvePosition()[0], fullPCA.getAvePosition()[1], fullPCA.getAvePosition()[2]);
1833  Eigen::Vector3f axisDirVec(fullPCA.getEigenVectors().row(2));
1834 
1835  using DocaToPCAPair = std::pair<float, const reco::PrincipalComponents*>;
1836  using DocaToPCAVec = std::vector<DocaToPCAPair>;
1837 
1838  DocaToPCAVec docaToPCAVec;
1839 
1840  // Outer loop over views
1841  for (const auto& idxToPCA : idxToPCAMap) {
1842  const reco::PrincipalComponents* pca = idxToPCA.second;
1843 
1844  // Now we need to calculate the doca and poca...
1845  // Start by getting this hits position
1846  Eigen::Vector3f pcaPos(
1847  pca->getAvePosition()[0], pca->getAvePosition()[1], pca->getAvePosition()[2]);
1848 
1849  // Form a TVector from this to the cluster average position
1850  Eigen::Vector3f pcaToAveVec = pcaPos - avePosition;
1851 
1852  // With this we can get the arclength to the doca point
1853  float arclenToPoca = pcaToAveVec.dot(axisDirVec);
1854 
1855  docaToPCAVec.emplace_back(arclenToPoca, pca);
1856  }
1857 
1858  std::sort(docaToPCAVec.begin(), docaToPCAVec.end(), [](const auto& left, const auto& right) {
1859  return left.first < right.first;
1860  });
1861 
1862  // Set up the space point creation
1863  // Right now error matrix is uniform...
1864  double spError[] = {1., 0., 1., 0., 0., 1.};
1865  double chisq = 1.;
1866 
1867  const reco::PrincipalComponents* lastPCA(nullptr);
1868 
1869  // Set up to loop through the clusters
1870  for (const auto& docaToPCAPair : docaToPCAVec) {
1871  // Recover the PCA for this cluster
1872  const reco::PrincipalComponents* curPCA = docaToPCAPair.second;
1873 
1874  if (lastPCA) {
1875  double lastPointPos[] = {
1876  lastPCA->getAvePosition()[0], lastPCA->getAvePosition()[1], lastPCA->getAvePosition()[2]};
1877  size_t lastPointBin = output.artVertexPointVector->size();
1878  double curPointPos[] = {
1879  curPCA->getAvePosition()[0], curPCA->getAvePosition()[1], curPCA->getAvePosition()[2]};
1880  size_t curPointBin = lastPointBin + 1;
1881 
1882  output.artVertexPointVector->emplace_back(lastPointPos, spError, chisq, lastPointBin);
1883  output.artVertexPointVector->emplace_back(curPointPos, spError, chisq, curPointBin);
1884 
1885  Eigen::Vector3f distVec(curPointPos[0] - lastPointPos[0],
1886  curPointPos[1] - lastPointPos[1],
1887  curPointPos[2] - lastPointPos[2]);
1888 
1889  output.artVertexEdgeVector->emplace_back(
1890  distVec.norm(), lastPointBin, curPointBin, output.artEdgeVector->size());
1891  }
1892 
1893  lastPCA = curPCA;
1894  }
1895  }
1896 
1897 } // namespace lar_cluster3d
process_name vertex
Definition: cheaterreco.fcl:51
reco::HitPairListPtr & getBestHitPairListPtr()
Definition: Cluster3D.h:479
std::unique_ptr< std::vector< recob::Edge > > artPathEdgeVector
float m_buildNeighborhoodTime
Keeps track of time to build epsilon neighborhood.
std::unique_ptr< art::Assns< recob::SpacePoint, recob::PFParticle > > artPFPartSPAssociations
double std(const std::vector< short > &wf, const double ped_mean, size_t start, size_t nsample)
Definition: UtilFunc.cxx:42
void PCAAnalysis_calc3DDocas(const reco::HitPairListPtr &hitPairVector, const reco::PrincipalComponents &pca) const
std::list< reco::ClusterHit3D > HitPairList
Definition: Cluster3D.h:338
std::unique_ptr< std::vector< recob::SpacePoint > > artVertexPointVector
std::unique_ptr< lar_cluster3d::IClusterModAlg > m_clusterPathAlg
Algorithm to do cluster path finding.
bool getSvdOK() const
Definition: Cluster3D.h:243
void makeEdgeSpacePointAssns(std::vector< recob::Edge > &edgeVector, RecobSpacePointVector &spacePointVector, art::Assns< recob::SpacePoint, recob::Edge > &edgeSPAssociations, const std::string &path="")
Utilities related to art service access.
An object to define a &quot;edge&quot; which is used to connect space points in a triangulation algorithm...
void PCAAnalysis_3D(const reco::HitPairListPtr &hitPairList, reco::PrincipalComponents &pca, bool skeletonOnly=false) const
Class managing the creation of a new recob::Cluster object.
std::unique_ptr< std::vector< recob::Seed > > artSeedVector
static constexpr Sample_t transform(Sample_t sample)
process_name cluster
Definition: cheaterreco.fcl:51
std::string m_vertexInstance
Special instance name for vertex points.
static constexpr size_t kPFParticlePrimary
Define index to signify primary particle.
Definition: PFParticle.h:61
std::list< HalfEdge > HalfEdgeList
Definition: DCEL.h:184
void makePFPartSpacePointAssns(std::vector< recob::SpacePoint > &spacePointVector, art::Assns< recob::SpacePoint, recob::PFParticle > &pfPartSPAssociations, size_t spacePointStart, const std::string &instance="")
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:157
std::unique_ptr< art::Assns< recob::Seed, recob::Hit > > artSeedHitAssociations
float m_parallelHitsTransWid
Cut on transverse width of cluster (PCA 2nd eigenvalue)
HalfEdge * getTwinHalfEdge() const
Definition: DCEL.h:163
Unknown view.
Definition: geo_types.h:136
Cluster3D class.
Definition: SkeletonAlg.h:27
float m_parallelHitsCosAng
Cut for PCA 3rd axis angle to X axis.
bool aParallelHitsCluster(const reco::PrincipalComponents &pca) const
There are several places we will want to know if a candidate cluster is a &quot;parallel hits&quot; type of clu...
Declaration of signal hit object.
walls no right
Definition: selectors.fcl:105
reco::PrincipalComponents & getSkeletonPCA()
Definition: Cluster3D.h:477
void clearStatusBits(unsigned bits) const
Definition: Cluster3D.h:180
void findTrackSeeds(art::Event &evt, reco::ClusterParameters &cluster, IHit3DBuilder::RecobHitToPtrMap &hitToPtrMap, std::vector< recob::Seed > &seedVec, art::Assns< recob::Seed, recob::Hit > &seedHitAssns) const
An interface to the seed finding algorithm.
void makeClusterHitAssns(RecobHitVector &recobHits)
SkeletonAlg m_skeletonAlg
Skeleton point finder.
void makePFPartEdgeAssns(std::vector< recob::Edge > &edgeVector, art::Assns< recob::Edge, recob::PFParticle > &pfPartEdgeAssociations, size_t edgeStart, const std::string &instance="")
std::unique_ptr< std::vector< recob::Edge > > artVertexEdgeVector
const std::string instance
int m_hits3D
Keeps track of the number of 3D hits made.
float getTotalCharge() const
Definition: Cluster3D.h:161
Hit has been rejected for any reason.
Definition: Cluster3D.h:98
reco::EdgeList & getBestEdgeList()
Definition: Cluster3D.h:480
int getNumHitsUsed() const
Definition: Cluster3D.h:244
HoughSeedFinderAlg class.
std::list< HitPairListPtr > HitPairListPtrList
Definition: Cluster3D.h:336
A utility class used in construction of 3D clusters.
Definition: Cluster3D.h:299
reco::HitPairListPtr & getHitPairListPtr()
Definition: Cluster3D.h:475
void ProduceArtClusters(util::GeometryUtilities const &gser, ArtOutputHandler &output, reco::HitPairList &hitPairList, reco::ClusterParametersList &clusterParametersList, IHit3DBuilder::RecobHitToPtrMap &hitToPtrMap) const
Top level output routine, allows checking cluster status.
bool m_onlyMakSpacePoints
If true we don&#39;t do the full cluster 3D processing.
void makePFPartSeedAssns(size_t numSeedsStart)
process_name hit
Definition: cheaterreco.fcl:51
reco::PlaneToClusterParamsMap & getClusterParams()
Definition: Cluster3D.h:473
HoughSeedFinderAlg m_seedFinderAlg
Seed finder.
std::unique_ptr< lar_cluster3d::IClusterAlg > m_clusterAlg
Algorithm to do 3D space point clustering.
std::unique_ptr< art::Assns< recob::PFParticle, recob::Seed > > artPFPartSeedAssociations
float m_makeHitsTime
Keeps track of time to build 3D hits.
std::unique_ptr< std::vector< recob::PFParticle > > artPFParticleVector
void makeSpacePointHitAssns(std::vector< recob::SpacePoint > &spacePointVector, RecobHitVector &recobHits, art::Assns< recob::Hit, recob::SpacePoint > &spHitAssns, const std::string &path="")
void GetSkeletonHits(const reco::HitPairListPtr &inputHitList, reco::HitPairListPtr &skeletonHitList) const
Return the skeleton hits from the input list.
art::Ptr< recob::Edge > makeEdgePtr(size_t index, const std::string &instance="")
bool findTrackSeeds(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, SeedHitPairListPairVec &seedHitMap) const override
Given the list of hits this will search for candidate Seed objects and return them.
This is an algorithm for finding recob::Seed objects in 3D clusters.
float getDocaToAxis() const
Definition: Cluster3D.h:168
static const SentryArgument_t Sentry
An instance of the sentry object.
Definition: Cluster.h:182
float m_dbscanTime
Keeps track of time to run DBScan.
Hit has been used in Cluster Splitting MST.
Definition: Cluster3D.h:111
void ImportHits(util::GeometryUtilities const &gser, Iter begin, Iter end)
Calls SetHits() with the hits in the sequence.
BEGIN_PROLOG triggeremu_data_config_icarus settings PMTADCthresholds sequence::icarus_stage0_multiTPC_TPC physics sequence::icarus_stage0_EastHits_TPC physics sequence::icarus_stage0_WestHits_TPC physics producers purityana0 caloskimCalorimetryCryoE physics caloskimCalorimetryCryoW physics path
ClusterParametersList & daughterList()
Definition: Cluster3D.h:449
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
std::unique_ptr< art::Assns< recob::Hit, recob::SpacePoint > > artPPHitAssociations
Overrides another ClusterParamsAlgBase class with selected constants.
This is an algorithm for finding recob::Seed objects in 3D clusters.
const EigenValues & getEigenValues() const
Definition: Cluster3D.h:245
art::PtrMaker< recob::Edge > fEdgePtrMaker
double distance(geo::Point_t const &point, CathodeDesc_t const &cathode)
Returns the distance of a point from the cathode.
std::string m_extremeInstance
Instance name for the extreme points.
reco::PrincipalComponents & getFullPCA()
Definition: Cluster3D.h:476
size_t countUltimateDaughters(reco::ClusterParameters &clusterParameters) const
Count number of end of line daughters.
float getChargeAsymmetry() const
Definition: Cluster3D.h:167
art::PtrVector< recob::Hit > RecobHitVector
art::PtrMaker< recob::SpacePoint > fSPPtrMakerPath
This is an algorithm for finding recob::Seed objects in 3D clusters.
void splitClustersWithMST(reco::ClusterParameters &clusterParameters, reco::ClusterParametersList &clusterParametersList) const
Attempt to split clusters by using a minimum spanning tree.
std::unique_ptr< std::vector< recob::SpacePoint > > artPathPointVector
void InitializeMonitoring()
Initialize the internal monitoring.
float m_finishTime
Keeps track of time to run output module.
Helper functions to create a cluster.
Wrapper for ClusterParamsAlgBase objects to accept diverse input.
std::unique_ptr< std::vector< recob::SpacePoint > > artSpacePointVector
float m_clusterMergeTime
Keeps track of the time to merge clusters.
const Eigen::Vector3f & getAvePosition() const
Definition: Cluster3D.h:247
std::unique_ptr< art::Assns< recob::PFParticle, recob::PCAxis > > artPFPartAxisAssociations
std::unique_ptr< art::Assns< recob::Cluster, recob::Hit > > artClusterAssociations
PCASeedFinderAlg m_pcaSeedFinderAlg
Use PCA axis to find seeds.
ArtOutputHandler(art::Event &evt, std::string &pathName, std::string &vertexName, std::string &extremeName)
std::list< ConvexHullKinkTuple > ConvexHullKinkTupleList
Definition: Cluster3D.h:354
void AverageSkeletonPositions(reco::HitPairListPtr &skeletonHitList) const
Modifies the position of input skeleton hits by averaging along the &quot;best&quot; wire direction.
Wrapper for ClusterParamsAlgBase objects to accept arbitrary input.
ParallelHitsSeedFinderAlg m_parallelHitsAlg
Deal with parallel hits clusters.
std::unique_ptr< std::vector< recob::Edge > > artEdgeVector
void MakeAndSavePCAPoints(ArtOutputHandler &, const reco::PrincipalComponents &, IdxToPCAMap &) const
std::list< const reco::ClusterHit3D * > HitPairListPtr
Definition: Cluster3D.h:334
walls no left
Definition: selectors.fcl:105
std::unique_ptr< art::Assns< recob::Hit, recob::SpacePoint > > artSPHitAssociations
This provides an art tool interface definition for 3D Cluster algorithms.
std::unique_ptr< art::Assns< recob::SpacePoint, recob::Edge > > artEdgePPAssociations
int m_hits
Keeps track of the number of hits seen.
void produce(art::Event &evt) override
Declaration of cluster object.
size_t FindAndStoreDaughters(util::GeometryUtilities const &gser, ArtOutputHandler &output, reco::ClusterParameters &clusterParameters, size_t pfParticleParent, IdxToPCAMap &idxToPCAMap, IHit3DBuilder::RecobHitToPtrMap &hitToPtrMap, Hit3DToSPPtrMap &hit3DToSPPtrMap, Hit3DToSPPtrMap &best3DToSPPtrMap) const
This will produce art output for daughters noting that it needs to be done recursively.
bool operator()(const std::pair< float, const reco::ClusterHit3D * > &left, const std::pair< float, const reco::ClusterHit3D * > &right)
void MakeAndSaveSpacePoints(ArtOutputHandler &output, std::vector< recob::SpacePoint > &spacePointVec, art::Assns< recob::Hit, recob::SpacePoint > &spHitAssns, reco::HitPairListPtr &clusHitPairVector, IHit3DBuilder::RecobHitToPtrMap &hitToPtrMap, Hit3DToSPPtrMap &hit3DToSPPtrMap, const std::string &path="") const
Special routine to handle creating and saving space points.
This provides an art tool interface definition for tools which construct 3D hits used in 3D clusterin...
This header file defines the interface to a principal components analysis designed to be used within ...
std::unique_ptr< art::Assns< recob::SpacePoint, recob::PFParticle > > artPFPartPPAssociations
std::unique_ptr< art::Assns< recob::PFParticle, recob::Cluster > > artPFPartClusAssociations
float m_pathFindingTime
Keeps track of the path finding time.
bool findTrackHits(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, reco::HitPairListPtrList &hitPairListPtrList) const
Given the list of hits this will return the sets of hits which belong on the same line...
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.
bool operator()(const reco::ClusterHit3D *hit3D) const
art::Ptr< recob::SpacePoint > makeSpacePointPtr(size_t index, const std::string &instance="")
std::unique_ptr< lar_cluster3d::IClusterParametersBuilder > m_clusterBuilder
Common cluster builder tool.
Cluster3D(fhicl::ParameterSet const &pset)
art::PtrMaker< recob::SpacePoint > fSPPtrMaker
std::unique_ptr< std::vector< recob::SpacePoint > > artExtremePointVector
Hierarchical representation of particle flow.
Definition: PFParticle.h:44
bool bitsAreSet(const unsigned int &bitsToCheck) const
Definition: Cluster3D.h:176
const float getAveHitDoca() const
Definition: Cluster3D.h:248
bool findTrackSeeds(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, SeedHitPairListPairVec &seedHitMap) const override
Given the list of hits this will search for candidate Seed objects and return them.
void makePFPartClusterAssns(size_t clusterStart)
Header file to define the interface to the SkeletonAlg.
Hit has been made into Space Point.
Definition: Cluster3D.h:102
std::vector< SeedHitPairListPair > SeedHitPairListPairVec
This provides an art tool interface definition for 3D Cluster algorithms.
void MakeAndSaveVertexPoints(ArtOutputHandler &, dcel2d::VertexList &, dcel2d::HalfEdgeList &) const
Special routine to handle creating and saving space points &amp; edges associated to voronoi diagrams...
std::unique_ptr< std::vector< recob::Cluster > > artClusterVector
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG simSlidingORM6O6 effSlidingORW output
std::unique_ptr< lar_cluster3d::IClusterModAlg > m_clusterMergeAlg
Algorithm to do cluster merging.
size_t ConvertToArtOutput(util::GeometryUtilities const &gser, ArtOutputHandler &output, reco::ClusterParameters &clusterParameters, size_t pfParticleParent, IHit3DBuilder::RecobHitToPtrMap &hitToPtrMap, Hit3DToSPPtrMap &hit3DToSPPtrMap, Hit3DToSPPtrMap &best3DToSPPtrMap) const
Produces the art output from all the work done in this producer module.
std::unique_ptr< lar_cluster3d::IHit3DBuilder > m_hit3DBuilderAlg
Builds the 3D hits to operate on.
T copy(T const &v)
PrincipalComponentsAlg m_pcaAlg
Principal Components algorithm.
std::unordered_map< const reco::ClusterHit2D *, ClusterToHitPairSetMap > Hit2DToClusterMap
Definition: Cluster3D.h:510
art::PtrMaker< recob::Edge > fEdgePtrMakerPath
ClusterHit2DVec m_hitVector
Definition: Cluster3D.h:326
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
art::ServiceHandle< art::TFileService > tfs
Algorithm collection class computing cluster parameters.
float getHitChiSquare() const
Definition: Cluster3D.h:165
Definition of the Cluster3D class.
float m_artHitsTime
Keeps track of time to recover hits.
Interface to class computing cluster parameters.
std::unique_ptr< std::vector< recob::PCAxis > > artPCAxisVector
TCEvent evt
Definition: DataStructs.cxx:8
Vertex * getTargetVertex() const
Definition: DCEL.h:161
Eigen::Vector3f Coords
Definition: DCEL.h:46
std::unique_ptr< art::Assns< recob::SpacePoint, recob::Edge > > artEdgeSPAssociations
std::unordered_map< const reco::ClusterHit3D *, art::Ptr< recob::SpacePoint >> Hit3DToSPPtrMap
void MakeAndSaveKinkPoints(ArtOutputHandler &output, reco::ConvexHullKinkTupleList &clusHitPairVector) const
Special routine to handle creating and saving space points.
PCASeedFinderAlg class.
std::vector< std::vector< double > > EigenVectors
Definition: PCAxis.h:29
bool empty(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:555
std::unordered_map< const recob::Hit *, art::Ptr< recob::Hit >> RecobHitToPtrMap
Defines a structure mapping art representation to internal.
Definition: IHit3DBuilder.h:62
void PrepareEvent(const art::Event &evt)
Event Preparation.
const Coords & getCoords() const
Definition: DCEL.h:69
std::list< Vertex > VertexList
Definition: DCEL.h:182
void splitClustersWithHough(reco::ClusterParameters &clusterParameters, reco::ClusterParametersList &clusterParametersList) const
Attempt to split clusters using the output of the Hough Filter.
std::string m_pathInstance
Special instance for path points.
art framework interface to geometry description
BEGIN_PROLOG could also be cout
auto const detProp
float m_totalTime
Keeps track of total execution time.
std::list< ClusterParameters > ClusterParametersList
Definition: Cluster3D.h:403
const EigenVectors & getEigenVectors() const
Definition: Cluster3D.h:246
bool findTrackSeeds(reco::HitPairListPtr &hitPairListPtr, reco::PrincipalComponents &inputPCA, SeedHitPairListPairVec &seedHitPairVec) const override
Given the list of hits this will search for candidate Seed objects and return them.
std::unique_ptr< art::Assns< recob::Edge, recob::PFParticle > > artPFPartEdgeAssociations
std::unique_ptr< art::Assns< recob::Edge, recob::PFParticle > > artPFPartPathEdgeAssociations
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:179
std::map< size_t, const reco::PrincipalComponents * > IdxToPCAMap
Special routine to handle creating and saving space points &amp; edges PCA points.
art::PtrVector< recob::SpacePoint > RecobSpacePointVector
bool m_enableMonitoring
Turn on monitoring of this algorithm.