ATLAS Offline Software
Loading...
Searching...
No Matches
Trk::FitParameters Class Reference

#include <FitParameters.h>

Collaboration diagram for Trk::FitParameters:

Public Member Functions

 FitParameters (const Perigee &perigee)
 FitParameters (double d0, double z0, double cosPhi, double sinPhi, double cotTheta, double ptInv0, const PerigeeSurface &surface)
 FitParameters (const FitParameters &parameters)=default
 FitParameters (FitParameters &&)=default
FitParametersoperator= (const FitParameters &)=default
FitParametersoperator= (FitParameters &&)=default
 ~FitParameters (void)=default
void addAlignment (bool constrained, double localAngle, double localOffset)
void addScatterer (double phi, double theta)
double alignmentAngle (int alignment) const
double alignmentAngleConstraint (int alignment) const
double alignmentOffset (int alignment) const
double alignmentOffsetConstraint (int alignment) const
const SurfaceassociatedSurface (void) const
double cosPhi (void) const
double cosTheta (void) const
double cotTheta (void) const
void covariance (Amg::MatrixX *finalCovariance, const Amg::MatrixX *fullCovariance)
double d0 (void) const
void d0 (double value)
double difference (int param) const
const Amg::VectorXdifferences (void) const
Amg::Vector3D direction (void) const
bool extremeMomentum (void) const
void extremeMomentum (bool value)
const Amg::MatrixXfinalCovariance (void) const
int firstAlignmentParameter (void) const
void firstAlignmentParameter (int value)
int firstScatteringParameter (void) const
void firstScatteringParameter (int value)
bool fitEnergyDeposit (void) const
void fitEnergyDeposit (double minEnergyDeposit)
bool fitMomentum (void) const
void fitMomentum (bool value)
const Amg::MatrixXfullCovariance (void) const
TrackSurfaceIntersection intersection (void) const
int numberAlignments (void) const
void numberAlignments (int numberAlignments)
int numberOscillations (void) const
int numberParameters (void) const
void numberParameters (int numberParameters)
int numberScatterers (void) const
void numberScatterers (int numberScatterers)
void performCutStep (double cutStep)
bool phiInstability (void) const
double ptInv0 (void) const
const Amg::MatrixX parameterDifference (const Amg::VectorX &parameters) const
Perigeeperigee (void) const
const Amg::Vector3Dposition (void) const
void print (MsgStream &log) const
void printCovariance (MsgStream &log) const
void printVerbose (MsgStream &log) const
double qOverP (void) const
void qOverP (double value)
double qOverP1 (void) const
void qOverP1 (double value)
void reset (const FitParameters &parameters)
void resetOscillations (void)
double scattererPhi (int scatterer) const
double scattererTheta (int scatterer) const
ScatteringAngles scatteringAngles (const FitMeasurement &fitMeasurement, int scatterer=-1) const
void setPhiInstability (void)
double sinPhi (void) const
double sinTheta (void) const
PerigeestartingPerigee (void) const
TrackParameterstrackParameters (MsgStream &log, const FitMeasurement &measurement, bool withCovariance=false)
void update (const Amg::VectorX &differences)
void update (Amg::Vector3D position, Amg::Vector3D direction, double qOverP, const Amg::MatrixX &leadingCovariance)
const Amg::Vector3Dvertex (void) const
double z0 (void) const

Private Attributes

std::vector< double > m_alignmentAngle
std::vector< double > m_alignmentAngleConstraint
std::vector< double > m_alignmentOffset
std::vector< double > m_alignmentOffsetConstraint
double m_cosPhi
double m_cosPhi1
double m_cosTheta
double m_cosTheta1
double m_cotTheta
double m_d0
Amg::VectorX m_differences {}
bool m_eigen {}
bool m_extremeMomentum
Amg::MatrixXm_finalCovariance
int m_firstAlignmentParameter
int m_firstScatteringParameter
bool m_fitEnergyDeposit
bool m_fitMomentum
const Amg::MatrixXm_fullCovariance
double m_minEnergyDeposit
int m_numberAlignments
int m_numberOscillations
int m_numberParameters
int m_numberScatterers
double m_oldDifference
const Perigeem_perigee
bool m_phiInstability
Amg::Vector3D m_position
double m_qOverP
double m_qOverP1
std::vector< double > m_scattererPhi
std::vector< double > m_scattererTheta
double m_sinPhi
double m_sinPhi1
double m_sinTheta
double m_sinTheta1
const Surfacem_surface
Amg::Vector3D m_vertex
double m_z0

Detailed Description

Definition at line 29 of file FitParameters.h.

Constructor & Destructor Documentation

◆ FitParameters() [1/4]

Trk::FitParameters::FitParameters ( const Perigee & perigee)

Definition at line 26 of file FitParameters.cxx.

27 : m_cosPhi1(0.),
28 m_cosTheta1(0.),
29 m_d0(perigee.parameters()[Trk::d0]),
31 m_extremeMomentum(false),
32 m_finalCovariance(nullptr),
35 m_fitEnergyDeposit(false),
36 m_fitMomentum(true),
37 m_fullCovariance(nullptr),
45 m_phiInstability(false),
46 m_qOverP1(0.),
47 m_sinPhi1(0.),
48 m_sinTheta1(0.),
49 m_surface(&perigee.associatedSurface()),
50 m_vertex(perigee.associatedSurface().center()),
51 m_z0(perigee.position().z()) {
52 Amg::Vector3D momentum = perigee.momentum();
53 const double ptInv0 = 1. / momentum.perp();
54 m_cosPhi = ptInv0 * momentum.x();
55 m_sinPhi = ptInv0 * momentum.y();
57 m_sinTheta = 1. / std::sqrt(1. + m_cotTheta * m_cotTheta);
59 m_qOverP = perigee.charge() * ptInv0 * m_sinTheta;
61 m_vertex.y() + m_d0 * m_cosPhi, m_z0);
62}
const Surface * m_surface
Amg::Vector3D m_position
Amg::Vector3D m_vertex
Perigee * perigee(void) const
Amg::MatrixX * m_finalCovariance
const Perigee * m_perigee
double ptInv0(void) const
const Amg::MatrixX * m_fullCovariance
Amg::VectorX m_differences
Eigen::Matrix< double, 3, 1 > Vector3D
@ d0
Definition ParamDefs.h:63

