SND@LHC Software
Loading...
Searching...
No Matches
RKTrackRep.h
Go to the documentation of this file.
1/* Copyright 2008-2010, Technische Universitaet Muenchen,
2 Authors: Christian Hoeppner & Sebastian Neubert & Johannes Rauch
3
4 This file is part of GENFIT.
5
6 GENFIT is free software: you can redistribute it and/or modify
7 it under the terms of the GNU Lesser General Public License as published
8 by the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 GENFIT is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU Lesser General Public License for more details.
15
16 You should have received a copy of the GNU Lesser General Public License
17 along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18*/
19
24#ifndef genfit_RKTrackRep_h
25#define genfit_RKTrackRep_h
26
27#include "AbsTrackRep.h"
28#include "StateOnPlane.h"
29#include "RKTools.h"
30#include "StepLimits.h"
31
32
33namespace genfit {
34
38struct RKStep {
39 MatStep matStep_; // material properties and stepsize
40 M1x7 state7_; // 7D state vector
42
44 memset(state7_, 0x00, 7*sizeof(double));
45 }
46};
47
48
52struct ExtrapStep {
53 M7x7 jac7_; // 5D jacobian of transport
54 M7x7 noise7_; // 5D noise matrix
55
57 memset(jac7_, 0, sizeof(M7x7));
58 memset(noise7_, 0, sizeof(M7x7));
59 }
60};
61
62
70class RKTrackRep : public AbsTrackRep {
71
72
73 public:
74
75 RKTrackRep();
76 RKTrackRep(int pdgCode, char propDir = 0);
77
78 virtual ~RKTrackRep();
79
80 virtual AbsTrackRep* clone() const {return new RKTrackRep(*this);}
81
82 virtual double extrapolateToPlane(StateOnPlane& state,
83 const SharedPlanePtr& plane,
84 bool stopAtBoundary = false,
85 bool calcJacobianNoise = false) const;
86
88
89 virtual double extrapolateToLine(StateOnPlane& state,
90 const TVector3& linePoint,
91 const TVector3& lineDirection,
92 bool stopAtBoundary = false,
93 bool calcJacobianNoise = false) const;
94
95 virtual double extrapolateToPoint(StateOnPlane& state,
96 const TVector3& point,
97 bool stopAtBoundary = false,
98 bool calcJacobianNoise = false) const {
99 return extrapToPoint(state, point, NULL, stopAtBoundary, calcJacobianNoise);
100 }
101
102 virtual double extrapolateToPoint(StateOnPlane& state,
103 const TVector3& point,
104 const TMatrixDSym& G, // weight matrix (metric)
105 bool stopAtBoundary = false,
106 bool calcJacobianNoise = false) const {
107 return extrapToPoint(state, point, &G, stopAtBoundary, calcJacobianNoise);
108 }
109
110 virtual double extrapolateToCylinder(StateOnPlane& state,
111 double radius,
112 const TVector3& linePoint = TVector3(0.,0.,0.),
113 const TVector3& lineDirection = TVector3(0.,0.,1.),
114 bool stopAtBoundary = false,
115 bool calcJacobianNoise = false) const;
116
117 virtual double extrapolateToSphere(StateOnPlane& state,
118 double radius,
119 const TVector3& point = TVector3(0.,0.,0.),
120 bool stopAtBoundary = false,
121 bool calcJacobianNoise = false) const;
122
123 virtual double extrapolateBy(StateOnPlane& state,
124 double step,
125 bool stopAtBoundary = false,
126 bool calcJacobianNoise = false) const;
127
128
129 unsigned int getDim() const {return 5;}
130
131 virtual TVector3 getPos(const StateOnPlane& state) const;
132
133 virtual TVector3 getMom(const StateOnPlane& state) const;
134 virtual void getPosMom(const StateOnPlane& state, TVector3& pos, TVector3& mom) const;
135
136 virtual double getMomMag(const StateOnPlane& state) const;
137 virtual double getMomVar(const MeasuredStateOnPlane& state) const;
138
139 virtual TMatrixDSym get6DCov(const MeasuredStateOnPlane& state) const;
140 virtual void getPosMomCov(const MeasuredStateOnPlane& state, TVector3& pos, TVector3& mom, TMatrixDSym& cov) const;
141 virtual double getCharge(const StateOnPlane& state) const;
142 virtual double getQop(const StateOnPlane& state) const {return state.getState()(0);}
143 double getSpu(const StateOnPlane& state) const;
144
145 virtual void getForwardJacobianAndNoise(TMatrixD& jacobian, TMatrixDSym& noise, TVectorD& deltaState) const;
146
147 virtual void getBackwardJacobianAndNoise(TMatrixD& jacobian, TMatrixDSym& noise, TVectorD& deltaState) const;
148
149 std::vector<genfit::MatStep> getSteps() const;
150
151 virtual double getRadiationLenght() const;
152
153 virtual double getTOF() const;
154
155
156 virtual void setPosMom(StateOnPlane& state, const TVector3& pos, const TVector3& mom) const;
157 virtual void setPosMom(StateOnPlane& state, const TVectorD& state6) const;
158 virtual void setPosMomErr(MeasuredStateOnPlane& state, const TVector3& pos, const TVector3& mom, const TVector3& posErr, const TVector3& momErr) const;
159 virtual void setPosMomCov(MeasuredStateOnPlane& state, const TVector3& pos, const TVector3& mom, const TMatrixDSym& cov6x6) const;
160 virtual void setPosMomCov(MeasuredStateOnPlane& state, const TVectorD& state6, const TMatrixDSym& cov6x6) const;
161
162 virtual void setChargeSign(StateOnPlane& state, double charge) const;
163 virtual void setQop(StateOnPlane& state, double qop) const {state.getState()(0) = qop;}
164
165 void setSpu(StateOnPlane& state, double spu) const;
166
168
175 double RKPropagate(M1x7& state7,
176 M7x7* jacobian,
177 M1x3& SA,
178 double S,
179 bool varField = true,
180 bool calcOnlyLastRowOfJ = false) const;
181
182 virtual bool isSameType(const AbsTrackRep* other);
183 virtual bool isSame(const AbsTrackRep* other);
184
185 private:
186
187 void initArrays() const;
188
189 virtual double extrapToPoint(StateOnPlane& state,
190 const TVector3& point,
191 const TMatrixDSym* G = NULL, // weight matrix (metric)
192 bool stopAtBoundary = false,
193 bool calcJacobianNoise = false) const;
194
195 void getState7(const StateOnPlane& state, M1x7& state7) const;
196 void getState5(StateOnPlane& state, const M1x7& state7) const; // state7 must already lie on plane of state!
197
198 void transformPM7(const MeasuredStateOnPlane& state,
199 M7x7& out7x7) const;
200
201 void calcJ_pM_5x7(M5x7& J_pM, const TVector3& U, const TVector3& V, const M1x3& pTilde, double spu) const;
202
203 void transformPM6(const MeasuredStateOnPlane& state,
204 M6x6& out6x6) const;
205
206 void transformM7P(const M7x7& in7x7,
207 const M1x7& state7,
208 MeasuredStateOnPlane& state) const; // plane must already be set!
209
210 void calcJ_Mp_7x5(M7x5& J_Mp, const TVector3& U, const TVector3& V, const TVector3& W, const M1x3& A) const;
211
212 void calcForwardJacobianAndNoise(const M1x7& startState7, const DetPlane& startPlane,
213 const M1x7& destState7, const DetPlane& destPlane) const;
214
215 void transformM6P(const M6x6& in6x6,
216 const M1x7& state7,
217 MeasuredStateOnPlane& state) const; // plane and charge must already be set!
218
220
229 bool RKutta(const M1x4& SU,
230 const DetPlane& plane,
231 double charge,
232 M1x7& state7,
233 M7x7* jacobianT,
234 double& coveredDistance, // signed
235 bool& checkJacProj,
236 M7x7& noiseProjection,
237 StepLimits& limits,
238 bool onlyOneStep = false,
239 bool calcOnlyLastRowOfJ = false) const;
240
241 double estimateStep(const M1x7& state7,
242 const M1x4& SU,
243 const DetPlane& plane,
244 const double& charge,
245 double& relMomLoss,
246 StepLimits& limits) const;
247
248 TVector3 pocaOnLine(const TVector3& linePoint,
249 const TVector3& lineDirection,
250 const TVector3& point) const;
251
253
261 double Extrap(const DetPlane& startPlane, // plane where Extrap starts
262 const DetPlane& destPlane, // plane where Extrap has to extrapolate to
263 double charge,
264 bool& isAtBoundary,
265 M1x7& state7,
266 bool fillExtrapSteps,
267 TMatrixDSym* cov = nullptr,
268 bool onlyOneStep = false,
269 bool stopAtBoundary = false,
270 double maxStep = 1.E99) const;
271
272 void checkCache(const StateOnPlane& state, const SharedPlanePtr* plane) const;
273
274 double momMag(const M1x7& state7) const;
275
276
279 mutable std::vector<RKStep> RKSteps_;
280 mutable int RKStepsFXStart_;
281 mutable int RKStepsFXStop_;
282 mutable std::vector<ExtrapStep> ExtrapSteps_;
283
284 mutable TMatrixD fJacobian_;
285 mutable TMatrixDSym fNoise_;
286
287 mutable bool useCache_;
288 mutable unsigned int cachePos_;
289
290 // auxiliary variables and arrays
291 // needed in Extrap()
295 mutable M7x7 J_MMT_;
296
297 public:
298
299 ClassDef(RKTrackRep, 1)
300
301};
302
303} /* End of namespace genfit */
306#endif // genfit_RKTrackRep_h
Abstract base class for a track representation.
Definition AbsTrackRep.h:66
virtual double extrapolateToLine(StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const =0
Extrapolates the state to the POCA to a line, and returns the extrapolation length and,...
Detector plane.
Definition DetPlane.h:61
StateOnPlane with additional covariance matrix.
AbsTrackRep with 5D track parameterization in plane coordinates: (q/p, u', v', u, v)
Definition RKTrackRep.h:70
virtual TVector3 getMom(const StateOnPlane &state) const
Get the cartesian momentum vector of a state.
void initArrays() const
std::vector< RKStep > RKSteps_
state where the last extrapolation has ended
Definition RKTrackRep.h:279
StepLimits limits_
Definition RKTrackRep.h:292
virtual void setPosMomErr(MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TVector3 &posErr, const TVector3 &momErr) const
Set position and momentum and error of state.
void calcJ_Mp_7x5(M7x5 &J_Mp, const TVector3 &U, const TVector3 &V, const TVector3 &W, const M1x3 &A) const
std::vector< ExtrapStep > ExtrapSteps_
Definition RKTrackRep.h:282
double getSpu(const StateOnPlane &state) const
virtual void getBackwardJacobianAndNoise(TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const
Get the jacobian and noise matrix of the last extrapolation if it would have been done in opposite di...
virtual void setQop(StateOnPlane &state, double qop) const
Set charge/momentum.
Definition RKTrackRep.h:163
M7x7 noiseProjection_
noise matrix of the last extrapolation
Definition RKTrackRep.h:294
virtual void setChargeSign(StateOnPlane &state, double charge) const
Set the sign of the charge according to charge.
int RKStepsFXStart_
RungeKutta steps made in the last extrapolation.
Definition RKTrackRep.h:280
double estimateStep(const M1x7 &state7, const M1x4 &SU, const DetPlane &plane, const double &charge, double &relMomLoss, StepLimits &limits) const
virtual double getQop(const StateOnPlane &state) const
Get charge over momentum.
Definition RKTrackRep.h:142
void setSpu(StateOnPlane &state, double spu) const
virtual double getCharge(const StateOnPlane &state) const
Get the (fitted) charge of a state. This is not always equal the pdg charge (e.g. if the charge sign ...
virtual void getPosMom(const StateOnPlane &state, TVector3 &pos, TVector3 &mom) const
Get cartesian position and momentum vector of a state.
void calcForwardJacobianAndNoise(const M1x7 &startState7, const DetPlane &startPlane, const M1x7 &destState7, const DetPlane &destPlane) const
double Extrap(const DetPlane &startPlane, const DetPlane &destPlane, double charge, bool &isAtBoundary, M1x7 &state7, bool fillExtrapSteps, TMatrixDSym *cov=nullptr, bool onlyOneStep=false, bool stopAtBoundary=false, double maxStep=1.E99) const
Handles propagation and material effects.
bool RKutta(const M1x4 &SU, const DetPlane &plane, double charge, M1x7 &state7, M7x7 *jacobianT, double &coveredDistance, bool &checkJacProj, M7x7 &noiseProjection, StepLimits &limits, bool onlyOneStep=false, bool calcOnlyLastRowOfJ=false) const
Propagates the particle through the magnetic field.
double momMag(const M1x7 &state7) const
virtual double extrapolateToSphere(StateOnPlane &state, double radius, const TVector3 &point=TVector3(0., 0., 0.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state to the sphere surface, and returns the extrapolation length and,...
void transformPM7(const MeasuredStateOnPlane &state, M7x7 &out7x7) const
virtual TMatrixDSym get6DCov(const MeasuredStateOnPlane &state) const
Get the 6D covariance.
void getState7(const StateOnPlane &state, M1x7 &state7) const
void transformM7P(const M7x7 &in7x7, const M1x7 &state7, MeasuredStateOnPlane &state) const
virtual double extrapolateBy(StateOnPlane &state, double step, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state by step (cm) and returns the extrapolation length and, via reference,...
std::vector< genfit::MatStep > getSteps() const
Get stepsizes and material properties of crossed materials of the last extrapolation.
virtual double extrapToPoint(StateOnPlane &state, const TVector3 &point, const TMatrixDSym *G=NULL, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
virtual void setPosMom(StateOnPlane &state, const TVector3 &pos, const TVector3 &mom) const
Set position and momentum of state.
virtual double getTOF() const
Get the time of flight [ns] of the last extrapolation.
virtual double extrapolateToCylinder(StateOnPlane &state, double radius, const TVector3 &linePoint=TVector3(0., 0., 0.), const TVector3 &lineDirection=TVector3(0., 0., 1.), bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state to the cylinder surface, and returns the extrapolation length and,...
void transformPM6(const MeasuredStateOnPlane &state, M6x6 &out6x6) const
void checkCache(const StateOnPlane &state, const SharedPlanePtr *plane) const
virtual void setPosMomCov(MeasuredStateOnPlane &state, const TVector3 &pos, const TVector3 &mom, const TMatrixDSym &cov6x6) const
Set position, momentum and covariance of state.
StateOnPlane lastStartState_
Definition RKTrackRep.h:277
unsigned int getDim() const
Get the dimension of the state vector used by the track representation.
Definition RKTrackRep.h:129
virtual ~RKTrackRep()
Definition RKTrackRep.cc:75
virtual bool isSame(const AbsTrackRep *other)
check if other is of same type (e.g. RKTrackRep) and has same pdg code.
virtual void getForwardJacobianAndNoise(TMatrixD &jacobian, TMatrixDSym &noise, TVectorD &deltaState) const
Get the jacobian and noise matrix of the last extrapolation.
virtual double getRadiationLenght() const
Get the accumulated X/X0 (path / radiation length) of the material crossed in the last extrapolation.
TMatrixD fJacobian_
steps made in Extrap during last extrapolation
Definition RKTrackRep.h:284
virtual double getMomVar(const MeasuredStateOnPlane &state) const
get the variance of the absolute value of the momentum .
unsigned int cachePos_
use cached RKSteps_ for extrapolation
Definition RKTrackRep.h:288
void calcJ_pM_5x7(M5x7 &J_pM, const TVector3 &U, const TVector3 &V, const M1x3 &pTilde, double spu) const
virtual double getMomMag(const StateOnPlane &state) const
get the magnitude of the momentum in GeV.
void getState5(StateOnPlane &state, const M1x7 &state7) const
virtual AbsTrackRep * clone() const
Clone the trackRep.
Definition RKTrackRep.h:80
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
TVector3 pocaOnLine(const TVector3 &linePoint, const TVector3 &lineDirection, const TVector3 &point) const
virtual double extrapolateToPoint(StateOnPlane &state, const TVector3 &point, const TMatrixDSym &G, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state to the POCA to a point in the metric of G, and returns the extrapolation lengt...
Definition RKTrackRep.h:102
virtual double extrapolateToPlane(StateOnPlane &state, const SharedPlanePtr &plane, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state to plane, and returns the extrapolation length and, via reference,...
Definition RKTrackRep.cc:97
virtual TVector3 getPos(const StateOnPlane &state) const
Get the cartesian position of a state.
TMatrixDSym fNoise_
Definition RKTrackRep.h:285
void transformM6P(const M6x6 &in6x6, const M1x7 &state7, MeasuredStateOnPlane &state) const
StateOnPlane lastEndState_
state where the last extrapolation has started
Definition RKTrackRep.h:278
virtual bool isSameType(const AbsTrackRep *other)
check if other is of same type (e.g. RKTrackRep).
virtual double extrapolateToLine(StateOnPlane &state, const TVector3 &linePoint, const TVector3 &lineDirection, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state to the POCA to a line, and returns the extrapolation length and,...
virtual double extrapolateToPoint(StateOnPlane &state, const TVector3 &point, bool stopAtBoundary=false, bool calcJacobianNoise=false) const
Extrapolates the state to the POCA to a point, and returns the extrapolation length and,...
Definition RKTrackRep.h:95
double RKPropagate(M1x7 &state7, M7x7 *jacobian, M1x3 &SA, double S, bool varField=true, bool calcOnlyLastRowOfJ=false) const
The actual Runge Kutta propagation.
A state with arbitrary dimension defined in a DetPlane.
const TVectorD & getState() const
Helper to store different limits on the stepsize for the RKTRackRep.
Definition StepLimits.h:54
Matrix inversion tools.
Definition AbsBField.h:29
double M6x6[6 *6]
Definition RKTools.h:38
double M5x7[5 *7]
Definition RKTools.h:44
double M1x7[1 *7]
Definition RKTools.h:36
double M1x4[1 *4]
Definition RKTools.h:34
double M7x5[7 *5]
Definition RKTools.h:42
double M1x3[1 *3]
Definition RKTools.h:33
boost::shared_ptr< genfit::DetPlane > SharedPlanePtr
Shared Pointer to a DetPlane.
double M7x7[7 *7]
Definition RKTools.h:39
Helper for RKTrackRep.
Definition RKTrackRep.h:52
Simple struct containing MaterialProperties and stepsize in the material.
Definition AbsTrackRep.h:42
Helper for RKTrackRep.
Definition RKTrackRep.h:38
StepLimits limits_
Definition RKTrackRep.h:41
MatStep matStep_
Definition RKTrackRep.h:39