All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Public Member Functions | Private Member Functions | Private Attributes | List of all members
mvapid::MVAAlg Class Reference

#include <MVAAlg.h>

Classes

struct  SortedObj
 
struct  SumDistance2
 

Public Member Functions

 MVAAlg (fhicl::ParameterSet const &pset)
 
void GetDetectorEdges ()
 
void GetWireNormals ()
 
void RunPID (art::Event &evt, std::vector< anab::MVAPIDResult > &result, art::Assns< recob::Track, anab::MVAPIDResult, void > &trackAssns, art::Assns< recob::Shower, anab::MVAPIDResult, void > &showerAssns)
 

Private Member Functions

int IsInActiveVol (const TVector3 &pos)
 
void PrepareEvent (const art::Event &event, const detinfo::DetectorClocksData &clockData)
 
void FitAndSortTrack (art::Ptr< recob::Track > track, int &isStoppingReco, SortedObj &sortedObj)
 
void SortShower (art::Ptr< recob::Shower > shower, int &isStoppingReco, mvapid::MVAAlg::SortedObj &sortedShower)
 
void RunPCA (std::vector< art::Ptr< recob::Hit >> &hits, std::vector< double > &eVals, std::vector< double > &eVecs)
 
void _Var_Shape (const SortedObj &track, double &coreHaloRatio, double &concentration, double &conicalness)
 