◆ FitParameters() [2/4]

Trk::FitParameters::FitParameters ( double d0,
double z0,
double cosPhi,
double sinPhi,
double cotTheta,
double ptInv0,
const PerigeeSurface & surface )

Definition at line 64 of file FitParameters.cxx.

68 m_cosPhi1(0.),
69 m_cosTheta1(0.),
71 m_d0(d0),
72 m_extremeMomentum(false),
73 m_finalCovariance(nullptr),
76 m_fitEnergyDeposit(false),
77 m_fitMomentum(true),
78 m_fullCovariance(nullptr),
85 m_perigee(nullptr),
86 m_phiInstability(false),
87 m_qOverP1(0.),
89 m_sinPhi1(0.),
90 m_sinTheta1(0.),
91 m_surface(&surface),
92 m_vertex(surface.center()),
93 m_z0(z0) {
94 m_sinTheta = 1. / std::sqrt(1. + m_cotTheta * m_cotTheta);
97 m_z0 = z0;
99 m_vertex.y() + m_d0 * m_cosPhi, m_z0);
100}
double z0(void) const
double sinPhi(void) const
double d0(void) const
double cotTheta(void) const
double cosPhi(void) const

◆ FitParameters() [3/4]

Trk::FitParameters::FitParameters ( const FitParameters & parameters)
default

◆ FitParameters() [4/4]

Trk::FitParameters::FitParameters ( FitParameters && )
default

◆ ~FitParameters()

Trk::FitParameters::~FitParameters ( void )
default

Member Function Documentation

◆ addAlignment()

void Trk::FitParameters::addAlignment ( bool constrained,
double localAngle,
double localOffset )

Definition at line 103 of file FitParameters.cxx.

104 {
107 if (constrained) {
110 }
112}
double angle(const GeoTrf::Vector2D &a, const GeoTrf::Vector2D &b)
std::vector< double > m_alignmentOffsetConstraint
std::vector< double > m_alignmentOffset
std::vector< double > m_alignmentAngleConstraint
std::vector< double > m_alignmentAngle

◆ addScatterer()

void Trk::FitParameters::addScatterer ( double phi,
double theta )

Definition at line 114 of file FitParameters.cxx.

114 {
118}
std::vector< double > m_scattererTheta
std::vector< double > m_scattererPhi
@ theta
Definition ParamDefs.h:66
@ phi
Definition ParamDefs.h:75

◆ alignmentAngle()

double Trk::FitParameters::alignmentAngle ( int alignment) const
inline

Definition at line 156 of file FitParameters.h.

156 {
158}

◆ alignmentAngleConstraint()

double Trk::FitParameters::alignmentAngleConstraint ( int alignment) const
inline

Definition at line 160 of file FitParameters.h.

160 {
162}

◆ alignmentOffset()

double Trk::FitParameters::alignmentOffset ( int alignment) const
inline

Definition at line 164 of file FitParameters.h.

164 {
166}

◆ alignmentOffsetConstraint()

double Trk::FitParameters::alignmentOffsetConstraint ( int alignment) const
inline

Definition at line 168 of file FitParameters.h.

168 {
170}

◆ associatedSurface()

const Surface * Trk::FitParameters::associatedSurface ( void ) const

Definition at line 120 of file FitParameters.cxx.

120 {
121 if (!m_perigee)
122 return nullptr;
123 return &m_perigee->associatedSurface();
124}

◆ cosPhi()

double Trk::FitParameters::cosPhi ( void ) const
inline

Definition at line 172 of file FitParameters.h.

172 {
173 return m_cosPhi;
174}

◆ cosTheta()

double Trk::FitParameters::cosTheta ( void ) const
inline

Definition at line 176 of file FitParameters.h.

176 {
177 return m_cosTheta;
178}

◆ cotTheta()

double Trk::FitParameters::cotTheta ( void ) const
inline

Definition at line 180 of file FitParameters.h.

180 {
181 return m_cotTheta;
182}

◆ covariance()

void Trk::FitParameters::covariance ( Amg::MatrixX * finalCovariance,
const Amg::MatrixX * fullCovariance )

Definition at line 126 of file FitParameters.cxx.

127 {
130}
const Amg::MatrixX * finalCovariance(void) const
const Amg::MatrixX * fullCovariance(void) const

◆ d0() [1/2]

void Trk::FitParameters::d0 ( double value)

Definition at line 132 of file FitParameters.cxx.

132 {
133 m_d0 = value;
135 m_vertex.y() + m_d0 * m_cosPhi, m_z0);
136}

◆ d0() [2/2]

double Trk::FitParameters::d0 ( void ) const
inline

Definition at line 184 of file FitParameters.h.

184 {
185 return m_d0;
186}

◆ difference()

double Trk::FitParameters::difference ( int param) const
inline

Use as_const to avoid compiler warning

Definition at line 188 of file FitParameters.h.

188 {
189 if (m_differences.size() == 0) {
190 return 0.;
191 }
193 return std::as_const(m_differences)(param);
194}

◆ differences()

const Amg::VectorX & Trk::FitParameters::differences ( void ) const
inline

Definition at line 196 of file FitParameters.h.

196 {
197 return m_differences;
198}

◆ direction()

Amg::Vector3D Trk::FitParameters::direction ( void ) const
inline

Definition at line 200 of file FitParameters.h.

200 {
202 m_cosTheta);
203}

◆ extremeMomentum() [1/2]

void Trk::FitParameters::extremeMomentum ( bool value)

Definition at line 138 of file FitParameters.cxx.

138 {
141}

◆ extremeMomentum() [2/2]

bool Trk::FitParameters::extremeMomentum ( void ) const
inline

Definition at line 205 of file FitParameters.h.

205 {
206 return m_extremeMomentum;
207}

◆ finalCovariance()

