PandaRoot
GFAbsTrackRep.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 /* Copyright 2008-2010, Technische Universitaet Muenchen,
14  Authors: Christian Hoeppner & Sebastian Neubert
15 
16  This file is part of GENFIT.
17 
18  GENFIT is free software: you can redistribute it and/or modify
19  it under the terms of the GNU Lesser General Public License as published
20  by the Free Software Foundation, either version 3 of the License, or
21  (at your option) any later version.
22 
23  GENFIT is distributed in the hope that it will be useful,
24  but WITHOUT ANY WARRANTY; without even the implied warranty of
25  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
26  GNU Lesser General Public License for more details.
27 
28  You should have received a copy of the GNU Lesser General Public License
29  along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
30 */
31 
36 #ifndef GFABSTRACKREP_H
37 #define GFABSTRACKREP_H
38 
39 #include <vector>
40 #include <list>
41 #include <iostream>
42 
43 #include "TMatrixT.h"
44 #include "TVector3.h"
45 
46 #include "GFDetPlane.h"
47 
48 class GFAbsRecoHit;
49 
92 class GFAbsTrackRep : public TObject {
93 
94  /*----- Data members -----*/
95  protected:
97  unsigned int fDimension;
98 
100  TMatrixT<double> fState;
101 
103  TMatrixT<double> fCov;
104 
106  double fChiSqu;
107  unsigned int fNdf;
108 
112  bool fInverted;
113 
115  TMatrixT<double> fFirstState;
116  TMatrixT<double> fFirstCov;
117 
118  TMatrixT<double> fLastState;
119  TMatrixT<double> fLastCov;
122 
123  // detector plane where the track parameters are given
125 
126  public:
127  virtual GFAbsTrackRep *clone() const = 0;
128 
129  virtual GFAbsTrackRep *prototype() const = 0;
130 
132 
137  virtual double extrapolate(const GFDetPlane &plane, TMatrixT<double> &statePred);
138 
139  public:
140  GFAbsTrackRep();
141  GFAbsTrackRep(int);
142  virtual ~GFAbsTrackRep();
143 
145 
155  virtual void extrapolateToPoint(const TVector3 &point, TVector3 &poca, TVector3 &normVec);
156 
158 
165  virtual void extrapolateToLine(const TVector3 &point1, const TVector3 &point2, TVector3 &poca, TVector3 &normVec, TVector3 &poca_onwire);
166 
168 
172  virtual double stepalong(double h, TVector3 &point, TVector3 &dir);
173 
175 
178  virtual double extrapolate(const GFDetPlane &plane, TMatrixT<double> &statePred, TMatrixT<double> &covPred) = 0;
179 
181 
184  double extrapolate(const GFDetPlane &plane);
185 
187  unsigned int getDim() const { return fDimension; }
188 
189  virtual void Print() const;
190 
191  inline TMatrixT<double> getState() const { return fState; }
192  inline TMatrixT<double> getCov() const { return fCov; }
193 
194  double getStateElem(int i) const { return fState(i, 0); }
195  double getCovElem(int i, int j) const { return fCov(i, j); }
196 
197  virtual TVector3 getPos(const GFDetPlane &pl) = 0;
198  virtual TVector3 getMom(const GFDetPlane &pl) = 0;
199  virtual void getPosMom(const GFDetPlane &pl, TVector3 &pos, TVector3 &mom) = 0;
200 
202 
206  virtual void getPosMomCov(const GFDetPlane &pl, TVector3 &pos, TVector3 &mom, TMatrixT<double> &cov);
207 
208  virtual double getCharge() const = 0;
209 
210  TVector3 getPos() { return getPos(fRefPlane); }
211  TVector3 getMom() { return getMom(fRefPlane); }
212  void getPosMomCov(TVector3 &pos, TVector3 &mom, TMatrixT<double> &c) { getPosMomCov(fRefPlane, pos, mom, c); }
213 
214  inline TMatrixT<double> getFirstState() const { return fFirstState; }
215  inline TMatrixT<double> getFirstCov() const { return fFirstCov; }
216  inline GFDetPlane getFirstPlane() const { return fFirstPlane; }
217  inline TMatrixT<double> getLastState() const { return fLastState; }
218  inline TMatrixT<double> getLastCov() const { return fLastCov; }
219  inline GFDetPlane getLastPlane() const { return fLastPlane; }
220  inline double getChiSqu() const { return fChiSqu; }
222  inline double getRedChiSqu() const
223  {
224  if (getNDF() > 0)
225  return getChiSqu() / getNDF();
226  return 0;
227  }
228  inline unsigned int getNDF() const
229  {
230  if (fNdf > getDim())
231  return fNdf - getDim();
232  return 0;
233  }
245  virtual void setData(const TMatrixT<double> &st, const GFDetPlane &pl, const TMatrixT<double> *cov = nullptr, const TMatrixT<double> *aux = nullptr)
246  {
247  fState = st;
248  fRefPlane = pl;
249  if (cov != nullptr)
250  fCov = *cov;
251  if (aux != nullptr) {
252  ;
253  }
254  }
255  inline void setCov(const TMatrixT<double> &aCov) { fCov = aCov; }
256  inline void setFirstState(const TMatrixT<double> &aState) { fFirstState = aState; }
257  inline void setFirstCov(const TMatrixT<double> &aCov) { fFirstCov = aCov; }
258  inline void setFirstPlane(const GFDetPlane &aPlane)
259  {
260  fFirstPlane = aPlane;
261  ;
262  }
263  inline void setLastState(const TMatrixT<double> &aState) { fLastState = aState; }
264  inline void setLastCov(const TMatrixT<double> &aCov) { fLastCov = aCov; }
265  inline void setLastPlane(const GFDetPlane &aPlane)
266  {
267  fLastPlane = aPlane;
268  ;
269  }
270 
271  const GFDetPlane &getReferencePlane() const { return fRefPlane; }
272 
273  inline void setChiSqu(double aChiSqu) { fChiSqu = aChiSqu; }
274  inline void setNDF(unsigned int n) { fNdf = n; }
275  inline void addChiSqu(double aChiSqu) { fChiSqu += aChiSqu; }
276  inline void addNDF(unsigned int n) { fNdf += n; }
277  inline void setStatusFlag(int _val) { fStatusFlag = _val; }
278 
279  virtual void switchDirection() = 0;
280 
282  bool setInverted(bool f = true)
283  {
284  fInverted = f;
285  return true;
286  }
287 
288  inline bool getStatusFlag() { return fStatusFlag; }
289 
290  virtual void reset();
291 
297  virtual bool hasAuxInfo() { return false; }
298 
308  virtual const TMatrixT<double> *getAuxInfo(const GFDetPlane &)
309  { //(const GFDetPlane& pl)
310  return nullptr;
311  }
312 
313  private:
314  void Abort(std::string method);
315 
316  ClassDef(GFAbsTrackRep, 3)
317 };
318 
319 #endif
320 
TVector3 getMom()
virtual void extrapolateToLine(const TVector3 &point1, const TVector3 &point2, TVector3 &poca, TVector3 &normVec, TVector3 &poca_onwire)
This method extrapolates to the point of closest approach to a line.
virtual double stepalong(double h, TVector3 &point, TVector3 &dir)
make step of h cm along the track
Base Class for genfit track representations. Defines interface for track parameterizations.
Definition: GFAbsTrackRep.h:92
unsigned int getDim() const
returns dimension of state vector
void addNDF(unsigned int n)
void setFirstState(const TMatrixT< double > &aState)
Detector plane genfit geometry class.
Definition: GFDetPlane.h:70
virtual void Print() const
virtual void getPosMom(const GFDetPlane &pl, TVector3 &pos, TVector3 &mom)=0
virtual ~GFAbsTrackRep()
TMatrixT< double > fLastState
void setFirstPlane(const GFDetPlane &aPlane)
TMatrixT< double > getCov() const
TMatrixT< double > getState() const
bool setInverted(bool f=true)
Deprecated. Should be removed soon.
TMatrixT< double > fFirstCov
TMatrixT< double > getLastCov() const
GFDetPlane getLastPlane() const
virtual void getPosMomCov(const GFDetPlane &pl, TVector3 &pos, TVector3 &mom, TMatrixT< double > &cov)
method which gets position, momentum and 6x6 covariance matrix
unsigned int i
Definition: P4_F32vec4.h:33
virtual void setData(const TMatrixT< double > &st, const GFDetPlane &pl, const TMatrixT< double > *cov=nullptr, const TMatrixT< double > *aux=nullptr)
Puts the track representation in a given state.
double getStateElem(int i) const
unsigned int fDimension
Dimensionality of track representation.
Definition: GFAbsTrackRep.h:97
unsigned int fNdf
TMatrixT< double > fFirstState
state, cov and plane for first and last point in fit
virtual GFAbsTrackRep * prototype() const =0
unsigned int getNDF() const
GFDetPlane getFirstPlane() const
void setCov(const TMatrixT< double > &aCov)
virtual const TMatrixT< double > * getAuxInfo(const GFDetPlane &)
Get auxillary information from the track representation.
bool getStatusFlag()
void setLastPlane(const GFDetPlane &aPlane)
virtual bool hasAuxInfo()
See if the track representation has auxillary information stored.
virtual void reset()
void setLastCov(const TMatrixT< double > &aCov)
TMatrixT< double > getLastState() const
bool fInverted
specifies the direction of flight of the particle
GFDetPlane fFirstPlane
void addChiSqu(double aChiSqu)
GFDetPlane fRefPlane
TVector3 getPos()
virtual double extrapolate(const GFDetPlane &plane, TMatrixT< double > &statePred)
returns the tracklength spanned in this extrapolation
TMatrixT< double > getFirstCov() const
Base Class for representing a Hit in GENFIT.
Definition: GFAbsRecoHit.h:83
GFDetPlane fLastPlane
double getRedChiSqu() const
returns chi2/ndf
double getChiSqu() const
void setNDF(unsigned int n)
void setChiSqu(double aChiSqu)
double getCovElem(int i, int j) const
int fStatusFlag
status of track representation: 0 means everything&#39;s OK
void setFirstCov(const TMatrixT< double > &aCov)
float f
Definition: P4_F32vec4.h:32
void setLastState(const TMatrixT< double > &aState)
void setStatusFlag(int _val)
const GFDetPlane & getReferencePlane() const
TMatrixT< double > fState
The vector of track parameters.
double fChiSqu
chiSqu of the track fit
virtual void extrapolateToPoint(const TVector3 &point, TVector3 &poca, TVector3 &normVec)
This method is to extrapolate the track to point of closest approach to a point in space...
virtual void switchDirection()=0
TMatrixT< double > fLastCov
virtual GFAbsTrackRep * clone() const =0
void getPosMomCov(TVector3 &pos, TVector3 &mom, TMatrixT< double > &c)
virtual double getCharge() const =0
TMatrixT< double > getFirstState() const
TMatrixT< double > fCov
The covariance matrix.