Class Reference for E1039 Core & Analysis Software
KalmanFitter.h
Go to the documentation of this file.
1 #ifndef _KALMANFITTER_H
2 #define _KALMANFITTER_H
3 
4 /*
5 KalmanFitter.h
6 
7 Definition of the class KalmanFitter, refine the measurement of KalmanTrack
8 
9 Author: Kun Liu, liuk@fnal.gov
10 Created: 10-18-2012
11 */
12 
13 #include <list>
14 #include <vector>
15 
16 #include <GlobalConsts.h>
17 
18 #include "KalmanUtil.h"
19 #include "KalmanTrack.h"
20 #include "KalmanFilter.h"
21 
22 class PHField;
23 
25 {
26 public:
27  KalmanFitter(const PHField* field, const TGeoManager *geom);
28  //~KalmanFitter();
29 
31  void setControlParameter(int nMaxIteration, double tolerance) { _max_iteration = nMaxIteration; _tolerance = tolerance; }
32 
36  int processOneTrack(KalmanTrack& _track);
37  void updateTrack(KalmanTrack& _track);
38 
40  void init();
41 
43  int initNodeList(KalmanTrack& _track);
44 
46  void updateAlignment();
47 
49  void setStartingParameter(KalmanTrack& _track);
50  void setStartingParameter(Node& _node);
51 
53  bool initSmoother(Node& _node);
54 
56  std::list<Node>& getNodeList() { return _nodes; }
57 
59  double findVertex();
60 
62  double getMomentumInitial(double& px, double& py, double& pz) { return _nodes.front().getSmoothed().get_mom(px, py, pz); }
63  double getMomentumFinal(double& px, double& py, double& pz) { return _nodes.back().getSmoothed().get_mom(px, py, pz); }
64 
65  double getPositionInitial(double& x, double& y, double& z) { return _nodes.front().getSmoothed().get_pos(x, y, z); }
66  double getPositionFinal(double& x, double& y, double& z) { return _nodes.back().getSmoothed().get_pos(x, y, z); }
67 
68  double getChisq() { return _chisq; }
69 
70  const TrkPar& getTrkParInitial() { return _nodes.front().getSmoothed(); }
71  const TrkPar& getTrkParFinal() { return _nodes.back().getSmoothed(); }
72 
73  //TrkPar getTrkPar(double z) { std::cout << "Will be implemented later" << std::endl; }
74 
75 private:
77  std::list<Node> _nodes;
78 
80  double _chisq;
81 
83  KalmanFilter *_kmfit;
84 
86  double rM_20[nChamberPlanes];
87  double rM_21[nChamberPlanes];
88  double rM_22[nChamberPlanes];
89  double z_planes[nChamberPlanes];
90 
92  int _max_iteration;
93  double _tolerance;
94 };
95 
96 #endif
#define nChamberPlanes
Definition: GlobalConsts.h:6
int processOneTrack(KalmanTrack &_track)
void setStartingParameter(KalmanTrack &_track)
Set the starting parameters.
std::list< Node > & getNodeList()
return the node list
Definition: KalmanFitter.h:56
bool initSmoother(Node &_node)
Initialize the smoother.
int initNodeList(KalmanTrack &_track)
Initialize the node list.
double getPositionInitial(double &x, double &y, double &z)
Definition: KalmanFitter.h:65
void updateTrack(KalmanTrack &_track)
double getPositionFinal(double &x, double &y, double &z)
Definition: KalmanFitter.h:66
double getMomentumFinal(double &px, double &py, double &pz)
Definition: KalmanFitter.h:63
double getMomentumInitial(double &px, double &py, double &pz)
Get the final results – temporary interfaces, just for debugging purposes.
Definition: KalmanFitter.h:62
const TrkPar & getTrkParInitial()
Definition: KalmanFitter.h:70
const TrkPar & getTrkParFinal()
Definition: KalmanFitter.h:71
KalmanFitter(const PHField *field, const TGeoManager *geom)
void setControlParameter(int nMaxIteration, double tolerance)
Set the convergence control parameters.
Definition: KalmanFitter.h:31
double getChisq()
Definition: KalmanFitter.h:68
void updateAlignment()
Update the actual z position of each node according to current fit results.
void init()
Initialize the kalman filter.
double findVertex()
Vertex finder.
transient DST object for field storage and access
Definition: PHField.h:14