LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
EMShower3D_module.cc
Go to the documentation of this file.
1 // Class: EMShower3D
3 // Module Type: producer
4 // File: EMShower3D_module.cc
5 // Authors: dorota.stefan@cern.ch robert.sulej@cern.ch
7 
15 #include "fhiclcpp/ParameterSet.h"
20 
23 
29 
31 
35 
36 #include <memory>
37 
40 
41 // ROOT includes
42 #include "TLorentzVector.h"
43 #include "TMathBase.h"
44 
45 #include <iostream>
46 #include <fstream>
47 
48 struct IniSeg
49 {
50  size_t idcl1;
51  size_t idcl2;
52  size_t idcl3;
53  size_t view1;
54  size_t view2;
55  size_t view3;
57  std::vector< art::Ptr< recob::Hit > > hits1;
58  std::vector< art::Ptr< recob::Hit > > hits2;
59  std::vector< art::Ptr< recob::Hit > > hits3;
60 };
61 
62 namespace ems
63 {
64  class EMShower3D;
65 }
66 
68 public:
69  explicit EMShower3D(fhicl::ParameterSet const & p);
70 
71  EMShower3D(EMShower3D const &) = delete;
72  EMShower3D(EMShower3D &&) = delete;
73  EMShower3D & operator = (EMShower3D const &) = delete;
74  EMShower3D & operator = (EMShower3D &&) = delete;
75 
76  void beginJob() override;
77 
78  void produce(art::Event & e) override;
79 
80  void reconfigure(fhicl::ParameterSet const& p);
81 
82 
83 private:
84  recob::Track ConvertFrom(pma::Track3D& src);
85  recob::Track ConvertFrom2(pma::Track3D& src);
86  recob::Cluster ConvertFrom(const std::vector< art::Ptr<recob::Hit> > & src);
87 
88  std::vector< ems::DirOfGamma* > CollectShower2D(art::Event const & e);
89 
90  void Link(art::Event const & e, std::vector< ems::DirOfGamma* > input);
91 
92  // Remove tracks which are too close to each other
93  void Reoptimize();
94 
95  void Make3DSeg(art::Event const & e, std::vector< ems::DirOfGamma* > pair);
96 
97  bool Validate(art::Event const & e, const pma::Track3D& src, size_t plane);
98  bool Validate(std::vector< ems::DirOfGamma* > input, size_t id1, size_t id2, size_t c1, size_t c2, size_t plane3);
99 
100  void FilterOutSmallParts(
101  double r2d,
102  const std::vector< art::Ptr<recob::Hit> >& hits_in,
103  std::vector< art::Ptr<recob::Hit> >& hits_out);
104 
105  bool GetCloseHits(
106  double r2d,
107  const std::vector< art::Ptr<recob::Hit> >& hits_in,
108  std::vector<size_t>& used,
109  std::vector< art::Ptr<recob::Hit> >& hits_out);
110 
111  bool Has(const std::vector<size_t>& v, size_t idx);
112 
113  size_t LinkCandidates(art::Event const & e, std::vector< ems::DirOfGamma* > input, size_t id);
114 
115  std::vector< IniSeg > fInisegs;
116  std::vector< IniSeg > fSeltracks;
117  std::vector< IniSeg > fPMA3D;
118 
119  std::vector< std::vector< art::Ptr<recob::Hit> > > fClusters;
120 
121  std::vector< size_t > fClustersNotUsed;
122  std::vector< size_t > fTracksNotUsed;
123 
124  unsigned int fTrkIndex; unsigned int fClIndex;
125  unsigned int fIniIndex;
126 
127  std::string fCluModuleLabel;
128  std::string fTrk3DModuleLabel;
129 
132 
134 };
135 
136 
138  : fProjectionMatchingAlg(p.get< fhicl::ParameterSet >("ProjectionMatchingAlg")),
139  fCalorimetryAlg(p.get< fhicl::ParameterSet >("CalorimetryAlg"))
140 {
141  reconfigure(p);
142 
143  produces< std::vector<recob::Track> >();
144  produces< std::vector<recob::Vertex> >();
145  produces< std::vector<recob::Cluster> >();
146  produces< std::vector<recob::SpacePoint> >();
147  produces< art::Assns<recob::Track, recob::Hit> >();
148  produces< art::Assns<recob::Track, recob::Vertex> >();
149  produces< art::Assns<recob::Cluster, recob::Hit> >();
150  produces< art::Assns<recob::Track, recob::SpacePoint> >();
151  produces< art::Assns<recob::SpacePoint, recob::Hit> >();
152  produces< art::Assns<recob::Track, recob::Cluster> >();
153 }
154 
156 {
157 }
158 
160 {
161  fCluModuleLabel = p.get< std::string >("ClustersModuleLabel");
162  //fProjectionMatchingAlg.reconfigure(p.get< fhicl::ParameterSet >("ProjectionMatchingAlg")); // use constructor-time config only
163  fTrk3DModuleLabel = p.get< std::string >("Trk3DModuleLabel");
164 
165  return;
166 }
167 
169 {
170 
171  return 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, src.size(), 0.0F, 0.0F, fClIndex, src[0]->View(), src[0]->WireID().planeID());
172 }
173 
175 {
176  auto const* geom = lar::providerFrom<geo::Geometry>();
177  auto const* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
178 
179  double avdrift = (src.front()->Point3D().X() + src.back()->Point3D().X()) * 0.5;
180  unsigned int nplanes = geom->Nplanes(src.front()->TPC(), src.front()->Cryo());
181  size_t nusedhitsmax = 0; int bestplane = -1;
182  for (unsigned int p = 0; p < nplanes; ++p)
183  {
184  unsigned int nusedP = 0;
185  fProjectionMatchingAlg.selectInitialHits(src, p, &nusedP);
186 
187  if (nusedP > nusedhitsmax)
188  {
189  nusedhitsmax = nusedP;
190  bestplane = int(p);
191  }
192  }
193 
194  std::vector< std::vector< double > > vdedx;
195  std::vector< double > dedx;
196 
197  for (unsigned int p = 0; p < nplanes; ++p)
198  {
199  unsigned int nusedP = 0;
200  double dqdxplane = fProjectionMatchingAlg.selectInitialHits(src, p, &nusedP);
201  double timeP = detprop->ConvertXToTicks(avdrift,
202  p,
203  src.front()->TPC(),
204  src.front()->Cryo());
205  double dEdxplane = fCalorimetryAlg.dEdx_AREA(dqdxplane, timeP, p);
206  dedx.push_back(dEdxplane);
207  if (int(p) == bestplane) dedx.push_back(1);
208  else dedx.push_back(0);
209  vdedx.push_back(dedx);
210  }
211 
212  std::vector< TVector3 > xyz, dircos;
213 
214  for (size_t i = 0; i < src.size(); i++)
215  {
216  xyz.push_back(src[i]->Point3D());
217 
218  if (i < src.size() - 1)
219  {
220  TVector3 dc(src[i + 1]->Point3D());
221  dc -= src[i]->Point3D();
222  dc *= 1.0 / dc.Mag();
223  dircos.push_back(dc);
224  }
225  else dircos.push_back(dircos.back());
226  }
227 
228  return recob::Track(xyz, dircos, vdedx, std::vector< double >(2, util::kBogusD), fIniIndex);
229 }
230 
232 {
233 
234  auto const* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
235  auto const* geom = lar::providerFrom<geo::Geometry>();
236 
237  double avdrift = (src.front()->Point3D().X() + src.back()->Point3D().X()) * 0.5;
238  unsigned int nplanes = geom->Nplanes(src.front()->TPC(), src.front()->Cryo());
239  size_t nusedhitsmax = 0; int bestplane = -1;
240  for (unsigned int p = 0; p < nplanes; ++p)
241  {
242  unsigned int nusedP = 0;
243  fProjectionMatchingAlg.selectInitialHits(src, p, &nusedP);
244 
245  if (nusedP > nusedhitsmax)
246  {
247  nusedhitsmax = nusedP;
248  bestplane = int(p);
249  }
250  }
251 
252  std::vector< std::vector< double > > vdedx;
253  std::vector< double > dedx;
254 
255  for (unsigned int p = 0; p < nplanes; ++p)
256  {
257  unsigned int nusedP = 0;
258  double dqdxplane = fProjectionMatchingAlg.selectInitialHits(src, p, &nusedP);
259  double timeP = detprop->ConvertXToTicks(avdrift,
260  p,
261  src.front()->TPC(),
262  src.front()->Cryo());
263  double dEdxplane = fCalorimetryAlg.dEdx_AREA(dqdxplane, timeP, p);
264  dedx.push_back(dEdxplane);
265  if (int(p) == bestplane) dedx.push_back(1);
266  else dedx.push_back(0);
267  vdedx.push_back(dedx);
268  }
269 
270  std::vector< TVector3 > xyz, dircos;
271 
272  for (size_t i = 0; i < src.size(); i++)
273  {
274  xyz.push_back(src[i]->Point3D());
275 
276  if (i < src.size() - 1)
277  {
278  TVector3 dc(src[i + 1]->Point3D());
279  dc -= src[i]->Point3D();
280  dc *= 1.0 / dc.Mag();
281  dircos.push_back(dc);
282  }
283  else dircos.push_back(dircos.back());
284  }
285 
286  return recob::Track(xyz, dircos, vdedx, std::vector< double >(2, util::kBogusD), fIniIndex);
287 }
288 
290 {
292  fSeltracks.clear();
293  fInisegs.clear();
294  fClusters.clear();
295  fPMA3D.clear();
296  fClustersNotUsed.clear();
297 
298  std::unique_ptr< std::vector< recob::Track > > tracks(new std::vector< recob::Track >);
299  std::unique_ptr< std::vector< recob::Vertex > > vertices(new std::vector< recob::Vertex >);
300  std::unique_ptr< std::vector< recob::Cluster > > clusters(new std::vector< recob::Cluster >);
301  std::unique_ptr< std::vector< recob::SpacePoint > > allsp(new std::vector< recob::SpacePoint >);
302 
303  std::unique_ptr< art::Assns< recob::Track, recob::Hit > > trk2hit(new art::Assns< recob::Track, recob::Hit >);
304  std::unique_ptr< art::Assns< recob::Track, recob::Vertex > > trk2vtx(new art::Assns< recob::Track, recob::Vertex >);
305  std::unique_ptr< art::Assns< recob::Cluster, recob::Hit > > cl2hit(new art::Assns< recob::Cluster, recob::Hit >);
306  std::unique_ptr< art::Assns< recob::Track, recob::Cluster > > trk2cl(new art::Assns< recob::Track, recob::Cluster >);
307  std::unique_ptr< art::Assns< recob::Track, recob::SpacePoint > > trk2sp(new art::Assns< recob::Track, recob::SpacePoint >);
308  std::unique_ptr< art::Assns< recob::SpacePoint, recob::Hit > > sp2hit(new art::Assns< recob::SpacePoint, recob::Hit >);
309 
310 
312  {
313  fClustersNotUsed.clear();
314  fInisegs.clear();
316 
317  for (size_t id = 0; id < fCluListHandle->size(); id++)
318  {
319  std::vector< art::Ptr<recob::Hit> > hitlist;
320  hitlist = fb.at(id);
321 
322  if (hitlist.size() > 5)
323  fClustersNotUsed.push_back(id);
324  }
325 
326  std::vector< ems::DirOfGamma* > showernviews = CollectShower2D(e);
327 
328  Link(e, showernviews);
329 
330  while (fInisegs.size())
331  {
332  fSeltracks.push_back(fInisegs[0]);
333  fInisegs.erase(fInisegs.begin() + 0);
334  }
335 
336  Reoptimize();
337 
338  // conversion from pma track to recob::track
339 
340  size_t spStart = 0, spEnd = 0;
341  double sp_pos[3], sp_err[6], vtx_pos[3];
342  for (size_t i = 0; i < 6; i++) sp_err[i] = 1.0;
343 
344  fTrkIndex = 0;
345 
346  for (auto const trk : fSeltracks)
347  {
348  tracks->push_back(ConvertFrom(*(trk.track)));
349 
350  vtx_pos[0] = trk.track->front()->Point3D().X();
351  vtx_pos[1] = trk.track->front()->Point3D().Y();
352  vtx_pos[2] = trk.track->front()->Point3D().Z();
353  vertices->push_back(recob::Vertex(vtx_pos, fTrkIndex));
354 
355  fTrkIndex++;
356 
357  std::vector< art::Ptr< recob::Cluster > > cl2d;
358  cl2d.push_back( art::Ptr< recob::Cluster >(fCluListHandle, trk.idcl1) );
359  cl2d.push_back( art::Ptr< recob::Cluster >(fCluListHandle, trk.idcl2) );
360 
361  std::vector< art::Ptr< recob::Hit > > hits2d;
363 
364  spStart = allsp->size();
365  for (int h = trk.track->size() - 1; h >= 0; h--)
366  {
367  pma::Hit3D* h3d = (*trk.track)[h];
368  hits2d.push_back(h3d->Hit2DPtr());
369 
370  if ((h == 0) ||
371  (sp_pos[0] != h3d->Point3D().X()) ||
372  (sp_pos[1] != h3d->Point3D().Y()) ||
373  (sp_pos[2] != h3d->Point3D().Z()))
374  {
375  if (sp_hits.size()) // hits assigned to the previous sp
376  {
377  util::CreateAssn(*this, e, *allsp, sp_hits, *sp2hit);
378  sp_hits.clear();
379  }
380  sp_pos[0] = h3d->Point3D().X();
381  sp_pos[1] = h3d->Point3D().Y();
382  sp_pos[2] = h3d->Point3D().Z();
383  allsp->push_back(recob::SpacePoint(sp_pos, sp_err, 1.0));
384  }
385  sp_hits.push_back(h3d->Hit2DPtr());
386  }
387  if (sp_hits.size()) // hits assigned to the last sp
388  {
389  util::CreateAssn(*this, e, *allsp, sp_hits, *sp2hit);
390  }
391  spEnd = allsp->size();
392 
393  if (vertices->size())
394  {
395  size_t vtx_idx = (size_t)(vertices->size() - 1);
396  util::CreateAssn(*this, e, *tracks, *vertices, *trk2vtx, vtx_idx, vtx_idx + 1);
397  }
398 
399  if (cl2d.size())
400  {
401  util::CreateAssn(*this, e, *tracks, cl2d, *trk2cl);
402  }
403 
404  if (hits2d.size())
405  {
406  util::CreateAssn(*this, e, *tracks, *allsp, *trk2sp, spStart, spEnd);
407  util::CreateAssn(*this, e, *tracks, hits2d, *trk2hit);
408  }
409  }
410 
411  fIniIndex = fTrkIndex + 1;
412  for (auto const trk : fPMA3D)
413  {
414  tracks->push_back(ConvertFrom2(*(trk.track)));
415 
416  fIniIndex++;
417 
418  std::vector< art::Ptr< recob::Cluster > > cl2d;
419  cl2d.push_back( art::Ptr< recob::Cluster >(fCluListHandle, trk.idcl1) );
420  cl2d.push_back( art::Ptr< recob::Cluster >(fCluListHandle, trk.idcl2) );
421 
422  std::vector< art::Ptr< recob::Hit > > hits2d;
424 
425  spStart = allsp->size();
426  for (int h = trk.track->size() - 1; h >= 0; h--)
427  {
428  pma::Hit3D* h3d = (*trk.track)[h];
429  hits2d.push_back(h3d->Hit2DPtr());
430 
431  if ((h == 0) ||
432  (sp_pos[0] != h3d->Point3D().X()) ||
433  (sp_pos[1] != h3d->Point3D().Y()) ||
434  (sp_pos[2] != h3d->Point3D().Z()))
435  {
436  if (sp_hits.size()) // hits assigned to the previous sp
437  {
438  util::CreateAssn(*this, e, *allsp, sp_hits, *sp2hit);
439  sp_hits.clear();
440  }
441  sp_pos[0] = h3d->Point3D().X();
442  sp_pos[1] = h3d->Point3D().Y();
443  sp_pos[2] = h3d->Point3D().Z();
444  allsp->push_back(recob::SpacePoint(sp_pos, sp_err, 1.0));
445  }
446  sp_hits.push_back(h3d->Hit2DPtr());
447  }
448  if (sp_hits.size()) // hits assigned to the last sp
449  {
450  util::CreateAssn(*this, e, *allsp, sp_hits, *sp2hit);
451  }
452  spEnd = allsp->size();
453 
454 
455  if (cl2d.size())
456  {
457  util::CreateAssn(*this, e, *tracks, cl2d, *trk2cl);
458  }
459 
460  if (hits2d.size())
461  {
462  util::CreateAssn(*this, e, *tracks, *allsp, *trk2sp, spStart, spEnd);
463  util::CreateAssn(*this, e, *tracks, hits2d, *trk2hit);
464  }
465  }
466 
467  // create cluster from hits, which were an input to find initial part of the cascade.
468  fClIndex = 0;
469  for (auto const& cl : fClusters)
470  if (cl.size())
471  {
472  clusters->push_back(ConvertFrom(cl));
473  fClIndex++;
474 
475  util::CreateAssn(*this, e, *clusters, cl, *cl2hit);
476  }
477 
478 
479  for (unsigned int i = 0; i < showernviews.size(); i++)
480  delete showernviews[i];
481 
482  for (unsigned int i = 0; i < fSeltracks.size(); i++)
483  delete fSeltracks[i].track;
484 
485  for (unsigned int i = 0; i < fInisegs.size(); i++)
486  delete fInisegs[i].track;
487 
488  for (unsigned int i = 0; i < fPMA3D.size(); i++)
489  delete fPMA3D[i].track;
490 
491 }
492 
493  e.put(std::move(tracks));
494  e.put(std::move(vertices));
495  e.put(std::move(clusters));
496  e.put(std::move(allsp));
497 
498  e.put(std::move(trk2hit));
499  e.put(std::move(trk2vtx));
500  e.put(std::move(cl2hit));
501  e.put(std::move(trk2cl));
502  e.put(std::move(trk2sp));
503  e.put(std::move(sp2hit));
504 
505 }
506 
508 {
509  if (!fSeltracks.size()) return;
510  const float min_dist = 3.0F;
511  size_t ta = 0;
512  while (ta < (fSeltracks.size() - 1))
513  {
514  size_t tb = ta + 1;
515  bool found = false;
516  while (tb < fSeltracks.size())
517  {
518  if (ta == tb) {tb++; continue;}
519 
520  TVector3 p1 = fSeltracks[ta].track->front()->Point3D();
521  TVector3 p2 = fSeltracks[tb].track->front()->Point3D();
522  float dist = std::sqrt(pma::Dist2(p1, p2));
523 
524  if (dist < min_dist)
525  if ((fSeltracks[ta].idcl1 == fSeltracks[tb].idcl1) || (fSeltracks[ta].idcl1 == fSeltracks[tb].idcl2) ||
526  (fSeltracks[ta].idcl2 == fSeltracks[tb].idcl1) || (fSeltracks[ta].idcl2 == fSeltracks[tb].idcl2))
527  {
528  found = true;
529  size_t view3 = fSeltracks[ta].view1; size_t idcl3 = fSeltracks[ta].idcl1;
530  std::vector< art::Ptr< recob::Hit > > hits3 = fSeltracks[ta].hits1;
531  std::vector< art::Ptr< recob::Hit > > hits = fSeltracks[ta].hits1;
532  for (size_t h = 0; h < fSeltracks[ta].hits2.size(); ++h)
533  hits.push_back(fSeltracks[ta].hits2[h]);
534 
535  if ((fSeltracks[tb].view1 != fSeltracks[ta].view1) && (fSeltracks[tb].view1 != fSeltracks[ta].view2))
536  {
537  view3 = fSeltracks[tb].view1;
538  for (size_t h = 0; h < fSeltracks[tb].hits1.size(); ++h)
539  hits.push_back(fSeltracks[tb].hits1[h]);
540  }
541  if ((fSeltracks[tb].view2 != fSeltracks[ta].view1) && (fSeltracks[tb].view2 != fSeltracks[ta].view2))
542  {
543  view3 = fSeltracks[tb].view2;
544  for (size_t h = 0; h < fSeltracks[tb].hits2.size(); ++h)
545  hits.push_back(fSeltracks[tb].hits2[h]);
546  }
547 
548 
549  if ((view3 == fSeltracks[ta].view1) || (view3 == fSeltracks[ta].view2))
550  {
551  delete fSeltracks[ta].track;
552  fSeltracks.erase(fSeltracks.begin() + ta);
553  }
554  else
555  {
557 
558  if (pma::Dist2(track->back()->Point3D(), fSeltracks[ta].track->front()->Point3D()) <
559  pma::Dist2(track->front()->Point3D(), fSeltracks[ta].track->front()->Point3D()))
560  track->Flip();
561 
562  IniSeg initrack;
563  initrack.idcl1 = fSeltracks[ta].idcl1; initrack.idcl3 = idcl3;
564  initrack.view1 = fSeltracks[ta].view1; initrack.view3 = view3;
565  initrack.hits1 = fSeltracks[ta].hits1; initrack.hits3 = hits3;
566  initrack.idcl2 = fSeltracks[ta].idcl2;
567  initrack.view2 = fSeltracks[ta].view2;
568  initrack.hits2 = fSeltracks[ta].hits2;
569  initrack.track = track;
570 
571 
572  delete fSeltracks[tb].track; delete fSeltracks[ta].track;
573  fSeltracks.erase(fSeltracks.begin() + tb); fSeltracks.erase(fSeltracks.begin() + ta);
574  fSeltracks.push_back(initrack);
575 
576  }
577  }
578 
579  if (found) break;
580  tb++;
581  }
582 
583  if (!found) ta++;
584  }
585 }
586 
587 std::vector< ems::DirOfGamma* > ems::EMShower3D::CollectShower2D(art::Event const & e)
588 {
589  std::vector< ems::DirOfGamma* > input;
590 
592  {
594  for (unsigned int c = 0; c < fCluListHandle->size(); c++)
595  {
596  std::vector< art::Ptr<recob::Hit> > hitlist;
597  hitlist = fb.at(c);
598 
599  if (hitlist.size() > 5)
600  {
601  std::vector< art::Ptr<recob::Hit> > hits_out;
602  FilterOutSmallParts(2.0, hitlist, hits_out);
603 
604  if (hits_out.size() > 5)
605  {
606  fClusters.push_back(hits_out);
607 
608  ems::DirOfGamma * sh = new ems::DirOfGamma(hits_out, 14, c);
609 
610  if (sh->GetHits2D().size())
611  input.push_back(sh);
612  }
613  }
614  }
615  }
616 
617  return input;
618 }
619 
620 void ems::EMShower3D::Link(art::Event const & e, std::vector< ems::DirOfGamma* > input)
621 {
622  auto const* detprop = lar::providerFrom<detinfo::DetectorPropertiesService>();
623 
624  std::vector< std::vector< size_t > > saveids;
625  std::vector< size_t > saveidsnotusedcls;
626  size_t i = 0;
627 
628  while (i < input.size())
629  {
630  if (!input[i]->GetCandidates().size()){i++; continue;}
631 
632  double mindist = 1.0; // cm
633  std::vector< ems::DirOfGamma* > pairs;
634 
635  size_t startview = input[i]->GetFirstHit()->WireID().Plane;
636  size_t tpc = input[i]->GetFirstHit()->WireID().TPC;
637  size_t cryo = input[i]->GetFirstHit()->WireID().Cryostat;
638 
639  float t1 = detprop->ConvertTicksToX(input[i]->GetFirstHit()->PeakTime(), startview, tpc, cryo);
640 
641  unsigned int idsave = 0;
642  for (unsigned int j = 0; j < input.size(); j++)
643  {
644  if (!input[j]->GetCandidates().size()) continue;
645 
646  size_t secondview = input[j]->GetFirstHit()->WireID().Plane;
647  size_t tpc_j = input[j]->GetFirstHit()->WireID().TPC;
648  size_t cryo_j = input[j]->GetFirstHit()->WireID().Cryostat;
649 
650  if ((i != j) && (secondview != startview) && (tpc == tpc_j) && (cryo == cryo_j))
651  {
652  float t2 = detprop->ConvertTicksToX(input[j]->GetFirstHit()->PeakTime(), secondview, tpc_j, cryo_j);
653  float dist = fabs(t2 - t1);
654 
655  if (dist < mindist)
656  {
657  mindist = dist;
658  pairs.clear();
659  pairs.push_back(input[i]); pairs.push_back(input[j]);
660  idsave = j;
661  }
662  }
663 
664  }
665 
666  bool exist = false;
667  for (unsigned int v = 0; v < saveids.size(); v++)
668  if ((saveids[v][0] == i) || (saveids[v][0] == idsave))
669  if ((saveids[v][1] == i) || (saveids[v][1] == idsave))
670  exist = true;
671 
672 
673  if (pairs.size())
674  {
675  if (!exist) Make3DSeg(e, pairs);
676 
677  std::vector< size_t > ids;
678  ids.push_back(i); ids.push_back(idsave);
679  saveids.push_back(ids);
680  }
681  else
682  {
683  saveidsnotusedcls.push_back(i);
684  }
685 
686  i++;
687  }
688 
689  i = 0;
690  while(i < saveidsnotusedcls.size())
691  {
692  LinkCandidates(e, input, i);
693  i++;
694  }
695 }
696 
697 size_t ems::EMShower3D::LinkCandidates(art::Event const & e, std::vector< ems::DirOfGamma* > input, size_t id)
698 {
700 
701  size_t index = id; bool found = false;
702 
703  if (input[id]->GetCandidates().size() < 2) { return index; }
704 
705  double mindist = 3.0; // cm
706  std::vector< ems::DirOfGamma* > pairs;
707 
708  size_t idcsave = 0; size_t idcjsave = 0;
709  size_t c = 0; size_t idsave = 0;
710  while (c < input[id]->GetCandidates().size())
711  {
712 
713  size_t startview = input[id]->GetCandidates()[c].GetPlane();
714  size_t tpc = input[id]->GetCandidates()[c].GetTPC();
715  size_t cryo = input[id]->GetCandidates()[c].GetCryo();
716 
717  float t1 = input[id]->GetCandidates()[c].GetPosition().Y(); // y --> drift in 2D space.
718 
719  // loop over 2D showers
720  for (size_t j = 0; j < input.size(); ++j)
721  {
722  if (!input[j]->GetCandidates().size()) continue;
723  if (j == id) continue;
724 
725  // loop over candidates
726  for (size_t cj = 0; cj < input[j]->GetCandidates().size(); ++cj)
727  {
728  size_t secondview = input[j]->GetCandidates()[cj].GetPlane();
729  size_t tpc_j = input[j]->GetCandidates()[cj].GetTPC();
730  size_t cryo_j = input[j]->GetCandidates()[cj].GetCryo();
731 
732  size_t thirdview = startview;
733 
734  const geo::CryostatGeo& cryostat = geom->Cryostat(cryo);
735  for (size_t p = 0; p < cryostat.MaxPlanes(); p++)
736  if ((p == startview) || (p == secondview)) {continue;}
737  else {thirdview = p; break;}
738 
739 
740  if ((startview != secondview) && (tpc == tpc_j) && (cryo == cryo_j))// && Validate(input, id, cj, thirdview))
741  {
742  float t2 = input[j]->GetCandidates()[cj].GetPosition().Y();
743  float dist = fabs(t2 - t1);
744 
745  if ((dist < mindist) && Validate(input, id, j, c, cj, thirdview))
746  {
747  mindist = dist;
748  pairs.clear();
749  pairs.push_back(input[id]); pairs.push_back(input[j]);
750  idsave = j; index = j;
751  idcsave = c; idcjsave = cj;
752  found = true;
753  }
754  }
755  }
756  }
757 
758  c++;
759  }
760 
761  if (found && pairs.size())
762  {
763  input[id]->SetIdCandidate(idcsave);
764  input[idsave]->SetIdCandidate(idcjsave);
765  Make3DSeg(e, pairs);
766  }
767 
768  return index;
769 }
770 
771 void ems::EMShower3D::Make3DSeg(art::Event const & e, std::vector< ems::DirOfGamma* > pair)
772 {
773  if (pair.size() < 2) return;
774 
775  // to build a track correctly 2d hits must belong to the same tpc
776  size_t tpc1 = pair[0]->GetFirstHit()->WireID().TPC;
777  size_t tpc2 = pair[1]->GetFirstHit()->WireID().TPC;
778 
779  std::vector< art::Ptr< recob::Hit > > vec1 = pair[0]->GetIniHits();
780  std::vector< art::Ptr< recob::Hit > > vec2 = pair[1]->GetIniHits();
781 
782  if ((vec1.size() < 3) && (vec2.size() < 3)) return;
783 
784  std::vector< art::Ptr<recob::Hit> > hitscl1uniquetpc;
785  std::vector< art::Ptr<recob::Hit> > hitscl2uniquetpc;
786 
787  if (tpc1 == tpc2)
788  for (size_t i = 0; i < vec1.size(); ++i)
789  for (size_t j = 0; j < vec2.size(); ++j)
790  if ((vec1[i]->WireID().TPC == vec2[j]->WireID().TPC) && (tpc1 == vec2[j]->WireID().TPC))
791  {
792  hitscl1uniquetpc.push_back(vec1[i]);
793  hitscl2uniquetpc.push_back(vec2[j]);
794  }
795 
796  if ((hitscl1uniquetpc.size() > 2) && (hitscl2uniquetpc.size() > 2))
797  {
798  pma::Track3D* trk = fProjectionMatchingAlg.buildSegment(hitscl1uniquetpc, hitscl2uniquetpc);
799 
800  //pma::Track3D* trk = fProjectionMatchingAlg.buildSegment(vec1, vec2);
801 
802  // turn the track that front is at vertex - easier to handle associations.
803  if ((trk->back()->Hit2DPtr() == pair[0]->GetFirstHit())
804  || (trk->back()->Hit2DPtr() == pair[1]->GetFirstHit())) trk->Flip();
805 
806 
807  IniSeg initrack;
808  initrack.idcl1 = pair[0]->GetIdCl();
809  initrack.view1 = pair[0]->GetFirstHit()->WireID().Plane;
810  initrack.hits1 = hitscl1uniquetpc;
811  initrack.idcl2 = pair[1]->GetIdCl();
812  initrack.view2 = pair[1]->GetFirstHit()->WireID().Plane;
813  initrack.hits2 = hitscl2uniquetpc;
814  initrack.track = trk;
815 
816  fInisegs.push_back(initrack);
817 
818 
819  }
820 }
821 
822 bool ems::EMShower3D::Validate(std::vector< ems::DirOfGamma* > input, size_t id1, size_t id2, size_t c1, size_t c2, size_t plane3)
823 {
824  bool result = false;
825  if (id1 == id2) return false;
826 
827  std::vector< art::Ptr< recob::Hit > > vec1 = input[id1]->GetCandidates()[c1].GetIniHits();
828  std::vector< art::Ptr< recob::Hit > > vec2 = input[id2]->GetCandidates()[c2].GetIniHits();
829 
830  if ((vec1.size() < 3) || (vec2.size() < 3)) return false;
831 
832  std::vector< art::Ptr<recob::Hit> > hitscl1uniquetpc;
833  std::vector< art::Ptr<recob::Hit> > hitscl2uniquetpc;
834 
835  size_t tpc = vec1[0]->WireID().TPC;
836  for (size_t i = 0; i < vec1.size(); ++i)
837  for (size_t j = 0; j < vec2.size(); ++j)
838  if ((vec1[i]->WireID().TPC == tpc) && (vec2[j]->WireID().TPC == tpc))
839  {
840  hitscl1uniquetpc.push_back(vec1[i]);
841  hitscl2uniquetpc.push_back(vec2[j]);
842  }
843 
844  if ((hitscl1uniquetpc.size() < 3) || (hitscl2uniquetpc.size() < 3)) return false;
845 
846  pma::Track3D* track = fProjectionMatchingAlg.buildSegment(hitscl1uniquetpc, hitscl2uniquetpc);
847  for (size_t i = 0; i < input.size(); ++i)
848  {
849  std::vector< Hit2D* > hits2dcl = input[i]->GetHits2D();
850  for (size_t h = 0; h < hits2dcl.size(); ++h)
851  {
852  TVector2 pfront = pma::GetProjectionToPlane(track->front()->Point3D(), plane3, track->FrontTPC(), track->FrontCryo());
853  TVector2 pback = pma::GetProjectionToPlane(track->back()->Point3D(), plane3, track->BackTPC(), track->BackCryo());
854  if ( (pma::Dist2(hits2dcl[h]->GetPointCm(), pfront) < 1.0F) &&
855  (pma::Dist2(hits2dcl[h]->GetPointCm(), pback) < 1.0F) )
856  {result = true; break;}
857  }
858  }
859  delete track;
860  return result;
861 }
862 
863 
864 
865 bool ems::EMShower3D::Validate(art::Event const & e, const pma::Track3D& src, size_t plane)
866 {
867  bool result = false;
868 
870  std::vector< art::Ptr<recob::Hit> > hitscl;
871  for (size_t id = 0; id < fClustersNotUsed.size(); id++)
872  {
873  std::vector< art::Ptr<recob::Hit> > hits = fbc.at(fClustersNotUsed[id]);
874  for (size_t i = 0; i < hits.size(); i++) hitscl.push_back(hits[i]);
875  }
876 
877 
878  if (fProjectionMatchingAlg.validate(src, hitscl) > 0.2) result = true;
879 
880  return result;
881 }
882 
883 bool ems::EMShower3D::Has(const std::vector<size_t>& v, size_t idx)
884 {
885  for (auto c : v) if (c == idx) return true;
886  return false;
887 }
888 
890  double r2d,
891  const std::vector< art::Ptr<recob::Hit> >& hits_in,
892  std::vector<size_t>& used,
893  std::vector< art::Ptr<recob::Hit> >& hits_out)
894 {
895 
896  hits_out.clear();
897 
898  const double gapMargin = 5.0; // can be changed to f(id_tpc1, id_tpc2)
899  size_t idx = 0;
900 
901  while ((idx < hits_in.size()) && Has(used, idx)) idx++;
902 
903  if (idx < hits_in.size())
904  {
905  hits_out.push_back(hits_in[idx]);
906  used.push_back(idx);
907 
908 
909  double r2d2 = r2d*r2d;
910  double gapMargin2 = sqrt(2 * gapMargin*gapMargin);
911  gapMargin2 = (gapMargin2 + r2d)*(gapMargin2 + r2d);
912 
913  bool collect = true;
914  while (collect)
915  {
916  collect = false;
917  for (size_t i = 0; i < hits_in.size(); i++)
918  if (!Has(used, i))
919  {
920  art::Ptr<recob::Hit> hi = hits_in[i];
921  TVector2 hi_cm = pma::WireDriftToCm(hi->WireID().Wire, hi->PeakTime(), hi->WireID().Plane, hi->WireID().TPC, hi->WireID().Cryostat);
922 
923  bool accept = false;
924  //for (auto const& ho : hits_out)
925  for (size_t idx_o = 0; idx_o < hits_out.size(); idx_o++)
926  {
927  art::Ptr<recob::Hit> ho = hits_out[idx_o];
928 
929  double d2 = pma::Dist2(
930  hi_cm, pma::WireDriftToCm(ho->WireID().Wire, ho->PeakTime(), ho->WireID().Plane, ho->WireID().TPC, ho->WireID().Cryostat));
931 
932  if (hi->WireID().TPC == ho->WireID().TPC)
933  {
934  if (d2 < r2d2) { accept = true; break; }
935  }
936  else
937  {
938  if (d2 < gapMargin2) { accept = true; break; }
939  }
940  }
941  if (accept)
942  {
943  collect = true;
944  hits_out.push_back(hi);
945  used.push_back(i);
946  }
947  }
948  }
949  return true;
950  }
951  else return false;
952 }
953 
954 
956  double r2d,
957  const std::vector< art::Ptr<recob::Hit> >& hits_in,
958  std::vector< art::Ptr<recob::Hit> >& hits_out)
959 {
960  size_t min_size = hits_in.size() / 5;
961  if (min_size < 3) min_size = 3;
962 
963  std::vector<size_t> used;
964  std::vector< art::Ptr<recob::Hit> > close_hits;
965 
966  while (GetCloseHits(r2d, hits_in, used, close_hits))
967  {
968  if (close_hits.size() > min_size)
969  for (auto h : close_hits) hits_out.push_back(h);
970  }
971 }
size_t LinkCandidates(art::Event const &e, std::vector< ems::DirOfGamma * > input, size_t id)
unsigned int TPC(void) const
Definition: PmaHit3D.h:57
void Link(art::Event const &e, std::vector< ems::DirOfGamma * > input)
TTree * t1
Definition: plottest35.C:26
unsigned int Cryo(void) const
Definition: PmaHit3D.h:56
unsigned int MaxPlanes() const
Returns the largest number of planes among the TPCs in this cryostat.
Implementation of the Projection Matching Algorithm.
double Dist2(const TVector2 &v1, const TVector2 &v2)
Definition: Utilities.cxx:19
unsigned int FrontTPC(void) const
Definition: PmaTrack3D.h:104
double dEdx_AREA(art::Ptr< recob::Hit > hit, double pitch, double T0=0) const
unsigned int BackTPC(void) const
Definition: PmaTrack3D.h:107
bool Flip(std::vector< pma::Track3D * > &allTracks)
Definition: PmaTrack3D.cxx:499
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
Declaration of signal hit object.
size_t idcl1
std::vector< IniSeg > fPMA3D
::fhicl::TupleAs< Point(::geo::Length_t,::geo::Length_t,::geo::Length_t)> Point3D
Atom object for reading a 3D point or vector (centimeters).
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:130
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:313
void reconfigure(fhicl::ParameterSet const &p)
Geometry information for a single cryostat.
Definition: CryostatGeo.h:36
art::Ptr< recob::Hit > const & Hit2DPtr(void) const
Definition: PmaHit3D.h:46
Definition of vertex object for LArSoft.
Definition: Vertex.h:35
std::vector< size_t > fClustersNotUsed
unsigned int fIniIndex
size_t view3
unsigned int fClIndex
size_t view2
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
art::Handle< std::vector< recob::Cluster > > fCluListHandle
void produce(art::Event &e) override
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
std::vector< size_t > fTracksNotUsed
double selectInitialHits(pma::Track3D &trk, unsigned int view=geo::kZ, unsigned int *nused=0) const
void hits()
Definition: readHits.C:15
TCanvas * c1
Definition: plotHisto.C:7
unsigned int fTrkIndex
EMShower3D(fhicl::ParameterSet const &p)
TCanvas * c2
Definition: plot_hist.C:75
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
void beginJob()
Definition: Breakpoints.cc:14
void push_back(Ptr< U > const &p)
Definition: PtrVector.h:441
std::string fTrk3DModuleLabel
std::vector< art::Ptr< recob::Hit > > hits3
size_t idcl3
recob::Track ConvertFrom2(pma::Track3D &src)
parameter set interface
unsigned int BackCryo(void) const
Definition: PmaTrack3D.h:108
T get(std::string const &key) const
Definition: ParameterSet.h:231
pma::Track3D * buildSegment(const std::vector< art::Ptr< recob::Hit > > &hits_1, const std::vector< art::Ptr< recob::Hit > > &hits_2=std::vector< art::Ptr< recob::Hit > >()) const
void FilterOutSmallParts(double r2d, const std::vector< art::Ptr< recob::Hit > > &hits_in, std::vector< art::Ptr< recob::Hit > > &hits_out)
TFile fb("Li6.root")
size_t idcl2
pma::Hit3D *& back()
Definition: PmaTrack3D.h:74
CryostatGeo const & Cryostat(geo::CryostatID const &cryoid) const
Returns the specified cryostat.
void beginJob() override
bool CreateAssn(PRODUCER const &prod, art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t indx=UINT_MAX)
Creates a single one-to-one association.
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
Declaration of cluster object.
Provides recob::Track data product.
calo::CalorimetryAlg fCalorimetryAlg
unsigned int FrontCryo(void) const
Definition: PmaTrack3D.h:105
bool Validate(art::Event const &e, const pma::Track3D &src, size_t plane)
Definition: DirOfGamma.h:18
size_type size() const
Definition: PtrVector.h:308
TTree * t2
Definition: plottest35.C:36
bool GetCloseHits(double r2d, const std::vector< art::Ptr< recob::Hit > > &hits_in, std::vector< size_t > &used, std::vector< art::Ptr< recob::Hit > > &hits_out)
TVector3 const & Point3D(void) const
Definition: PmaHit3D.h:48
std::vector< std::vector< art::Ptr< recob::Hit > > > fClusters
Implementation of the Projection Matching Algorithm.
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:219
Utility object to perform functions of association.
void Make3DSeg(art::Event const &e, std::vector< ems::DirOfGamma * > pair)
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
pma::ProjectionMatchingAlg fProjectionMatchingAlg
TVector2 WireDriftToCm(unsigned int wire, float drift, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:307
TVector2 GetProjectionToPlane(const TVector3 &p, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:291
bool Has(const std::vector< size_t > &v, size_t idx)
size_t view1
TrackCollectionProxyElement< TrackCollProxy > Track
Proxy to an element of a proxy collection of recob::Track objects.
Definition: Track.h:1035
std::vector< art::Ptr< recob::Hit > > hits2
constexpr double kBogusD
obviously bogus double value
std::vector< art::Ptr< recob::Hit > > hits1
std::vector< IniSeg > fSeltracks
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
pma::Track3D * track
Float_t e
Definition: plot.C:34
recob::Track ConvertFrom(pma::Track3D &src)
size_t size() const
Definition: PmaTrack3D.h:76
std::string fCluModuleLabel
std::vector< Hit2D * > const & GetHits2D(void) const
Definition: DirOfGamma.h:126
Float_t track
Definition: plot.C:34
void clear()
Definition: PtrVector.h:537
pma::Hit3D *& front()
Definition: PmaTrack3D.h:72
double validate(const pma::Track3D &trk, const std::vector< art::Ptr< recob::Hit > > &hits) const
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:51
std::vector< IniSeg > fInisegs
std::vector< ems::DirOfGamma * > CollectShower2D(art::Event const &e)