LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
MVAAlg.h
Go to the documentation of this file.
1 // \fileMVAAlg.h
3 // m.haigh@warwick.ac.uk
5 #ifndef MVAAlg_H
6 #define MVAAlg_H
7 
8 #include <map>
9 #include <string>
10 #include <vector>
11 
12 namespace art {
13  class Event;
14 }
17 namespace fhicl {
18  class ParameterSet;
19 }
20 
22 
24 
25 namespace detinfo {
26  class DetectorClocksData;
27  class DetectorPropertiesData;
28 }
29 namespace recob {
30  class Hit;
31  class Shower;
32  class SpacePoint;
33  class Track;
34 }
35 #include "Math/Vector3Dfwd.h"
36 #include "TGraph2D.h"
37 #include "TLorentzVector.h"
38 #include "TMVA/Reader.h"
39 #include "TVector3.h"
40 
41 namespace mvapid {
42 
43  //---------------------------------------------------------------
44  class MVAAlg {
45  public:
46  struct SortedObj {
47  TVector3 start, end, dir;
48  double length;
49  std::map<double, const art::Ptr<recob::Hit>> hitMap;
50  };
51 
52  struct SumDistance2 {
53  // the TGraph is a data member of the object
54  TGraph2D* fGraph;
55 
56  SumDistance2(TGraph2D* g) : fGraph(g) {}
57 
58  // implementation of the function to be minimized
59  double operator()(const double* p)
60  {
61 
62  ROOT::Math::XYZVector x0(p[0], p[2], p[4]);
63  ROOT::Math::XYZVector u(p[1], p[3], p[5]);
64 
65  u = u.Unit();
66  double* x = fGraph->GetX();
67  double* y = fGraph->GetY();
68  double* z = fGraph->GetZ();
69  int npoints = fGraph->GetN();
70  double sum = 0;
71  for (int i = 0; i < npoints; ++i) {
72  ROOT::Math::XYZVector xp(x[i], y[i], z[i]);
73  sum += ((xp - x0).Cross(u)).Mag2();
74  }
75  return sum;
76  }
77  };
78 
79  MVAAlg(fhicl::ParameterSet const& pset);
80 
81  void GetDetectorEdges();
82 
83  void GetWireNormals();
84 
85  void RunPID(art::Event& evt,
86  std::vector<anab::MVAPIDResult>& result,
89 
90  private:
91  int IsInActiveVol(const TVector3& pos);
92 
93  void PrepareEvent(const art::Event& event, const detinfo::DetectorClocksData& clockData);
94 
95  void FitAndSortTrack(art::Ptr<recob::Track> track, int& isStoppingReco, SortedObj& sortedObj);
96 
97  //void SortShower(art::Ptr<recob::Shower> shower,TVector3 dir,int& isStoppingReco,
98  // mvapid::MVAAlg::SortedObj& sortedShower);
99  void SortShower(art::Ptr<recob::Shower> shower,
100  int& isStoppingReco,
101  mvapid::MVAAlg::SortedObj& sortedShower);
102 
103  void RunPCA(std::vector<art::Ptr<recob::Hit>>& hits,
104  std::vector<double>& eVals,
105  std::vector<double>& eVecs);
106 
107  void _Var_Shape(const SortedObj& track,
108  double& coreHaloRatio,
109  double& concentration,
110  double& conicalness);
111 
112  double CalcSegmentdEdxFrac(const detinfo::DetectorClocksData& clock_data,
113  const detinfo::DetectorPropertiesData& det_prop,
114  const SortedObj& track,
115  double start,
116  double end);
117 
118  double CalcSegmentdEdxDist(const detinfo::DetectorClocksData& clock_data,
119  const detinfo::DetectorPropertiesData& det_prop,
120  const SortedObj& track,
121  double start,
122  double end);
123 
124  double CalcSegmentdEdxDistAtEnd(const detinfo::DetectorClocksData& clock_data,
125  const detinfo::DetectorPropertiesData& det_prop,
126  const mvapid::MVAAlg::SortedObj& track,
127  double distAtEnd);
128 
129  int LinFit(const art::Ptr<recob::Track> track, TVector3& trackPoint, TVector3& trackDir);
130 
131  int LinFitShower(const art::Ptr<recob::Shower> shower,
132  TVector3& showerPoint,
133  TVector3& showerDir);
134 
136 
137  double fEventT0;
138 
139  double fDetMinX, fDetMaxX, fDetMinY, fDetMaxY, fDetMinZ, fDetMaxZ;
140 
141  std::map<int, double> fNormToWiresY;
142  std::map<int, double> fNormToWiresZ;
143 
144  std::string fTrackLabel;
145  std::string fShowerLabel;
146  std::string fHitLabel;
147  std::string fSpacePointLabel;
148  std::string fTrackingLabel;
149 
150  std::vector<art::Ptr<recob::Track>> fTracks;
151  std::vector<art::Ptr<recob::Shower>> fShowers;
152  std::vector<art::Ptr<recob::SpacePoint>> fSpacePoints;
153  std::vector<art::Ptr<recob::Hit>> fHits;
154 
155  std::map<art::Ptr<recob::Track>, std::vector<art::Ptr<recob::Hit>>> fTracksToHits;
156  std::map<art::Ptr<recob::Track>, std::vector<art::Ptr<recob::SpacePoint>>> fTracksToSpacePoints;
157  std::map<art::Ptr<recob::Shower>, std::vector<art::Ptr<recob::Hit>>> fShowersToHits;
158  std::map<art::Ptr<recob::Shower>, std::vector<art::Ptr<recob::SpacePoint>>>
160  std::map<art::Ptr<recob::Hit>, art::Ptr<recob::SpacePoint>> fHitsToSpacePoints;
161  std::map<art::Ptr<recob::SpacePoint>, art::Ptr<recob::Hit>> fSpacePointsToHits;
162 
164 
165  TMVA::Reader fReader;
166 
167  std::vector<std::string> fMVAMethods;
168  std::vector<std::string> fWeightFiles;
169 
171 
172  TLorentzVector fVertex4Vect;
173 
174  }; // class MVAAlg
175 
176 } // namespace mvapid
177 
178 #endif // ifndef MVAAlg_H
Float_t x
Definition: compare.C:6
std::vector< art::Ptr< recob::Shower > > fShowers
Definition: MVAAlg.h:151
anab::MVAPIDResult fResHolder
Definition: MVAAlg.h:163
const calo::CalorimetryAlg fCaloAlg
Definition: MVAAlg.h:135
Reconstruction base classes.
SumDistance2(TGraph2D *g)
Definition: MVAAlg.h:56
double operator()(const double *p)
Definition: MVAAlg.h:59
std::map< art::Ptr< recob::Hit >, art::Ptr< recob::SpacePoint > > fHitsToSpacePoints
Definition: MVAAlg.h:160
std::vector< std::string > fMVAMethods
Definition: MVAAlg.h:167
Float_t y
Definition: compare.C:6
Double_t z
Definition: plot.C:276
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::Hit > > > fTracksToHits
Definition: MVAAlg.h:155
std::string fSpacePointLabel
Definition: MVAAlg.h:147
TMVA::Reader fReader
Definition: MVAAlg.h:165
bool fCheatVertex
Definition: MVAAlg.h:170
std::vector< art::Ptr< recob::Hit > > fHits
Definition: MVAAlg.h:153
decltype(auto) constexpr end(T &&obj)
ADL-aware version of std::end.
Definition: StdUtils.h:77
Definition: MVAAlg.h:41
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:289
void hits()
Definition: readHits.C:15
parameter set interface
std::map< double, const art::Ptr< recob::Hit > > hitMap
Definition: MVAAlg.h:49
std::map< int, double > fNormToWiresY
Definition: MVAAlg.h:141
double fDetMinZ
Definition: MVAAlg.h:139
std::vector< art::Ptr< recob::SpacePoint > > fSpacePoints
Definition: MVAAlg.h:152
TLorentzVector fVertex4Vect
Definition: MVAAlg.h:172
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::Hit > > > fShowersToHits
Definition: MVAAlg.h:157
General LArSoft Utilities.
std::string fTrackingLabel
Definition: MVAAlg.h:148
std::string fTrackLabel
Definition: MVAAlg.h:144
std::string fHitLabel
Definition: MVAAlg.h:146
Contains all timing reference information for the detector.
TDirectory * dir
Definition: macro.C:5
std::map< art::Ptr< recob::Track >, std::vector< art::Ptr< recob::SpacePoint > > > fTracksToSpacePoints
Definition: MVAAlg.h:156
std::vector< art::Ptr< recob::Track > > fTracks
Definition: MVAAlg.h:150
Definition: MVAAlg.h:12
std::string fShowerLabel
Definition: MVAAlg.h:145
std::map< art::Ptr< recob::Shower >, std::vector< art::Ptr< recob::SpacePoint > > > fShowersToSpacePoints
Definition: MVAAlg.h:159
TCEvent evt
Definition: DataStructs.cxx:8
std::map< art::Ptr< recob::SpacePoint >, art::Ptr< recob::Hit > > fSpacePointsToHits
Definition: MVAAlg.h:161
Double_t sum
Definition: plot.C:31
Float_t track
Definition: plot.C:35
std::map< int, double > fNormToWiresZ
Definition: MVAAlg.h:142
Definition: Track.hh:42
double fEventT0
Definition: MVAAlg.h:137
std::vector< std::string > fWeightFiles
Definition: MVAAlg.h:168
Event finding and building.