double CalcSegmentdEdxFrac (const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
 
double CalcSegmentdEdxDist (const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
 
double CalcSegmentdEdxDistAtEnd (const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const mvapid::MVAAlg::SortedObj &track, double distAtEnd)
 
int LinFit (const art::Ptr< recob::Track > track, TVector3 &trackPoint, TVector3 &trackDir)
 
int LinFitShower (const art::Ptr< recob::Shower > shower, TVector3 &showerPoint, TVector3 &showerDir)
 

Private Attributes

const calo::CalorimetryAlg fCaloAlg
 
double fEventT0
 
double fDetMinX
 
double fDetMaxX
 
double fDetMinY
 
double fDetMaxY
 
double fDetMinZ
 
double fDetMaxZ
 
std::map< int, double > fNormToWiresY
 
std::map< int, double > fNormToWiresZ
 
std::string fTrackLabel
 
std::string fShowerLabel
 
std::string fHitLabel
 
std::string fSpacePointLabel
 
std::string fTrackingLabel
 
std::vector< art::Ptr
< recob::Track > > 
fTracks
 
std::vector< art::Ptr
< recob::Shower > > 
fShowers
 
std::vector< art::Ptr
< recob::SpacePoint > > 
fSpacePoints
 
std::vector< art::Ptr
< recob::Hit > > 
fHits
 
std::map< art::Ptr
< recob::Track >, std::vector
< art::Ptr< recob::Hit > > > 
fTracksToHits
 
std::map< art::Ptr
< recob::Track >, std::vector
< art::Ptr< recob::SpacePoint > > > 
fTracksToSpacePoints
 
std::map< art::Ptr
< recob::Shower >, std::vector
< art::Ptr< recob::Hit > > > 
fShowersToHits
 
std::map< art::Ptr
< recob::Shower >, std::vector
< art::Ptr< recob::SpacePoint > > > 
fShowersToSpacePoints
 
std::map< art::Ptr< recob::Hit >
, art::Ptr< recob::SpacePoint > > 
fHitsToSpacePoints
 
std::map< art::Ptr
< recob::SpacePoint >
, art::Ptr< recob::Hit > > 
fSpacePointsToHits
 
anab::MVAPIDResult fResHolder
 
TMVA::Reader fReader
 
std::vector< std::string > fMVAMethods
 
std::vector< std::string > fWeightFiles
 
bool fCheatVertex
 
TLorentzVector fVertex4Vect
 

Detailed Description

Definition at line 40 of file MVAAlg.h.

Constructor & Destructor Documentation

mvapid::MVAAlg::MVAAlg ( fhicl::ParameterSet const &  pset)

Definition at line 33 of file MVAAlg.cxx.

34  : fCaloAlg(pset.get<fhicl::ParameterSet>("CalorimetryAlg")), fReader("")
35 {
36  fHitLabel = pset.get<std::string>("HitLabel");
37  fTrackLabel = pset.get<std::string>("TrackLabel");
38  fShowerLabel = pset.get<std::string>("ShowerLabel");
39  fSpacePointLabel = pset.get<std::string>("SpacePointLabel");
40  fTrackingLabel = pset.get<std::string>("TrackingLabel", "");
41 
42  fCheatVertex = pset.get<bool>("CheatVertex", false);
43 
44  fReader.AddVariable("evalRatio", &fResHolder.evalRatio);
45  fReader.AddVariable("coreHaloRatio", &fResHolder.coreHaloRatio);
46  fReader.AddVariable("concentration", &fResHolder.concentration);
47  fReader.AddVariable("conicalness", &fResHolder.conicalness);
48  fReader.AddVariable("dEdxStart", &fResHolder.dEdxStart);
49  fReader.AddVariable("dEdxEnd", &fResHolder.dEdxEnd);
50  fReader.AddVariable("dEdxEndRatio", &fResHolder.dEdxEndRatio);
51 
52  fMVAMethods = pset.get<std::vector<std::string>>("MVAMethods");
53  std::vector<std::string> weightFileBnames = pset.get<std::vector<std::string>>("WeightFiles");
54 
55  cet::search_path searchPath("FW_SEARCH_PATH");
56  for (auto fileIter = weightFileBnames.begin(); fileIter != weightFileBnames.end(); ++fileIter) {
57  std::string fileWithPath;
58  if (!searchPath.find_file(*fileIter, fileWithPath)) {
59  fWeightFiles.clear();
60  fMVAMethods.clear();
61  throw cet::exception("MVAPID") << "Unable to find weight file " << *fileIter
62  << " in search path " << searchPath.to_string() << std::endl;
63  }
64  fWeightFiles.push_back(fileWithPath);
65  }
66 
67  if (fMVAMethods.size() != fWeightFiles.size()) {
68  std::cerr << "Mismatch in number of MVA methods and weight files!" << std::endl;
69  exit(1);
70  }
71 
72  for (unsigned int iMethod = 0; iMethod != fMVAMethods.size(); ++iMethod) {
73  fReader.BookMVA(fMVAMethods[iMethod], fWeightFiles[iMethod]);
74  }
75 }
anab::MVAPIDResult fResHolder
Definition: MVAAlg.h:160
const calo::CalorimetryAlg fCaloAlg
Definition: MVAAlg.h:132
BEGIN_PROLOG could also be cerr
std::vector< std::string > fMVAMethods
Definition: MVAAlg.h:164
std::string fSpacePointLabel
Definition: MVAAlg.h:144
TMVA::Reader fReader
Definition: MVAAlg.h:162
bool fCheatVertex
Definition: MVAAlg.h:167
std::string fTrackingLabel
Definition: MVAAlg.h:145
std::string fTrackLabel
Definition: MVAAlg.h:141
std::string fHitLabel
Definition: MVAAlg.h:143
std::string fShowerLabel
Definition: MVAAlg.h:142
std::vector< std::string > fWeightFiles
Definition: MVAAlg.h:165

Member Function Documentation

void mvapid::MVAAlg::_Var_Shape ( const SortedObj track,
double &  coreHaloRatio,
double &  concentration,
double &  conicalness 
)
private

Definition at line 547 of file MVAAlg.cxx.

551 {
552 
553  static const unsigned int conMinHits = 3;
554  static const double minCharge = 0.1;
555  static const double conFracRange = 0.2;
556  static const double MoliereRadius = 10.1;
557  static const double MoliereRadiusFraction = 0.2;
558 
559  double totalCharge = 0;
560  double totalChargeStart = 0;
561  double totalChargeEnd = 0;
562 
563  double chargeCore = 0;
564  double chargeHalo = 0;
565  double chargeCon = 0;
566  unsigned int nHits = 0;
567 
568  //stuff for conicalness
569  double chargeConStart = 0;
570  double chargeConEnd = 0;
571  unsigned int nHitsConStart = 0;
572  unsigned int nHitsConEnd = 0;
573 
574  for (auto hitIter = track.hitMap.begin(); hitIter != track.hitMap.end(); ++hitIter) {
575  if (fHitsToSpacePoints.count(hitIter->second)) {
576  art::Ptr<recob::SpacePoint> sp = fHitsToSpacePoints.at(hitIter->second);
577 
578  double distFromTrackFit = ((TVector3(sp->XYZ()) - track.start).Cross(track.dir)).Mag();
579 
580  ++nHits;
581 
582  if (distFromTrackFit < MoliereRadiusFraction * MoliereRadius)
583  chargeCore += hitIter->second->Integral();
584  else
585  chargeHalo += hitIter->second->Integral();
586 
587  totalCharge += hitIter->second->Integral();
588 
589  chargeCon += hitIter->second->Integral() / std::max(1.E-2, distFromTrackFit);
590  if (hitIter->first / track.length < conFracRange) {
591  chargeConStart += distFromTrackFit * distFromTrackFit * hitIter->second->Integral();
592  ++nHitsConStart;
593  totalChargeStart += hitIter->second->Integral();
594  }
595  else if (1. - hitIter->first / track.length < conFracRange) {
596  chargeConEnd += distFromTrackFit * distFromTrackFit * hitIter->second->Integral();
597  ++nHitsConEnd;
598  totalChargeEnd += hitIter->second->Integral();
599  }
600  }
601  }
602 
603  coreHaloRatio = chargeHalo / TMath::Max(1.0E-3, chargeCore);
604  coreHaloRatio = TMath::Min(100.0, coreHaloRatio);
605  concentration = chargeCon / totalCharge;
606  if (nHitsConStart >= conMinHits && nHitsConEnd >= conMinHits && totalChargeEnd > minCharge &&
607  sqrt(chargeConStart) > minCharge && totalChargeStart > minCharge) {
608  conicalness = (sqrt(chargeConEnd) / totalChargeEnd) / (sqrt(chargeConStart) / totalChargeStart);
609  }
610  else {
611  conicalness = 1.;
612  }
613 }
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:157
process_name E
process_name use argoneut_mc_hitfinder track
double mvapid::MVAAlg::CalcSegmentdEdxDist ( const detinfo::DetectorClocksData clock_data,
const detinfo::DetectorPropertiesData det_prop,
const SortedObj track,
double  start,
double  end 
)
private

Definition at line 639 of file MVAAlg.cxx.

644 {
645  art::ServiceHandle<geo::Geometry const> geom;
646 
647  double totaldEdx = 0;
648  unsigned int nHits = 0;
649 
650  //Loop over hits again to calculate average dE/dx and shape variables
651  for (auto hitIter = track.hitMap.begin(); hitIter != track.hitMap.end(); ++hitIter) {
652 
653  if (hitIter->first < start) continue;
654  if (hitIter->first >= end) break;
655 
656  art::Ptr<recob::Hit> hit = hitIter->second;
657 
658  //Pitch to use in dEdx calculation
659  double yzPitch =
660  geom->WirePitch(hit->WireID().Plane,
661  hit->WireID().TPC); //pitch not taking into account angle of track or shower
662  double xComponent, pitch3D;
663 
664  TVector3 dir = track.dir;
665 
666  //This assumes equal numbers of TPCs in each cryostat and equal numbers of planes in each TPC
667  int planeKey = hit->WireID().Cryostat * geom->NTPC(0) * geom->Nplanes(0, 0) +
668  hit->WireID().TPC * geom->Nplanes(0, 0) + hit->WireID().Plane;
669 
670  if (fNormToWiresY.count(planeKey) && fNormToWiresZ.count(planeKey)) {
671  TVector3 normToWires(0.0, fNormToWiresY.at(planeKey), fNormToWiresZ.at(planeKey));
672  yzPitch =
673  geom->WirePitch(hit->WireID().Plane, hit->WireID().TPC) / fabs(dir.Dot(normToWires));
674  }
675 
676  xComponent = yzPitch * dir[0] / sqrt(dir[1] * dir[1] + dir[2] * dir[2]);
677  pitch3D = sqrt(xComponent * xComponent + yzPitch * yzPitch);
678 
679  double dEdx = fCaloAlg.dEdx_AREA(clock_data, det_prop, *hit, pitch3D, fEventT0);
680  if (dEdx < 50.) {
681  ++nHits;
682  totaldEdx += dEdx;
683  }
684  }
685 
686  return nHits ? totaldEdx / nHits : 0;
687 }
const calo::CalorimetryAlg fCaloAlg
Definition: MVAAlg.h:132
process_name use argoneut_mc_hitfinder track
process_name hit
Definition: cheaterreco.fcl:51
std::map< int, double > fNormToWiresY
Definition: MVAAlg.h:138
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
float dEdx(detinfo::DetectorClocksData const &clockData, detinfo::DetectorPropertiesData const &detProp, const TCSlice &slc, TP3D &tp3d)
Definition: PFPUtils.cxx:2687
tuple dir
Definition: dropbox.py:28
double dEdx_AREA(detinfo::DetectorClocksData const &clock_data, detinfo::DetectorPropertiesData const &det_prop, recob::Hit const &hit, double pitch, double T0=0) const
std::map< int, double > fNormToWiresZ
Definition: MVAAlg.h:139
double fEventT0
Definition: MVAAlg.h:134
double mvapid::MVAAlg::CalcSegmentdEdxDistAtEnd ( const detinfo::DetectorClocksData clock_data,
const detinfo::DetectorPropertiesData det_prop,
const mvapid::MVAAlg::SortedObj track,
double  distAtEnd 
)
private

Definition at line 628 of file MVAAlg.cxx.

632 {
633 
634  double trackLength = (track.end - track.start).Mag();
635  return CalcSegmentdEdxDist(clock_data, det_prop, track, trackLength - distAtEnd, trackLength);
636 }
double CalcSegmentdEdxDist(const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
Definition: MVAAlg.cxx:639
double mvapid::MVAAlg::CalcSegmentdEdxFrac ( const detinfo::DetectorClocksData clock_data,
const detinfo::DetectorPropertiesData det_prop,
const SortedObj track,
double  start,
double  end 
)
private

Definition at line 616 of file MVAAlg.cxx.

621 {
622 
623  double trackLength = (track.end - track.start).Mag();
624  return CalcSegmentdEdxDist(clock_data, det_prop, track, start * trackLength, end * trackLength);
625 }
process_name use argoneut_mc_hitfinder track
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
double CalcSegmentdEdxDist(const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
Definition: MVAAlg.cxx:639
void mvapid::MVAAlg::FitAndSortTrack ( art::Ptr< recob::Track track,
int &  isStoppingReco,
SortedObj sortedObj 
)
private

Definition at line 348 of file MVAAlg.cxx.

351 {
352 
353  sortedTrack.hitMap.clear();
354  TVector3 trackPoint, trackDir;
355  this->LinFit(track, trackPoint, trackDir);
356 
357  TVector3 nearestPointStart, nearestPointEnd;
358 
359  //For single-particle events can opt to cheat vertex from start of primary trajectory.
360  //Ok since in real events it should be possible to identify the true vertex.
361  if (fCheatVertex) {
362  if ((track->End<TVector3>() - fVertex4Vect.Vect()).Mag() >
363  (track->Vertex<TVector3>() - fVertex4Vect.Vect()).Mag()) {
364  nearestPointStart =
365  trackPoint +
366  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
367  nearestPointEnd = trackPoint + trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) /
368  trackDir.Mag2());
369  isStoppingReco = this->IsInActiveVol(track->End<TVector3>());
370  }
371  else {
372  nearestPointStart =
373  trackPoint +
374  trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) / trackDir.Mag2());
375  nearestPointEnd =
376  trackPoint +
377  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
378  isStoppingReco = this->IsInActiveVol(track->Vertex<TVector3>());
379  trackDir *= -1.;
380  }
381  }
382  else {
383  if (track->End<TVector3>().Z() >=
384  track->Vertex<TVector3>().Z()) { //Otherwise assume particle is forward-going for now...
385  nearestPointStart =
386  trackPoint +
387  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
388  nearestPointEnd = trackPoint + trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) /
389  trackDir.Mag2());
390  isStoppingReco = this->IsInActiveVol(track->End<TVector3>());
391  }
392  else {
393  nearestPointStart =
394  trackPoint +
395  trackDir * (trackDir.Dot(track->End<TVector3>() - trackPoint) / trackDir.Mag2());
396  nearestPointEnd =
397  trackPoint +
398  trackDir * (trackDir.Dot(track->Vertex<TVector3>() - trackPoint) / trackDir.Mag2());
399  isStoppingReco = this->IsInActiveVol(track->Vertex<TVector3>());
400  }
401 
402  if (trackDir.Z() <= 0) {
403  trackDir.SetX(-trackDir.X());
404  trackDir.SetY(-trackDir.Y());
405  trackDir.SetZ(-trackDir.Z());
406  }
407  }
408 
409  sortedTrack.start = nearestPointStart;
410  sortedTrack.end = nearestPointEnd;
411  sortedTrack.dir = trackDir;
412  sortedTrack.length = (nearestPointEnd - nearestPointStart).Mag();
413 
414  std::vector<art::Ptr<recob::Hit>> hits = fTracksToHits[track];
415 
416  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
417 
418  if (!fHitsToSpacePoints.count(*hitIter)) continue;
419  art::Ptr<recob::SpacePoint> sp = fHitsToSpacePoints.at(*hitIter);
420 
421  TVector3 nearestPoint =
422  trackPoint + trackDir * (trackDir.Dot(TVector3(sp->XYZ()) - trackPoint) / trackDir.Mag2());
423  double lengthAlongTrack = (nearestPointStart - nearestPoint).Mag();
424  sortedTrack.hitMap.insert(std::pair<double, art::Ptr<recob::Hit>>(lengthAlongTrack, *hitIter));
425  }
426 }
int LinFit(const art::Ptr< recob::Track > track, TVector3 &trackPoint, TVector3 &trackDir)
Definition: MVAAlg.cxx:690
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:157
int IsInActiveVol(const TVector3 &pos)
Definition: MVAAlg.cxx:78
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::Hit > > > fTracksToHits
Definition: MVAAlg.h:152
process_name use argoneut_mc_hitfinder track
bool fCheatVertex
Definition: MVAAlg.h:167
TLorentzVector fVertex4Vect
Definition: MVAAlg.h:169
void mvapid::MVAAlg::GetDetectorEdges ( )

