MyMarlinTPC  170316
RowTripletBasedTrackFinderProcessor.h
Go to the documentation of this file.
1 #ifndef ROWTRIPLETBASEDTRACKFINDERPROCESSOR_H
2 #define ROWTRIPLETBASEDTRACKFINDERPROCESSOR_H
3 
4 #include <marlin/Processor.h>
5 #include <marlin/Global.h>
6 #include <lcio.h>
7 #include <string>
8 
9 #include <gear/TPCModule.h>
10 
11 // LCIO:
12 #include "EVENT/TrackerHit.h"
13 
14 #include "TFile.h"
15 #include "TTree.h"
16 #include "TMath.h"
17 #include "TVectorD.h"
18 #include "TMatrixD.h"
19 #include "TMatrixDSym.h"
20 
28 namespace marlintpc {
30 
77 class RowTripletBasedTrackFinderProcessor: public marlin::Processor {
78  public:
79 
80  virtual Processor* newProcessor() {
82  }
83 
85 
86  virtual void init();
87 
88  virtual void processRunHeader(EVENT::LCRunHeader* run);
89 
90  virtual void processEvent(EVENT::LCEvent* evt);
91 
92  virtual void check(EVENT::LCEvent* evt);
93 
94  virtual void end();
95 
96  protected:
97  /* the place for protected and private member data and functions */
98  std::string _inputColName;
99  std::string _outputColName;
101  double _distCut;
102  double _trpCut;
103  double _posCut;
104  double _dirCut;
105  double _segCut;
106  double _unusedCut;
107  int _maxGap;
109  bool _refAtPCA;
110 
111  private:
112 
113  double _Bzc;
114  std::map<int, std::vector<int> > _modNeighbours;
115  double _Xcenter;
116  double _Ycenter;
117 
118  bool areNeighbourModules(gear::TPCModule*, gear::TPCModule*, bool);
119 };
120 
122 
125 class rb_Hit {
126  public:
127  rb_Hit(const int iHit, const int mod, const int row, const EVENT::TrackerHit& aHit);
128  void print() const;
129  int getHitNum() const;
130  int getMod() const;
131  int getRow() const;
132  double getX() const;
133  double getY() const;
134  double getZ() const;
135  void getPos(double*) const;
136  double getVarXY(const double = 0.) const;
137  double getVarZ(const double = 0.) const;
138  double getPhiMeas() const;
139  bool getUsed() const;
140  void setUsed(const double);
141  double getCosBeta() const;
142  double getDistXY(const double, const double, const double, const double) const;
143  double getDistZ(const double) const;
144 
145  private:
146  const int _hitNum;
147  const int _hitModule;
148  const int _hitRow;
149  const double _hitX;
150  const double _hitY;
151  const double _hitZ;
152  const double _hitPhiMeas;
153  const double _hitVarXY;
154  const double _hitVarR;
155  const double _hitVarZ;
156  bool _used;
157  double _cosb;
158 
159  double _calcPhiMeas(const EVENT::TrackerHit& aHit);
160  double _getVarXY(FloatVec);
161  double _getVarR(FloatVec);
162  double _getVarZ(FloatVec);
163 
164 };
165 
166 typedef std::vector<rb_Hit*> hitListType;
167 typedef std::map<int, hitListType> rowHitMapType;
168 typedef std::map<int, rowHitMapType> modHitMapType;
169 
171 
176 class rb_Doublet {
177  public:
179  bool match(rb_Hit*, const double, const double) const;
180  rb_Hit* getHit(const int) const;
181  void getParameters(double&, double&, double&, double&, double&, double&, double&, double&, double&, double&) const;
182 
183  private:
184  rb_Hit* _hit[2];
185  double _phiMeas;
186  double _cosPhi;
187  double _sinPhi;
188  double _xav;
189  double _yav;
190  double _zav;
191  double _phi;
192  double _tanl;
193  double _length;
194  double _cosb;
195  double _der2XY;
196  double _der2Z;
197  double _varXY;
198  double _varZ;
199 };
200 
202 
207 class rb_Triplet {
208  public:
210  int getRow() const;
211  double getX() const;
212  double getY() const;
213  double getZ() const;
214  void getPos(double*) const;
215  double getVarXYPos(const double = 0.) const;
216  double getVarXYDir() const;
217  double getVarZPos(const double = 0.) const;
218  double getVarZDir() const;
219  double getPhiMeas() const;
220  double getPhi() const;
221  double getTanl() const;
222  void getPosCloseTo(const double, const double, double*, double&) const;
223  rb_Hit* getHit(const int) const;
224  bool match(rb_Triplet*, const double, const double) const;
225  double getCosBeta() const;
226 
227  private:
228  rb_Hit* _hit[3];
229  double _phiMeas;
230  double _posX;
231  double _posY;
232  double _posZ;
233  double _varXY[2];
234  double _phi;
235  double _tanl;
236  double _varZ[2];
237  double _length;
238  double _cosb;
239 };
240 
241 typedef std::vector<rb_Triplet*> trpListType;
242 
244 
249 class rb_Segment {
250  public:
251  rb_Segment(double, trpListType);
252  rb_Segment(double, hitListType, unsigned int = 1);
253  rb_Segment(std::vector<rb_Segment*> segments);
254  double getBzc() const;
255  int getFirstRow() const;
256  int getLastRow() const;
257  void getRefPoint(double*) const;
258  int getNdf() const;
259  double getChi2() const;
260  void getPar(double *) const;
261  const hitListType& getHitList() const;
262  const TVectorD& getPar() const;
263  const TMatrixDSym& getCov() const;
264  bool match(rb_Segment*, const double, const int) const;
265  bool match(rb_Hit*, const double) const;
266  void match(hitListType&, hitListType&, const double) const;
267  void getLCIOStateAtRefPoint(const double*, TVectorD&, TMatrixDSym&) const;
268  void fillRowMap(std::map<int, int>&) const;
269  int getId() const;
270  void setId(int);
271  void addHit(rb_Hit *);
272 
273  private:
274  int _segId;
275  hitListType _hitList;
276  double _bzc;
277  double _refX;
278  double _refY;
279  double _refZ;
280  int _npar;
281  int _ndf;
282  double _chi2;
283  TVectorD _parameters;
284  TMatrixDSym _covariance;
285  void calcRefPoint();
286  void fitSegment(double);
287 };
288 
289 typedef std::vector<rb_Segment*> segListType;
290 typedef std::map<int, segListType> modSegMapType;
291 
293 
299  public:
301  void addIndex(int);
302  void addMatch(std::pair<int, int>);
303  std::map<int, std::vector<int> > getClasses();
304 
305  private:
306  std::map<int, int> _index;
307  std::vector<std::pair<int, int> > _matches;
308 };
309 
310 } // namespace marlintpc
311 #endif // ROWTRIPLETBASEDTRACKFINDERPROCESSOR_H
double _segCut
Chi2/Ndf cut for segment matching (Ndf=4 for B off, 5 for B on)
const double _hitVarR
variance in radial direction
std::map< int, std::vector< int > > _modNeighbours
neighbour modules
int _npar
number of parameters (5: helix, 4: line)
std::vector< rb_Triplet * > trpListType
std::map< int, segListType > modSegMapType
double _chi2
chi2 from segment fit
double _phiMeas
measurement direction in XY (tangential to row)
double _refY
Y of reference point.
double _posCut
Chi2 cut for triplet position matching in XY and Z.
double _trpCut
Chi2 cut for triplet definition on XY and Z residuals.
bool areNeighbourModules(gear::TPCModule *, gear::TPCModule *, bool)
Check if modules are neighbours.
double _bzc
magnetic field strength (Bz*c)
std::string _outputColName
Name of the output collection.
TMatrixDSym _covariance
covariance matrix
double _phiMeas
(average) measurement direction in XY
double _unusedCut
Chi2 cut for matching unused hits in XY and Z.
double _distCut
Coarse cut on XY and Z residuals for triplet preselection.
double _dirCut
Chi2 cut for triplet direction matching in XY and Z.
int _ndf
number of degrees of freedom of segment fit
bool _refAtPCA
Use Pca as reference point (else 1. hit)
const int _hitNum
input hit collection index
std::vector< rb_Segment * > segListType
std::vector< rb_Hit * > hitListType
std::vector< std::pair< int, int > > _matches
list of matches
int _maxGap
Cut for (row) distance to segment for matching unused hits in XY and Z.
const double _hitPhiMeas
measurement direction in XY (tangential to row)
double _refZ
Z of reference point.
std::map< int, rowHitMapType > modHitMapType
double _bfieldScaleFactor
scale factor for magnetic field (default: 1.0)
hitListType _hitList
pointers to hits
double _refX
X of reference point.
virtual void processEvent(EVENT::LCEvent *evt)
Process event.
std::map< int, hitListType > rowHitMapType
const double _hitVarZ
variance for Z measurement
const double _hitVarXY
variance for XY measurement