LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
MergeEMShower3D_module.cc
Go to the documentation of this file.
1 // Class: MergeEMShower3D
3 // Module Type: producer
4 // File: MergeEMShower3D_module.cc
5 // Authors: dorota.stefan@cern.ch robert.sulej@cern.ch
7 
15 #include "fhiclcpp/ParameterSet.h"
20 
22 
30 
35 
36 #include <memory>
37 
39 
40 // ROOT includes
41 #include "TTree.h"
42 #include "TLorentzVector.h"
43 #include "TMathBase.h"
44 
45 namespace ems
46 {
47  class MergeEMShower3D;
48  class ShowerInfo;
49  class bDistLess;
50  class ShowersCollection;
51 }
52 
54 {
55  public:
56  ShowerInfo(int key, int gid, bool hasvtx, std::vector<double> vensum, recob::Track const& trk);
57  double Pointsto(ems::ShowerInfo const& s1) const;
58  double Angleto(TVector3 const& pmapoint) const;
59 
60  bool HasConPoint(void) const {return fHasVtx;}
61  int GetKey(void) const {return fKey;}
62  int GetGid(void) const {return fGId;}
63  int GetPlaneId(void) const {return fPlaneId;}
64  double GetAdcSum(void) const {return fAdcSum;}
65  std::vector<double> GetEnTot(void) const {return fVentot; }
66 
67  TVector3 GetFront(void) const {return fFront;}
68  TVector3 GetEnd(void) const {return fEnd;}
69  TVector3 GetDir(void) const {return fDir;}
70  std::vector<double> GetDedx(void) const { return fVdqdx; }
71 
72  double GetP0Dist(void) const { return fP0Dist; }
73  void SetP0Dist(const TVector3 p)
74  {
75  if (fHasVtx) fP0Dist = std::sqrt( pma::Dist2(p, fFront) );
76  else fP0Dist = std::sqrt( pma::Dist2(p, 0.5 * (fFront + fEnd)) );
77  }
78 
79  private:
80  int fKey;
81  int fGId;
82  int fPlaneId;
83  bool fHasVtx;
84  double fAdcSum;
85  std::vector<double> fVentot;
86  double fP0Dist;
87 
88  TVector3 fFront;
89  TVector3 fEnd;
90  TVector3 fDir;
91 
92  std::vector<double> fVdqdx;
93 
94 };
95 
96 ems::ShowerInfo::ShowerInfo(int key, int gid, bool hasvtx, std::vector<double> vensum, recob::Track const& trk) :
97 fKey(key),
98 fGId(gid),
99 fVentot(vensum),
100 fP0Dist(0)
101 {
102  fDir = trk.VertexDirection();
103  fFront = trk.Vertex();
104  fEnd = trk.End();
105  fHasVtx = hasvtx;
106 
108  /*fVdqdx.push_back(trk.DQdxAtPoint(0, geo::kU));
109  fVdqdx.push_back(trk.DQdxAtPoint(0, geo::kV));
110  fVdqdx.push_back(trk.DQdxAtPoint(0, geo::kZ));*/
111 
112 
113  bool isplane = false;
114  unsigned int nplanes = geom->Nplanes();
115  for (unsigned int p = 0; p < nplanes; ++p)
116  {
117  geo::View_t view = geom->View(p);
118  fVdqdx.push_back(trk.DQdxAtPoint(0, view));
119  if (trk.DQdxAtPoint(1, view) > 0)
120  {fPlaneId = int(p); isplane = true;}
121  }
122  if (!isplane) fPlaneId = -1;
123 
124  /*if (trk.DQdxAtPoint(1, geo::kU) > 0) fPlaneId = geo::kU;
125  else if (trk.DQdxAtPoint(1, geo::kV) > 0) fPlaneId = geo::kV;
126  else if (trk.DQdxAtPoint(1, geo::kZ) > 0) fPlaneId = geo::kZ;
127  else fPlaneId = -1;*/
128 }
129 
131 {
132  double cosine0 = (this->fFront - s1.fFront) * (1.0 / (this->fFront - s1.fFront).Mag()) * s1.fDir;
133  double th0 = 180.0F * (std::acos(cosine0)) / TMath::Pi();
134 
135  std::vector< std::pair<TVector3, TVector3> > lines;
136  lines.push_back(std::pair<TVector3, TVector3>(this->fFront, this->fFront + this->fDir));
137  lines.push_back(std::pair<TVector3, TVector3>(s1.fFront, s1.fFront + s1.fDir));
138 
139  TVector3 isect, pm;
140  pma::SolveLeastSquares3D(lines, isect);
141  pm = pma::GetProjectionToSegment(isect, this->fFront, this->fFront + this->fDir);
142  double cosine1 = (pm - s1.fFront) * (1.0/(pm - s1.fFront).Mag()) * s1.fDir;
143  double th1 = 180.0F * (std::acos(cosine1)) / TMath::Pi();
144 
145  double thmin = th0;
146  if (th1 < th0) thmin = th1;
147 
148  return thmin;
149 }
150 
151 double ems::ShowerInfo::Angleto(TVector3 const& pmapoint) const
152 {
153  double cosine = (pmapoint - this->fFront) * (1.0 / (pmapoint - this->fFront).Mag()) * this->fDir;
154  double th = 180.0F * (std::acos(cosine)) / TMath::Pi();
155 
156  return th;
157 }
158 
160  public std::binary_function<const ems::ShowerInfo&, const ems::ShowerInfo&, bool>
161 {
162  public:
163  bool operator() (const ems::ShowerInfo& s1, const ems::ShowerInfo& s2)
164  {
165  return (s1.GetP0Dist() < s2.GetP0Dist());
166  }
167 };
168 
170 {
171  public:
173 
174  double MinAngle(ShowersCollection const& col2);
175  void Merge(ShowersCollection const& col2);
176  void Merge(ShowerInfo const& src);
177 
178  double Angle(TVector3 p0, TVector3 test);
179 
180  double Angle(TVector3 test);
181 
182  size_t Size() {return fParts.size();}
183 
184  std::vector< ShowerInfo > const & GetParts(void) const { return fParts; }
185 
186  bool IsClean() {return Clean;}
187 
188  int PlaneId();
189  TVector3 Front();
190  TVector3 Dir();
191  std::vector<double> DeDx(void);
192  std::vector<double> TotEn(void);
193 
195 
196  private:
197  bool Clean;
198  int GId;
199  std::vector< ShowerInfo > fParts;
200 };
201 
203 first(part)
204 {
205  Clean = true;
206  GId = part.GetGid();
207  if (GId == 9999)
208  Clean = false;
209 
210  fParts.push_back(part);
211 }
212 
214 {
215  if (fParts.size()) return fParts.front().GetPlaneId();
216  else return 0;
217 }
218 
220 {
221  if (fParts.size()) return fParts.front().GetFront();
222  else return TVector3(0, 0, 0);
223 }
224 
226 {
227  if (fParts.size()) return fParts.front().GetDir();
228  else return TVector3(0, 0, 0);
229 }
230 
231 std::vector<double> ems::ShowersCollection::DeDx()
232 {
233  if (fParts.size()) return fParts.front().GetDedx();
234  else return std::vector<double>(0.0);
235 }
236 
237 std::vector<double> ems::ShowersCollection::TotEn()
238 {
239  std::vector<double> en;
240 
241  double enkU = 0.0; double enkV = 0.0; double enkZ = 0.0;
242  for (size_t i = 0; i < fParts.size(); ++i)
243  {
244  std::vector<double> showeren = fParts[i].GetEnTot();
245  enkU += showeren[0];
246  enkV += showeren[1];
247  enkZ += showeren[2];
248  }
249 
250  en.push_back(enkU);
251  en.push_back(enkV);
252  en.push_back(enkZ);
253 
254  return en;
255 }
256 
257 double ems::ShowersCollection::Angle(TVector3 p0, TVector3 test)
258 {
259  p0 = Front() - p0;
260  p0 *= 1.0 / p0.Mag();
261  test -= Front();
262  test *= 1.0 / test.Mag();
263  double c = p0 * test;
264  if (c > 1.0) c = 1.0;
265  return 180.0 * acos(c) / TMath::Pi();
266 }
267 
269 {
270  test -= Front();
271  test *= 1.0 / test.Mag();
272  double c = Dir() * test;
273  if (c > 1.0) c = 1.0;
274  return 180.0 * acos(c) / TMath::Pi();
275 }
276 
277 
279 {
280  double a, a_min = 360.0;
281  for (size_t i = 0; i < fParts.size(); ++i)
282  {
283  const ShowerInfo& temp = fParts[i];
284  a = col2.first.Pointsto(temp);
285  if (a < a_min) a_min = a;
286  }
287  return a_min;
288 }
289 
291 {
292  for (size_t i = 0; i < col2.fParts.size(); ++i)
293  for (size_t j = 0; j < fParts.size(); ++j)
294  if (fParts[j].GetGid() != col2.fParts[i].GetGid())
295  {
296  Clean = false;
297  break;
298  }
299 
300  for (size_t i = 0; i < col2.fParts.size(); ++i)
301  fParts.push_back(col2.fParts[i]);
302 
303 }
304 
306 {
307  if (fParts.size())
308  {
309  for (auto const & p : fParts)
310  if (p.GetGid() != src.GetGid())
311  {
312  Clean = false;
313  break;
314  }
315  }
316  else
317  {
318  Clean = true;
319  GId = src.GetGid();
320  }
321  fParts.push_back(src);
322 }
323 
325 public:
326  explicit MergeEMShower3D(fhicl::ParameterSet const & p);
327 
328  MergeEMShower3D(MergeEMShower3D const &) = delete;
329  MergeEMShower3D(MergeEMShower3D &&) = delete;
330  MergeEMShower3D & operator = (MergeEMShower3D const &) = delete;
331  MergeEMShower3D & operator = (MergeEMShower3D &&) = delete;
332 
333  void beginJob() override;
334 
335  void produce(art::Event & e) override;
336 
337  void reconfigure(fhicl::ParameterSet const& p);
338 
339  void mcinfo(art::Event & evt);
340 
341  int getClusterBestId(const std::vector< art::Ptr<recob::Hit> >& v);
342 
343  int getGammaId(art::Event & evt, const size_t t);
344 
345  double getCos3D(const TVector3& p0, const recob::Track& trk);
346 
347  double getClusterAdcSum(const std::vector< art::Ptr<recob::Hit> >& v);
348 
349  std::vector< ems::ShowersCollection > collectshowers(art::Event & evt, std::vector< ems::ShowerInfo > & showers, bool refpoint);
350 
351  TVector3 getBestPoint(
352  const std::vector<recob::Track>& tracks,
353  const std::vector< ShowerInfo >& showers,
354  const TVector3& p0, const TVector3& p1,
355  double step);
356 
357 private:
358 
359  TTree* fEvTree;
360 // TTree* fClTree;
361 
363 
364  int fNParts[2];
365  int fNPMA; int fNConv; int fNTot;
366 // int fHasConvPt;
367  int fWhat, fEvWhat;
368  int fNMerged, fNCleanMerged;
369 
371  double fMcMom /*, fTrkLen *//*, fAdcSum */;
373 
374  TVector3 fPi0vtx;
375  TVector3 fRefPoint;
376 
377  std::string fHitsModuleLabel;
378  std::string fCluModuleLabel;
379  std::string fTrk3DModuleLabel;
380  std::string fVtxModuleLabel;
381 
384 
386 
387 };
388 
390  : fShowerEnergyAlg(p.get<fhicl::ParameterSet>("ShowerEnergyAlg"))
391 {
392  reconfigure(p);
393 
394  produces< std::vector<recob::Shower> >();
395  produces< std::vector<recob::Vertex> >();
396 
397  produces< art::Assns<recob::Shower, recob::Vertex> >();
398  produces< art::Assns<recob::Shower, recob::Cluster> >();
399  produces< art::Assns<recob::Shower, recob::Hit> >();
400 }
401 
403 {
405 
406  fEvTree = tfs->make<TTree>("ShowerTestEv", "pi0 tests");
407  fEvTree->Branch("fEvNumber", &fEvNumber, "fEvNumber/I");
408  fEvTree->Branch("fNPartsA", &fNParts[0], "fNPartsA/I");
409  fEvTree->Branch("fNPartsB", &fNParts[1], "fNPartsB/I");
410  fEvTree->Branch("fNTot", &fNTot, "fNTot/I");
411  fEvTree->Branch("fNPMA", &fNPMA, "fNPMA/I");
412  fEvTree->Branch("fNConv", &fNConv, "fNConv/I");
413  fEvTree->Branch("fEvWhat", &fEvWhat, "fEvWhat/I");
414  fEvTree->Branch("fMcMom", &fMcMom, "fMcMom/D");
415 
416  fEvTree->Branch("fDistvtxmcreco", &fDistvtxmcreco, "fDistvtxmcreco/D");
417  fEvTree->Branch("fNMerged", &fNMerged, "fNMerged/I");
418  fEvTree->Branch("fNCleanMerged", &fNCleanMerged, "fNCleanMerged/I");
419 }
420 
422 {
423  fHitsModuleLabel = p.get< std::string >("HitsModuleLabel");
424  fCluModuleLabel = p.get< std::string >("ClustersModuleLabel");
425  fTrk3DModuleLabel = p.get< std::string >("Trk3DModuleLabel");
426  fVtxModuleLabel = p.get< std::string >("VtxModuleLabel");
427 
428  fNarrowConeAngle = p.get< double >("NarrowConeAngle");
429  fWideConeAngle = p.get< double >("WideConeAngle");
430 
431  return;
432 }
433 
435 {
436  std::unique_ptr< std::vector< recob::Shower > > cascades(new std::vector< recob::Shower >);
437  std::unique_ptr< std::vector< recob::Vertex > > vertices(new std::vector< recob::Vertex >);
438 
439  std::unique_ptr< art::Assns< recob::Shower, recob::Vertex > > shs2vtx(new art::Assns< recob::Shower, recob::Vertex >);
440  std::unique_ptr< art::Assns< recob::Shower, recob::Cluster > > shs2cl(new art::Assns< recob::Shower, recob::Cluster >);
441  std::unique_ptr< art::Assns< recob::Shower, recob::Hit > > shs2hit(new art::Assns< recob::Shower, recob::Hit >);
442 
443  std::vector< ShowersCollection > gammawithconv;
444 
448  art::Handle< std::vector<recob::Hit> > hitListHandle;
449  if (evt.getByLabel(fTrk3DModuleLabel, trkListHandle) &&
450  evt.getByLabel(fVtxModuleLabel, vtxListHandle) &&
451  evt.getByLabel(fCluModuleLabel, cluListHandle) &&
452  evt.getByLabel(fHitsModuleLabel, hitListHandle))
453  {
454 
455  fEvNumber = evt.id().event();
456  fNMerged = 0; fNCleanMerged = 0; fNParts[0] = 0; fNParts[1] = 0;
457  fNPMA = 0; fNConv = 0; fNTot = 0; fMcMom = 0; fWhat = 0;
458  fDistvtxmcreco = -10.0;
459 
460  art::FindManyP< recob::Cluster > cluFromTrk(trkListHandle, evt, fTrk3DModuleLabel);
461  art::FindManyP< recob::Vertex > vtxFromTrk(trkListHandle, evt, fVtxModuleLabel);
462  art::FindManyP< recob::Hit > hitFromClu(cluListHandle, evt, fCluModuleLabel);
463 
464  // all reconstructed shower fragments
465  std::vector< ShowerInfo > showers;
466  for (size_t t = 0; t < trkListHandle->size(); ++t)
467  {
468  const recob::Track& trk = (*trkListHandle)[t];
469 
470  auto src_clu_list = cluFromTrk.at(t);
471  double adcsum = 0.0;
472  std::vector<double> ensum(3, 0.0);
473 
474  for (size_t c = 0; c < src_clu_list.size(); ++c)
475  {
476  std::vector< art::Ptr<recob::Hit> > v = hitFromClu.at(src_clu_list[c].key());
477  adcsum += getClusterAdcSum(v);
478 
479  if (src_clu_list[c]->View() == geo::kU) { ensum[0] += adcsum; }
480  else if (src_clu_list[c]->View() == geo::kV) { ensum[1] += adcsum; }
481  else if (src_clu_list[c]->View() == geo::kZ) { ensum[2] += adcsum; }
482  }
483 
484  auto cnv = vtxFromTrk.at(t);
485  if (cnv.size())
486  {
487  int gid = getGammaId(evt, t);
488  ShowerInfo si(t, gid, true, ensum, trk);
489  showers.push_back(si);
490  }
491  else
492  {
493  int gid = getGammaId(evt, t);
494  ShowerInfo si(t, gid, false, ensum, trk);
495  showers.push_back(si);
496  }
497  }
498 
499  fNTot = showers.size();
500 
502  gammawithconv = collectshowers(evt, showers, true); // true switches on refpoint
503 
504  // with pma segments reconstructed
505  // procceed with pma segments when two conversions
506  std::vector< ShowerInfo > pmaseg;
507  for (size_t i = 0; i < showers.size(); ++i)
508  if (!showers[i].HasConPoint())
509  pmaseg.push_back(showers[i]);
510 
511  fNPMA = pmaseg.size();
512 
513  const double bigcone = fWideConeAngle; // degree.
514 
515  if (gammawithconv.size())
516  for (size_t i = 0; i < pmaseg.size(); i++)
517  {
518  double a_min = bigcone; size_t best = 0;
519  for (size_t j = 0; j < gammawithconv.size(); j++)
520  {
521  TVector3 halfpoint = (pmaseg[i].GetFront() + pmaseg[i].GetEnd()) * 0.5;
522  double a = gammawithconv[j].first.Angleto(halfpoint);
523  if (a < a_min)
524  {
525  a_min = a; best = j;
526  }
527  }
528  if (a_min < bigcone)
529  gammawithconv[best].Merge(pmaseg[i]);
530  }
531 
532  mcinfo(evt);
533 
534 
535  for (size_t i = 0; i < gammawithconv.size(); ++i)
536  if (gammawithconv[i].IsClean()) fNCleanMerged++;
537 
538  fNMerged = gammawithconv.size();
539  if (gammawithconv.size() == 2)
540  {
541  fNParts[0] = gammawithconv[0].Size();
542  fNParts[1] = gammawithconv[1].Size();
543  }
544 
545  fDistvtxmcreco = std::sqrt(pma::Dist2(fRefPoint, fPi0vtx));
546  fEvTree->Fill();
547 
548  fVtxIndex = 0;
549  for (size_t i = 0; i < gammawithconv.size(); ++i)
550  {
551  int id = i;
552  TVector3 v0(0., 0., 0.);
553  TVector3 dir = gammawithconv[i].Dir();
554  TVector3 front = gammawithconv[i].Front();
555 
556  std::vector<double> dedx = gammawithconv[i].DeDx();
557  std::vector< double > vd;
558  recob::Shower cas(
559  dir, v0, front, v0,
560  vd, vd, dedx, vd, gammawithconv[i].PlaneId(), id);
561 
562 
563  std::vector< art::Ptr<recob::Cluster> > cls;
564  std::vector< art::Ptr<recob::Hit> > hits;
565  for (size_t p = 0; p < gammawithconv[i].Size(); p++)
566  {
567  int trkKey = gammawithconv[i].GetParts()[p].GetKey();
568  auto src_clu_list = cluFromTrk.at(trkKey);
569 
570  for (size_t c = 0; c < src_clu_list.size(); c++)
571  {
572  cls.push_back(src_clu_list[c]);
573 
574  auto v = hitFromClu.at(src_clu_list[c].key());
575  for (size_t h = 0; h < v.size(); h++) hits.push_back(v[h]);
576  }
577 
578  auto ver_list = vtxFromTrk.at(trkKey);
579  }
580 
581  std::vector<double> totE;
582  std::vector<double> totEerr;
583  for (int i = 0; i<3; ++i){
584  totE.push_back(fShowerEnergyAlg.ShowerEnergy(hits,i));
585  totEerr.push_back(0);
586  }
587  cas.set_total_energy(totE);
588  cas.set_total_energy_err(totEerr);
589  cascades->push_back(cas);
590 
591  double vtx_pos[3] = {front.X(), front.Y(), front.Z()};
592  vertices->push_back(recob::Vertex(vtx_pos, fVtxIndex));
593  fVtxIndex++;
594 
595  if (vertices->size())
596  {
597  size_t vtx_idx = (size_t)(vertices->size() - 1);
598  util::CreateAssn(*this, evt, *cascades, *vertices, *shs2vtx, vtx_idx, vtx_idx + 1);
599  }
600 
601  util::CreateAssn(*this, evt, *cascades, cls, *shs2cl);
602  util::CreateAssn(*this, evt, *cascades, hits, *shs2hit);
603 
604  }
605 
606  }
607 
608  evt.put(std::move(cascades));
609  evt.put(std::move(vertices));
610 
611  evt.put(std::move(shs2vtx));
612  evt.put(std::move(shs2cl));
613  evt.put(std::move(shs2hit));
614 }
615 
616 std::vector< ems::ShowersCollection > ems::MergeEMShower3D::collectshowers(art::Event & evt, std::vector< ems::ShowerInfo > & showers, bool refpoint)
617 {
618  std::vector< ems::ShowersCollection > gammawithconv;
619 
620  const double smallcone = fNarrowConeAngle; // degree.
621  bool merge = true;
622  if (refpoint)
623  {
626  if (evt.getByLabel(fTrk3DModuleLabel, trkListHandle))
627  {
628  double dsize[6];
629  geom->CryostatBoundaries(dsize, 0);
630  TVector3 geoP0(dsize[0], dsize[2], dsize[4]), geoP1(dsize[1], dsize[3], dsize[5]);
631  TVector3 pov = getBestPoint(*trkListHandle, showers, geoP0, geoP1, 5.0);
632  fRefPoint = pov;
633 
634  for (size_t is = 0; is < showers.size(); is++)
635  showers[is].SetP0Dist(pov);
636  std::sort(showers.begin(), showers.end(), ems::bDistLess());
637 
638  for (size_t is = 0; is < showers.size(); is++)
639  {
640  size_t m_idx = 0;
641  double a, a_min = 360.0;
642  for (size_t m = 0; m < gammawithconv.size(); ++m)
643  {
644  a = gammawithconv[m].Angle(showers[is].GetFront());
645  if ((a < fNarrowConeAngle) && (a < a_min))
646  {
647  m_idx = m; a_min = a;
648  }
649 
650  std::vector< std::pair<TVector3, TVector3> > lines;
651  lines.push_back(std::pair<TVector3, TVector3>(
652  showers[is].GetFront(), showers[is].GetFront() + showers[is].GetDir()));
653  lines.push_back(std::pair<TVector3, TVector3>(
654  gammawithconv[m].Front(), gammawithconv[m].Front() + gammawithconv[m].Dir()));
655 
656  TVector3 isect, pm;
657  pma::SolveLeastSquares3D(lines, isect);
658  pm = pma::GetProjectionToSegment(isect,
659  gammawithconv[m].Front(), gammawithconv[m].Front() + gammawithconv[m].Dir());
660 
661  if (pma::Dist2(pov, gammawithconv[m].Front()) < pma::Dist2(pov, pm))
662  {
663  a = gammawithconv[m].Angle(isect);
664  if ((a < fNarrowConeAngle) && (a < a_min))
665  {
666  m_idx = m; a_min = a;
667  }
668  }
669  }
670  if (a_min < fNarrowConeAngle)
671  {
672  gammawithconv[m_idx].Merge(showers[is]);
673  }
674  else if (showers[is].HasConPoint())
675  {
676  ems::ShowersCollection sc(showers[is]);
677  gammawithconv.push_back(sc);
678  }
679 
680  }
681  }
682  }
683  else
684  {
685  for (size_t i = 0; i < showers.size(); ++i)
686  if (showers[i].HasConPoint())
687  {
688  ems::ShowersCollection sc(showers[i]);
689  gammawithconv.push_back(sc);
690  }
691 
692  fNConv = gammawithconv.size();
693  while (merge)
694  {
695  merge = false;
696  size_t i = 0;
697  while (i < gammawithconv.size())
698  {
699  size_t best = 0; double a_min = smallcone;
700  for (size_t j = 0; j < gammawithconv.size(); j++)
701  if (i != j)
702  {
703  double a = gammawithconv[i].MinAngle(gammawithconv[j]);
704  if (a < a_min)
705  {
706  a_min = a; best = j; merge = true;
707  }
708  }
709  if (merge)
710  {
711  gammawithconv[i].Merge(gammawithconv[best]);
712  gammawithconv.erase(gammawithconv.begin() + best);
713  break;
714  }
715  i++;
716  }
717  }
718  }
719  return gammawithconv;
720 }
721 
722 
724 {
726 
727  std::map< int, size_t > ids;
728  for (auto ptr : v)
729  {
730  auto hid = bt_serv->HitToTrackIDEs(ptr);
731  if (hid.size()) ids[hid.front().trackID]++;
732  }
733 
734  int best_id = 9999;
735  size_t max = 0;
736  for (auto p : ids)
737  {
738  if ((p.second > (int)(0.8 * v.size())) && (p.second > max))
739  {
740  max = p.second; best_id = p.first;
741  }
742  }
743 
744  return best_id;
745 }
746 
748 {
749  int pi0_idx = -1;
750 
751  art::Handle< std::vector<simb::MCTruth> > mctruthListHandle;
752  std::vector< art::Ptr<simb::MCTruth> > mclist;
753  if (evt.getByLabel("generator", mctruthListHandle))
754  {
755  art::fill_ptr_vector(mclist, mctruthListHandle);
756  if (mclist.size())
757  {
758  size_t np = mclist[0]->NParticles();
759  for (size_t i = 0; i < np; ++i)
760  if (mclist[0]->GetParticle(i).PdgCode() == 111)
761  {pi0_idx = i; break;}
762  }
763  }
764 
765  if (pi0_idx < 0) return;
766 
767  const simb::MCParticle & pi0 = mclist[0]->GetParticle(pi0_idx);
768  TVector3 pi0_vtx(pi0.Vx(), pi0.Vy(), pi0.Vz());
769  fPi0vtx = pi0_vtx;
770  fMcMom = pi0.P();
771 }
772 
774 {
775  fWhat = 0;
776  int cid = 0, gid = 0;
780  art::Handle< std::vector<recob::Hit> > hitListHandle;
781  if (evt.getByLabel(fTrk3DModuleLabel, trkListHandle) &&
782  evt.getByLabel(fVtxModuleLabel, vtxListHandle) &&
783  evt.getByLabel(fCluModuleLabel, cluListHandle) &&
784  evt.getByLabel(fHitsModuleLabel, hitListHandle))
785  {
786  art::FindManyP< recob::Cluster > cluFromTrk(trkListHandle, evt, fTrk3DModuleLabel);
787  art::FindManyP< recob::Vertex > vtxFromTrk(trkListHandle, evt, fVtxModuleLabel);
788  art::FindManyP< recob::Hit > hitFromClu(cluListHandle, evt, fCluModuleLabel);
789 
790  auto src_clu_list = cluFromTrk.at(t);
791  for (size_t c = 0; c < src_clu_list.size(); ++c)
792  {
793  std::vector< art::Ptr<recob::Hit> > v = hitFromClu.at(src_clu_list[c].key());
794  cid = getClusterBestId(v);
795 
796  if ((fWhat == 0) && (cid == 9999)) fWhat = 1; // confused 2D
797  if ((fWhat == 0) && (gid == 0)) gid = cid;
798  if ((fWhat == 0) && (cid != gid)) fWhat = 2; // confused 3D
799  }
800  }
801 
802  if (fWhat == 0) return gid;
803  else return 9999;
804 }
805 
807  const std::vector< recob::Track >& tracks,
808  const std::vector< ShowerInfo >& showers,
809  const TVector3& p0, const TVector3& p1, double step)
810 {
811  TVector3 best, p;
812 
813  double f, fmin = 1.0e10;
814 
815  double dx = step;
816  double dy = step;
817  double dz = step;
818 
819  double x0 = p0.X();
820  while (x0 < p1.X())
821  {
822  double y0 = p0.Y();
823  while (y0 < p1.Y())
824  {
825  double z0 = p0.Z();
826  while (z0 < p1.Z())
827  {
828  f = 0.0;
829  //size_t nOK = 0;
830  TVector3 mid(0., 0., 0.);
831 
832  p.SetXYZ(x0, y0, z0);
833  for (size_t t = 0; t < tracks.size(); t++)
834  // if (showers[t].OK)
835  {
836  auto const & trk = tracks[t];
837 
838  double cos = -acos( getCos3D(p, trk) ) / (0.5*TMath::Pi());
839  //double cos = getCos3D(p, trk);
840 
841  if (showers[t].HasConPoint()) cos *= 3.0;
842  cos *= sqrt( showers[t].GetAdcSum() );
843 
844  mid += 0.5 * (trk.Vertex() + trk.End());
845 
846  f += cos;
847  // nOK++;
848  }
849  //if (!nOK) return best;
850 
851  // f /= nOK;
852  // mid *= 1.0 / nOK;
853 
854  f = -f + 0.0001 * sqrt(pma::Dist2(p, mid));
855  if (f < fmin)
856  {
857  fmin = f; best = p;
858  }
859 
860  z0 += dz;
861  }
862  y0 += dy;
863  }
864  x0 += dx;
865  }
866  return best;
867 }
868 
869 double ems::MergeEMShower3D::getCos3D(const TVector3& p0, const recob::Track& trk)
870 {
871  TVector3 p1, dir;
872 
873  if (pma::Dist2(p0, trk.Vertex()) < pma::Dist2(p0, trk.End()))
874  p1 = trk.Vertex();
875  else
876  p1 = trk.End();
877 
878  p1 -= p0;
879  double m = p1.Mag();
880  if (m > 0.0)
881  {
882  p1 *= 1.0 / m;
883  double c = fabs(p1 * trk.VertexDirection());
884  if (c > 1.0) c = 1.0;
885  return c;
886  }
887  else return 0.0;
888 }
889 
891 {
892  double sum = 0.0;
893  for (auto ptr : v)
894  {
895  sum += ptr->SummedADC();
896  }
897  return sum;
898 }
899 
double ShowerEnergy(std::vector< art::Ptr< recob::Hit > > const &hits, int plane)
Finds the total energy deposited by the shower in this view.
void reconfigure(fhicl::ParameterSet const &p)
TVector3 GetDir(void) const
Implementation of the Projection Matching Algorithm.
const double & DQdxAtPoint(unsigned int p, geo::View_t view=geo::kUnknown) const
Covariance matrices are either set or not.
Definition: Track.cxx:80
double Dist2(const TVector2 &v1, const TVector2 &v2)
Definition: Utilities.cxx:19
TVector3 VertexDirection() const
Covariance matrices are either set or not.
Definition: Track.h:247
const std::vector< sim::TrackIDE > HitToTrackIDEs(recob::Hit const &hit)
TVector3 GetEnd(void) const
enum geo::_plane_proj View_t
Enumerate the possible plane projections.
double Pointsto(ems::ShowerInfo const &s1) const
ShowersCollection(ShowerInfo const &part)
Planes which measure V.
Definition: geo_types.h:77
Declaration of signal hit object.
void set_total_energy(const std::vector< double > &q)
Definition: Shower.h:129
Planes which measure Z direction.
Definition: geo_types.h:79
double GetP0Dist(void) const
shower::ShowerEnergyAlg fShowerEnergyAlg
std::vector< double > fVentot
std::vector< double > GetDedx(void) const
ShowerInfo(int key, int gid, bool hasvtx, std::vector< double > vensum, recob::Track const &trk)
Definition of vertex object for LArSoft.
Definition: Vertex.h:35
void set_total_energy_err(const std::vector< double > &q)
Definition: Shower.h:130
double getCos3D(const TVector3 &p0, const recob::Track &trk)
void CryostatBoundaries(double *boundaries, geo::CryostatID const &cid) const
Returns the boundaries of the specified cryostat.
int getGammaId(art::Event &evt, const size_t t)
TFile f
Definition: plotHisto.C:6
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265
int getClusterBestId(const std::vector< art::Ptr< recob::Hit > > &v)
double Angleto(TVector3 const &pmapoint) const
Planes which measure U.
Definition: geo_types.h:76
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
Int_t max
Definition: plot.C:27
void hits()
Definition: readHits.C:15
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
double MinAngle(ShowersCollection const &col2)
void beginJob()
Definition: Breakpoints.cc:14
TVector2 GetProjectionToSegment(const TVector2 &p, const TVector2 &p0, const TVector2 &p1)
Definition: Utilities.cxx:167
parameter set interface
double P(const int i=0) const
Definition: MCParticle.h:238
T get(std::string const &key) const
Definition: ParameterSet.h:231
TString part[npart]
Definition: Style.C:32
std::vector< ShowerInfo > const & GetParts(void) const
std::vector< ems::ShowersCollection > collectshowers(art::Event &evt, std::vector< ems::ShowerInfo > &showers, bool refpoint)
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.
double SolveLeastSquares3D(const std::vector< std::pair< TVector3, TVector3 > > &lines, TVector3 &result)
Definition: Utilities.cxx:187
std::vector< double > DeDx(void)
void mcinfo(art::Event &evt)
double Angle(TVector3 p0, TVector3 test)
Declaration of cluster object.
Provides recob::Track data product.
View_t View(geo::PlaneID const &pid) const
Returns the view (wire orientation) on the channels of specified TPC plane.
std::vector< double > fVdqdx
double GetAdcSum(void) const
TVector3 getBestPoint(const std::vector< recob::Track > &tracks, const std::vector< ShowerInfo > &showers, const TVector3 &p0, const TVector3 &p1, double step)
Definition: DirOfGamma.h:18
MergeEMShower3D(fhicl::ParameterSet const &p)
int GetKey(void) const
double Vx(const int i=0) const
Definition: MCParticle.h:225
int GetGid(void) const
Implementation of the Projection Matching Algorithm.
T * make(ARGS...args) const
Utility object to perform functions of association.
int GetPlaneId(void) const
TDirectory * dir
Definition: macro.C:5
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
void Merge(ShowersCollection const &col2)
TVector3 Vertex() const
Covariance matrices are either set or not.
Definition: Track.h:245
double Vz(const int i=0) const
Definition: MCParticle.h:227
EventNumber_t event() const
Definition: EventID.h:117
std::vector< ShowerInfo > fParts
double getClusterAdcSum(const std::vector< art::Ptr< recob::Hit > > &v)
void SetP0Dist(const TVector3 p)
void fill_ptr_vector(std::vector< Ptr< T >> &ptrs, H const &h)
Definition: Ptr.h:464
Float_t e
Definition: plot.C:34
TVector3 GetFront(void) const
std::vector< double > TotEn(void)
TVector3 End() const
Covariance matrices are either set or not.
Definition: Track.h:246
bool HasConPoint(void) const
EventID id() const
Definition: Event.h:56
double Vy(const int i=0) const
Definition: MCParticle.h:226
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< double > GetEnTot(void) const
void produce(art::Event &e) override