Definition at line 91 of file MVAAlg.cxx.

92 {
93 
94  art::ServiceHandle<geo::Geometry const> geom;
95 
96  fDetMinX = 999999.9;
97  fDetMaxX = -999999.9;
98  fDetMinY = 999999.9;
99  fDetMaxY = -999999.9;
100  fDetMinZ = 999999.9;
101  fDetMaxZ = -999999.9;
102 
103  for (unsigned int t = 0; t < geom->TotalNTPC(); t++) {
104  if (geom->TPC(t).MinX() < fDetMinX) fDetMinX = geom->TPC(t).MinX();
105  if (geom->TPC(t).MaxX() > fDetMaxX) fDetMaxX = geom->TPC(t).MaxX();
106  if (geom->TPC(t).MinY() < fDetMinY) fDetMinY = geom->TPC(t).MinY();
107  if (geom->TPC(t).MaxY() > fDetMaxY) fDetMaxY = geom->TPC(t).MaxY();
108  if (geom->TPC(t).MinZ() < fDetMinZ) fDetMinZ = geom->TPC(t).MinZ();
109  if (geom->TPC(t).MaxZ() > fDetMaxZ) fDetMaxZ = geom->TPC(t).MaxZ();
110  }
111 }
double fDetMaxX
Definition: MVAAlg.h:136
double fDetMinY
Definition: MVAAlg.h:136
double fDetMinZ
Definition: MVAAlg.h:136
double fDetMinX
Definition: MVAAlg.h:136
double fDetMaxZ
Definition: MVAAlg.h:136
double fDetMaxY
Definition: MVAAlg.h:136
void mvapid::MVAAlg::GetWireNormals ( )