const Amg::MatrixX * Trk::FitParameters::finalCovariance ( void ) const
inline

Definition at line 209 of file FitParameters.h.

209 {
210 return m_finalCovariance;
211}

◆ firstAlignmentParameter() [1/2]

void Trk::FitParameters::firstAlignmentParameter ( int value)

◆ firstAlignmentParameter() [2/2]

int Trk::FitParameters::firstAlignmentParameter ( void ) const
inline

Definition at line 213 of file FitParameters.h.

213 {
215}

◆ firstScatteringParameter() [1/2]

void Trk::FitParameters::firstScatteringParameter ( int value)

◆ firstScatteringParameter() [2/2]

int Trk::FitParameters::firstScatteringParameter ( void ) const
inline

Definition at line 217 of file FitParameters.h.

217 {
219}

◆ fitEnergyDeposit() [1/2]

void Trk::FitParameters::fitEnergyDeposit ( double minEnergyDeposit)

Definition at line 143 of file FitParameters.cxx.

143 {
144 m_fitEnergyDeposit = true;
145 m_minEnergyDeposit = minEnergyDeposit;
146}

◆ fitEnergyDeposit() [2/2]

bool Trk::FitParameters::fitEnergyDeposit ( void ) const
inline

Definition at line 221 of file FitParameters.h.

221 {
222 return m_fitEnergyDeposit;
223}

◆ fitMomentum() [1/2]

void Trk::FitParameters::fitMomentum ( bool value)

Definition at line 148 of file FitParameters.cxx.

148 {
150}

◆ fitMomentum() [2/2]

bool Trk::FitParameters::fitMomentum ( void ) const
inline

Definition at line 225 of file FitParameters.h.

225 {
226 return m_fitMomentum;
227}

◆ fullCovariance()

const Amg::MatrixX * Trk::FitParameters::fullCovariance ( void ) const
inline

Definition at line 229 of file FitParameters.h.

229 {
230 return m_fullCovariance;
231}

◆ intersection()

TrackSurfaceIntersection Trk::FitParameters::intersection ( void ) const

Definition at line 152 of file FitParameters.cxx.

152 {
153 return TrackSurfaceIntersection(
156 0.);
157}

◆ numberAlignments() [1/2]

void Trk::FitParameters::numberAlignments ( int numberAlignments)

Definition at line 159 of file FitParameters.cxx.

159 {
161 if (!numberAlignments)
162 return;
163 m_alignmentAngle = std::vector<double>(numberAlignments, 0.);
164 m_alignmentAngleConstraint = std::vector<double>(numberAlignments, 0.);
165 m_alignmentOffset = std::vector<double>(numberAlignments, 0.);
166 m_alignmentOffsetConstraint = std::vector<double>(numberAlignments, 0.);
167}
int numberAlignments(void) const

◆ numberAlignments() [2/2]

int Trk::FitParameters::numberAlignments ( void ) const
inline

Definition at line 233 of file FitParameters.h.

233 {
234 return m_numberAlignments;
235}

◆ numberOscillations()

int Trk::FitParameters::numberOscillations ( void ) const
inline

Definition at line 237 of file FitParameters.h.

237 {
239}

◆ numberParameters() [1/2]

void Trk::FitParameters::numberParameters ( int numberParameters)

◆ numberParameters() [2/2]

int Trk::FitParameters::numberParameters ( void ) const
inline

Definition at line 241 of file FitParameters.h.

241 {
242 return m_numberParameters;
243}

◆ numberScatterers() [1/2]

void Trk::FitParameters::numberScatterers ( int numberScatterers)

Definition at line 176 of file FitParameters.cxx.

176 {
178 if (!numberScatterers)
179 return;
180 m_scattererPhi = std::vector<double>(numberScatterers, 0.);
181 m_scattererTheta = std::vector<double>(numberScatterers, 0.);
182}
int numberScatterers(void) const

◆ numberScatterers() [2/2]

int Trk::FitParameters::numberScatterers ( void ) const
inline

Definition at line 245 of file FitParameters.h.

245 {
246 return m_numberScatterers;
247}

◆ operator=() [1/2]

FitParameters & Trk::FitParameters::operator= ( const FitParameters & )
default

◆ operator=() [2/2]

FitParameters & Trk::FitParameters::operator= ( FitParameters && )
default

◆ parameterDifference()

const Amg::MatrixX Trk::FitParameters::parameterDifference ( const Amg::VectorX & parameters) const

Definition at line 184 of file FitParameters.cxx.

185 {
187 difference(0, 0) = parameters(0) - m_d0;
188 difference(0, 1) = parameters(1) - m_z0;
190
191 // FIXME: diff is now delta(theta) : check sign
192 // difference[0][3] = parameters[4] - m_cotTheta;
193 difference(0, 3) = std::sin(parameters(4)) * m_cosTheta -
194 std::cos(parameters(4)) * m_sinTheta;
195 difference(0, 4) = (parameters(5) - m_qOverP) * Gaudi::Units::TeV;
196 return difference;
197}
double difference(int param) const
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.

◆ performCutStep()

void Trk::FitParameters::performCutStep ( double cutStep)

Definition at line 199 of file FitParameters.cxx.

199 {
200 // revert parameters to previous parameter change with cutStep*value
201 // i.e. 0 < cutstep < 1 such that cutStep = 0 gives complete reversion
202
203 Amg::VectorX cutDifferences(m_differences);
204 cutDifferences *= (cutStep - 1.);
205 Amg::VectorX oldDifferences(m_differences);
206 oldDifferences *= cutStep;
207
208 // leave phi alone when unstable
209 if (m_phiInstability) {
210 cutDifferences(2) = 0.;
211 oldDifferences(2) = (m_differences)(2);
212 }
213
214 // apply cut
215 update(cutDifferences);
216 m_differences = Amg::VectorX(oldDifferences);
217
219 m_oldDifference = 0.;
220 // std::cout << " after cutstep " << std::endl;
221}
void update(const Amg::VectorX &differences)
Eigen::Matrix< double, Eigen::Dynamic, 1 > VectorX
Dynamic Vector - dynamic allocation.

