All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
larreco/larreco/RecoAlg/TrackMomentumCalculator.h
Go to the documentation of this file.
1 /// \file TrackMomentumCalculator.h
2 // \author sowjanyag@phys.ksu.edu
3 
4 #ifndef TrackMomentumCalculator_H
5 #define TrackMomentumCalculator_H
6 
7 #include "canvas/Persistency/Common/Ptr.h"
9 
10 #include "TGraph.h"
11 #include "TVector3.h"
12 
13 #include <optional>
14 #include <vector>
15 #include <tuple>
16 
17 class TPolyLine3D;
18 
19 namespace trkf {
20 
22  public:
23  TrackMomentumCalculator(double minLength = 100.0,
24  double maxLength = 1350.0);
25 
26  double GetTrackMomentum(double trkrange, int pdg) const;
27  double GetMomentumMultiScatterChi2(art::Ptr<recob::Track> const& trk, const bool checkValidPoints = false);
28  double GetMomentumMultiScatterLLHD(art::Ptr<recob::Track> const& trk);
29  double GetMuMultiScatterLLHD3(art::Ptr<recob::Track> const& trk, bool dir);
30  TVector3 GetMultiScatterStartingPoint(art::Ptr<recob::Track> const& trk);
31 
32  private:
33  bool plotRecoTracks_(std::vector<float> const& xxx,
34  std::vector<float> const& yyy,
35  std::vector<float> const& zzz);
36 
37  struct Segments {
38  std::vector<float> x, nx;
39  std::vector<float> y, ny;
40  std::vector<float> z, nz;
41  std::vector<float> L;
42  };
43 
44  std::optional<Segments> getSegTracks_(std::vector<float> const& xxx,
45  std::vector<float> const& yyy,
46  std::vector<float> const& zzz,
47  double seg_size);
48 
49  std::tuple<double, double, double> getDeltaThetaRMS_(Segments const& segments,
50  double thick) const;
51 
52  int getDeltaThetaij_(std::vector<float>& ei,
53  std::vector<float>& ej,
54  std::vector<float>& th,
55  std::vector<float>& ind,
56  Segments const& segments,
57  double thick) const;
58 
59  double my_g(double xx, double Q, double s) const;
60 
61  double my_mcs_llhd(std::vector<float> const& dEi,
62  std::vector<float> const& dEj,
63  std::vector<float> const& dthij,
64  std::vector<float> const& ind,
65  double x0, double x1) const;
66 
67  float seg_stop{-1.};
68  int n_seg{};
69 
70  float x_seg[50000];
71  float y_seg[50000];
72  float z_seg[50000];
73 
74  double find_angle(double vz, double vy) const;
75 
76  float steps_size{10.};
77  int n_steps{6};
78  std::vector<float> steps;
79 
80  double minLength;
81  double maxLength;
82 
83  // The following are objects that are created but not drawn or
84  // saved. This class should consider accepting a "debug"
85  // parameter where if it is specified, then the graphs will be
86  // created; otherwise, their creation is unnecessary and impedes
87  // efficiency.
88  //
89  // N.B. TPolyLine3D objects are owned by ROOT, and we thus refer
90  // to them by pointer. It is important that 'delete' is not
91  // called on the TPolyLine3D pointers during destruction of a
92  // TrackMomentumCalculator object.
93  TPolyLine3D* gr_reco_xyz{nullptr};
94  TGraph gr_reco_xy{};
95  TGraph gr_reco_yz{};
96  TGraph gr_reco_xz{};
97 
98  TPolyLine3D* gr_seg_xyz{nullptr};
99  TGraph gr_seg_xy{};
100  TGraph gr_seg_yz{};
101  TGraph gr_seg_xz{};
102 
103  };
104 
105 } // namespace trkf
106 
107 #endif // TrackMomentumCalculator_H
var pdg
Definition: selectors.fcl:14
double GetMomentumMultiScatterChi2(art::Ptr< recob::Track > const &trk, const bool checkValidPoints=false)
double my_g(double xx, double Q, double s) const
int getDeltaThetaij_(std::vector< float > &ei, std::vector< float > &ej, std::vector< float > &th, std::vector< float > &ind, Segments const &segments, double thick) const
std::optional< Segments > getSegTracks_(std::vector< float > const &xxx, std::vector< float > const &yyy, std::vector< float > const &zzz, double seg_size)
std::tuple< double, double, double > getDeltaThetaRMS_(Segments const &segments, double thick) const
Provides recob::Track data product.
tuple dir
Definition: dropbox.py:28
double GetMuMultiScatterLLHD3(art::Ptr< recob::Track > const &trk, bool dir)
TrackMomentumCalculator(double minLength=100.0, double maxLength=1350.0)
double my_mcs_llhd(std::vector< float > const &dEi, std::vector< float > const &dEj, std::vector< float > const &dthij, std::vector< float > const &ind, double x0, double x1) const
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
TVector3 GetMultiScatterStartingPoint(art::Ptr< recob::Track > const &trk)
double GetMomentumMultiScatterLLHD(art::Ptr< recob::Track > const &trk)
double GetTrackMomentum(double trkrange, int pdg) const
bool plotRecoTracks_(std::vector< float > const &xxx, std::vector< float > const &yyy, std::vector< float > const &zzz)