LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
CBAlgoTrackSeparate.cxx
Go to the documentation of this file.
1 #ifndef RECOTOOL_CBALGOTRACKSEPARATE_CXX
2 #define RECOTOOL_CBALGOTRACKSEPARATE_CXX
3 
4 #include "CBAlgoTrackSeparate.h"
5 
6 namespace cmtool {
7 
8  //----------------------------------------
10  //----------------------------------------
11  {
12 
13  //this just sets default values
14  SetVerbose(true);
15  SetDebug(true);
16 
17  //1e9 is huge; everything will be merged
18  SetMinNumHits(30);
19  SetMinAngleDiff(15.); //in degrees
20  SetMaxOpeningAngle(12.0); //in deg (parameter in rad!!)
21  SetMinLength(10.);
23  SetMaxWidth(10.);
24 
25 
26  //NOTE! Using this flag means all of the other crap
27  //(minNumHits, anglediff, openingangle, blah blah)
28  //is totally irrelevant. if we stick with this flag as the algo,
29  //we probably want to delete all of the old method
30  SetUseEP(true);
31  SetEPCutoff(0.99000);
33  _wire_2_cm = geou.WireToCm();
34  _time_2_cm = geou.TimeToCm();
35 
36  }
37 
38  //--------------------------------------------------------
39  bool CBAlgoTrackSeparate::Bool(const ::cluster::ClusterParamsAlg &cluster1,
40  const ::cluster::ClusterParamsAlg &cluster2)
41  //--------------------------------------------------------
42  {
43  //if you are using EP method for this algo:
44  if(_use_EP){
45  if(cluster1.GetParams().eigenvalue_principal > _ep_cut &&
46  cluster2.GetParams().eigenvalue_principal > _ep_cut)
47  return true;
48  }
49  //if you are using the original method for this algo:
50  else{
51 
52  //two clusters are considered un-mergable if:
53  //1) both have more than _MinNumHits
54  //2) opening angle for both < _MAxOpeningAngle
55  //3) diff. in direction of both < _MinAngleDiff
56 
57  size_t N_Hits1 = cluster1.GetHitVector().size();
58  size_t N_Hits2 = cluster2.GetHitVector().size();
59  auto start_point1 = cluster1.GetParams().start_point;
60  auto start_point2 = cluster2.GetParams().start_point;
61  double angle_2d1 = cluster1.GetParams().angle_2d;
62  double angle_2d2 = cluster2.GetParams().angle_2d;
63  double opening_angle1 = cluster1.GetParams().opening_angle;
64  double opening_angle2 = cluster2.GetParams().opening_angle;
65  Polygon2D PolyObject1 = cluster1.GetParams().PolyObject;
66  Polygon2D PolyObject2 = cluster2.GetParams().PolyObject;
67  double length1 = cluster1.GetParams().length;
68  double length2 = cluster2.GetParams().length;
69  double width1 = cluster1.GetParams().width;
70  double width2 = cluster2.GetParams().width;
71 
72  //first filter out low hits clusters
73  if ( (N_Hits1 > _MinNumHits) and
74  (N_Hits2 > _MinNumHits) ) {
75  if (_debug) {
76  std::cout << "Cluster1 Num Hits: " << N_Hits1 << std::endl;
77  std::cout << "\t Start: (" << start_point1.w << " " << start_point1.t << " )" << std::endl;
78  std::cout << "\t Opening ANgle " << opening_angle1*(360/(2*3.14)) << std::endl;
79  std::cout << "\t Angle2D: " << angle_2d1 << std::endl;
80  std::cout << "\t Length: " << length1 << std::endl;
81  std::cout << "\t Width: " << width1 << std::endl;
82  std::cout << "Cluster2 Num Hits: " << N_Hits2 << std::endl;
83  std::cout << "\t Start: (" << start_point2.w << " " << start_point2.t << " )" << std::endl;
84  std::cout << "\t Opening ANgle " << opening_angle2*(360/(2*3.14)) << std::endl;
85  std::cout << "\t Angle2D: " << angle_2d2 << std::endl;
86  std::cout << "\t Length: " << length2 << std::endl;
87  std::cout << "\t Width: " << width2 << std::endl;
88  }
89  if ( (N_Hits1 > _MinNumHits) and
90  (N_Hits2 > _MinNumHits) and
91  ( abs(angle_2d1 - angle_2d2) > _MinAngleDiff ) and
92  //( PolyObject1.Area()/N_Hits1 > _MinDensity ) and
93  //( PolyObject2.Area()/N_Hits2 > _MinDensity ) and
94  (opening_angle1 < _MaxOpeningAngle/(360/(2*3.14))) and
95  (opening_angle2 < _MaxOpeningAngle/(360/(2*3.14))) and
96  (width1 < _MaxWidth) and
97  (width2 < _MaxWidth) and
98  (length1 > _MinLength) and
99  (length2 > _MinLength) ){
100  if (_verbose) { std::cout << "*****************************************Separate with TrackSeparate!" << std::endl; }
101  return true;
102  }
103  }
104  }
105 
106  return false;
107 
108  }
109 
110  //-----------------------
112  //-----------------------
113  {
114  }
115 
116 }
117 
118 #endif
Double_t TimeToCm() const
Double_t WireToCm() const
void SetMaxWidth(float maxwidth)
virtual void Report()
Function to report what&#39;s going on per merging.
CBAlgoTrackSeparate()
Default constructor.
Class def header for a class CBAlgoTrackSeparate.
virtual bool Bool(const ::cluster::ClusterParamsAlg &cluster1, const ::cluster::ClusterParamsAlg &cluster2)
void SetVerbose(bool on)
Setter function for verbosity.
void SetMinPolyHitDensity(float mindensity)
void SetMaxOpeningAngle(float maxangle)
void SetMinLength(float minlength)
void SetMinAngleDiff(float anglesep)