Definition at line 114 of file MVAAlg.cxx.

115 {
116 
117  art::ServiceHandle<geo::Geometry const> geom;
118 
119  fNormToWiresY.clear();
120  fNormToWiresZ.clear();
121 
122  int planeKey;
123 
124  //Get normals to wires for each plane in the detector
125  //This assumes equal numbers of TPCs in each cryostat and equal numbers of planes in each TPC
126  for (geo::PlaneGeo const& plane : geom->IteratePlanes()) {
127  std::string id = std::string(plane.ID());
128  int pcryo = id.find("C");
129  int ptpc = id.find("T");
130  int pplane = id.find("P");
131  std::string scryo = id.substr(pcryo + 2, 2);
132  std::string stpc = id.substr(ptpc + 2, 2);
133  std::string splane = id.substr(pplane + 2, 2);
134  int icryo = std::stoi(scryo);
135  int itpc = std::stoi(stpc);
136  int iplane = std::stoi(splane);
137  planeKey = icryo * geom->NTPC(0) * geom->Nplanes(0, 0) + itpc * geom->Nplanes(0, 0) +
138  iplane; //single index for all planes in detector
139  fNormToWiresY.insert(
140  std::make_pair(planeKey, -plane.Wire(0).Direction().Z())); //y component of normal
141  fNormToWiresZ.insert(
142  std::make_pair(planeKey, plane.Wire(0).Direction().Y())); //z component of normal
143  }
144 }
std::map< int, double > fNormToWiresY
Definition: MVAAlg.h:138
Geometry information for a single wire plane.The plane is represented in the geometry by a solid whic...
Definition: PlaneGeo.h:82
std::map< int, double > fNormToWiresZ
Definition: MVAAlg.h:139
int mvapid::MVAAlg::IsInActiveVol ( const TVector3 &  pos)
private

Definition at line 78 of file MVAAlg.cxx.

79 {
80  const double fiducialDist = 5.0;
81 
82  if (pos.X() > (fDetMinX + fiducialDist) && pos.X() < (fDetMaxX - fiducialDist) &&
83  pos.Y() > (fDetMinY + fiducialDist) && pos.Y() < (fDetMaxY - fiducialDist) &&
84  pos.Z() > (fDetMinZ + fiducialDist) && pos.Z() < (fDetMaxZ - fiducialDist))
85  return 1;
86  else
87  return 0;
88 }
double fDetMaxX
Definition: MVAAlg.h:136
double fDetMinY
Definition: MVAAlg.h:136
double fDetMinZ
Definition: MVAAlg.h:136
double fDetMinX
Definition: MVAAlg.h:136
double fDetMaxZ
Definition: MVAAlg.h:136
double fDetMaxY
Definition: MVAAlg.h:136
int mvapid::MVAAlg::LinFit ( const art::Ptr< recob::Track track,
TVector3 &  trackPoint,
TVector3 &  trackDir 
)
private

Definition at line 690 of file MVAAlg.cxx.

691 {
692 
693  const std::vector<art::Ptr<recob::SpacePoint>>& sp = fTracksToSpacePoints.at(track);
694 
695  TGraph2D grFit(1);
696  unsigned int iPt = 0;
697  for (auto spIter = sp.begin(); spIter != sp.end(); ++spIter) {
698  TVector3 point = (*spIter)->XYZ();
699  grFit.SetPoint(iPt++, point.X(), point.Y(), point.Z());
700  }
701 
702  //Lift from the ROOT line3Dfit.C tutorial
703  ROOT::Fit::Fitter fitter;
704  // make the functor object
705  mvapid::MVAAlg::SumDistance2 sdist(&grFit);
706 
707  ROOT::Math::Functor fcn(sdist, 6);
708 
709  //Initial fit parameters from track start and end...
710  TVector3 trackStart = track->Vertex<TVector3>();
711  TVector3 trackEnd = track->End<TVector3>();
712  trackDir = (trackEnd - trackStart).Unit();
713 
714  TVector3 x0 = trackStart - trackDir;
715  TVector3 u = trackDir;
716 
717  double pStart[6] = {x0.X(), u.X(), x0.Y(), u.Y(), x0.Z(), u.Z()};
718 
719  fitter.SetFCN(fcn, pStart);
720 
721  bool ok = fitter.FitFCN();
722  if (!ok) {
723  trackPoint.SetXYZ(x0.X(), x0.Y(), x0.Z());
724  trackDir.SetXYZ(u.X(), u.Y(), u.Z());
725  trackDir = trackDir.Unit();
726  return 1;
727  }
728  else {
729  const ROOT::Fit::FitResult& result = fitter.Result();
730  const double* parFit = result.GetParams();
731  trackPoint.SetXYZ(parFit[0], parFit[2], parFit[4]);
732  trackDir.SetXYZ(parFit[1], parFit[3], parFit[5]);
733  trackDir = trackDir.Unit();
734  return 0;
735  }
736 }
process_name use argoneut_mc_hitfinder track
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::SpacePoint > > > fTracksToSpacePoints
Definition: MVAAlg.h:153
int mvapid::MVAAlg::LinFitShower ( const art::Ptr< recob::Shower shower,
TVector3 &  showerPoint,
TVector3 &  showerDir 
)
private