◆ perigee()

Perigee * Trk::FitParameters::perigee ( void ) const

Definition at line 223 of file FitParameters.cxx.

223 {
224 // copy 'final' covariance
225 AmgSymMatrix(5) covMatrix = AmgSymMatrix(5)(*m_finalCovariance);
226 const double pT = std::abs(m_sinTheta / m_qOverP);
227 double charge = 1.;
228 if (m_qOverP < 0.)
229 charge = -1.;
230 const Amg::Vector3D momentum(pT * m_cosPhi, pT * m_sinPhi, pT * m_cotTheta);
231
232 if (m_surface) {
233 return new Perigee(m_position, momentum, charge,
234 dynamic_cast<const Trk::PerigeeSurface&>(*m_surface),
235 std::move(covMatrix));
236 } else {
237 return new Perigee(m_position, momentum, charge, m_vertex,
238 std::move(covMatrix));
239 }
240}
double charge(const T &p)
Definition AtlasPID.h:997
#define AmgSymMatrix(dim)
if(febId1==febId2)
Eigen::Matrix< double, 3, 1 > Vector3D
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee

◆ phiInstability()

bool Trk::FitParameters::phiInstability ( void ) const

Definition at line 242 of file FitParameters.cxx.

242 {
243 return m_phiInstability;
244}

◆ position()

const Amg::Vector3D & Trk::FitParameters::position ( void ) const
inline

Definition at line 253 of file FitParameters.h.

253 {
254 return m_position;
255}

◆ print()

void Trk::FitParameters::print ( MsgStream & log) const

Definition at line 246 of file FitParameters.cxx.

246 {
247 log << std::setiosflags(std::ios::fixed | std::ios::right) << std::setw(16)
248 << std::setprecision(1) << m_position.perp() << " perigee radius"
249 << std::setw(10) << std::setprecision(3) << m_d0 << " d0" << std::setw(11)
250 << std::setprecision(3) << m_position.z() << " z0" << std::setw(11)
251 << std::setprecision(6) << std::atan2(m_sinPhi, m_cosPhi) << " phi"
252 << std::setw(11) << std::setprecision(6) << std::acos(m_cosTheta)
253 << " theta" << std::setw(13) << std::setprecision(3)
254 << m_sinTheta / (m_qOverP * Gaudi::Units::GeV) << " pT (GeV)";
255}

◆ printCovariance()

void Trk::FitParameters::printCovariance ( MsgStream & log) const

Definition at line 257 of file FitParameters.cxx.

257 {
258 double error00 = 0.;
260 if (cov(0, 0) > 0.)
261 error00 = std::sqrt(cov(0, 0));
262 double error11 = 0.;
263 if (cov(1, 1) > 0.)
264 error11 = std::sqrt(cov(1, 1));
265 double error22 = 0.;
266 if (cov(2, 2) > 0.)
267 error22 = std::sqrt(cov(2, 2));
268 double error33 = 0.;
269 if (cov(3, 3) > 0.)
270 error33 = std::sqrt(cov(3, 3));
271 double error44 = 0.;
272 if (cov(4, 4) > 0.)
273 error44 = std::sqrt(cov(4, 4));
274 double correl02 = 0.;
275 const double denom02 = cov(0, 0) * cov(2, 2);
276 if (denom02 > 0.)
277 correl02 = cov(0, 2) / std::sqrt(denom02);
278 double correl13 = 0.;
279 const double denom13 = cov(1, 1) * cov(3, 3);
280 if (denom13 > 0.)
281 correl13 = cov(1, 3) / std::sqrt(denom13);
282 log << std::setiosflags(std::ios::fixed | std::ios::right) << std::setw(10)
283 << std::setprecision(3) << error00 << std::setw(14)
284 << std::setprecision(3) << error11 << std::setw(14)
285 << std::setprecision(6) << error22 << std::setw(15)
286 << std::setprecision(6) << error33 << std::setw(19)
287 << std::setprecision(3)
288 << (error44 * m_sinTheta) / (m_qOverP * m_qOverP * Gaudi::Units::GeV)
289 << " (covariance)" << std::setw(9) << std::setprecision(3) << correl02
290 << " Cd0phi" << std::setw(8) << std::setprecision(3) << correl13
291 << " Cz0theta" << std::endl;
292}

◆ printVerbose()

void Trk::FitParameters::printVerbose ( MsgStream & log) const

Definition at line 294 of file FitParameters.cxx.

