LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
MultiEMShowers_module.cc
Go to the documentation of this file.
1 // Class: MultiEMShowers
3 // Module Type: analyzer
4 // File: MultiEMShowers_module.cc
5 // Author: dorota.stefan@cern.ch robert.sulej@cern.ch
7 
13 #include "art_root_io/TFileService.h"
16 #include "fhiclcpp/ParameterSet.h"
18 
28 
29 #include <algorithm>
30 #include <limits>
31 #include <memory>
32 
33 // ROOT includes
34 #include "TLorentzVector.h"
35 #include "TTree.h"
36 
37 namespace ems {
38  class MCinfo;
39  class MultiEMShowers;
40 }
41 
42 class ems::MCinfo {
43 public:
44  MCinfo();
45  void Info();
46  void Findtpcborders();
47 
48  int GetNgammas() const { return fNgammas; }
49 
50  double GetMompi0() const { return fMompi0; }
51  double GetMomGamma1() const { return fGammamom1; }
52  double GetMomGamma2() const { return fGammamom2; }
53 
54  double GetCosine() { return fCosine; }
55 
56  TVector3 const& GetPrimary() const { return fPrimary; }
57  TVector3 const& GetPospi0() const { return fPi0pos; }
58  TVector3 const& GetPosgamma1() const { return fConvgamma1; }
59  TVector3 const& GetPosgamma2() const { return fConvgamma2; }
60 
61  TVector3 const& GetDirgamma1() const { return fDirgamma1; }
62  TVector3 const& GetDirgamma2() const { return fDirgamma2; }
63 
64  bool IsInside1() const { return fInside1; }
65  bool IsInside2() const { return fInside2; }
66 
67  bool IsCompton() const { return fCompton; }
68 
69 private:
70  bool insideFidVol(const TLorentzVector& pvtx) const;
71  double fMinx{std::numeric_limits<double>::max()};
72  double fMaxx{std::numeric_limits<double>::min()};
73  double fMiny{std::numeric_limits<double>::max()};
74  double fMaxy{std::numeric_limits<double>::min()};
75  double fMinz{std::numeric_limits<double>::max()};
76  double fMaxz{std::numeric_limits<double>::min()};
77 
78  double fFidVolCut;
79 
80  int fNgammas;
81 
82  double fMompi0;
83  double fGammamom1;
84  bool fInside1;
85  double fGammamom2;
86  bool fInside2;
87 
88  double fCosine;
89 
90  bool fCompton;
91 
92  TVector3 fPrimary;
93  TVector3 fPi0pos;
94  TVector3 fConvgamma1;
95  TVector3 fConvgamma2;
96  TVector3 fDirgamma1;
97  TVector3 fDirgamma2;
98 };
99 
101 {
102  Info();
103  Findtpcborders();
104 }
105 
107 {
109  for (const geo::TPCGeo& tpcg : geom->Iterate<geo::TPCGeo>()) {
110  fMinx = std::min(fMinx, tpcg.MinX());
111  fMaxx = std::max(fMaxx, tpcg.MaxX());
112  fMiny = std::min(fMiny, tpcg.MinY());
113  fMaxy = std::max(fMaxy, tpcg.MaxY());
114  fMinz = std::min(fMinz, tpcg.MinZ());
115  fMaxz = std::max(fMaxz, tpcg.MaxZ());
116  }
117 }
118 
120 {
121  fMompi0 = 0.0;
122  fPi0pos.SetXYZ(0, 0, 0);
123  fNgammas = 0;
124  fCosine = 0.0;
125  fInside1 = false;
126  fInside2 = false;
127  fCompton = false;
128 
129  fGammamom1 = 0.0;
130  fGammamom2 = 0.0;
131  fConvgamma1.SetXYZ(0, 0, 0);
132  fConvgamma2.SetXYZ(0, 0, 0);
133  fDirgamma1.SetXYZ(0, 0, 0);
134  fDirgamma2.SetXYZ(0, 0, 0);
135 
137  const sim::ParticleList& plist = pi_serv->ParticleList();
138  for (sim::ParticleList::const_iterator ipar = plist.begin(); ipar != plist.end(); ++ipar) {
139  const simb::MCParticle* particle = ipar->second;
140 
141  if (particle->Process() != "primary") continue;
142 
143  TLorentzVector posvec = particle->Position();
144  TVector3 pose(posvec.X(), posvec.Y(), posvec.Z());
145  fPrimary = pose;
146 
147  if (particle->PdgCode() == 111) {
148  fMompi0 = particle->P();
149 
150  TLorentzVector posvec3 = particle->Position();
151  TVector3 pospi0(posvec3.X(), posvec3.Y(), posvec3.Z());
152  fPi0pos = pospi0;
153 
154  if (particle->NumberDaughters() != 2) continue;
155 
156  const simb::MCParticle* daughter1 = pi_serv->TrackIdToParticle_P(particle->Daughter(0));
157  if (daughter1->PdgCode() != 22) continue;
158 
159  const simb::MCParticle* daughter2 = pi_serv->TrackIdToParticle_P(particle->Daughter(1));
160  if (daughter2->PdgCode() != 22) continue;
161 
162  fNgammas = particle->NumberDaughters();
163  TLorentzVector mom1 = pi_serv->TrackIdToParticle_P(particle->Daughter(0))->Momentum();
164  TLorentzVector mom2 = pi_serv->TrackIdToParticle_P(particle->Daughter(1))->Momentum();
165 
166  // compton process
167  if (daughter1->EndProcess() == "phot") fCompton = true;
168  if (daughter2->EndProcess() == "phot") fCompton = true;
169 
170  TVector3 mom1vec3(mom1.Px(), mom1.Py(), mom1.Pz());
171  fGammamom1 = pi_serv->TrackIdToParticle_P(particle->Daughter(0))->P();
172  TVector3 mom2vec3(mom2.Px(), mom2.Py(), mom2.Pz());
173  fGammamom2 = pi_serv->TrackIdToParticle_P(particle->Daughter(1))->P();
174 
175  TLorentzVector pos1 = pi_serv->TrackIdToParticle_P(particle->Daughter(0))->EndPosition();
176  TLorentzVector pos2 = pi_serv->TrackIdToParticle_P(particle->Daughter(1))->EndPosition();
177 
178  if (insideFidVol(pos1)) fInside1 = true;
179  if (insideFidVol(pos2)) fInside2 = true;
180 
181  fConvgamma1.SetXYZ(pos1.X(), pos1.Y(), pos1.Z());
182  fConvgamma2.SetXYZ(pos2.X(), pos2.Y(), pos2.Z());
183 
184  TVector3 vecnorm1 = mom1vec3.Unit();
185  fDirgamma1 = vecnorm1;
186  TVector3 vecnorm2 = mom2vec3.Unit();
187  fDirgamma2 = vecnorm2;
188 
190  }
191  else {
192  fNgammas = particle->NumberDaughters();
193  }
194  }
195 }
196 
197 bool ems::MCinfo::insideFidVol(const TLorentzVector& pvtx) const
198 {
199 
200  bool inside = false;
201  //x
202  double dista = fabs(fMinx - pvtx.X());
203  double distb = fabs(pvtx.X() - fMaxx);
204  if ((pvtx.X() > fMinx) && (pvtx.X() < fMaxx) && (dista > fFidVolCut) && (distb > fFidVolCut))
205  inside = true;
206  //y
207  dista = fabs(fMaxy - pvtx.Y());
208  distb = fabs(pvtx.Y() - fMiny);
209  if (inside && (pvtx.Y() > fMiny) && (pvtx.Y() < fMaxy) && (dista > fFidVolCut) &&
210  (distb > fFidVolCut))
211  inside = true;
212  else
213  inside = false;
214 
215  //z
216  dista = fabs(fMaxz - pvtx.Z());
217  distb = fabs(pvtx.Z() - fMinz);
218  if (inside && (pvtx.Z() > fMinz) && (pvtx.Z() < fMaxz) && (dista > fFidVolCut) &&
219  (distb > fFidVolCut))
220  inside = true;
221  else
222  inside = false;
223 
224  return inside;
225 }
226 
228 public:
229  explicit MultiEMShowers(fhicl::ParameterSet const& p);
230 
231  MultiEMShowers(MultiEMShowers const&) = delete;
232  MultiEMShowers(MultiEMShowers&&) = delete;
233  MultiEMShowers& operator=(MultiEMShowers const&) = delete;
234  MultiEMShowers& operator=(MultiEMShowers&&) = delete;
235 
236 private:
237  void beginJob() override;
238  void endJob() override;
239 
240  void analyze(art::Event const& e) override;
241 
242  bool convCluster(art::Event const& evt);
243  double getMinDist(detinfo::DetectorPropertiesData const& detProp,
245  TVector3 const& convmc,
246  size_t view,
247  size_t tpc,
248  size_t cryo);
253 
254  // ROOT
255  TTree* fEvTree;
257  int fNGroups;
258 
259  // mc
260  double fPi0mom;
261  double fGmom1;
262  double fGmom2;
263  double fMcth;
264  int fNgammas;
266  int fEvComp;
268  int fEvInput;
269  TVector3 fGdir1;
270  TVector3 fGdir2;
271  TVector3 fPrimary;
272 
273  //reco
274  int fEvReco;
276  int fEv2Good;
277  int fCountph;
279 
280  TTree* fShTree;
281  TTree* fRecoTree;
282  double fStartX;
283  double fStartY;
284  double fStartZ;
285  double fDedxZ;
286  double fDedxV;
287  double fDedxU;
288  double fMCrecovtx;
289  double fMCrecoTh;
292  double fRecth;
293  double fRecthgood;
296  double fGdirmcreco1;
297  double fGdirmcreco2;
300 
306 };
307 
309 {
310  fConvGood = 0;
311  fConvWrong = 0;
312  fConvBothGood = 0;
313  fEvFidVol = 0;
314  fEvComp = 0;
315  fEvGMomCut = 0;
316  fEvReco = 0;
317  fEvInput = 0;
318  fEv2Groups = 0;
319  fEv2Good = 0;
320  fHitsModuleLabel = p.get<art::InputTag>("HitsModuleLabel");
321  fCluModuleLabel = p.get<art::InputTag>("ClustersModuleLabel");
322  fTrk3DModuleLabel = p.get<art::InputTag>("Trk3DModuleLabel");
323  fVtxModuleLabel = p.get<art::InputTag>("VtxModuleLabel");
324  fShsModuleLabel = p.get<art::InputTag>("ShsModuleLabel");
325 }
326 
328 {
330 
331  fEvTree = tfs->make<TTree>("MultiShowers", "showers3d");
332  fEvTree->Branch("fEvNumber", &fEvNumber, "fEvNumber/I");
333  fEvTree->Branch("fNGroups", &fNGroups, "fNGroups/I");
334  fEvTree->Branch("fPi0mom", &fPi0mom, "fPi0mom/D");
335  fEvTree->Branch("fNgammas", &fNgammas, "fNgammas/I");
336  fEvTree->Branch("fGmom1", &fGmom1, "fGmom1/D");
337  fEvTree->Branch("fGmom2", &fGmom2, "fGmom2/D");
338  fEvTree->Branch("fMcth", &fMcth, "fMcth/D");
339  fEvTree->Branch("fRecth", &fRecth, "fRecth/D");
340  fEvTree->Branch("fMCrecovtx", &fMCrecovtx, "fMCrecovtx/D");
341  fEvTree->Branch("fMCrecoTh", &fMCrecoTh, "fMCrecoTh/D");
342  fEvTree->Branch("fConvGood", &fConvGood, "fConvGood/I");
343  fEvTree->Branch("fConvWrong", &fConvWrong, "fConvWrong/I");
344  fEvTree->Branch("fConvBothGood", &fConvBothGood, "fConvBothGood/I");
345  fEvTree->Branch("fGammasInside", &fGammasInside, "fGammasInside/I");
346  fEvTree->Branch("fCountph", &fCountph, "fCountph/I");
347  fEvTree->Branch("fCountreco", &fCountreco, "fCountreco/I");
348 
349  fRecoTree = tfs->make<TTree>("Cascades", "conv points");
350  fRecoTree->Branch("fRecth", &fRecth, "fRecth/D");
351  fRecoTree->Branch("fMCrecovtx", &fMCrecovtx, "fMCrecovtx/D");
352  fRecoTree->Branch("fMCrecoTh", &fMCrecoTh, "fMCrecoTh/D");
353  fRecoTree->Branch("fRecthgood", &fRecthgood, "fRecthgood/D");
354  fRecoTree->Branch("fMCrecovtxgood", &fMCrecovtxgood, "fMCrecovtxgood/D");
355  fRecoTree->Branch("fMCrecoThgood", &fMCrecoThgood, "fMCrecoThgood/D");
356  fRecoTree->Branch("fGdirmcreco1", &fGdirmcreco1, "fGdirmcreco1/D");
357  fRecoTree->Branch("fGdirmcreco2", &fGdirmcreco2, "fGdirmcreco2/D");
358  fRecoTree->Branch("fGdirmcreco1good", &fGdirmcreco1good, "fGdirmcreco1good/D");
359  fRecoTree->Branch("fGdirmcreco2good", &fGdirmcreco2good, "fGdirmcreco2good/D");
360 
361  fShTree = tfs->make<TTree>("Shower", "conv point");
362 
363  fShTree->Branch("fStartX", &fStartX, "fStartX/D");
364  fShTree->Branch("fStartY", &fStartY, "fStartY/D");
365  fShTree->Branch("fStartZ", &fStartZ, "fStartZ/D");
366  fShTree->Branch("fDedxZ", &fDedxZ, "fDedxZ/D");
367  fShTree->Branch("fDedxV", &fDedxV, "fDedxV/D");
368  fShTree->Branch("fDedxU", &fDedxU, "fDedxU/D");
369 }
370 
372 {
373  mf::LogInfo log("MultiEMShower");
374  log << "******************** fEvFidVol = " << fEvFidVol << "\n";
375  log << "******************** fEvGMomCut = " << fEvGMomCut << "\n";
376  log << "******************** fEvComp = " << fEvComp << "\n";
377  log << "******************** fEvReco = " << fEvReco << "\n";
378  log << "******************** fEvInput = " << fEvInput << "\n";
379  log << "******************** fEv2Groups = " << fEv2Groups << "\n";
380  log << "******************** fEv2Good = " << fEv2Good << "\n";
381  if (fEvInput)
382  log << "******************** reco % = " << double(fEvReco) / double(fEvInput) << "\n";
383 }
384 
386 {
387  fEvNumber = e.id().event();
388  fNGroups = 0;
389  fStartX = 0.0;
390  fStartY = 0.0;
391  fStartZ = 0.0;
392  fPi0mom = 0.0;
393  fNgammas = 0;
394  fDistConvrecomc1 = 0.0;
395  fDistConvrecomc2 = 0.0;
396  fMCrecovtx = -400.0;
397  fMCrecovtxgood = -400.0;
398  fRecth = -400.0;
399  fRecthgood = -400.0;
400  fMCrecoTh = -400.0;
401  fMCrecoThgood = -400.0;
402  fGammasInside = 0;
403  fCountph = 0;
404  fCountreco = 0;
405  fGdirmcreco1 = 0.0;
406  fGdirmcreco2 = 0.0;
407  fGdirmcreco1good = 0.0;
408  fGdirmcreco2good = 0.0;
409  fDedxZ = 0.0;
410  fDedxV = 0.0;
411  fDedxU = 0.0;
412 
413  ems::MCinfo mc;
414  fPrimary = mc.GetPrimary();
415  fPi0mom = mc.GetMompi0();
416  fGmom1 = mc.GetMomGamma1();
417  fGmom2 = mc.GetMomGamma2();
418  fGdir1 = mc.GetDirgamma1();
419  fGdir2 = mc.GetDirgamma2();
420  fNgammas = mc.GetNgammas();
421  TVector3 pospi0 = mc.GetPospi0();
422 
423  double cosinemc = mc.GetCosine();
424  fMcth = 180.0F * (std::acos(cosinemc)) / TMath::Pi();
425  TVector3 convp[2];
426  convp[0] = mc.GetPosgamma1();
427  convp[1] = mc.GetPosgamma2();
428  const double maxdist = 2.0; //cm
429 
430  // check whether two photons are inside fid vol
431  if (mc.IsInside1() && mc.IsInside2()) {
432  fGammasInside = 1;
433  fEvFidVol++;
434  }
435 
436  if ((fGmom1 > 0.1) && (fGmom2 > 0.1)) fEvGMomCut++;
437 
438  if (mc.IsCompton()) fEvComp++;
439 
445  if (e.getByLabel(fShsModuleLabel, shsListHandle) &&
446  e.getByLabel(fTrk3DModuleLabel, trkListHandle) &&
447  e.getByLabel(fVtxModuleLabel, vtxListHandle) &&
448  e.getByLabel(fCluModuleLabel, cluListHandle) &&
449  e.getByLabel(fHitsModuleLabel, hitListHandle)) {
450  art::FindManyP<recob::Cluster> cluFromShs(shsListHandle, e, fShsModuleLabel);
451  art::FindManyP<recob::Cluster> cluFromTrk(trkListHandle, e, fTrk3DModuleLabel);
452  art::FindManyP<recob::Vertex> vtxFromTrk(trkListHandle, e, fVtxModuleLabel);
453  art::FindManyP<recob::Hit> hitFromClu(cluListHandle, e, fCluModuleLabel);
454 
455  fNGroups = shsListHandle->size();
456 
457  fCountph = 0;
458  if (fNgammas == 2) // pi0
459  {
460  int idph = -1;
461  for (size_t s = 0; s < shsListHandle->size(); ++s) {
462  const recob::Shower& sh = (*shsListHandle)[s];
463  double mindist = maxdist;
464  bool found = false;
465 
466  for (int i = 0; i < fNgammas; ++i) {
467  double dist = sqrt(pma::Dist2(sh.ShowerStart(), convp[i]));
468  if ((dist < mindist) && (idph != i)) {
469  mindist = dist;
470  idph = i;
471  found = true;
472  }
473  }
474  if (found) {
475  fConvGood++;
476  fCountph++;
477  }
478  else {
479  fConvWrong++;
480  }
481  }
482  if (fCountph == 2) fConvBothGood++;
483 
484  // plot a few variables if there are 2 showers
485  if (fCountph == 2)
486  for (size_t s = 0; s < shsListHandle->size(); ++s) {
487  const recob::Shower& sh = (*shsListHandle)[s];
488  TVector3 pos = sh.ShowerStart();
489  fStartX = pos.X();
490  fStartY = pos.Y();
491  fStartZ = pos.Z();
492  std::vector<double> const& vecdedx = sh.dEdx();
493 
494  if (vecdedx.size() == 3) {
495  fDedxZ = vecdedx[0];
496  fDedxV = vecdedx[1];
497  fDedxU = vecdedx[2];
498  }
499 
500  fShTree->Fill();
501  }
502  }
503  else // other than pi0
504  {
505  for (size_t s = 0; s < shsListHandle->size(); ++s) {
506  const recob::Shower& sh = (*shsListHandle)[s];
507  double mindist = maxdist;
508 
509  double dist = sqrt(pma::Dist2(sh.ShowerStart(), fPrimary));
510  if (dist < mindist) {
511  TVector3 pos = sh.ShowerStart();
512  fStartX = pos.X();
513  fStartY = pos.Y();
514  fStartZ = pos.Z();
515  std::vector<double> vecdedx = sh.dEdx();
516  if (vecdedx.size() == 3) {
517  fDedxZ = vecdedx[0];
518  fDedxV = vecdedx[1];
519  fDedxU = vecdedx[2];
520  }
521  }
522 
523  fShTree->Fill();
524  }
525  }
526  // compute the crossing point
527 
528  //cut from mc and clusters
529 
530  if (mc.IsInside1() && mc.IsInside2() && (fGmom1 > 0.1) && (fGmom2 > 0.1) && (!mc.IsCompton()) &&
531  convCluster(e)) {
532  fCountreco = 1;
533  if (fNGroups == 2) fEv2Groups++;
534  if ((fNGroups == 2) && (fCountph == 2)) fEv2Good++;
535  // cut from reco
536  //if (countph == 2)
537  if (fNGroups == 2) {
538  std::vector<std::pair<TVector3, TVector3>> lines;
539  const recob::Shower& sh1 = (*shsListHandle)[0];
540  const recob::Shower& sh2 = (*shsListHandle)[1];
541 
542  std::pair<TVector3, TVector3> frontback1(sh1.ShowerStart(),
543  sh1.ShowerStart() + sh1.Direction());
544  std::pair<TVector3, TVector3> frontback2(sh2.ShowerStart(),
545  sh2.ShowerStart() + sh2.Direction());
546  lines.push_back(frontback1);
547  lines.push_back(frontback2);
548 
549  TVector3 result;
550  pma::SolveLeastSquares3D(lines, result); // mse.
551 
552  double dist1_0 = pma::Dist2(result, sh1.ShowerStart());
553  double dist2_0 = pma::Dist2(result, sh1.ShowerStart() + sh1.Direction());
554  double dist1_1 = pma::Dist2(result, sh2.ShowerStart());
555  double dist2_1 = pma::Dist2(result, sh2.ShowerStart() + sh2.Direction());
556  if ((dist1_0 > dist2_0) || (dist1_1 > dist2_1)) {}
557  else {
558  fMCrecovtx = std::sqrt(pma::Dist2(pospi0, result));
559 
560  if (fCountph == 2) fMCrecovtxgood = fMCrecovtx;
561 
562  double cosine_reco = sh1.Direction() * sh2.Direction();
563  fRecth = 180.0F * (std::acos(cosine_reco)) / TMath::Pi();
564 
565  fGdirmcreco1 = fGdir1 * sh1.Direction();
566  fGdirmcreco2 = fGdir2 * sh2.Direction();
567  if (fCountph == 2) {
570  }
571 
572  if (fCountph == 2) fRecthgood = fRecth;
573 
574  fMCrecoTh = fRecth - fMcth;
575 
576  if (fCountph == 2) fMCrecoThgood = fMCrecoTh;
577 
578  fEvReco++;
579  fRecoTree->Fill();
580  }
581  }
582  fEvInput++;
583  //fRecoTree->Fill();
584  }
585  }
586 
587  fEvTree->Fill();
588 }
589 
590 // true if there are clusters corresponding to mc conversion points
592 {
593  ems::MCinfo mc;
594  TVector3 const convp[2]{mc.GetPosgamma1(), mc.GetPosgamma2()};
595 
596  geo::Point_t const vtx{convp[0].X(), convp[0].Y(), convp[0].Z()};
597 
599  geo::TPCID idtpc = geom->FindTPCAtPosition(vtx);
600  size_t cryoid = geom->PositionToCryostatID(vtx).Cryostat;
601 
604 
605  //map: conversion point, vec of id clusters in each view
606  std::map<size_t, std::vector<size_t>> used;
607 
608  art::FindManyP<recob::Hit> fbc(cluListHandle, evt, fCluModuleLabel);
609 
610  auto const detProp = art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt);
611  double maxdist = 1.0; // 1 cm
612  if (geom->HasTPC(idtpc)) {
613  const geo::CryostatGeo& cryostat = geom->Cryostat(geo::CryostatID(cryoid));
614  if (evt.getByLabel(fHitsModuleLabel, hitListHandle) &&
615  evt.getByLabel(fCluModuleLabel, cluListHandle)) {
616  size_t conv = 0;
617  while (conv < 2) {
618  for (size_t view = 0; view < cryostat.MaxPlanes(); view++) {
619 
620  double mindist = maxdist;
621  int clid = 0;
622  for (size_t c = 0; c < cluListHandle->size(); ++c) {
623 
624  bool exist = false;
625  for (auto const& ids : used)
626  for (auto i : ids.second)
627  if (i == c) exist = true;
628  if (exist) continue;
629 
630  std::vector<art::Ptr<recob::Hit>> hits = fbc.at(c);
631  if (hits.size() < 20) continue;
632  if (hits[0]->WireID().Plane != view) continue;
633 
634  double dist = getMinDist(detProp, hits, convp[conv], view, idtpc.TPC, cryoid);
635  if (dist < mindist) {
636  mindist = dist;
637  clid = c;
638  }
639  }
640  if (mindist < maxdist) used[conv].push_back(clid);
641  }
642  conv++;
643  }
644  }
645  }
646  bool result = false;
647 
648  if (used.size() > 1)
649  for (auto const& ids : used) {
650  if (ids.second.size() > 1)
651  result = true;
652  else {
653  result = false;
654  break;
655  }
656  }
657 
658  return result;
659 }
660 
663  TVector3 const& convmc,
664  size_t view,
665  size_t tpc,
666  size_t cryo)
667 {
668  double mindist = 9999;
669  // MC vertex projected to view
670  TVector2 proj = pma::GetProjectionToPlane(convmc, view, tpc, cryo);
671 
672  // loop over hits to find the closest to MC 2d vtx
673  for (size_t h = 0; h < v.size(); ++h) {
674  if ((v[h]->WireID().Plane == view) && (v[h]->WireID().TPC == tpc)) {
675  TVector2 hpoint =
676  pma::WireDriftToCm(detProp, v[h]->WireID().Wire, v[h]->PeakTime(), view, tpc, cryo);
677 
678  double dist = pma::Dist2(proj, hpoint);
679  if (dist < mindist) { mindist = dist; }
680  }
681  }
682 
683  return mindist;
684 }
685 
details::range_type< T > Iterate() const
Initializes the specified ID with the ID of the first cryostat.
Definition: GeometryCore.h:541
const TVector3 & ShowerStart() const
Definition: Shower.h:197
TVector3 const & GetDirgamma2() const
const TLorentzVector & Position(const int i=0) const
Definition: MCParticle.h:220
int PdgCode() const
Definition: MCParticle.h:213
double SolveLeastSquares3D(const std::vector< std::pair< TVector3, TVector3 >> &lines, TVector3 &result)
Definition: Utilities.cxx:158
unsigned int MaxPlanes() const
Returns the largest number of planes among the TPCs in this cryostat.
double Dist2(const TVector2 &v1, const TVector2 &v2)
Definition: Utilities.cxx:39
const simb::MCParticle * TrackIdToParticle_P(int id) const
MultiEMShowers(fhicl::ParameterSet const &p)
Declaration of signal hit object.
list_type::const_iterator const_iterator
Definition: ParticleList.h:132
Geometry information for a single TPC.
Definition: TPCGeo.h:36
void analyze(art::Event const &e) override
bool convCluster(art::Event const &evt)
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
TVector2 WireDriftToCm(detinfo::DetectorPropertiesData const &detProp, unsigned int wire, float drift, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:286
Geometry information for a single cryostat.
Definition: CryostatGeo.h:43
std::string Process() const
Definition: MCParticle.h:216
CryostatGeo const & Cryostat(CryostatID const &cryoid=cryostat_zero) const
Returns the specified cryostat.
double GetMomGamma2() const
int NumberDaughters() const
Definition: MCParticle.h:218
double GetMompi0() const
int Daughter(const int i) const
Definition: MCParticle.cxx:118
double getMinDist(detinfo::DetectorPropertiesData const &detProp, std::vector< art::Ptr< recob::Hit >> const &v, TVector3 const &convmc, size_t view, size_t tpc, size_t cryo)
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
void hits()
Definition: readHits.C:15
TPCID FindTPCAtPosition(Point_t const &point) const
Returns the ID of the TPC at specified location.
const std::vector< double > & dEdx() const
Definition: Shower.h:235
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:65
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
std::string EndProcess() const
Definition: MCParticle.h:217
void beginJob()
Definition: Breakpoints.cc:14
TVector3 const & GetPosgamma1() const
double P(const int i=0) const
Definition: MCParticle.h:235
T get(std::string const &key) const
Definition: ParameterSet.h:314
int GetNgammas() const
iterator begin()
Definition: ParticleList.h:305
Provides recob::Track data product.
const TVector3 & Direction() const
Definition: Shower.h:188
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
Declaration of cluster object.
TVector3 const & GetPrimary() const
TVector3 const & GetDirgamma1() const
Definition: DirOfGamma.h:22
const sim::ParticleList & ParticleList() const
bool IsInside2() const
CryostatID PositionToCryostatID(Point_t const &point) const
Returns the ID of the cryostat at specified location.
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double >, ROOT::Math::GlobalCoordinateSystemTag > Point_t
Type for representation of position in physical 3D space.
Definition: geo_vectors.h:180
bool IsCompton() const
Implementation of the Projection Matching Algorithm.
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
TVector3 const & GetPospi0() const
double GetMomGamma1() const
constexpr double dist(const TReal *x, const TReal *y, const unsigned int dimension)
TVector2 GetProjectionToPlane(const TVector3 &p, unsigned int plane, unsigned int tpc, unsigned int cryo)
Definition: Utilities.cxx:263
Float_t proj
Definition: plot.C:35
EventNumber_t event() const
Definition: EventID.h:116
TCEvent evt
Definition: DataStructs.cxx:8
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
Float_t e
Definition: plot.C:35
TVector3 const & GetPosgamma2() const
recob::tracking::Plane Plane
Definition: TrackState.h:17
bool IsInside1() const
bool HasTPC(TPCID const &tpcid) const
Returns whether we have the specified TPC.
Definition: GeometryCore.h:699
EventID id() const
Definition: Event.cc:23
art framework interface to geometry description
The data type to uniquely identify a cryostat.
Definition: geo_types.h:192
bool insideFidVol(const TLorentzVector &pvtx) const