Definition at line 739 of file MVAAlg.cxx.

742 {
743 
744  const std::vector<art::Ptr<recob::SpacePoint>>& sp = fShowersToSpacePoints.at(shower);
745 
746  TGraph2D grFit(1);
747  unsigned int iPt = 0;
748  for (auto spIter = sp.begin(); spIter != sp.end(); ++spIter) {
749  TVector3 point = (*spIter)->XYZ();
750  grFit.SetPoint(iPt++, point.X(), point.Y(), point.Z());
751  }
752 
753  //Lift from the ROOT line3Dfit.C tutorial
754  ROOT::Fit::Fitter fitter;
755  // make the functor object
756  mvapid::MVAAlg::SumDistance2 sdist(&grFit);
757 
758  ROOT::Math::Functor fcn(sdist, 6);
759 
760  //Initial fit parameters from shower start and end...
761  TVector3 showerStart = shower->ShowerStart();
762  showerDir = shower->Direction().Unit();
763 
764  TVector3 x0 = showerStart - showerDir;
765  TVector3 u = showerDir;
766 
767  double pStart[6] = {x0.X(), u.X(), x0.Y(), u.Y(), x0.Z(), u.Z()};
768 
769  fitter.SetFCN(fcn, pStart);
770 
771  bool ok = fitter.FitFCN();
772  if (!ok) {
773  showerPoint.SetXYZ(x0.X(), x0.Y(), x0.Z());
774  showerDir.SetXYZ(u.X(), u.Y(), u.Z());
775  showerDir = showerDir.Unit();
776  return 1;
777  }
778  else {
779  const ROOT::Fit::FitResult& result = fitter.Result();
780  const double* parFit = result.GetParams();
781  showerPoint.SetXYZ(parFit[0], parFit[2], parFit[4]);
782  showerDir.SetXYZ(parFit[1], parFit[3], parFit[5]);
783  showerDir = showerDir.Unit();
784  return 0;
785  }
786 }
process_name shower
Definition: cheaterreco.fcl:51
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::SpacePoint > > > fShowersToSpacePoints
Definition: MVAAlg.h:156
void mvapid::MVAAlg::PrepareEvent ( const art::Event &  event,
const detinfo::DetectorClocksData clockData 
)
private

Definition at line 248 of file MVAAlg.cxx.

249 {
250 
251  fHits.clear();
252  fSpacePoints.clear();
253  fTracks.clear();
254  fShowers.clear();
255  fSpacePointsToHits.clear();
256  fHitsToSpacePoints.clear();
257  fTracksToHits.clear();
258  fTracksToSpacePoints.clear();
259  fShowersToHits.clear();
260  fShowersToSpacePoints.clear();
261 
262  fEventT0 = trigger_offset(clockData);
263 
264  art::Handle<std::vector<recob::Hit>> hitsHandle;
265  evt.getByLabel(fHitLabel, hitsHandle);
266 
267  for (unsigned int iHit = 0; iHit < hitsHandle->size(); ++iHit) {
268  const art::Ptr<recob::Hit> hit(hitsHandle, iHit);
269  fHits.push_back(hit);
270  }
271 
272  art::Handle<std::vector<recob::Track>> tracksHandle;
273  evt.getByLabel(fTrackLabel, tracksHandle);
274 
275  for (unsigned int iTrack = 0; iTrack < tracksHandle->size(); ++iTrack) {
276  const art::Ptr<recob::Track> track(tracksHandle, iTrack);
277  fTracks.push_back(track);
278  }
279 
280  art::Handle<std::vector<recob::Shower>> showersHandle;
281  evt.getByLabel(fShowerLabel, showersHandle);
282 
283  for (unsigned int iShower = 0; iShower < showersHandle->size(); ++iShower) {
284  const art::Ptr<recob::Shower> shower(showersHandle, iShower);
285  fShowers.push_back(shower);
286  }
287 
288  art::Handle<std::vector<recob::SpacePoint>> spHandle;
289  evt.getByLabel(fSpacePointLabel, spHandle);
290 
291  for (unsigned int iSP = 0; iSP < spHandle->size(); ++iSP) {
292  const art::Ptr<recob::SpacePoint> spacePoint(spHandle, iSP);
293  fSpacePoints.push_back(spacePoint);
294  }
295 
296  art::FindManyP<recob::Hit> findTracksToHits(fTracks, evt, fTrackLabel);
297  art::FindManyP<recob::Hit> findShowersToHits(fShowers, evt, fShowerLabel);
298  art::FindOneP<recob::Hit> findSPToHits(fSpacePoints, evt, fSpacePointLabel);
299 
300  for (unsigned int iSP = 0; iSP < fSpacePoints.size(); ++iSP) {
301  const art::Ptr<recob::SpacePoint> spacePoint = fSpacePoints.at(iSP);
302 
303  const art::Ptr<recob::Hit> hit = findSPToHits.at(iSP);
304  fSpacePointsToHits[spacePoint] = hit;
305  fHitsToSpacePoints[hit] = spacePoint;
306  }
307 
308  for (unsigned int iTrack = 0; iTrack < fTracks.size(); ++iTrack) {
309  const art::Ptr<recob::Track> track = fTracks.at(iTrack);
310 
311  const std::vector<art::Ptr<recob::Hit>> trackHits = findTracksToHits.at(iTrack);
312 
313  for (unsigned int iHit = 0; iHit < trackHits.size(); ++iHit) {
314  const art::Ptr<recob::Hit> hit = trackHits.at(iHit);
315  fTracksToHits[track].push_back(hit);
316  if (fHitsToSpacePoints.count(hit)) {
317  fTracksToSpacePoints[track].push_back(fHitsToSpacePoints.at(hit));
318  }
319  }
320  }
321 
322  for (unsigned int iShower = 0; iShower < fShowers.size(); ++iShower) {
323  const art::Ptr<recob::Shower> shower = fShowers.at(iShower);
324  const std::vector<art::Ptr<recob::Hit>> showerHits = findShowersToHits.at(iShower);
325 
326  for (unsigned int iHit = 0; iHit < showerHits.size(); ++iHit) {
327  const art::Ptr<recob::Hit> hit = showerHits.at(iHit);
328  fShowersToHits[shower].push_back(hit);
329  if (fHitsToSpacePoints.count(hit)) {
330  fShowersToSpacePoints[shower].push_back(fHitsToSpacePoints.at(hit));
331  }
332  }
333  }
334 
335  if (fCheatVertex) {
336  art::Handle<std::vector<simb::MCParticle>> partHandle;
337  evt.getByLabel(fTrackingLabel, partHandle);
338 
339  if (partHandle->size() == 0 || partHandle->at(0).TrackId() != 1) {
340  std::cout << "Error, ID of first track in largeant list is not 0" << std::endl;
341  exit(1);
342  }
343  fVertex4Vect = partHandle->at(0).Position();
344  }
345 }
std::vector< art::Ptr< recob::Shower > > fShowers
Definition: MVAAlg.h:148
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:157
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::Hit > > > fTracksToHits
Definition: MVAAlg.h:152
process_name use argoneut_mc_hitfinder track
std::string fSpacePointLabel
Definition: MVAAlg.h:144
process_name hit
Definition: cheaterreco.fcl:51
bool fCheatVertex
Definition: MVAAlg.h:167
process_name shower
Definition: cheaterreco.fcl:51
std::vector< art::Ptr< recob::Hit > > fHits
Definition: MVAAlg.h:150
std::vector< art::Ptr< recob::SpacePoint > > fSpacePoints
Definition: MVAAlg.h:149
TLorentzVector fVertex4Vect
Definition: MVAAlg.h:169
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::Hit > > > fShowersToHits
Definition: MVAAlg.h:154
std::string fTrackingLabel
Definition: MVAAlg.h:145
std::string fTrackLabel
Definition: MVAAlg.h:141
std::string fHitLabel
Definition: MVAAlg.h:143
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::SpacePoint > > > fTracksToSpacePoints
Definition: MVAAlg.h:153
std::vector< art::Ptr< recob::Track > > fTracks
Definition: MVAAlg.h:147
int trigger_offset(DetectorClocksData const &data)
std::string fShowerLabel
Definition: MVAAlg.h:142
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::SpacePoint > > > fShowersToSpacePoints
Definition: MVAAlg.h:156
TCEvent evt
Definition: DataStructs.cxx:8
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > fSpacePointsToHits
Definition: MVAAlg.h:158
double fEventT0
Definition: MVAAlg.h:134
BEGIN_PROLOG could also be cout
void mvapid::MVAAlg::RunPCA ( std::vector< art::Ptr< recob::Hit >> &  hits,
std::vector< double > &  eVals,
std::vector< double > &  eVecs 
)
private

