LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
ClusterMergeHelper.cxx
Go to the documentation of this file.
1 //
3 // ClusterMergeHelper source
4 //
6 
7 #include "ClusterMergeHelper.h"
8 
13 
19 
20 namespace cluster {
21 
22  //####################################################################################################
24  util::GeometryUtilities const& gser,
26  //####################################################################################################
27  {
28  fInputClusters.clear();
29  fOutputClusters.clear();
30 
31  std::vector<std::vector<util::PxHit>> px_clusters(clusters.size(), std::vector<util::PxHit>());
32 
33  fInputClusters.resize(clusters.size(), std::vector<art::Ptr<recob::Hit>>());
34 
35  for (size_t cluster_index = 0; cluster_index < clusters.size(); ++cluster_index) {
36 
37  px_clusters.at(cluster_index).resize(clusters.at(cluster_index).size(), util::PxHit());
38 
39  fInputClusters.at(cluster_index).resize(clusters.at(cluster_index).size());
40 
41  for (size_t hit_index = 0; hit_index < clusters.at(cluster_index).size(); ++hit_index) {
42 
43  px_clusters.at(cluster_index).at(hit_index).plane =
44  clusters.at(cluster_index).at(hit_index)->WireID().Plane;
45  px_clusters.at(cluster_index).at(hit_index).w =
46  clusters.at(cluster_index).at(hit_index)->WireID().Wire * fGeoU.WireToCm();
47  px_clusters.at(cluster_index).at(hit_index).t =
48  clusters.at(cluster_index).at(hit_index)->PeakTime() * fGeoU.TimeToCm();
49  px_clusters.at(cluster_index).at(hit_index).charge =
50  clusters.at(cluster_index).at(hit_index)->Integral();
51 
52  fInputClusters.at(cluster_index).at(hit_index) = clusters.at(cluster_index).at(hit_index);
53  }
54  }
55 
56  SetClusters(gser, px_clusters);
57  }
58 
59  //##################################################################################################
61  const art::Event& evt,
62  const std::string& cluster_module_label)
63  //##################################################################################################
64  {
65 
67  evt.getByLabel(cluster_module_label, clusters_h);
68 
69  if (!(clusters_h.isValid()))
70 
71  throw cet::exception(__FUNCTION__)
72  << "\033[93m"
73  << " Failed to retrieve recob::Cluster with label: " << cluster_module_label.c_str()
74  << "\033[00m" << std::endl;
75 
76  std::vector<std::vector<art::Ptr<recob::Hit>>> cluster_hits_v;
77  cluster_hits_v.reserve(clusters_h->size());
78 
79  art::FindManyP<recob::Hit> hit_m(clusters_h, evt, cluster_module_label);
80 
81  for (size_t i = 0; i < clusters_h->size(); ++i)
82  cluster_hits_v.push_back(hit_m.at(i));
83 
84  SetClusters(gser, cluster_hits_v);
85  }
86 
87  //################################
89  //################################
90  {
91  if (fMgr.GetClusters().size())
92  throw cet::exception(__PRETTY_FUNCTION__)
93  << "\033[93m"
94  << "Merged cluster set not empty... Called Process() twice?"
95  << "\033[00m" << std::endl;
96 
97  fMgr.Process(gser);
98 
99  // Now create output clusters
100  auto res = fMgr.GetBookKeeper();
101 
102  std::vector<std::vector<unsigned short>> out_clusters = res.GetResult();
103  fOutputClusters.clear();
104 
105  fOutputClusters.reserve(out_clusters.size());
106 
107  for (auto const& cluster_index_v : out_clusters) {
108 
109  std::vector<art::Ptr<recob::Hit>> out_cluster;
110 
111  for (auto const& cluster_index : cluster_index_v) {
112  out_cluster.reserve(out_cluster.size() + fInputClusters.at(cluster_index).size());
113  for (auto const& hit_ptr : fInputClusters.at(cluster_index))
114  out_cluster.push_back(hit_ptr);
115  }
116  fOutputClusters.push_back(out_cluster);
117  }
118  }
119 
120  //######################################################################################################
121  const std::vector<std::vector<art::Ptr<recob::Hit>>>& ClusterMergeHelper::GetMergedClusterHits()
122  const
123  //######################################################################################################
124  {
125  if (!fOutputClusters.size())
126  throw cet::exception(__FUNCTION__)
127  << "\033[93m"
128  << "You must call Process() before calling " << __FUNCTION__ << " to retrieve result."
129  << "\033[00m" << std::endl;
130 
131  return fOutputClusters;
132  }
133 
134  //#####################################################################################
135  const std::vector<cluster::ClusterParamsAlg>& ClusterMergeHelper::GetMergedCPAN() const
136  //#####################################################################################
137  {
138  if (!fOutputClusters.size())
139  throw cet::exception(__FUNCTION__)
140  << "\033[93m"
141  << "You must call Process() before calling " << __FUNCTION__ << " to retrieve result."
142  << "\033[00m" << std::endl;
143 
144  return fMgr.GetClusters();
145  }
146 
148  art::Event& ev,
149  std::vector<recob::Cluster>& out_clusters,
151  {
152 
153  if (!fOutputClusters.size())
154 
155  throw cet::exception(__FUNCTION__)
156  << "\033[93m"
157  << "You must call Process() before calling " << __FUNCTION__ << " to retrieve result."
158  << "\033[00m" << std::endl;
159 
161 
162  // Store output
163  for (size_t out_index = 0; out_index < GetMergedCPAN().size(); ++out_index) {
164 
165  // To save typing let's just retrieve const cluster_params instance
166  const cluster_params& res = GetMergedCPAN()[out_index].GetParams();
167 
168  // this "algo" is actually parroting its cluster_params
169  LazyClusterParamsAlg algo(res);
170 
171  std::vector<art::Ptr<recob::Hit>> const& hits = GetMergedClusterHits().at(out_index);
172 
173  // the full plane needed but not a part of cluster_params...
174  // get the one from the first hit
175  geo::PlaneID plane; // invalid by default
176  if (!hits.empty()) plane = hits.front()->WireID().planeID();
177 
178  // View_t needed but not a part of cluster_params, so retrieve it here
179  geo::View_t view_id = geo->Plane(plane).View();
180 
181  // Push back a new cluster data product with parameters copied from cluster_params
182  out_clusters.emplace_back(res.start_point.w / fGeoU.WireToCm(), // start_wire
183  0., // sigma_start_wire
184  res.start_point.t / fGeoU.TimeToCm(), // start_tick
185  0., // sigma_start_tick
186  algo.StartCharge(gser).value(),
187  algo.StartAngle().value(),
188  algo.StartOpeningAngle().value(),
189  res.end_point.w / fGeoU.WireToCm(), // end_wire
190  0., // sigma_end_wire
191  res.end_point.t / fGeoU.TimeToCm(), // end_tick
192  0., // sigma_end_tick
193  algo.EndCharge(gser).value(),
194  algo.EndAngle().value(),
195  algo.EndOpeningAngle().value(),
196  algo.Integral().value(),
197  algo.IntegralStdDev().value(),
198  algo.SummedADC().value(),
199  algo.SummedADCStdDev().value(),
200  algo.NHits(),
201  algo.MultipleHitDensity(),
202  algo.Width(gser),
203  out_clusters.size(), // ID
204  view_id,
205  plane,
207 
208  util::CreateAssn(ev, out_clusters, GetMergedClusterHits().at(out_index), assns);
209  }
210  }
211 
212 } // namespace cluster
Algorithm class inheriting pre-computed results.
Algorithm class inheriting cluster parameters.
std::vector< std::vector< art::Ptr< recob::Hit > > > fInputClusters
Input clusters in terms of a vector of art::Ptr<recob::Hit> collection.
Measure_t EndAngle() override
Computes the angle of the cluster.
size_t NHits() override
Returns the number of hits in the cluster.
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
std::vector< std::vector< unsigned short > > GetResult() const
Measure_t SummedADC() override
Computes the total charge of the cluster from Hit::SummedADC()
const std::vector< cluster::ClusterParamsAlg > & GetClusters() const
A method to obtain output clusters.
Definition: CMergeManager.h:49
Cluster finding and building.
static const SentryArgument_t Sentry
An instance of the sentry object.
Definition: Cluster.h:174
const std::vector< std::vector< art::Ptr< recob::Hit > > > & GetMergedClusterHits() const
Utility method to retrieve merged clusters in terms of a vector of art::Ptr<recob::Hit> ...
Measure_t EndOpeningAngle() override
Computes the opening angle at the start or end of the cluster.
Measure_t Integral() override
Computes the total charge of the cluster from Hit::Integral()
bool isValid() const noexcept
Definition: Handle.h:203
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
PlaneGeo const & Plane(PlaneID const &planeid) const
Returns the specified wire.
View_t View() const
Which coordinate does this plane measure.
Definition: PlaneGeo.h:166
void hits()
Definition: readHits.C:15
double t
Definition: PxUtils.h:10
void AppendResult(util::GeometryUtilities const &gser, art::Event &ev, std::vector< recob::Cluster > &out_clusters, art::Assns< recob::Cluster, recob::Hit > &assns) const
Utility method to append result set to user&#39;s data product storage.
Measure_t SummedADCStdDev() override
Computes the standard deviation on the charge of the cluster hits.
std::vector< std::vector< art::Ptr< recob::Hit > > > fOutputClusters
Output clusters in terms of a vector of art::Ptr<recob::Hit> collection.
float MultipleHitDensity() override
Fraction of wires in the cluster with more than one hit.
void Process(util::GeometryUtilities const &gser)
Function to execute CMergeManager::Process()
util::PxPoint start_point
start point
Definition: ClusterParams.h:23
::util::GeometryUtilities fGeoU
GeometryUtilities.
Measure_t IntegralStdDev() override
Computes the standard deviation on the charge of the cluster hits.
double w
Definition: PxUtils.h:9
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 getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Utility object to perform functions of association.
util::PxPoint end_point
end point
Definition: ClusterParams.h:24
Measure_t StartOpeningAngle() override
Computes the opening angle at the start or end of the cluster.
void Process(util::GeometryUtilities const &gser)
A method to execute the main action, to be called per event.
Measure_t StartAngle() override
Computes the angle of the cluster.
const std::vector< cluster::ClusterParamsAlg > & GetMergedCPAN() const
Utility method to retrieve merged clusters in terms of a vector of CPAN.
Measure_t StartCharge(util::GeometryUtilities const &gser) override
Computes the charge on the first and last wire of the track.
void SetClusters(util::GeometryUtilities const &gser, const std::vector< std::vector< art::Ptr< recob::Hit >>> &clusters)
Utility method to set cluster input information to CMergeManager from LArSoft data product (vector of...
const CMergeBookKeeper & GetBookKeeper() const
A method to obtain book keeper.
Definition: CMergeManager.h:52
TCEvent evt
Definition: DataStructs.cxx:8
Measure_t EndCharge(util::GeometryUtilities const &gser) override
Computes the charge on the first and last wire of the track.
::cmtool::CMergeManager fMgr
CMergeManager instance.
Namespace collecting geometry-related classes utilities.
float Width(util::GeometryUtilities const &) override
Computes the width of the cluster.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33