EmTrackClusterId3out_module.cc
Go to the documentation of this file.
1 /////////////////////////////////////////////////////////////////////////////////
2 // Class: EmTrackClusterId
3 // Module Type: producer
4 // File: EmTrackClusterId_module.cc
5 // Authors: dorota.stefan@cern.ch pplonski86@gmail.com robert.sulej@cern.ch
6 //
7 // Module applies CNN to 2D image made of deconvoluted wire waveforms in order
8 // to distinguish EM-like activity from track-like objects. New clusters of
9 // hits are produced to include also unclustered hits and tag everything in
10 // a common way.
11 // NOTE: This module uses 3-output CNN models, see EmTrackMichelClusterId for
12 // usage of 4-output models and EmTrackClusterId2out_module.cc for 2-output
13 // models.
14 //
15 /////////////////////////////////////////////////////////////////////////////////
16 
24 #include "fhiclcpp/types/Atom.h"
26 #include "fhiclcpp/types/Table.h"
28 
34 
37 
38 #include <memory>
39 
40 namespace nnet {
41 
43 public:
44 
45  // these types to be replaced with use of feature proposed in redmine #12602
46  typedef std::unordered_map< unsigned int, std::vector< size_t > > view_keymap;
47  typedef std::unordered_map< unsigned int, view_keymap > tpc_view_keymap;
48  typedef std::unordered_map< unsigned int, tpc_view_keymap > cryo_tpc_view_keymap;
49 
50  struct Config {
51  using Name = fhicl::Name;
53 
55  Name("PointIdAlg")
56  };
58  Name("BatchSize"),
59  Comment("number of samples processed in one batch")
60  };
61 
63  Name("WireLabel"),
64  Comment("tag of deconvoluted ADC on wires (recob::Wire)")
65  };
66 
68  Name("HitModuleLabel"),
69  Comment("tag of hits to be EM/track tagged")
70  };
71 
73  Name("ClusterModuleLabel"),
74  Comment("tag of clusters to be used as a source of EM/track tagged new clusters (incl. single-hit clusters ) using accumulated results from hits")
75  };
76 
78  Name("TrackModuleLabel"),
79  Comment("tag of 3D tracks to be EM/track tagged using accumulated results from hits in the best 2D projection")
80  };
81 
83  Name("Views"),
84  Comment("tag clusters in selected views only, or in all views if empty list")
85  };
86  };
88  explicit EmTrackClusterId(Parameters const & p);
89 
90  EmTrackClusterId(EmTrackClusterId const &) = delete;
92  EmTrackClusterId & operator = (EmTrackClusterId const &) = delete;
94 
95 private:
96  void produce(art::Event & e) override;
97 
98  bool isViewSelected(int view) const;
99 
100  size_t fBatchSize;
102  anab::MVAWriter<3> fMVAWriter; // <-------------- using 3-output CNN model
103 
109 
110  std::vector< int > fViews;
111 
112  art::InputTag fNewClustersTag; // input tag for the clusters produced by this module
113 };
114 // ------------------------------------------------------
115 
117  EDProducer{config},
118  fBatchSize(config().BatchSize()),
120  fMVAWriter(producesCollector(), "emtrack"),
121  fWireProducerLabel(config().WireLabel()),
122  fHitModuleLabel(config().HitModuleLabel()),
123  fClusterModuleLabel(config().ClusterModuleLabel()),
124  fTrackModuleLabel(config().TrackModuleLabel()),
125  fViews(config().Views()),
126 
128  config.get_PSet().get<std::string>("module_label"), "",
130 {
132 
133  if (!fClusterModuleLabel.label().empty())
134  {
135  produces< std::vector<recob::Cluster> >();
136  produces< art::Assns<recob::Cluster, recob::Hit> >();
137 
139  fDoClusters = true;
140  }
141  else { fDoClusters = false; }
142 
143  if (!fTrackModuleLabel.label().empty())
144  {
146  fDoTracks = true;
147  }
148  else { fDoTracks = false; }
149 }
150 // ------------------------------------------------------
151 
153 {
154  mf::LogVerbatim("EmTrackClusterId") << "next event: " << evt.run() << " / " << evt.id().event();
155 
156  auto wireHandle = evt.getValidHandle< std::vector<recob::Wire> >(fWireProducerLabel);
157 
158  unsigned int cryo, tpc, view;
159 
160  // ******************* get and sort hits ********************
161  auto hitListHandle = evt.getValidHandle< std::vector<recob::Hit> >(fHitModuleLabel);
162  std::vector< art::Ptr<recob::Hit> > hitPtrList;
163  art::fill_ptr_vector(hitPtrList, hitListHandle);
164 
166  for (auto const& h : hitPtrList)
167  {
168  view = h->WireID().Plane;
169  if (!isViewSelected(view)) continue;
170 
171  cryo = h->WireID().Cryostat;
172  tpc = h->WireID().TPC;
173 
174  hitMap[cryo][tpc][view].push_back(h.key());
175  }
176 
177  // ********************* classify hits **********************
178  auto hitID = fMVAWriter.initOutputs<recob::Hit>(fHitModuleLabel, hitPtrList.size(), fPointIdAlg.outputLabels());
179 
180  std::vector< char > hitInFA(hitPtrList.size(), 0); // tag hits in fid. area as 1, use 0 for hits close to the projectrion edges
181  for (auto const & pcryo : hitMap)
182  {
183  cryo = pcryo.first;
184  for (auto const & ptpc : pcryo.second)
185  {
186  tpc = ptpc.first;
187  for (auto const & pview : ptpc.second)
188  {
189  view = pview.first;
190  if (!isViewSelected(view)) continue; // should not happen, hits were selected
191 
192  fPointIdAlg.setWireDriftData(*wireHandle, view, tpc, cryo);
193 
194  // (1) do all hits in this plane ------------------------------------------------
195  for (size_t idx = 0; idx < pview.second.size(); idx += fBatchSize)
196  {
197  std::vector< std::pair<unsigned int, float> > points;
198  std::vector< size_t > keys;
199  for (size_t k = 0; k < fBatchSize; ++k)
200  {
201  if (idx + k >= pview.second.size()) { break; } // careful about the tail
202 
203  size_t h = pview.second[idx+k]; // h is the Ptr< recob::Hit >::key()
204  const recob::Hit & hit = *(hitPtrList[h]);
205  points.emplace_back(hit.WireID().Wire, hit.PeakTime());
206  keys.push_back(h);
207  }
208 
209  auto batch_out = fPointIdAlg.predictIdVectors(points);
210  if (points.size() != batch_out.size())
211  {
212  throw cet::exception("EmTrackClusterId") << "hits processing failed" << std::endl;
213  }
214 
215  for (size_t k = 0; k < points.size(); ++k)
216  {
217  size_t h = keys[k];
218  fMVAWriter.setOutput(hitID, h, batch_out[k]);
219  if (fPointIdAlg.isInsideFiducialRegion(points[k].first, points[k].second))
220  { hitInFA[h] = 1; }
221  }
222  } // hits done ------------------------------------------------------------------
223  }
224  }
225  }
226 
227  // (2) do clusters when hits are ready in all planes ----------------------------------------
228  if (fDoClusters)
229  {
230  // **************** prepare for new clusters ****************
231  auto clusters = std::make_unique< std::vector< recob::Cluster > >();
232  auto clu2hit = std::make_unique< art::Assns< recob::Cluster, recob::Hit > >();
233 
234  // ************** get and sort input clusters ***************
235  auto cluListHandle = evt.getValidHandle< std::vector<recob::Cluster> >(fClusterModuleLabel);
236  std::vector< art::Ptr<recob::Cluster> > cluPtrList;
237  art::fill_ptr_vector(cluPtrList, cluListHandle);
238 
240  for (auto const& c : cluPtrList)
241  {
242  view = c->Plane().Plane;
243  if (!isViewSelected(view)) continue;
244 
245  cryo = c->Plane().Cryostat;
246  tpc = c->Plane().TPC;
247 
248  cluMap[cryo][tpc][view].push_back(c.key());
249  }
250 
252 
253  unsigned int cidx = 0; // new clusters index
254  art::FindManyP< recob::Hit > hitsFromClusters(cluListHandle, evt, fClusterModuleLabel);
255  std::vector< bool > hitUsed(hitPtrList.size(), false); // tag hits used in clusters
256  for (auto const & pcryo : cluMap)
257  {
258  cryo = pcryo.first;
259  for (auto const & ptpc : pcryo.second)
260  {
261  tpc = ptpc.first;
262  for (auto const & pview : ptpc.second)
263  {
264  view = pview.first;
265  if (!isViewSelected(view)) continue; // should not happen, clusters were pre-selected
266 
267  for (size_t c : pview.second) // c is the Ptr< recob::Cluster >::key()
268  {
269  auto v = hitsFromClusters.at(c);
270  if (v.empty()) continue;
271 
272  for (auto const & hit : v)
273  {
274  if (hitUsed[hit.key()]) { mf::LogWarning("EmTrackClusterId") << "hit already used in another cluster"; }
275  hitUsed[hit.key()] = true;
276  }
277 
278  auto vout = fMVAWriter.getOutput<recob::Hit>(v,
279  [&](art::Ptr<recob::Hit> const & ptr) { return (float)hitInFA[ptr.key()]; });
280 
281  float pvalue = vout[0] / (vout[0] + vout[1]);
282  mf::LogVerbatim("EmTrackClusterId") << "cluster in tpc:" << tpc << " view:" << view
283  << " size:" << v.size() << " p:" << pvalue;
284 
285  clusters->emplace_back(
286  recob::Cluster(0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F,
287  v.size(), 0.0F, 0.0F, cidx, (geo::View_t)view, v.front()->WireID().planeID()));
288  util::CreateAssn(*this, evt, *clusters, v, *clu2hit);
289  cidx++;
290 
291  fMVAWriter.addOutput(cluID, vout); // add copy of the input cluster
292  }
293 
294  // (2b) make single-hit clusters --------------------------------------------
295  for (size_t h : hitMap[cryo][tpc][view]) // h is the Ptr< recob::Hit >::key()
296  {
297  if (hitUsed[h]) continue;
298 
299  auto vout = fMVAWriter.getOutput<recob::Hit>(h);
300  float pvalue = vout[0] / (vout[0] + vout[1]);
301 
302  mf::LogVerbatim("EmTrackClusterId") << "single hit in tpc:" << tpc << " view:" << view
303  << " wire:" << hitPtrList[h]->WireID().Wire << " drift:" << hitPtrList[h]->PeakTime() << " p:" << pvalue;
304 
305  art::PtrVector< recob::Hit > cluster_hits;
306  cluster_hits.push_back(hitPtrList[h]);
307  clusters->emplace_back(
308  recob::Cluster(0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F, 0.0F,
309  1, 0.0F, 0.0F, cidx, (geo::View_t)view, hitPtrList[h]->WireID().planeID()));
310  util::CreateAssn(*this, evt, *clusters, cluster_hits, *clu2hit);
311  cidx++;
312 
313  fMVAWriter.addOutput(cluID, vout); // add single-hit cluster tagging unclutered hit
314  }
315  mf::LogVerbatim("EmTrackClusterId") << "...produced " << cidx - pview.second.size() << " single-hit clusters.";
316  }
317  }
318  }
319 
320  evt.put(std::move(clusters));
321  evt.put(std::move(clu2hit));
322  } // all clusters done ----------------------------------------------------------------------
323 
324  // (3) do tracks when all hits in all cryo/tpc/plane are done -------------------------------
325  if (fDoTracks)
326  {
327  auto trkListHandle = evt.getValidHandle< std::vector<recob::Track> >(fTrackModuleLabel);
328  art::FindManyP< recob::Hit > hitsFromTracks(trkListHandle, evt, fTrackModuleLabel);
329  std::vector< std::vector< art::Ptr<recob::Hit> > > trkHitPtrList(trkListHandle->size());
330  for (size_t t = 0; t < trkListHandle->size(); ++t)
331  {
332  auto v = hitsFromTracks.at(t);
333  size_t nh[3] = { 0, 0, 0 };
334  for (auto const & hptr : v) { ++nh[hptr->View()]; }
335  size_t best_view = 2; // collection
336  if ((nh[0] >= nh[1]) && (nh[0] > 2 * nh[2])) best_view = 0; // ind1
337  if ((nh[1] >= nh[0]) && (nh[1] > 2 * nh[2])) best_view = 1; // ind2
338 
339  size_t k = 0;
340  while (!isViewSelected(best_view))
341  {
342  best_view = (best_view + 1) % 3;
343  if (++k > 3) { throw cet::exception("EmTrackClusterId") << "No views selected at all?" << std::endl; }
344  }
345 
346  for (auto const & hptr : v)
347  {
348  if (hptr->View() == best_view) trkHitPtrList[t].emplace_back(hptr);
349  }
350  }
351 
352  auto trkID = fMVAWriter.initOutputs<recob::Track>(fTrackModuleLabel, trkHitPtrList.size(), fPointIdAlg.outputLabels());
353  for (size_t t = 0; t < trkHitPtrList.size(); ++t) // t is the Ptr< recob::Track >::key()
354  {
355  auto vout = fMVAWriter.getOutput<recob::Hit>(trkHitPtrList[t],
356  [&](art::Ptr<recob::Hit> const & ptr) { return (float)hitInFA[ptr.key()]; });
357  fMVAWriter.setOutput(trkID, t, vout);
358  }
359  }
360  // tracks done ------------------------------------------------------------------------------
361 
362  fMVAWriter.saveOutputs(evt);
363 }
364 // ------------------------------------------------------
365 
367 {
368  if (fViews.empty()) return true;
369  else
370  {
371  for (auto k : fViews) if (k == view) { return true; }
372  return false;
373  }
374 }
375 // ------------------------------------------------------
376 
378 
379 }
bool isInsideFiducialRegion(unsigned int wire, float drift) const
Definition: PointIdAlg.cxx:328
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
std::unordered_map< unsigned int, view_keymap > tpc_view_keymap
std::unordered_map< unsigned int, std::vector< size_t > > view_keymap
fhicl::Atom< art::InputTag > TrackModuleLabel
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
std::string string
Definition: nybbler.cc:12
geo::WireID WireID() const
Definition: Hit.h:233
void setOutput(FVector_ID id, size_t key, std::array< float, N > const &values)
Definition: MVAWriter.h:175
Declaration of signal hit object.
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
void produce(art::Event &e) override
Set of hits with a 2D structure.
Definition: Cluster.h:71
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:576
bool setWireDriftData(const std::vector< recob::Wire > &wires, unsigned int plane, unsigned int tpc, unsigned int cryo)
bool isViewSelected(int view) const
std::unordered_map< unsigned int, tpc_view_keymap > cryo_tpc_view_keymap
fhicl::Atom< art::InputTag > HitModuleLabel
std::string const & label() const noexcept
Definition: InputTag.cc:76
fhicl::Atom< art::InputTag > WireLabel
std::vector< std::string > const & outputLabels(void) const
network output labels
Definition: PointIdAlg.h:122
FVector_ID initOutputs(std::string const &dataTag, size_t dataSize, std::vector< std::string > const &names=std::vector< std::string >(N,""))
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:69
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:435
static Config * config
Definition: config.cpp:1054
def move(depos, offset)
Definition: depos.py:107
fhicl::Atom< art::InputTag > ClusterModuleLabel
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
Definition: DataViewImpl.h:491
RunNumber_t run() const
Definition: DataViewImpl.cc:82
std::array< float, N > getOutput(std::vector< art::Ptr< T > > const &items) const
Definition: MVAWriter.h:189
Declaration of cluster object.
void saveOutputs(art::Event &evt)
Check consistency and save all the results in the event.
Definition: MVAWriter.h:325
Provides recob::Track data product.
Detector simulation of raw signals on wires.
ProducesCollector & producesCollector() noexcept
void addOutput(FVector_ID id, std::array< float, N > const &values)
Definition: MVAWriter.h:180
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.
p
Definition: test.py:223
AdcCodeMitigator::Name Name
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:218
#define Comment
Utility object to perform functions of association.
std::vector< std::vector< float > > predictIdVectors(std::vector< std::pair< unsigned int, float > > points) const
Definition: PointIdAlg.cxx:246
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
EventNumber_t event() const
Definition: EventID.h:116
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:48
constexpr PlaneID const & planeID() const
Definition: geo_types.h:634
TCEvent evt
Definition: DataStructs.cxx:7
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:291
ProductID put(std::unique_ptr< PROD > &&edp, FullSemantic< Level::Run > const semantic)
Definition: DataViewImpl.h:730
EventID id() const
Definition: Event.cc:37
EmTrackClusterId & operator=(EmTrackClusterId const &)=delete
Track from a non-cascading particle.A recob::Track consists of a recob::TrackTrajectory, plus additional members relevant for a "fitted" track:
Definition: Track.h:49
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
QTextStream & endl(QTextStream &s)
h
training ###############################
Definition: train_cnn.py:186