Definition at line 522 of file MVAAlg.cxx.

525 {
526  TPrincipal* principal = new TPrincipal(3, "D");
527 
528  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
529 
530  if (fHitsToSpacePoints.count(*hitIter)) {
531  principal->AddRow(fHitsToSpacePoints.at(*hitIter)->XYZ());
532  }
533  }
534 
535  // PERFORM PCA
536  principal->MakePrincipals();
537  // GET EIGENVALUES AND EIGENVECTORS
538  for (unsigned int i = 0; i < 3; ++i) {
539  eVals.push_back(principal->GetEigenValues()->GetMatrixArray()[i]);
540  }
541 
542  for (unsigned int i = 0; i < 9; ++i) {
543  eVecs.push_back(principal->GetEigenVectors()->GetMatrixArray()[i]);
544  }
545 }
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:157
void mvapid::MVAAlg::RunPID ( art::Event &  evt,
std::vector< anab::MVAPIDResult > &  result,
art::Assns< recob::Track, anab::MVAPIDResult, void > &  trackAssns,
art::Assns< recob::Shower, anab::MVAPIDResult, void > &  showerAssns 
)

Definition at line 147 of file MVAAlg.cxx.

151 {
152  auto const clockData =
153  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
154  auto const detProp =
155  art::ServiceHandle<detinfo::DetectorPropertiesService const>()->DataFor(evt, clockData);
156  this->PrepareEvent(evt, clockData);
157 
158  for (auto trackIter = fTracks.begin(); trackIter != fTracks.end(); ++trackIter) {
159  mvapid::MVAAlg::SortedObj sortedObj;
160 
161  std::vector<double> eVals, eVecs;
162  int isStoppingReco;
163  this->RunPCA(fTracksToHits[*trackIter], eVals, eVecs);
164  double evalRatio;
165  if (eVals[0] < 0.0001)
166  evalRatio = 0.0;
167  else
168  evalRatio = std::sqrt(eVals[1] * eVals[1] + eVals[2] * eVals[2]) / eVals[0];
169  this->FitAndSortTrack(*trackIter, isStoppingReco, sortedObj);
170  double coreHaloRatio, concentration, conicalness;
171  this->_Var_Shape(sortedObj, coreHaloRatio, concentration, conicalness);
172  double dEdxStart = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0., 0.05);
173  double dEdxEnd = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.9, 1.0);
174  double dEdxPenultimate = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.8, 0.9);
175 
176  fResHolder.isTrack = 1;
177  fResHolder.isStoppingReco = isStoppingReco;
178  fResHolder.nSpacePoints = sortedObj.hitMap.size();
179  fResHolder.trackID = (*trackIter)->ID();
180  fResHolder.evalRatio = evalRatio;
181  fResHolder.concentration = concentration;
182  fResHolder.coreHaloRatio = coreHaloRatio;
183  fResHolder.conicalness = conicalness;
184  fResHolder.dEdxStart = dEdxStart;
185  fResHolder.dEdxEnd = dEdxEnd;
186  if (dEdxPenultimate < 0.1)
187  fResHolder.dEdxEndRatio = 1.0;
188  else
189  fResHolder.dEdxEndRatio = dEdxEnd / dEdxPenultimate;
190  fResHolder.length = sortedObj.length;
191 
192  for (auto methodIter = fMVAMethods.begin(); methodIter != fMVAMethods.end(); ++methodIter) {
193  fResHolder.mvaOutput[*methodIter] = fReader.EvaluateMVA(*methodIter);
194  }
195  result.push_back(fResHolder);
196  util::CreateAssn(evt, result, *trackIter, trackAssns);
197  }
198 
199  for (auto showerIter = fShowers.begin(); showerIter != fShowers.end(); ++showerIter) {
200  mvapid::MVAAlg::SortedObj sortedObj;
201 
202  std::vector<double> eVals, eVecs;
203  int isStoppingReco;
204 
205  this->RunPCA(fShowersToHits[*showerIter], eVals, eVecs);
206 
207  double evalRatio;
208  if (eVals[0] < 0.0001)
209  evalRatio = 0.0;
210  else
211  evalRatio = std::sqrt(eVals[1] * eVals[1] + eVals[2] * eVals[2]) / eVals[0];
212 
213  this->SortShower(*showerIter, isStoppingReco, sortedObj);
214 
215  double coreHaloRatio, concentration, conicalness;
216  this->_Var_Shape(sortedObj, coreHaloRatio, concentration, conicalness);
217  double dEdxStart = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0., 0.05);
218  double dEdxEnd = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.9, 1.0);
219  double dEdxPenultimate = CalcSegmentdEdxFrac(clockData, detProp, sortedObj, 0.8, 0.9);
220 
221  fResHolder.isTrack = 0;
222  fResHolder.isStoppingReco = isStoppingReco;
223  fResHolder.nSpacePoints = sortedObj.hitMap.size();
225  (*showerIter)->ID() + 1000; //For the moment label showers by adding 1000 to ID
226 
227  fResHolder.evalRatio = evalRatio;
228  fResHolder.concentration = concentration;
229  fResHolder.coreHaloRatio = coreHaloRatio;
230  fResHolder.conicalness = conicalness;
231  fResHolder.dEdxStart = dEdxStart;
232  fResHolder.dEdxEnd = dEdxEnd;
233  if (dEdxPenultimate < 0.1)
234  fResHolder.dEdxEndRatio = 1.0;
235  else
236  fResHolder.dEdxEndRatio = dEdxEnd / dEdxPenultimate;
237  fResHolder.length = sortedObj.length;
238 
239  for (auto methodIter = fMVAMethods.begin(); methodIter != fMVAMethods.end(); ++methodIter) {
240  fResHolder.mvaOutput[*methodIter] = fReader.EvaluateMVA(*methodIter);
241  }
242  result.push_back(fResHolder);
243  util::CreateAssn(evt, result, *showerIter, showerAssns);
244  }
245 }
std::vector< art::Ptr< recob::Shower > > fShowers
Definition: MVAAlg.h:148
anab::MVAPIDResult fResHolder
Definition: MVAAlg.h:160
void RunPCA(std::vector< art::Ptr< recob::Hit >> &hits, std::vector< double > &eVals, std::vector< double > &eVecs)
Definition: MVAAlg.cxx:522
std::vector< std::string > fMVAMethods
Definition: MVAAlg.h:164
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::Hit > > > fTracksToHits
Definition: MVAAlg.h:152
void PrepareEvent(const art::Event &event, const detinfo::DetectorClocksData &clockData)
Definition: MVAAlg.cxx:248
TMVA::Reader fReader
Definition: MVAAlg.h:162
void SortShower(art::Ptr< recob::Shower > shower, int &isStoppingReco, mvapid::MVAAlg::SortedObj &sortedShower)
Definition: MVAAlg.cxx:431
std::map< double, const art::Ptr< recob::Hit > > hitMap
Definition: MVAAlg.h:45
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::Hit > > > fShowersToHits
Definition: MVAAlg.h:154
double CalcSegmentdEdxFrac(const detinfo::DetectorClocksData &clock_data, const detinfo::DetectorPropertiesData &det_prop, const SortedObj &track, double start, double end)
Definition: MVAAlg.cxx:616
bool CreateAssn(art::Event &evt, std::vector< T > const &a, art::Ptr< U > const &b, art::Assns< U, T > &assn, std::string a_instance, size_t index=UINT_MAX)
Creates a single one-to-one association.
void _Var_Shape(const SortedObj &track, double &coreHaloRatio, double &concentration, double &conicalness)
Definition: MVAAlg.cxx:547
std::vector< art::Ptr< recob::Track > > fTracks
Definition: MVAAlg.h:147
TCEvent evt
Definition: DataStructs.cxx:8
std::map< std::string, double > mvaOutput
Definition: MVAPIDResult.h:27
void FitAndSortTrack(art::Ptr< recob::Track > track, int &isStoppingReco, SortedObj &sortedObj)
Definition: MVAAlg.cxx:348
unsigned int trackID
Definition: MVAPIDResult.h:21
auto const detProp
void mvapid::MVAAlg::SortShower ( art::Ptr< recob::Shower shower,
int &  isStoppingReco,
mvapid::MVAAlg::SortedObj sortedShower 
)
private