294 {
295 log << std::endl;
296
297 if (m_differences.size() != 0) {
299 log << " dParams ====" << std::setiosflags(std::ios::fixed)
300 << std::setw(10) << std::setprecision(4) << differences(0) << " (0) "
301 << std::setw(10) << std::setprecision(4) << differences(1) << " (1) "
302 << std::setw(10) << std::setprecision(5) << differences(2) << " (2) "
303 << std::setw(10) << std::setprecision(5) << differences(3) << " (3) "
304 << std::setw(13) << std::setprecision(9)
305 << differences(4) * Gaudi::Units::GeV / Gaudi::Units::TeV << " (4) ";
307 log << std::setiosflags(std::ios::fixed) << std::setw(13)
308 << std::setprecision(9)
309 << differences(5) * Gaudi::Units::GeV / Gaudi::Units::TeV << " (5) ";
310 log << std::endl;
311
312 if (m_numberAlignments) {
313 log << " dAlign ==== ";
314 unsigned param = m_firstAlignmentParameter;
315 for (int scat = 0; scat < m_numberAlignments; ++scat) {
316 ++param;
317 if (scat % 5 == 0 && scat > 0)
318 log << std::endl << " ";
319 log << std::setiosflags(std::ios::fixed) << std::setw(10)
320 << std::setprecision(6) << differences(param);
321 ++param;
322 log << std::setiosflags(std::ios::fixed) << std::setw(10)
323 << std::setprecision(6) << differences(param) << " ("
324 << std::setw(2) << scat << "A) ";
325 }
326 log << std::endl;
327 }
328
329 if (m_numberScatterers) {
330 log << " dScat ==== ";
331 unsigned param = m_firstScatteringParameter;
332 for (int scat = 0; scat < m_numberScatterers; ++scat) {
333 ++param;
334 if (scat % 5 == 0 && scat > 0)
335 log << std::endl << " ";
336 log << std::setiosflags(std::ios::fixed) << std::setw(10)
337 << std::setprecision(6) << differences(param);
338 ++param;
339 log << std::setiosflags(std::ios::fixed) << std::setw(10)
340 << std::setprecision(6) << differences(param) << " ("
341 << std::setw(2) << scat << "S) ";
342 }
343 log << std::endl;
344 }
345 }
346
347 log << std::setiosflags(std::ios::fixed | std::ios::right)
348 << " parameters: " << std::setw(12) << std::setprecision(4) << m_d0
349 << " transverse impact " << std::setw(10) << std::setprecision(4)
350 << m_position.z() << " z0 " << std::setw(10) << std::setprecision(6)
351 << std::atan2(m_sinPhi, m_cosPhi) << std::setw(10) << std::setprecision(6)
352 << m_cotTheta << " phi,cotTheta " << std::setw(13)
353 << std::setprecision(9) << m_qOverP / m_sinTheta << " inverse pT "
354 << std::setw(12) << std::setprecision(6)
355 << m_sinTheta / (m_qOverP * Gaudi::Units::GeV) << " signed pT ";
356 if (m_fitEnergyDeposit) {
357 // TODO: should give fitted energy loss
358 log << std::endl
359 << " E before/after energy deposit" << std::setw(12)
360 << std::setprecision(3) << 1. / std::abs(m_qOverP * Gaudi::Units::GeV)
361 << std::setw(12) << std::setprecision(3)
362 << 1. / std::abs(m_qOverP1 * Gaudi::Units::GeV);
363 }
364 if (m_numberAlignments) {
365 log << std::endl << " alignment number, angle, offset: ";
366 for (int align = 0; align < m_numberAlignments; ++align) {
367 log << std::setiosflags(std::ios::fixed) << std::setw(6) << align
368 << std::setw(10) << std::setprecision(6) << m_alignmentAngle[align]
369 << std::setw(10) << std::setprecision(6) << m_alignmentOffset[align];
370 if ((align + 1) % 5 == 0)
371 log << std::endl << " ";
372 }
373 }
374 if (m_numberScatterers) {
375 log << std::endl << " scatterer number, delta(phi), delta(theta): ";
376 for (int scat = 0; scat < m_numberScatterers; ++scat) {
377 log << std::setiosflags(std::ios::fixed) << std::setw(6) << scat
378 << std::setw(10) << std::setprecision(6) << m_scattererPhi[scat]
379 << std::setw(10) << std::setprecision(6) << m_scattererTheta[scat];
380 if ((scat + 1) % 5 == 0)
381 log << std::endl << " ";
382 }
383 }
384 log << std::endl;
385}
const Amg::VectorX & differences(void) const

◆ ptInv0()

double Trk::FitParameters::ptInv0 ( void ) const
inline

Definition at line 249 of file FitParameters.h.

249 {
250 return m_qOverP / m_sinTheta;
251}

◆ qOverP() [1/2]

void Trk::FitParameters::qOverP ( double value)

Definition at line 387 of file FitParameters.cxx.

387 {
388 m_qOverP = value;
389}

◆ qOverP() [2/2]

double Trk::FitParameters::qOverP ( void ) const
inline

Definition at line 257 of file FitParameters.h.

257 {
258 return m_qOverP;
259}

◆ qOverP1() [1/2]

void Trk::FitParameters::qOverP1 ( double value)

Definition at line 391 of file FitParameters.cxx.

391 {
393}

◆ qOverP1() [2/2]

double Trk::FitParameters::qOverP1 ( void ) const
inline

Definition at line 261 of file FitParameters.h.

261 {
262 return m_qOverP1;
263}

◆ reset()

void Trk::FitParameters::reset ( const FitParameters & parameters)

Definition at line 395 of file FitParameters.cxx.

395 {
396 // method is needed to complement copy in places where design uses
397 // a reference to a FitParameter pointer
398 // essentially a copy, with history of previous iteration removed
399 m_cosPhi = parameters.m_cosPhi;
400 m_cosPhi1 = parameters.m_cosPhi1;
401 m_cosTheta = parameters.m_cosTheta;
402 m_cosTheta1 = parameters.m_cosTheta1;
403 m_cotTheta = parameters.m_cotTheta;
404 m_d0 = parameters.m_d0;
405 m_extremeMomentum = parameters.m_extremeMomentum;
406 m_finalCovariance = parameters.m_finalCovariance;
407 m_firstScatteringParameter = parameters.m_firstScatteringParameter;
408 m_firstScatteringParameter = parameters.m_firstScatteringParameter;
409 m_fitEnergyDeposit = parameters.m_fitEnergyDeposit;
410 m_fitMomentum = parameters.m_fitMomentum;
411 m_fullCovariance = parameters.m_fullCovariance;
412 m_minEnergyDeposit = parameters.m_minEnergyDeposit;
413 m_numberAlignments = parameters.m_numberAlignments;
414 m_numberOscillations = parameters.m_numberOscillations;
415 m_numberParameters = parameters.m_numberParameters;
416 m_numberScatterers = parameters.m_numberScatterers;
417 m_oldDifference = parameters.m_oldDifference;
418 m_perigee = parameters.m_perigee;
419 m_phiInstability = parameters.m_phiInstability;
420 m_position = parameters.m_position;
421 m_qOverP = parameters.m_qOverP;
422 m_qOverP1 = parameters.m_qOverP1;
423 m_sinPhi = parameters.m_sinPhi;
424 m_sinPhi1 = parameters.m_sinPhi1;
425 m_sinTheta = parameters.m_sinTheta;
426 m_sinTheta1 = parameters.m_sinTheta1;
427 m_vertex = parameters.vertex();
428 m_z0 = parameters.m_z0;
429 for (int s = 0; s != m_numberAlignments; ++s) {
430 m_alignmentAngle[s] = parameters.m_alignmentAngle[s];
431 m_alignmentAngleConstraint[s] = parameters.m_alignmentAngleConstraint[s];
432 m_alignmentOffset[s] = parameters.m_alignmentOffset[s];
433 m_alignmentOffsetConstraint[s] = parameters.m_alignmentOffsetConstraint[s];
434 }
435 for (int s = 0; s != m_numberScatterers; ++s) {
436 m_scattererPhi[s] = parameters.m_scattererPhi[s];
437 m_scattererTheta[s] = parameters.m_scattererTheta[s];
438 }
439
440 // restore difference history
441 if (parameters.m_differences.size() != 0) {
442 m_differences = Amg::VectorX(parameters.m_differences);
443 } else {
445 m_differences.setZero();
446 }
447
449 m_oldDifference = 0.;
450}

◆ resetOscillations()

void Trk::FitParameters::resetOscillations ( void )
inline

Definition at line 281 of file FitParameters.h.

281 {
283}

◆ scattererPhi()

double Trk::FitParameters::scattererPhi ( int scatterer) const
inline

Definition at line 265 of file FitParameters.h.

265 {
266 return m_scattererPhi[scatterer];
267}

◆ scattererTheta()

double Trk::FitParameters::scattererTheta ( int scatterer) const
inline

Definition at line 269 of file FitParameters.h.

269 {
270 return m_scattererTheta[scatterer];
271}

◆ scatteringAngles()

ScatteringAngles Trk::FitParameters::scatteringAngles ( const FitMeasurement & fitMeasurement,
int scatterer = -1 ) const

Definition at line 452 of file FitParameters.cxx.

453 {
454 // scattering sigma used in chi2 computation
455 const double scattererSigmaTheta = 1. / fitMeasurement.weight();
456 const double scattererSigmaPhi =
457 scattererSigmaTheta /
458 fitMeasurement.intersection(FittedTrajectory).direction().perp();
459 if (scatterer < 0) {
460 return {0., 0., scattererSigmaPhi, scattererSigmaTheta};
461 } else {
462 return {m_scattererPhi[scatterer], m_scattererTheta[scatterer],
463 scattererSigmaPhi, scattererSigmaTheta};
464 }
465}
@ FittedTrajectory

◆ setPhiInstability()

void Trk::FitParameters::setPhiInstability ( void )

Definition at line 467 of file FitParameters.cxx.

467 {
468 m_phiInstability = true;
469}

◆ sinPhi()

double Trk::FitParameters::sinPhi ( void ) const
inline

Definition at line 273 of file FitParameters.h.

273 {
274 return m_sinPhi;
275}

◆ sinTheta()

double Trk::FitParameters::sinTheta ( void ) const
inline

Definition at line 277 of file FitParameters.h.

277 {
278 return m_sinTheta;
279}

◆ startingPerigee()

Perigee * Trk::FitParameters::startingPerigee ( void ) const

Definition at line 471 of file FitParameters.cxx.

471 {
472 // create momentum
473 const double pT = std::abs(m_sinTheta / m_qOverP);
474 double charge = 1.;
475 if (m_qOverP < 0.)
476 charge = -1.;
477 const Amg::Vector3D momentum(pT * m_cosPhi, pT * m_sinPhi, pT * m_cotTheta);
478
479 return new Perigee(m_position, momentum, charge, m_vertex);
480}

◆ trackParameters()

TrackParameters * Trk::FitParameters::trackParameters ( MsgStream & log,
const FitMeasurement & measurement,
bool withCovariance = false )

Definition at line 482 of file FitParameters.cxx.

483 {
484 // make checks necessary for the TrackParameters to be computed
485 // 1) a Surface is required
486 if (!measurement.surface()) {
487 log << MSG::WARNING
488 << "FitParameters::trackParameters - measurement lacks Surface"
489 << endmsg;
490 return nullptr;
491 }
492
493 // 2) a SurfaceIntersection is required
494 if (!measurement.hasIntersection(FittedTrajectory)) {
495 log << MSG::WARNING
496 << "FitParameters::trackParameters - invalid measurement" << endmsg;
497 return nullptr;
498 }
499
500 // 3) the intersection position has to lie sufficiently close to the Surface
501 const TrackSurfaceIntersection& intersection =
502 measurement.intersection(FittedTrajectory);
503 Amg::Vector2D localPos;
504 if (!measurement.surface()->globalToLocal(
505 intersection.position(), intersection.direction(), localPos)) {
506 log << MSG::WARNING
507 << "FitParameters::trackParameters - globalToLocal failure" << endmsg;
508 return nullptr;
509 }
510
511 // cache parameters at EnergyDeposit
512 if (measurement.isEnergyDeposit()) {
513 m_sinTheta1 = intersection.direction().perp();
514 m_cosPhi1 = intersection.direction().x() / m_sinTheta1;
515 m_sinPhi1 = intersection.direction().y() / m_sinTheta1;
516 m_cosTheta1 = intersection.direction().z();
517 m_qOverP1 = measurement.qOverP();
518 }
519
520 // propagate full covariance to form localCovariance
521 std::optional<AmgSymMatrix(5)> covMatrix = std::nullopt;
522 if (withCovariance && (measurement.isDrift() || measurement.isCluster() ||
523 measurement.isPerigee())) {
525 const double sigma = 1. / measurement.weight();
526 const double sigma2 = 1. / measurement.weight2();
527 int const lastParameter = measurement.lastParameter();
528 Amg::MatrixX jacobian = Amg::MatrixX::Zero(5, lastParameter);
529 for (int i = 0; i != lastParameter; ++i) {
530 jacobian(0, i) = sigma * measurement.derivative(i);
531 jacobian(1, i) = sigma2 * measurement.derivative2(i);
532 }
533
534 // phi,theta derivs into jac
535 jacobian(2, 2) = 1.;
536 jacobian(3, 3) = 1.;
538 while (++i < lastParameter) {
539 jacobian(2, i) = 1.;
540 ++i;
541 jacobian(3, i) = 1.;
542 }
543
544 // only if fit to curvature
545 if (m_fitMomentum && m_qOverP) {
546 const double sinTheta = direction.perp();
547 if (m_fitEnergyDeposit && measurement.afterCalo()) {
548 const double deltaPhi =
550 const double deltaTheta =
552 jacobian(0, 5) *= Gaudi::Units::TeV;
553 jacobian(1, 5) *= Gaudi::Units::TeV;
554 jacobian(2, 5) = deltaPhi / measurement.qOverP();
555 jacobian(3, 5) = deltaTheta / measurement.qOverP();
556 jacobian(4, 5) = measurement.qOverP() / m_qOverP1;
557 } else {
558 const double deltaPhi =
559 (direction.y() * m_cosPhi - direction.x() * m_sinPhi) / sinTheta;
560 const double deltaTheta =
562 jacobian(0, 4) *= Gaudi::Units::TeV;
563 jacobian(1, 4) *= Gaudi::Units::TeV;
564 jacobian(2, 4) = deltaPhi / measurement.qOverP();
565 jacobian(3, 4) = deltaTheta / measurement.qOverP();
566 jacobian(4, 4) = measurement.qOverP() / m_qOverP;
567 }
568 } else {
569 jacobian(4, 4) = 1.;
570 }
571
572 // similarity transform
574 jacobian * m_fullCovariance->block(0, 0, lastParameter, lastParameter) *
575 jacobian.transpose());
576 }
577
578 double phi = intersection.direction().phi();
579 double theta = intersection.direction().theta();
580 if (measurement.isFlipped()) {
581 if (phi > 0.) {
582 phi -= M_PI;
583 } else {
584 phi += M_PI;
585 }
586 theta = M_PI - theta;
587 }
588
589 if (!measurement.surface()) {
590 log << MSG::WARNING
591 << "FitParameters::trackParameters - unrecognized surface" << endmsg;
592 return nullptr;
593 }
594 // finally can create the appropriate TrackParameters based on the surface
595 return measurement.surface()
596 ->createUniqueTrackParameters(localPos[locR], localPos[locZ], phi, theta,
597 measurement.qOverP(), std::move(covMatrix))
598 .release();
599}
#define M_PI
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
#define endmsg
Amg::Vector3D direction(void) const
TrackSurfaceIntersection intersection(void) const
double sinTheta(void) const
Eigen::Matrix< double, 2, 1 > Vector2D
@ locR
Definition ParamDefs.h:44
@ locZ
local cylindrical
Definition ParamDefs.h:42

◆ update() [1/2]

void Trk::FitParameters::update ( Amg::Vector3D position,
Amg::Vector3D direction,
double qOverP,
const Amg::MatrixX & leadingCovariance )

Definition at line 670 of file FitParameters.cxx.

672 {
673 // update parameters after leading material corrections
675 m_sinTheta = direction.perp();
676 const double sinThetaInv = 1. / m_sinTheta;
677 m_cotTheta = sinThetaInv * m_cosTheta;
678 m_cosPhi = sinThetaInv * direction.x();
679 m_sinPhi = sinThetaInv * direction.y();
680 m_cotTheta = sinThetaInv * direction.z();
681 m_cosTheta = direction.z();
683 m_z0 = position.z();
684 m_d0 = (m_vertex.x() - position.x()) * m_sinPhi -
685 (m_vertex.y() - position.y()) * m_cosPhi;
686
687 // update covariance
688 // (*m_finalCovariance) += leadingCovariance; // ??
689 for (int i = 0; i != 5; ++i) {
690 for (int j = 0; j != 5; ++j) {
691 (*m_finalCovariance)(i, j) += leadingCovariance(i, j);
692 }
693 }
694}
double qOverP(void) const
const Amg::Vector3D & position(void) const

◆ update() [2/2]

void Trk::FitParameters::update ( const Amg::VectorX & differences)

Definition at line 601 of file FitParameters.cxx.

601 {
602 // keep update values in case of cutStep procedure
605 } else {
607 }
610
611 // misalignment parameters
612 std::vector<double>::iterator a = m_alignmentAngle.begin();
613 std::vector<double>::iterator o = m_alignmentOffset.begin();
614 int align = m_firstAlignmentParameter;
615 for (int i = 0; i != m_numberAlignments; ++i) {
616 (*a++) += differences(++align);
617 (*o++) += differences(++align);
618 }
619
620 // scattering angles
621 std::vector<double>::iterator p = m_scattererPhi.begin();
622 std::vector<double>::iterator t = m_scattererTheta.begin();
624 for (int i = 0; i != m_numberScatterers; ++i) {
625 (*p++) += differences(++scat);
626 (*t++) += differences(++scat);
627 }
628
629 // qOverP, cotTheta
630 if (m_fitMomentum)
631 m_qOverP += differences(4) / Gaudi::Units::TeV;
633
634 // impose charge conservation and decreasing energy
635 if (m_fitEnergyDeposit) {
636 m_qOverP1 += differences(5) / Gaudi::Units::TeV;
637 const double deposit = 1. / std::abs(m_qOverP) - 1. / std::abs(m_qOverP1);
638 if (std::abs(deposit) < std::abs(m_minEnergyDeposit) ||
639 deposit * m_minEnergyDeposit < 0. || m_qOverP * m_qOverP1 < 0.) {
640 m_qOverP = 1. / (1. / std::abs(m_qOverP1) + m_minEnergyDeposit);
641 if (m_qOverP1 < 0.)
643 }
644 }
645
646 // protect phi against some rounding instabilities
647 double sinDPhi = differences(2);
648 double cosDPhi = 0.;
649 if (std::abs(sinDPhi) < 1.0) {
650 cosDPhi = std::sqrt(1. - sinDPhi * sinDPhi);
651 } else {
652 if (sinDPhi > 0.) {
653 sinDPhi = 1.0;
654 } else {
655 sinDPhi = -1.0;
656 }
657 }
658
659 const double cosPhi = m_cosPhi * cosDPhi - m_sinPhi * sinDPhi;
660 m_sinPhi = m_sinPhi * cosDPhi + m_cosPhi * sinDPhi;
662 m_z0 += differences(1);
663 m_d0 += differences(0);
664 m_sinTheta = 1. / std::sqrt(1. + m_cotTheta * m_cotTheta);
667 m_vertex.y() + m_d0 * m_cosPhi, m_z0);
668}
static Double_t a

◆ vertex()

const Amg::Vector3D & Trk::FitParameters::vertex ( void ) const
inline

Definition at line 285 of file FitParameters.h.

285 {
286 return m_vertex;
287}

◆ z0()

double Trk::FitParameters::z0 ( void ) const
inline

Definition at line 289 of file FitParameters.h.

289 {
290 return m_position.z();
291}

Member Data Documentation

◆ m_alignmentAngle

std::vector<double> Trk::FitParameters::m_alignmentAngle
private

Definition at line 113 of file FitParameters.h.

◆ m_alignmentAngleConstraint

std::vector<double> Trk::FitParameters::m_alignmentAngleConstraint
private

Definition at line 114 of file FitParameters.h.

◆ m_alignmentOffset

std::vector<double> Trk::FitParameters::m_alignmentOffset
private

Definition at line 115 of file FitParameters.h.

◆ m_alignmentOffsetConstraint

std::vector<double> Trk::FitParameters::m_alignmentOffsetConstraint
private

Definition at line 116 of file FitParameters.h.

◆ m_cosPhi

double Trk::FitParameters::m_cosPhi
private

Definition at line 117 of file FitParameters.h.

◆ m_cosPhi1

double Trk::FitParameters::m_cosPhi1
private

Definition at line 118 of file FitParameters.h.

◆ m_cosTheta

double Trk::FitParameters::m_cosTheta
private

Definition at line 119 of file FitParameters.h.

◆ m_cosTheta1

double Trk::FitParameters::m_cosTheta1
private

Definition at line 120 of file FitParameters.h.

◆ m_cotTheta

double Trk::FitParameters::m_cotTheta
private

Definition at line 121 of file FitParameters.h.

◆ m_d0

double Trk::FitParameters::m_d0
private

Definition at line 122 of file FitParameters.h.

◆ m_differences

Amg::VectorX Trk::FitParameters::m_differences {}
private

Definition at line 123 of file FitParameters.h.

123{};

◆ m_eigen

bool Trk::FitParameters::m_eigen {}
private

Definition at line 124 of file FitParameters.h.

124{};

◆ m_extremeMomentum

bool Trk::FitParameters::m_extremeMomentum
private

Definition at line 125 of file FitParameters.h.

◆ m_finalCovariance

Amg::MatrixX* Trk::FitParameters::m_finalCovariance
private

Definition at line 126 of file FitParameters.h.

◆ m_firstAlignmentParameter

int Trk::FitParameters::m_firstAlignmentParameter
private

Definition at line 127 of file FitParameters.h.

◆ m_firstScatteringParameter

int Trk::FitParameters::m_firstScatteringParameter
private

Definition at line 128 of file FitParameters.h.

◆ m_fitEnergyDeposit

bool Trk::FitParameters::m_fitEnergyDeposit
private

Definition at line 129 of file FitParameters.h.

◆ m_fitMomentum

bool Trk::FitParameters::m_fitMomentum
private

Definition at line 130 of file FitParameters.h.

◆ m_fullCovariance

const Amg::MatrixX* Trk::FitParameters::m_fullCovariance
private

Definition at line 131 of file FitParameters.h.

◆ m_minEnergyDeposit

double Trk::FitParameters::m_minEnergyDeposit
private

Definition at line 132 of file FitParameters.h.

◆ m_numberAlignments

int Trk::FitParameters::m_numberAlignments
private

Definition at line 133 of file FitParameters.h.

◆ m_numberOscillations

int Trk::FitParameters::m_numberOscillations
private

Definition at line 134 of file FitParameters.h.

◆ m_numberParameters

int Trk::FitParameters::m_numberParameters
private

Definition at line 135 of file FitParameters.h.

◆ m_numberScatterers

int Trk::FitParameters::m_numberScatterers
private

Definition at line 136 of file FitParameters.h.

◆ m_oldDifference

double Trk::FitParameters::m_oldDifference
private

Definition at line 137 of file FitParameters.h.

◆ m_perigee

const Perigee* Trk::FitParameters::m_perigee
private

Definition at line 138 of file FitParameters.h.

◆ m_phiInstability

bool Trk::FitParameters::m_phiInstability
private

Definition at line 139 of file FitParameters.h.

◆ m_position

Amg::Vector3D Trk::FitParameters::m_position
private

Definition at line 140 of file FitParameters.h.

◆ m_qOverP

double Trk::FitParameters::m_qOverP
private

Definition at line 141 of file FitParameters.h.

◆ m_qOverP1

double Trk::FitParameters::m_qOverP1
private

Definition at line 142 of file FitParameters.h.

◆ m_scattererPhi

std::vector<double> Trk::FitParameters::m_scattererPhi
private

Definition at line 143 of file FitParameters.h.

◆ m_scattererTheta

std::vector<double> Trk::FitParameters::m_scattererTheta
private

Definition at line 144 of file FitParameters.h.

◆ m_sinPhi

double Trk::FitParameters::m_sinPhi
private

Definition at line 145 of file FitParameters.h.

◆ m_sinPhi1

double Trk::FitParameters::m_sinPhi1
private

Definition at line 146 of file FitParameters.h.

◆ m_sinTheta

double Trk::FitParameters::m_sinTheta
private

Definition at line 147 of file FitParameters.h.

◆ m_sinTheta1

double Trk::FitParameters::m_sinTheta1
private

Definition at line 148 of file FitParameters.h.

◆ m_surface

const Surface* Trk::FitParameters::m_surface
private

Definition at line 149 of file FitParameters.h.

◆ m_vertex

Amg::Vector3D Trk::FitParameters::m_vertex
private

Definition at line 150 of file FitParameters.h.

◆ m_z0

double Trk::FitParameters::m_z0
private

Definition at line 151 of file FitParameters.h.


The documentation for this class was generated from the following files: