LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
TrackMomentumCalculator.h
Go to the documentation of this file.
1 // \author sowjanyag@phys.ksu.edu
3 
4 #ifndef TrackMomentumCalculator_H
5 #define TrackMomentumCalculator_H
6 
10 
11 #include "TAxis.h"
12 #include "TGraph.h"
13 #include "TGraphErrors.h"
14 #include "TMath.h"
15 #include "TMatrixDSym.h"
16 #include "TMatrixDSymEigen.h"
17 #include "TPolyLine3D.h"
18 #include "TSpline.h"
19 #include "TVector3.h"
20 
21 #include "Math/Factory.h"
22 #include "Math/Minimizer.h"
23 #include "Minuit2/FCNBase.h"
24 #include "Minuit2/FunctionMinimum.h"
25 #include "Minuit2/Minuit2Minimizer.h"
26 #include "Minuit2/MnMigrad.h"
27 #include "Minuit2/MnPrint.h"
28 #include "Minuit2/MnUserParameterState.h"
29 #include "Minuit2/MnUserParameters.h"
30 
31 #include <cmath>
32 #include <iostream>
33 #include <optional>
34 #include <vector>
35 #include <tuple>
36 
37 namespace trkf {
38 
40  public:
41  TrackMomentumCalculator(double minLength = 100.0,
42  double maxLength = 1350.0);
43 
44  double GetTrackMomentum(double trkrange, int pdg) const;
47  double GetMuMultiScatterLLHD3(art::Ptr<recob::Track> const& trk, bool dir);
49 
50  private:
51  bool plotRecoTracks_(std::vector<float> const& xxx,
52  std::vector<float> const& yyy,
53  std::vector<float> const& zzz);
54 
55  struct Segments {
56  std::vector<float> x, nx;
57  std::vector<float> y, ny;
58  std::vector<float> z, nz;
59  std::vector<float> L;
60  };
61 
62  std::optional<Segments> getSegTracks_(std::vector<float> const& xxx,
63  std::vector<float> const& yyy,
64  std::vector<float> const& zzz,
65  double seg_size);
66 
67  std::tuple<double, double, double> getDeltaThetaRMS_(Segments const& segments,
68  double thick) const;
69 
70  int getDeltaThetaij_(std::vector<float>& ei,
71  std::vector<float>& ej,
72  std::vector<float>& th,
73  std::vector<float>& ind,
74  Segments const& segments,
75  double thick) const;
76 
77  double my_g(double xx, double Q, double s) const;
78 
79  double my_mcs_llhd(std::vector<float> const& dEi,
80  std::vector<float> const& dEj,
81  std::vector<float> const& dthij,
82  std::vector<float> const& ind,
83  double x0, double x1) const;
84 
85  float seg_stop{-1.};
86  int n_seg{};
87 
88  float x_seg[50000];
89  float y_seg[50000];
90  float z_seg[50000];
91 
92  double find_angle(double vz, double vy) const;
93 
94  float steps_size{10.};
95  int n_steps{6};
96  std::vector<float> steps;
97 
98  double minLength;
99  double maxLength;
100 
101  // The following are objects that are created but not drawn or
102  // saved. This class should consider accepting a "debug"
103  // parameter where if it is specified, then the graphs will be
104  // created; otherwise, their creation is unnecessary and impedes
105  // efficiency.
106  //
107  // N.B. TPolyLine3D objects are owned by ROOT, and we thus refer
108  // to them by pointer. It is important that 'delete' is not
109  // called on the TPolyLine3D pointers during destruction of a
110  // TrackMomentumCalculator object.
111  TPolyLine3D* gr_reco_xyz{nullptr};
112  TGraph gr_reco_xy{};
113  TGraph gr_reco_yz{};
114  TGraph gr_reco_xz{};
115 
116  TPolyLine3D* gr_seg_xyz{nullptr};
117  TGraph gr_seg_xy{};
118  TGraph gr_seg_yz{};
119  TGraph gr_seg_xz{};
120 
121  };
122 
123 } // namespace trkf
124 
125 #endif // TrackMomentumCalculator_H
Float_t s
Definition: plot.C:23
Double_t xx
Definition: macro.C:12
double my_g(double xx, double Q, double s) const
double GetMomentumMultiScatterChi2(art::Ptr< recob::Track > const &trk)
int getDeltaThetaij_(std::vector< float > &ei, std::vector< float > &ej, std::vector< float > &th, std::vector< float > &ind, Segments const &segments, double thick) const
Float_t x1[n_points_granero]
Definition: compare.C:5
std::optional< Segments > getSegTracks_(std::vector< float > const &xxx, std::vector< float > const &yyy, std::vector< float > const &zzz, double seg_size)
Provides recob::Track data product.
double find_angle(double vz, double vy) const
std::tuple< double, double, double > getDeltaThetaRMS_(Segments const &segments, double thick) const
double GetMuMultiScatterLLHD3(art::Ptr< recob::Track > const &trk, bool dir)
TDirectory * dir
Definition: macro.C:5
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
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)