Definition at line 431 of file MVAAlg.cxx.

434 {
435  sortedShower.hitMap.clear();
436 
437  std::vector<art::Ptr<recob::Hit>> hits = fShowersToHits[shower];
438 
439  TVector3 showerEnd(0, 0, 0);
440  double furthestHitFromStart = -999.9;
441  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
442 
443  if (!fHitsToSpacePoints.count(*hitIter)) continue;
444  art::Ptr<recob::SpacePoint> sp = fHitsToSpacePoints.at(*hitIter);
445  if ((TVector3(sp->XYZ()) - shower->ShowerStart()).Mag() > furthestHitFromStart) {
446  showerEnd = TVector3(sp->XYZ());
447  furthestHitFromStart = (TVector3(sp->XYZ()) - shower->ShowerStart()).Mag();
448  }
449  }
450 
451  TVector3 showerPoint, showerDir;
452  this->LinFitShower(shower, showerPoint, showerDir);
453 
454  TVector3 nearestPointStart, nearestPointEnd;
455 
456  //Ensure that shower is fitted in correct direction (assuming for now that particle moves in +z direction)
457 
458  if (fCheatVertex) {
459  if ((showerEnd - fVertex4Vect.Vect()).Mag() >
460  (shower->ShowerStart() - fVertex4Vect.Vect()).Mag()) {
461  nearestPointStart =
462  showerPoint +
463  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
464  nearestPointEnd =
465  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
466  isStoppingReco = this->IsInActiveVol(showerEnd);
467  }
468  else {
469  nearestPointStart =
470  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
471  nearestPointEnd =
472  showerPoint +
473  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
474  isStoppingReco = this->IsInActiveVol(shower->ShowerStart());
475  showerDir *= -1.;
476  }
477  }
478  else {
479  if (showerEnd.Z() >= shower->ShowerStart().Z()) {
480  nearestPointStart =
481  showerPoint +
482  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
483  nearestPointEnd =
484  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
485  isStoppingReco = this->IsInActiveVol(showerEnd);
486  }
487  else {
488  nearestPointStart =
489  showerPoint + showerDir * (showerDir.Dot(showerEnd - showerPoint) / showerDir.Mag2());
490  nearestPointEnd =
491  showerPoint +
492  showerDir * (showerDir.Dot(shower->ShowerStart() - showerPoint) / showerDir.Mag2());
493  isStoppingReco = this->IsInActiveVol(shower->ShowerStart());
494  }
495 
496  if (showerDir.Z() <= 0) {
497  showerDir.SetX(-showerDir.X());
498  showerDir.SetY(-showerDir.Y());
499  showerDir.SetZ(-showerDir.Z());
500  }
501  }
502 
503  sortedShower.start = nearestPointStart;
504  sortedShower.end = nearestPointEnd;
505  sortedShower.dir = showerDir;
506  sortedShower.length = (nearestPointEnd - nearestPointStart).Mag();
507 
508  for (auto hitIter = hits.begin(); hitIter != hits.end(); ++hitIter) {
509 
510  if (!fHitsToSpacePoints.count(*hitIter)) continue;
511  art::Ptr<recob::SpacePoint> sp = fHitsToSpacePoints.at(*hitIter);
512 
513  TVector3 nearestPoint =
514  showerPoint +
515  showerDir * (showerDir.Dot(TVector3(sp->XYZ()) - showerPoint) / showerDir.Mag2());
516  double lengthAlongShower = (nearestPointStart - nearestPoint).Mag();
517  sortedShower.hitMap.insert(
518  std::pair<double, art::Ptr<recob::Hit>>(lengthAlongShower, *hitIter));
519  }
520 }
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:157
int IsInActiveVol(const TVector3 &pos)
Definition: MVAAlg.cxx:78
bool fCheatVertex
Definition: MVAAlg.h:167
process_name shower
Definition: cheaterreco.fcl:51
int LinFitShower(const art::Ptr< recob::Shower > shower, TVector3 &showerPoint, TVector3 &showerDir)
Definition: MVAAlg.cxx:739
std::map< double, const art::Ptr< recob::Hit > > hitMap
Definition: MVAAlg.h:45
TLorentzVector fVertex4Vect
Definition: MVAAlg.h:169
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::Hit > > > fShowersToHits
Definition: MVAAlg.h:154

Member Data Documentation

const calo::CalorimetryAlg mvapid::MVAAlg::fCaloAlg
private

Definition at line 132 of file MVAAlg.h.

bool mvapid::MVAAlg::fCheatVertex
private

Definition at line 167 of file MVAAlg.h.

double mvapid::MVAAlg::fDetMaxX
private

Definition at line 136 of file MVAAlg.h.

double mvapid::MVAAlg::fDetMaxY
private

Definition at line 136 of file MVAAlg.h.

double mvapid::MVAAlg::fDetMaxZ
private

Definition at line 136 of file MVAAlg.h.

double mvapid::MVAAlg::fDetMinX
private

Definition at line 136 of file MVAAlg.h.

double mvapid::MVAAlg::fDetMinY
private

Definition at line 136 of file MVAAlg.h.

double mvapid::MVAAlg::fDetMinZ
private

Definition at line 136 of file MVAAlg.h.

double mvapid::MVAAlg::fEventT0
private

Definition at line 134 of file MVAAlg.h.

std::string mvapid::MVAAlg::fHitLabel
private

Definition at line 143 of file MVAAlg.h.

std::vector<art::Ptr<recob::Hit> > mvapid::MVAAlg::fHits
private

Definition at line 150 of file MVAAlg.h.

std::map<art::Ptr<recob::Hit>, art::Ptr<recob::SpacePoint> > mvapid::MVAAlg::fHitsToSpacePoints
private

Definition at line 157 of file MVAAlg.h.

std::vector<std::string> mvapid::MVAAlg::fMVAMethods
private

Definition at line 164 of file MVAAlg.h.

std::map<int, double> mvapid::MVAAlg::fNormToWiresY
private

Definition at line 138 of file MVAAlg.h.

std::map<int, double> mvapid::MVAAlg::fNormToWiresZ
private

Definition at line 139 of file MVAAlg.h.

TMVA::Reader mvapid::MVAAlg::fReader
private

Definition at line 162 of file MVAAlg.h.

anab::MVAPIDResult mvapid::MVAAlg::fResHolder
private

Definition at line 160 of file MVAAlg.h.

std::string mvapid::MVAAlg::fShowerLabel
private

Definition at line 142 of file MVAAlg.h.

std::vector<art::Ptr<recob::Shower> > mvapid::MVAAlg::fShowers
private

Definition at line 148 of file MVAAlg.h.

std::map<art::Ptr<recob::Shower>, std::vector<art::Ptr<recob::Hit> > > mvapid::MVAAlg::fShowersToHits
private

Definition at line 154 of file MVAAlg.h.

std::map<art::Ptr<recob::Shower>, std::vector<art::Ptr<recob::SpacePoint> > > mvapid::MVAAlg::fShowersToSpacePoints
private

Definition at line 156 of file MVAAlg.h.

std::string mvapid::MVAAlg::fSpacePointLabel
private

Definition at line 144 of file MVAAlg.h.

std::vector<art::Ptr<recob::SpacePoint> > mvapid::MVAAlg::fSpacePoints
private

Definition at line 149 of file MVAAlg.h.

std::map<art::Ptr<recob::SpacePoint>, art::Ptr<recob::Hit> > mvapid::MVAAlg::fSpacePointsToHits
private

Definition at line 158 of file MVAAlg.h.

std::string mvapid::MVAAlg::fTrackingLabel
private

Definition at line 145 of file MVAAlg.h.

std::string mvapid::MVAAlg::fTrackLabel
private

Definition at line 141 of file MVAAlg.h.

std::vector<art::Ptr<recob::Track> > mvapid::MVAAlg::fTracks
private

Definition at line 147 of file MVAAlg.h.

std::map<art::Ptr<recob::Track>, std::vector<art::Ptr<recob::Hit> > > mvapid::MVAAlg::fTracksToHits
private

Definition at line 152 of file MVAAlg.h.

std::map<art::Ptr<recob::Track>, std::vector<art::Ptr<recob::SpacePoint> > > mvapid::MVAAlg::fTracksToSpacePoints
private

Definition at line 153 of file MVAAlg.h.

TLorentzVector mvapid::MVAAlg::fVertex4Vect
private

Definition at line 169 of file MVAAlg.h.

std::vector<std::string> mvapid::MVAAlg::fWeightFiles
private

Definition at line 165 of file MVAAlg.h.


The documentation for this class was generated from the following files: