LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
DisambigAlg.cxx
Go to the documentation of this file.
1 //
3 // DisambigAlg.cxx
4 //
5 // tylerdalion@gmail.com
6 //
7 // description
8 //
9 //
10 // Since we are letting the module make the final column of hits, a disambiguated
11 // hit in this algorithm is a pair< art::Ptr<recob::Hit>, geo::WireID >
12 //
13 //
15 
16 //Framework includes:
18 #include "fhiclcpp/ParameterSet.h"
20 
27 
28 #include <cmath>
29 #include <cstdlib>
30 #include <map>
31 
32 #include "range/v3/view.hpp"
33 
34 using lar::to_element;
35 using ranges::views::filter;
36 using ranges::views::transform;
37 
38 namespace apa {
39 
41  : fAPAGeo(p.get<fhicl::ParameterSet>("APAGeometryAlg"))
42  {
43  fCrawl = p.get<bool>("Crawl");
44  fUseEndP = p.get<bool>("UseEndP");
45  fCompareViews = p.get<bool>("CompareViews");
46  fCloseHitsRadius = p.get<double>("CloseHitsRadius");
47  fMaxEndPDegRange = p.get<double>("MaxEndPDegRange");
48  fNChanJumps = p.get<unsigned int>("NChanJumps");
49  }
50 
51  //----------------------------------------------------------
52  //----------------------------------------------------------
54  detinfo::DetectorPropertiesData const& detProp,
55  art::Handle<std::vector<recob::Hit>> ChannelHits)
56  {
57  fUeffSoFar.clear();
58  fVeffSoFar.clear();
59  fnUSoFar.clear();
60  fnVSoFar.clear();
61  fnDUSoFar.clear();
62  fnDVSoFar.clear();
63  fChannelToHits.clear();
64  fAPAToUVHits.clear();
65  fAPAToZHits.clear();
66  fAPAToHits.clear();
67  fAPAToEndPHits.clear();
68  fAPAToDHits.clear();
69  fDisambigHits.clear();
70  fChanTimeToWid.clear();
71  fHasBeenDisambiged.clear();
72 
73  std::vector<art::Ptr<recob::Hit>> ChHits;
74  art::fill_ptr_vector(ChHits, ChannelHits);
75 
76  fHasBeenDisambiged.clear();
77  unsigned int skipNoise(0);
78  // Map hits by channel/APA, initialize the disambiguation status map
79  for (size_t h = 0; h < ChHits.size(); h++) {
80  art::Ptr<recob::Hit> const& hit = ChHits[h];
81 
82  // **temporary** option to skip noise hits
83  try {
84  bt_serv->HitToXYZ(clockData, hit);
85  }
86  catch (...) {
87  skipNoise++;
88  continue;
89  }
90 
91  geo::View_t view = hit->View();
92  unsigned int apa(0), cryo(0);
93  fAPAGeo.ChannelToAPA(hit->Channel(), apa, cryo);
94  fAPAToHits[apa].push_back(hit);
95  if (view == geo::kZ) {
96  fAPAToZHits[apa].push_back(hit);
97  continue;
98  }
99  else if (view == geo::kU || view == geo::kV) {
100  std::pair<double, double> ChanTime(hit->Channel() * 1., hit->PeakTime() * 1.);
101  this->fHasBeenDisambiged[apa][ChanTime] = false;
102  fChannelToHits[hit->Channel()].push_back(hit);
103  fAPAToUVHits[apa].push_back(hit);
104  }
105  }
106 
107  if (skipNoise > 0)
108  mf::LogWarning("DisambigAlg")
109  << "\nSkipped " << skipNoise << " induction noise hits using the BackTrackerService.\n"
110  << "This is only to temporarily deal with the excessive amount of noise due to the bad "
111  "deconvolution.\n";
112 
113  mf::LogVerbatim("RunDisambig") << "\n~~~~~~~~~~~ Running Disambiguation ~~~~~~~~~~~\n";
114 
115  std::map<unsigned int, std::vector<art::Ptr<recob::Hit>>>::iterator APA_it;
116  for (APA_it = fAPAToUVHits.begin(); APA_it != fAPAToUVHits.end(); APA_it++) {
117  unsigned int apa = APA_it->first;
118 
119  mf::LogVerbatim("RunDisambig") << "APA " << apa << ":";
120 
121  fUeffSoFar[apa] = 0.;
122  fVeffSoFar[apa] = 0.;
123  fnUSoFar[apa] = 0;
124  fnVSoFar[apa] = 0;
125  fnDUSoFar[apa] = 0;
126  fnDVSoFar[apa] = 0;
127 
128  // Always run this...
129  this->TrivialDisambig(clockData, detProp, apa);
130  this->AssessDisambigSoFar(apa);
131  mf::LogVerbatim("RunDisambig")
132  << " Trivial Disambig --> " << fnDUSoFar[apa] << " / " << fnUSoFar[apa] << " U, "
133  << fnDVSoFar[apa] << " / " << fnVSoFar[apa] << " V";
134 
135  // ... and pick the rest with the configurations.
136  if (fCrawl) {
137  this->Crawl(apa);
138  this->AssessDisambigSoFar(apa);
139  mf::LogVerbatim("RunDisambig")
140  << " Crawl --> " << fnDUSoFar[apa] << " / " << fnUSoFar[apa] << " U, "
141  << fnDVSoFar[apa] << " / " << fnVSoFar[apa] << " V";
142  }
143 
144  if (fUseEndP) {
145  this->FindChanTimeEndPts(detProp, apa);
146  this->UseEndPts(detProp, apa); // does the crawl from inside
147  this->AssessDisambigSoFar(apa);
148  mf::LogVerbatim("RunDisambig")
149  << " Endpoint Crawl --> " << fnDUSoFar[apa] << " / " << fnUSoFar[apa] << " U, "
150  << fnDVSoFar[apa] << " / " << fnVSoFar[apa] << " V";
151  }
152 
153  if (fCompareViews) {
154  unsigned int nDisambig(1);
155  while (nDisambig > 0) {
156  nDisambig = this->CompareViews(detProp, apa);
157  this->Crawl(apa);
158  }
159  this->AssessDisambigSoFar(apa);
160  mf::LogVerbatim("RunDisambig")
161  << " Compare Views --> " << fnDUSoFar[apa] << " / " << fnUSoFar[apa] << " U, "
162  << fnDVSoFar[apa] << " / " << fnVSoFar[apa] << " V";
163  }
164 
165  // For now just buld a simple list to get from the module
166  for (size_t i = 0; i < fAPAToDHits[apa].size(); i++)
167  fDisambigHits.push_back(fAPAToDHits[apa][i]);
168 
169  } // end loop through APA
170  }
171 
172  //-------------------------------------------------
173  //-------------------------------------------------
175  geo::WireID wid,
176  unsigned int apa)
177  {
178  std::pair<double, double> ChanTime(hit->Channel() * 1., hit->PeakTime() * 1.);
179  if (fHasBeenDisambiged[apa][ChanTime]) return;
180 
181  if (!wid.isValid) {
182  mf::LogWarning("InvalidWireID") << "wid is invalid, hit not being made\n";
183  return;
184  }
185 
186  fAPAToDHits[apa].emplace_back(hit, wid);
187  fHasBeenDisambiged[apa][ChanTime] = true;
188  fChanTimeToWid[ChanTime] = wid;
189  }
190 
191  //----------------------------------------------------------
192  //----------------------------------------------------------
194  recob::Hit const& hitA,
195  recob::Hit const& hitB)
196  {
197  double AsT = hitA.PeakTimeMinusRMS();
198  double AeT = hitA.PeakTimePlusRMS();
199  double BsT = hitB.PeakTimeMinusRMS();
200  double BeT = hitB.PeakTimePlusRMS();
201 
202  if (hitA.View() == geo::kU) {
203  AsT -= detProp.TimeOffsetU();
204  AeT -= detProp.TimeOffsetU();
205  }
206  else if (hitA.View() == geo::kV) {
207  AsT -= detProp.TimeOffsetV();
208  AeT -= detProp.TimeOffsetV();
209  }
210  else if (hitA.View() == geo::kZ) {
211  AsT -= detProp.TimeOffsetZ();
212  AeT -= detProp.TimeOffsetZ();
213  }
214 
215  if (hitB.View() == geo::kU) {
216  BsT += detProp.TimeOffsetU();
217  BeT -= detProp.TimeOffsetU();
218  }
219  else if (hitB.View() == geo::kV) {
220  BsT -= detProp.TimeOffsetV();
221  BeT -= detProp.TimeOffsetV();
222  }
223  else if (hitA.View() == geo::kZ) { // FIXME: Shouldn't this be hitB, BsT, and BeT?
224  AsT -= detProp.TimeOffsetZ();
225  AeT -= detProp.TimeOffsetZ();
226  }
227 
228  return (AsT <= BsT && BsT <= AeT) || (AsT <= BeT && BeT <= AeT) || (BsT <= AsT && AsT <= BeT) ||
229  (BsT <= AeT && AeT <= BeT);
230  }
231 
232  //----------------------------------------------------------
233  //----------------------------------------------------------
235  detinfo::DetectorPropertiesData const& detProp,
236  unsigned int apa)
237  {
238  // Loop through ambiguous hits (U/V) in this APA
239  for (auto const& hitPtr : fAPAToUVHits[apa]) {
240  auto const& hit = *hitPtr;
241  raw::ChannelID_t chan = hit.Channel();
242  unsigned int peakT = hit.PeakTime();
243 
244  std::vector<geo::WireID> hitwids = geom->ChannelToWire(chan);
245  std::vector<bool> IsReasonableWid(hitwids.size(), false);
246  unsigned short nPossibleWids(0);
247  for (size_t w = 0; w < hitwids.size(); w++) {
248  geo::WireID wid = hitwids[w];
249 
250  double xyzStart[3] = {0.};
251  double xyzEnd[3] = {0.};
252  geom->WireEndPoints(wid, xyzStart, xyzEnd);
253  unsigned int side(wid.TPC % 2), cryo(wid.Cryostat);
254  double zminPos(xyzStart[2]), zmaxPos(xyzEnd[2]);
255 
256  // get appropriate x and y with tpc center
257  unsigned int tpc =
258  2 * apa + side - cryo * geom->NTPC(); // apa number does not reset per cryo
259  auto const tpcCenter = geom->TPC(geo::TPCID(cryo, tpc)).GetCenter();
260 
261  // get channel range
262  auto Min = tpcCenter;
263  Min.SetZ(zminPos);
264  auto Max = tpcCenter;
265  Max.SetZ(zmaxPos);
266  raw::ChannelID_t ZminChan = geom->NearestChannel(Min, geo::PlaneID{cryo, tpc, 2});
267  raw::ChannelID_t ZmaxChan = geom->NearestChannel(Max, geo::PlaneID{cryo, tpc, 2});
268 
269  for (auto const& zhit : fAPAToZHits[apa] | transform(to_element)) {
270  raw::ChannelID_t chan = zhit.Channel();
271  if (chan <= ZminChan || ZmaxChan <= chan) continue;
272 
273  if (this->HitsOverlapInTime(detProp, hit, zhit)) {
274  IsReasonableWid[w] = true;
275  nPossibleWids++;
276  break;
277  }
278  }
279 
280  } // end hit chan-wid loop
281 
282  if (nPossibleWids == 0) {
283  std::vector<double> xyz;
284  try {
285  xyz = bt_serv->HitToXYZ(clockData, hit);
286  } // TEMPORARY
287  catch (...) {
288  continue;
289  }
291  mf::LogWarning("UniqueTimeSeg")
292  << "U/V hit inconsistent with Z info; peak time is " << peakT << " in APA " << apa
293  << " on channel " << hit.Channel();
294  }
295  else if (nPossibleWids == 1) {
296  for (size_t d = 0; d < hitwids.size(); d++)
297  if (IsReasonableWid[d]) this->MakeDisambigHit(hitPtr, hitwids[d], apa);
298  }
299  else if (nPossibleWids == 2) {
301  }
302 
303  } // end ambig hits loop
304  }
305 
306  //----------------------------------------------------------
307  //----------------------------------------------------------
308  unsigned int DisambigAlg::MakeCloseHits(int ext, geo::WireID Dwid, double Dmin, double Dmax)
309  {
310  // Function to look, on a channel *ext* channels away from a
311  // disambiguated hit channel, for hits with time windows touching
312  // range *Dmin to Dmax*. If found, make such a hit to have a
313  // wireID adjacent to supplied *wid*. Returns number of NEW hits
314  // made.
315 
317  geo::View_t view = geom->View(Dchan);
318  if (view == geo::kZ)
319  throw cet::exception("MakeCloseHits") << "Function not meant for non-wrapped channels.\n";
320 
321  // Account for wrapping
322  raw::ChannelID_t firstChan = fAPAGeo.FirstChannelInView(view, Dchan);
323  unsigned int ChanPerView = fAPAGeo.ChannelsInView(view);
324  int tempchan = Dchan + ext; // need sign for the set of channels starting with channel 0
325  if (tempchan < (int)firstChan) tempchan += ChanPerView;
326  if (tempchan > (int)(firstChan + ChanPerView - 1)) tempchan -= ChanPerView;
327  raw::ChannelID_t chan = (raw::ChannelID_t)(tempchan);
328 
329  // There may just be no hits
330  if (fChannelToHits.count(chan) == 0) return 0;
331 
332  // There are close channel hits, so for each
333  unsigned int apa(0), cryo(0);
334  fAPAGeo.ChannelToAPA(chan, apa, cryo);
335  unsigned int MakeCount(0);
336  for (size_t i = 0; i < fChannelToHits[chan].size(); i++) {
337  art::Ptr<recob::Hit> closeHit = fChannelToHits[chan][i];
338  double st = closeHit->PeakTimeMinusRMS();
339  double et = closeHit->PeakTimePlusRMS();
340  std::vector<geo::WireID> wids = geom->ChannelToWire(chan);
341 
342  if (!(Dmin <= st && st <= Dmax) && !(Dmin <= et && et <= Dmax)) continue;
343 
344  // Found hit with window overlapping given range,
345  // now find the only reasonable wireID.
346  for (size_t w = 0; w < wids.size(); w++) {
347  if (wids[w].TPC != Dwid.TPC) continue;
348  if ((int)(wids[w].Wire) - (int)(Dwid.Wire) != ext) continue;
349 
350  // In this case, we have a unique wireID.
351  // Check to see if it has already been made - if so, do not incriment count
352  std::pair<double, double> ChanTime(closeHit->Channel() * 1., closeHit->PeakTime() * 1.);
353  if (!fHasBeenDisambiged[apa][ChanTime]) {
354  this->MakeDisambigHit(closeHit, wids[w], apa);
355  MakeCount++;
356  //std::cout << " Close hit found on channel " << chan << ", time " << st<<"-"<<et << "... \n";
357  //std::cout << " ... giving it wireID ("<< Dwid.Cryostat <<"," << Dwid.TPC
358  // << "," << Dwid.Plane << "," << Dwid.Wire << ")\n";
359  }
360  break;
361  } // end find right wireID
362 
363  } // end loop through all hits on chan
364 
365  return MakeCount;
366  }
367 
368  //----------------------------------------------------------
369  //----------------------------------------------------------
370  void DisambigAlg::Crawl(unsigned int apa)
371  {
372 
373  std::vector<art::Ptr<recob::Hit>> hits = fAPAToUVHits[apa];
374 
375  // repeat this method until stable
376  unsigned int nExtended(1);
377  while (nExtended > 0) {
378  nExtended = 0;
379 
380  // Look for any disambiguated hit ...
381  for (size_t h = 0; h < hits.size(); h++) {
382  std::pair<double, double> ChanTime(hits[h]->Channel() * 1., hits[h]->PeakTime() * 1.);
383  if (!fHasBeenDisambiged[apa][ChanTime]) continue;
384  double stD = hits[h]->PeakTimePlusRMS(-1.);
385  double etD = hits[h]->PeakTimePlusRMS(+1.);
386  double hitWindow = etD - stD;
387  geo::WireID Dwid = fChanTimeToWid[ChanTime];
388 
389  // ... and if any neighboring-channel hits are close enough in time,
390  // extend the disambiguation to the neighboring wire.
391  unsigned int extensions = 0;
392  for (unsigned int ext = 1; ext < fNChanJumps + 1; ext++) {
394  unsigned int N(0);
395  double timeExt = hitWindow * ext;
396  N += this->MakeCloseHits((int)(-ext), Dwid, stD - 5 - timeExt, etD + 5 + timeExt);
397  N += this->MakeCloseHits((int)(ext), Dwid, stD - 5 - timeExt, etD + 5 + timeExt);
398  extensions += N;
399  }
400  nExtended += extensions;
401 
402  } // end UV hits loop
403 
404  } // end while still disambiguating***
405 
406  // *** nested while loops allow disambigauation to hop the channel-wrap-boundary
407  }
408 
409  //----------------------------------------------------------
410  //----------------------------------------------------------
412  unsigned int apa)
413  {
416 
417  double pi = 3.14159265;
418  double fMaxEndPRadRange = fMaxEndPDegRange / 180. * (2 * pi);
419 
420  for (size_t h = 0; h < fAPAToHits[apa].size(); h++) {
421  art::Ptr<recob::Hit> centhit = fAPAToHits[apa][h];
422  geo::View_t view = centhit->View();
423  unsigned int plane = 0;
424  if (view == geo::kV) { plane = 1; }
425  else if (view == geo::kZ)
426  plane = 2;
427  std::vector<double> ChanTimeCenter(2, 0.);
428  unsigned int relchan = centhit->Channel() - fAPAGeo.FirstChannelInView(centhit->Channel());
429  ChanTimeCenter[0] = relchan * geom->WirePitch(view);
430  ChanTimeCenter[1] = detProp.ConvertTicksToX(centhit->PeakTime(),
431  plane,
432  apa * 2, // tpc doesnt matter
433  centhit->WireID().Cryostat);
434  //std::vector< art::Ptr<recob::Hit> > CloseHits;
435  std::vector<std::vector<double>> CloseHitsChanTime;
436  std::vector<double> FurthestCloseChanTime(2, 0.); //double maxDist = 0.;
437  std::vector<double> ClosestChanTime(2, 0.);
438  double minDist = fCloseHitsRadius + 1.;
439  double ChanDistRange = fAPAGeo.ChannelsInView(view) * geom->WirePitch(view);
440 
441  for (size_t c = 0; c < fAPAToHits[apa].size(); c++) {
442  art::Ptr<recob::Hit> closehit = fAPAToHits[apa][c];
443  if (view != closehit->View()) continue;
444  if (view == geo::kZ && centhit->WireID().TPC != closehit->WireID().TPC) continue;
445  unsigned int plane = 0;
446  if (view == geo::kV) { plane = 1; }
447  else if (view == geo::kZ)
448  plane = 2;
449  std::vector<double> ChanTimeClose(2, 0.);
450  unsigned int relchanclose =
451  closehit->Channel() - fAPAGeo.FirstChannelInView(closehit->Channel());
452  ChanTimeClose[0] = relchanclose * geom->WirePitch(view);
453  ChanTimeClose[1] = detProp.ConvertTicksToX(closehit->PeakTime(),
454  plane,
455  apa * 2, // tpc doesnt matter
456  closehit->WireID().Cryostat);
457  if (ChanTimeClose == ChanTimeCenter) continue; // move on if the same one
458 
459  double ChanDist = ChanTimeClose[0] - ChanTimeCenter[0];
460  if (ChanDist > ChanDistRange / 2) ChanDist = ChanDistRange - ChanDist;
461 
462  double distance = std::hypot(ChanDist, ChanTimeClose[1] - ChanTimeCenter[1]);
463 
464  if (distance <= fCloseHitsRadius) CloseHitsChanTime.push_back(ChanTimeClose);
465 
466  if (distance < minDist) {
467  ClosestChanTime = ChanTimeClose;
468  minDist = distance;
469  }
470 
471  } // end close-by hit loop
472 
473  if (CloseHitsChanTime.size() < 5) continue; // quick fix, to-be improved
474 
475  double minRad(2 * pi + 1.), maxRad(0.);
476  bool CloseToNegPi(false), CloseToPosPi(false);
477  for (size_t i = 0; i < CloseHitsChanTime.size(); i++) {
478  std::vector<double> ThisChanTime(CloseHitsChanTime[i]);
479  //if(ThisChanTime==ClosestChanTime) continue;
480  double ChanDist = ThisChanTime[0] - ChanTimeCenter[0];
481  if (ChanDist > ChanDistRange / 2) ChanDist = ChanDistRange - ChanDist;
482  double hitrad = std::atan2(ThisChanTime[1] - ChanTimeCenter[1], ChanDist);
483  if (hitrad > maxRad) maxRad = hitrad;
484  if (hitrad < minRad) minRad = hitrad;
485  if (hitrad + fMaxEndPRadRange > pi)
486  CloseToPosPi = true;
487  else if (hitrad - fMaxEndPRadRange < -pi)
488  CloseToNegPi = true;
489  }
490 
491  // activity at this boundary automatically kills the test, move boundary and redo
492  if (CloseToPosPi && CloseToNegPi) {
493  for (size_t i = 0; i < CloseHitsChanTime.size(); i++) {
494  std::vector<double> ThisChanTime(CloseHitsChanTime[i]);
495  //if(ThisChanTime==ClosestChanTime) continue;
496  double ChanDist = ThisChanTime[0] - ChanTimeCenter[0];
497  if (ChanDist > ChanDistRange / 2) ChanDist = ChanDistRange - ChanDist;
498  double hitrad = std::atan2(ThisChanTime[1] - ChanTimeCenter[1], ChanDist);
499  if (hitrad > 0) hitrad = pi - hitrad; // reflec across pi/2 line
500  if (hitrad < 0) hitrad = -pi - hitrad;
501  if (hitrad > maxRad) maxRad = hitrad;
502  if (hitrad < minRad) minRad = hitrad;
503  }
504  }
505 
506  if (maxRad - minRad < fMaxEndPRadRange) fAPAToEndPHits[apa].push_back(centhit);
507 
508  } // end UV hit loop
509 
510  if (fAPAToEndPHits[apa].size() == 0) return 0;
511  mf::LogVerbatim("FindChanTimeEndPts") << " Found " << fAPAToEndPHits[apa].size()
512  << " endpoint hits in apa " << apa << std::endl;
513  for (size_t ep = 0; ep < fAPAToEndPHits[apa].size(); ep++) {
514  art::Ptr<recob::Hit> epHit = fAPAToEndPHits[apa][ep];
515  mf::LogVerbatim("FindChanTimeEndPts") << " endP on channel " << epHit->Channel()
516  << " at time " << epHit->PeakTime() << std::endl;
517  }
518 
519  return fAPAToEndPHits[apa].size();
520  }
521 
522  //----------------------------------------------------------
523  //----------------------------------------------------------
524  void DisambigAlg::UseEndPts(detinfo::DetectorPropertiesData const& detProp, unsigned int apa)
525  {
526 
528 
529  if (fAPAToEndPHits[apa].size() == 0) {
530  mf::LogVerbatim("UseEndPts") << " APA " << apa << " has no endpoints.";
531  return;
532  }
533  std::vector<art::Ptr<recob::Hit>> const& endPts = fAPAToEndPHits[apa];
534 
535  std::vector<std::vector<art::Ptr<recob::Hit>>> EndPMatch;
536  unsigned short nZendPts(0);
537 
538  auto on_z_plane = [](art::Ptr<recob::Hit> const& hit) { return hit->View() == geo::kZ; };
539  auto not_on_z_plane = [](art::Ptr<recob::Hit> const& hit) { return hit->View() != geo::kZ; };
540  for (auto const& ZHitPtr : endPts | filter(on_z_plane)) {
541  auto const& ZHit = *ZHitPtr;
542  art::Ptr<recob::Hit> Uhit = ZHitPtr;
543  art::Ptr<recob::Hit> Vhit = ZHitPtr;
544  unsigned short Umatch(0), Vmatch(0);
545  ++nZendPts;
546 
547  // look for U and V hits overlapping in time
548  for (auto const& hitPtr : endPts | filter(not_on_z_plane)) {
549  auto const& hit = *hitPtr;
550  if (not HitsOverlapInTime(detProp, ZHit, hit)) continue;
551 
552  if (hit.View() == geo::kU) {
553  Uhit = hitPtr;
554  Umatch++;
555  }
556  else if (hit.View() == geo::kV) {
557  Vhit = hitPtr;
558  Vmatch++;
559  }
560  }
561 
562  unsigned int tpc(ZHit.WireID().TPC), cryo(ZHit.WireID().Cryostat);
563  auto const tpcCenter = geom->TPC(geo::TPCID{cryo, tpc}).GetCenter();
564 
565  if (Umatch == 1 && Vmatch == 1) {
566 
567  std::vector<double> yzEndPt =
568  fAPAGeo.ThreeChanPos(Uhit->Channel(), Vhit->Channel(), ZHit.Channel());
569  double intersect[3] = {tpcCenter.X(), yzEndPt[0], yzEndPt[1]};
570 
571  geo::WireID Uwid = fAPAGeo.NearestWireIDOnChan(intersect, Uhit->Channel(), 0, tpc, cryo);
572  geo::WireID Vwid = fAPAGeo.NearestWireIDOnChan(intersect, Vhit->Channel(), 1, tpc, cryo);
573  this->MakeDisambigHit(Uhit, Uwid, apa);
574  this->MakeDisambigHit(Vhit, Vwid, apa);
575  }
576  else if (Umatch == 1 && Vmatch != 1) {
577 
578  std::vector<geo::WireIDIntersection> widIntersects;
579  fAPAGeo.APAChannelsIntersect(Uhit->Channel(), ZHit.Channel(), widIntersects);
580  if (widIntersects.size() == 0)
581  continue;
582  else if (widIntersects.size() == 1) {
583  double intersect[3] = {tpcCenter.X(), widIntersects[0].y, widIntersects[0].z};
584  geo::WireID Uwid = fAPAGeo.NearestWireIDOnChan(intersect, Uhit->Channel(), 0, tpc, cryo);
585  this->MakeDisambigHit(Uhit, Uwid, apa);
586  }
587  else {
588  for (size_t i = 0; i < widIntersects.size(); i++) {
589  // compare to V hit times, see if only one makes sense
590  }
591  }
592  }
593  else if (Umatch == 1 && Vmatch != 1) {
594 
595  std::vector<geo::WireIDIntersection> widIntersects;
596  fAPAGeo.APAChannelsIntersect(Vhit->Channel(), ZHit.Channel(), widIntersects);
597  if (widIntersects.size() == 0)
598  continue;
599  else if (widIntersects.size() == 1) {
600  double intersect[3] = {tpcCenter.X(), widIntersects[0].y, widIntersects[0].z};
601  geo::WireID Vwid = fAPAGeo.NearestWireIDOnChan(intersect, Vhit->Channel(), 0, tpc, cryo);
602  this->MakeDisambigHit(Vhit, Vwid, apa);
603  }
604  }
605  }
606 
607  if (nZendPts == 0 && endPts.size() == 2 &&
608  this->HitsOverlapInTime(detProp, *endPts[0], *endPts[1])) {
609  std::vector<geo::WireIDIntersection> widIntersects;
610  fAPAGeo.APAChannelsIntersect(endPts[0]->Channel(), endPts[1]->Channel(), widIntersects);
611  if (widIntersects.size() == 1) {
612  unsigned int cryo = endPts[0]->WireID().Cryostat;
613  unsigned int tpc = widIntersects[0].TPC;
614  auto const tpcCenter = geom->TPC(geo::TPCID(cryo, tpc)).GetCenter();
615  double intersect[3] = {tpcCenter.X(), widIntersects[0].y, widIntersects[0].z};
616  unsigned int plane0(0), plane1(0);
617  if (endPts[0]->View() == geo::kV) plane0 = 1;
618  if (endPts[1]->View() == geo::kV) plane1 = 1;
619  geo::WireID wid0 =
620  fAPAGeo.NearestWireIDOnChan(intersect, endPts[0]->Channel(), plane0, tpc, cryo);
621  this->MakeDisambigHit(endPts[0], wid0, apa);
622  geo::WireID wid1 =
623  fAPAGeo.NearestWireIDOnChan(intersect, endPts[1]->Channel(), plane1, tpc, cryo);
624  this->MakeDisambigHit(endPts[1], wid1, apa);
625  }
626  }
627 
628  this->Crawl(apa);
629  }
630 
631  //----------------------------------------------------------
632  //----------------------------------------------------------
634  {
635  unsigned int nU(0), nV(0);
636  for (size_t h = 0; h < fAPAToUVHits[apa].size(); h++) {
638  if (hit->View() == geo::kU)
639  nU++;
640  else if (hit->View() == geo::kV)
641  nV++;
642  }
643 
644  unsigned int nDU(0), nDV(0);
645  for (size_t h = 0; h < fAPAToDHits[apa].size(); h++) {
646  art::Ptr<recob::Hit> hit = fAPAToDHits[apa][h].first;
647  if (hit->View() == geo::kU)
648  nDU++;
649  else if (hit->View() == geo::kV)
650  nDV++;
651  }
652 
653  fUeffSoFar[apa] = (nDU * 1.) / (nU * 1.);
654  fVeffSoFar[apa] = (nDV * 1.) / (nV * 1.);
655  fnUSoFar[apa] = nU;
656  fnVSoFar[apa] = nV;
657  fnDUSoFar[apa] = nDU;
658  fnDVSoFar[apa] = nDV;
659  }
660 
661  //----------------------------------------------------------
662  //----------------------------------------------------------
664  unsigned int apa)
665  {
666  unsigned int nDisambiguations(0);
667 
668  // loop through all hits that are still ambiguous
669  for (auto const& ambighitPtr : fAPAToUVHits[apa]) {
670  auto const& ambighit = *ambighitPtr;
671  raw::ChannelID_t ambigchan = ambighit.Channel();
672  std::pair<double, double> ambigChanTime(ambigchan * 1., ambighit.PeakTime());
673  if (fHasBeenDisambiged[apa][ambigChanTime]) continue;
674  geo::View_t view = ambighit.View();
675  std::vector<geo::WireID> ambigwids = geom->ChannelToWire(ambigchan);
676  std::vector<unsigned int> widDcounts(ambigwids.size(), 0);
677  std::vector<unsigned int> widAcounts(ambigwids.size(), 0);
678 
679  // loop through hits in the other view which are close in time
680  for (auto const& hit : fAPAToUVHits[apa] | transform(to_element)) {
681  if (hit.View() == view || !this->HitsOverlapInTime(detProp, ambighit, hit)) continue;
682 
683  // An other-view-hit overlaps in time, see what
684  // wids of the ambiguous hit's channels it overlaps
685  raw::ChannelID_t chan = hit.Channel();
686  std::vector<geo::WireID> wids = geom->ChannelToWire(chan);
687  std::pair<double, double> ChanTime(chan * 1., hit.PeakTime());
688  geo::WireIDIntersection widIntersect; // only so we can use the function
689  if (fHasBeenDisambiged[apa][ChanTime]) {
690  for (size_t a = 0; a < ambigwids.size(); a++)
691  if (ambigwids[a].TPC == fChanTimeToWid[ChanTime].TPC &&
692  geom->WireIDsIntersect(ambigwids[a], fChanTimeToWid[ChanTime], widIntersect))
693  widDcounts[a]++;
694  }
695  else {
696  // still might be able to glean disambiguation
697  // from the ambiguous hits at this time
698  for (size_t a = 0; a < ambigwids.size(); a++)
699  for (size_t w = 0; w < wids.size(); w++)
700  if (ambigwids[a].TPC == wids[w].TPC &&
701  geom->WireIDsIntersect(ambigwids[a], wids[w], widIntersect))
702  widAcounts[a]++;
703  }
704  } // end loop through close-time hits
705 
706  // For now, just make a hit if either ambig or disambig hits
707  // unanimously intersect a single wireID
708  unsigned int Dcount(0), Acount(0);
709  for (size_t d = 0; d < widDcounts.size(); d++)
710  Dcount += widDcounts[d];
711  for (size_t a = 0; a < widAcounts.size(); a++)
712  Acount += widAcounts[a];
713  for (size_t d = 0; d < widDcounts.size(); d++) {
714  if (Dcount == widDcounts[d] && Dcount > 0 && Acount == 0) {
715  this->MakeDisambigHit(ambighitPtr, ambigwids[d], apa);
716  nDisambiguations++;
717  }
718  }
719  } // end loop through still ambiguous hits
720 
721  return nDisambiguations;
722  }
723 
724 } //end namespace apa
std::map< unsigned int, std::map< std::pair< double, double >, bool > > fHasBeenDisambiged
Convenient way to keep track of disambiguation so far.
Definition: DisambigAlg.h:86
MaybeLogger_< ELseverityLevel::ELsev_info, true > LogVerbatim
std::map< raw::ChannelID_t, std::vector< art::Ptr< recob::Hit > > > fChannelToHits
Definition: DisambigAlg.h:75
unsigned int NTPC(CryostatID const &cryoid=cryostat_zero) const
Returns the total number of TPCs in the specified cryostat.
Definition: GeometryCore.h:686
constexpr to_element_t to_element
Definition: ToElement.h:25
std::vector< WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
Encapsulate the construction of a single cyostat.
double fCloseHitsRadius
Definition: DisambigAlg.h:105
bool HitsOverlapInTime(detinfo::DetectorPropertiesData const &detProp, recob::Hit const &hitA, recob::Hit const &hitB)
std::map< unsigned int, std::vector< art::Ptr< recob::Hit > > > fAPAToZHits
Definition: DisambigAlg.h:76
std::map< unsigned int, unsigned int > fnDVSoFar
Definition: DisambigAlg.h:62
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
Planes which measure V.
Definition: geo_types.h:136
apa::APAGeometryAlg fAPAGeo
Definition: DisambigAlg.h:69
Declaration of signal hit object.
std::map< unsigned int, unsigned int > fnVSoFar
Definition: DisambigAlg.h:60
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
bool isValid
Whether this ID points to a valid element.
Definition: geo_types.h:210
bool WireIDsIntersect(WireID const &wid1, WireID const &wid2, Point_t &intersection) const
Computes the intersection between two wires.
uint32_t FirstChannelInView(geo::View_t geoview, unsigned int apa, unsigned int cryo) const
void Crawl(unsigned int apa)
Extend what we disambiguation we do have in apa.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
Planes which measure Z direction.
Definition: geo_types.h:138
geo::View_t View() const
View for the plane of the hit.
Definition: Hit.h:276
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
void RunDisambig(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, art::Handle< std::vector< recob::Hit >> GausHits)
Definition: DisambigAlg.cxx:53
std::map< unsigned int, std::vector< art::Ptr< recob::Hit > > > fAPAToEndPHits
\ todo: Channel/APA to hits can be done in a unified way
Definition: DisambigAlg.h:79
std::vector< double > HitToXYZ(detinfo::DetectorClocksData const &clockData, const recob::Hit &hit) const
bool APAChannelsIntersect(uint32_t chan1, uint32_t chan2, std::vector< geo::WireIDIntersection > &IntersectVector) const
If the channels intersect, get all intersections.
art::ServiceHandle< geo::Geometry const > geom
Definition: DisambigAlg.h:70
Framework includes.
TPCGeo const & TPC(TPCID const &tpcid=tpc_zero) const
Returns the specified TPC.
Definition: GeometryCore.h:722
void TrivialDisambig(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, unsigned int apa)
Make the easiest and safest disambiguations in apa.
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:101
Planes which measure U.
Definition: geo_types.h:135
Point_t GetCenter() const
Returns the center of the TPC volume in world coordinates [cm].
Definition: TPCGeo.h:247
unsigned int ChannelToAPA(uint32_t chan) const
Get number of the APA containing the given channel.
void hits()
Definition: readHits.C:15
void AssessDisambigSoFar(unsigned int apa)
See how much disambiguation has been done in this apa so far.
parameter set interface
std::map< unsigned int, std::vector< art::Ptr< recob::Hit > > > fAPAToUVHits
Definition: DisambigAlg.h:76
T get(std::string const &key) const
Definition: ParameterSet.h:314
std::map< unsigned int, std::vector< art::Ptr< recob::Hit > > > fAPAToHits
Definition: DisambigAlg.h:77
DisambigAlg(fhicl::ParameterSet const &pset)
Definition: DisambigAlg.cxx:40
Float_t d
Definition: plot.C:235
unsigned int ChannelsInView(geo::View_t geoview) const
std::map< unsigned int, double > fVeffSoFar
Definition: DisambigAlg.h:58
void MakeDisambigHit(art::Ptr< recob::Hit > const &hit, geo::WireID, unsigned int apa)
Makes a disambiguated hit while keeping track of what has already been disambiguated.
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
unsigned int fNChanJumps
Number of channels the crawl can jump over.
Definition: DisambigAlg.h:104
void UseEndPts(detinfo::DetectorPropertiesData const &detProp, unsigned int apa)
float PeakTimeMinusRMS(float sigmas=+1.) const
Returns a time sigmas RMS away from the peak time.
Definition: Hit.h:290
Detector simulation of raw signals on wires.
double ConvertTicksToX(double ticks, int p, int t, int c) const
raw::ChannelID_t PlaneWireToChannel(WireID const &wireid) const
Returns the ID of the TPC channel connected to the specified wire.
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
decltype(auto) get(T &&obj)
ADL-aware version of std::to_string.
Definition: StdUtils.h:120
geo::WireID NearestWireIDOnChan(const double WorldLoc[3], uint32_t chan, unsigned int const plane, unsigned int const tpc=0, unsigned int const cstat=0) const
Contains all timing reference information for the detector.
bool fCrawl
\ todo: Write function that compares hits more detailedly
Definition: DisambigAlg.h:101
View_t View(PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
constexpr T pi()
Returns the constant pi (up to 35 decimal digits of precision)
art::ServiceHandle< cheat::BackTrackerService const > bt_serv
Definition: DisambigAlg.h:72
raw::ChannelID_t NearestChannel(Point_t const &worldLoc, PlaneID const &planeid) const
Returns the ID of the channel nearest to the specified position.
std::vector< std::pair< art::Ptr< recob::Hit >, geo::WireID > > fDisambigHits
The final list of hits to pass back to be made.
Definition: DisambigAlg.h:64
double fMaxEndPDegRange
Definition: DisambigAlg.h:107
unsigned int CompareViews(detinfo::DetectorPropertiesData const &detProp, unsigned int apa)
Compare U and V to see if one says something about the other.
std::map< unsigned int, double > fUeffSoFar
Definition: DisambigAlg.h:57
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:46
float PeakTimePlusRMS(float sigmas=+1.) const
Returns a time sigmas RMS away from the peak time.
Definition: Hit.h:285
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:306
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
std::vector< double > ThreeChanPos(uint32_t u, uint32_t v, uint32_t z) const
Find the center of the 3 intersections, choose best if multiple.
std::map< unsigned int, unsigned int > fnUSoFar
Definition: DisambigAlg.h:59
Float_t w
Definition: plot.C:20
unsigned int FindChanTimeEndPts(detinfo::DetectorPropertiesData const &detProp, unsigned int apa)
Basic endpoint-hit finder per apa.
raw::ChannelID_t Channel() const
ID of the readout channel the hit was extracted from.
Definition: Hit.h:268
std::map< unsigned int, std::vector< std::pair< art::Ptr< recob::Hit >, geo::WireID > > > fAPAToDHits
Hold the disambiguations per APA.
Definition: DisambigAlg.h:80
Length_t WirePitch(PlaneID const &planeid=plane_zero) const
Returns the distance between two consecutive wires.
unsigned int MakeCloseHits(int ext, geo::WireID wid, double Dmin, double Dmax)
Having disambiguated a time range on a wireID, extend to neighboring channels.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void WireEndPoints(WireID const &wireid, double *xyzStart, double *xyzEnd) const
Fills two arrays with the coordinates of the wire end points.
Encapsulate the construction of a single detector plane.
std::map< std::pair< double, double >, geo::WireID > fChanTimeToWid
If a hit is disambiguated, map its chan and peak time to the chosen wireID.
Definition: DisambigAlg.h:84
std::map< unsigned int, unsigned int > fnDUSoFar
Definition: DisambigAlg.h:61