ATLAS Offline Software
VKalExtPropagator.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 //
5 // The VKalExtPropagator object is created if ATHENA propagator exists
6 // and is supplied via jobOptions. A pointer to it is supplied to VKalVrtCore
7 // for every vertex fit via VKalVrtControlBase object.
8 // VKalVrtCore uses it to extrapolate tracks for a vertex fit.
9 //
10 // If ATHENA propagator doesn't exist, the VKalVrtCore uses its internal
11 // propagator without material.
12 //
13 //-------------------------------------------------------------------------
14 
15 
16 // Header include
18 
23 //-------------------------------------------------
24 #include<iostream>
25 
26 
27 //&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
28 // External propagator access for VKalVrt
29 
30 namespace Trk {
31 
32  // Constructor
34  {
35  m_extrapolator = nullptr;
36  m_vkalFitSvc = pnt;
37  }
39 
40 
42  {
43  m_extrapolator = Pnt;
44  }
45 
46 //Protection against exit outside ID volume
47 //
48  double VKalExtPropagator::Protection(const double *RefEnd,
49  const IVKalState& istate) const
50  {
51  const TrkVKalVrtFitter::State& state = static_cast<const TrkVKalVrtFitter::State&> (istate);
52 
53  double Xend=RefEnd[0] + state.m_refFrameX;
54  double Yend=RefEnd[1] + state.m_refFrameY;
55  double Zend=RefEnd[2] + state.m_refFrameZ;
56  double Rlim=sqrt(Xend*Xend+Yend*Yend) / m_vkalFitSvc->m_IDsizeR;
57  double Zlim=fabs(Zend) / m_vkalFitSvc->m_IDsizeZ;
58  double Scale = Rlim; if(Zlim>Rlim) Scale=Zlim;
59 //std::cout<<"relative TARG="<<RefEnd[0]<<","<<RefEnd[1]<<","<<RefEnd[2]
60 //<<" global ref.="<<m_vkalFitSvc->state.m_refFrameX<<","<<m_vkalFitSvc->state.m_refFrameY<<","<<m_vkalFitSvc->state.m_refFrameZ
61 //<<" Limits="<<m_vkalFitSvc->m_IDsizeR<<","<<m_vkalFitSvc->m_IDsizeZ<<" scale="<<Scale<<'\n';
62  return Scale;
63  }
64 
65  bool VKalExtPropagator::checkTarget( double *RefEnd,
66  const IVKalState& istate) const
67  {
68  //double targV[3]={ RefEnd[0], RefEnd[1], RefEnd[2]};
69  return Protection(RefEnd, istate) <= 1.;
70  }
71 /*----------------------------------------------------------------------------------*/
72 // Addres of this fuction is supplied to TrkVKalVrtCore for ATLAS InDet extrapolation
73 // with material. Each time VKalVrt needs propagation - it calls it.
74 //
75 // In case of usage with first measured point
76 // it always use this point as starting point,
77 // so ParOld,CovOld,RefStart are irrelevant
78 //
79 // VKalVrtCore works in relative coordinates wrt (state.m_refFrameX,state.m_refFrameY,state.m_refFrameZ)
80 // For ATLAS propagator the Core coordinates must be moved back to global ref.frame
81 /*------------------------------------------------------------------------------------*/
82  void VKalExtPropagator::Propagate( long int trkID, long int Charge,
83  double *ParOld, double *CovOld, double *RefStart,
84  double *RefEnd, double *ParNew, double *CovNew,
85  IVKalState& istate) const
86  {
87  TrkVKalVrtFitter::State& state = static_cast<TrkVKalVrtFitter::State&> (istate);
88 
89  int trkID_loc=trkID; if(trkID_loc<0)trkID_loc=0;
90 //std::cout<<__func__<<" Ext.Propagator TrkID="<<trkID<<"to (local!!!)="<<RefEnd[0]<<", "<<RefEnd[1]<<", "<<RefEnd[2]<<'\n';
91 //-----------
92  double vX=RefEnd[0]; double vY=RefEnd[1]; double vZ=RefEnd[2]; //relative coords
93  // Propagation target in GLOBAL frame
94  Amg::Vector3D endPointG( vX + state.m_refFrameX, vY + state.m_refFrameY, vZ + state.m_refFrameZ);
95 //
96 // ---- Make MeasuredPerigee from input. Mag.field at start point is used here
97 //
98  std::vector<double> PerigeeIni( ParOld, ParOld+5 );
99  std::vector<double> CovPerigeeIni( 15, 0. );
100  if( CovOld != nullptr) {
101 // for(int i=0; i<15;i++) CovPerigeeIni.push_back( CovOld[i] );
102  std::copy(CovOld,CovOld+15,CovPerigeeIni.begin() );
103  }else{
104 // for(int i=0; i<15;i++) CovPerigeeIni.push_back(0.);
105  CovPerigeeIni[0]=1.e6;CovPerigeeIni[2]=1.e6;CovPerigeeIni[5]=1.;CovPerigeeIni[9]=1.;CovPerigeeIni[14]=fabs(PerigeeIni[4]);
106  }
107  //--- This creates Perigee in GLOBAL frame from input in realtive coordinates
108  const Perigee* inpPer =
109  m_vkalFitSvc->CreatePerigee( RefStart[0], RefStart[1], RefStart[2], PerigeeIni, CovPerigeeIni, state).release();
110  const TrackParameters * inpPar= (const TrackParameters*) inpPer;
111 //
112 // ----- Magnetic field is taken at target point (GLOBAL calculated from relative frame input)
113 //
114  double fx,fy,BMAG_FIXED;
115  state.m_fitField.getMagFld(vX,vY,vZ,fx,fy,BMAG_FIXED);
116 //
117 //-------------------- Extrapolation itself
118 //
119  const Trk::TrackParameters* endPer = nullptr;
120  if(trkID<0){
121  endPer = myExtrapWithMatUpdate( trkID, inpPar, &endPointG, state);
122  }else{
123  endPer = myExtrapWithMatUpdate( trkID, inpPar, &endPointG, state);
124  }
125 //-----------------------------------
126  if( endPer == nullptr ) { // No extrapolation done!!!
127  ParNew[0]=0.; ParNew[1]=0.;ParNew[2]=0.;ParNew[3]=0.;ParNew[4]=0.;
128  delete inpPer; return;
129  }
130  const Perigee* mPer = dynamic_cast<const Perigee*>(endPer);
131  const AtaStraightLine* Line = dynamic_cast<const AtaStraightLine*>(endPer);
132  AmgVector(5) VectPerig; VectPerig.setZero();
133  const AmgSymMatrix(5) *CovMtx=nullptr;
134  if( mPer ){
135  VectPerig = mPer->parameters();
136  CovMtx = mPer->covariance();
137  }
138  if( Line ){
139  VectPerig = Line->parameters();
140  CovMtx = Line->covariance();
141  }
142  if( (Line==nullptr && mPer==nullptr) || CovMtx==nullptr ){
143  ParNew[0]=0.; ParNew[1]=0.;ParNew[2]=0.;ParNew[3]=0.;ParNew[4]=0.;
144  delete inpPer; return;
145  }
146 
147  if((*CovMtx)(0,0)<=0. || (*CovMtx)(1,1)<=0.){ //protection against bad error matrix
148  ParNew[0]=0.; ParNew[1]=0.;ParNew[2]=0.;ParNew[3]=0.;ParNew[4]=0.;
149  delete inpPer; delete endPer;
150  return;
151  }
152  double CovVertTrk[15];
153  long int locCharge=Charge;
154  CovVertTrk[ 0] =(*CovMtx)(0,0);
155  CovVertTrk[ 1] =(*CovMtx)(1,0);
156  CovVertTrk[ 2] =(*CovMtx)(1,1);
157  CovVertTrk[ 3] =(*CovMtx)(2,0);
158  CovVertTrk[ 4] =(*CovMtx)(2,1);
159  CovVertTrk[ 5] =(*CovMtx)(2,2);
160  CovVertTrk[ 6] =(*CovMtx)(3,0);
161  CovVertTrk[ 7] =(*CovMtx)(3,1);
162  CovVertTrk[ 8] =(*CovMtx)(3,2);
163  CovVertTrk[ 9] =(*CovMtx)(3,3);
164  CovVertTrk[10] =(*CovMtx)(4,0);
165  CovVertTrk[11] =(*CovMtx)(4,1);
166  CovVertTrk[12] =(*CovMtx)(4,2);
167  CovVertTrk[13] =(*CovMtx)(4,3);
168  CovVertTrk[14] =(*CovMtx)(4,4);
169 //std::cout<<" extrapPoint="<<endPer->position().x()<<", "<<endPer->position().y()<<", "<<endPer->position().y()<<'\n';
170 //std::cout<<" extrapCov="<<(*CovMtx)(0,0)<<", "<<(*CovMtx)(1,1)<<", "<<(*CovMtx)(2,2)<<
171 // ", "<<(*CovMtx)(3,3)<<", "<<(*CovMtx)(4,4)<<'\n';
172 
173  if(CovNew != nullptr) {
174  m_vkalFitSvc->VKalTransform( BMAG_FIXED, VectPerig(0), VectPerig(1),
175  VectPerig(2), VectPerig(3), VectPerig(4), CovVertTrk,
176  locCharge, &ParNew[0] , &CovNew[0]);
177  }else{
178  double CovVertTrkTmp[15];
179  m_vkalFitSvc->VKalTransform( BMAG_FIXED, VectPerig(0), VectPerig(1),
180  VectPerig(2), VectPerig(3), VectPerig(4), CovVertTrk,
181  locCharge, &ParNew[0] , CovVertTrkTmp);
182  }
183  delete inpPer; delete endPer;
184  }
185 
186 
187 /*--------------------------------------------------------------------------------------*/
188 /* Logic of material subtruction or addition to track error matrix is implemented here. */
189 /* What to do depends on direction of extrapolation and relative position of initial, */
190 /* final and track reference points */
191 /* */
192 /* All points are in GLOBAL Atlas frame here! */
193 /* */
194 /* Real track (TrkID>=0) extrapolation always starts from INITIAL track object */
195 /* and NOT from parameters provided by VKalVrtCore. More CPU time but better precision*/
196 /* Only combined (TrkID<0) tracks use VKalVrtCore parameters for start */
197 /*--------------------------------------------------------------------------------------*/
199  const TrackParameters *inpPer,
200  Amg::Vector3D * endPoint,
201  const IVKalState& istate) const
202  {
203  const TrkVKalVrtFitter::State& state = static_cast<const TrkVKalVrtFitter::State&> (istate);
204  const EventContext& ctx = (state.m_eventContext)
205  ? *(state.m_eventContext)
206  : Gaudi::Hive::currentContext();
207 
208  const Trk::TrackParameters* endPer=nullptr;
209  const Trk::TrackParameters* tmpPer=nullptr;
210  ParticleHypothesis prtType = pion;
211  //End surface
212  PerigeeSurface surfEnd( *endPoint );
213  //Initial point (global frame)
214  Amg::Vector3D iniPoint = inpPer->position();
215  Amg::Vector3D pmom=inpPer->momentum();
216  //Track reference point ( some point on track provided initially, global frame )
217  Amg::Vector3D refPoint(0.,0.,0.);
218  if(TrkID>=0)refPoint = state.m_trkControl.at(TrkID).trkRefGlobPos;
219  //
220  Amg::Vector3D step = (*endPoint) - iniPoint;
221  //
222  int Strategy = 0; if(TrkID>=0) Strategy = state.m_trkControl[TrkID].extrapolationType;
223  //
224  // Extrapolation for new track - no material at all
225  //
226  const TrackParameters *pntOnTrk=nullptr;
227  if (TrkID < 0) {
228  prtType = undefined;
229  if (pmom.dot(step) > 0.) {
230  endPer = m_extrapolator->extrapolateDirectly(ctx,
231  *inpPer, surfEnd, alongMomentum, true, pion).release();
232  } else {
233  endPer = m_extrapolator->extrapolateDirectly(ctx,
234  *inpPer, surfEnd, oppositeMomentum, true, pion).release();
235  }
236  return endPer;
237  }
238  pntOnTrk = dynamic_cast<const TrackParameters*>(
239  state.m_trkControl.at(TrkID).TrkPnt);
240  if (pntOnTrk == nullptr){
241  return endPer;
242  }
243  prtType = pion; // Pion hypothesis is always used for extrapolation
244  // Redefinition of starting point for extrapolation
245  iniPoint = pntOnTrk->position();
246  step = (*endPoint) - iniPoint;
247 
248  //
249  // Extrapolation for first measured point strategy. Start from it and
250  // always add material
251  //
252  if( Strategy == 0) {
254  if (pmom.dot(step) < 0) {
256  }
257  endPer = m_extrapolator->extrapolate(
258  ctx, *pntOnTrk, surfEnd, dir, true, prtType, addNoise).release();
259  return endPer;
260  }
261  //
262  // Extrapolation for any measured point
263  //
264  if (Strategy == 1) {
267  if (pmom.dot(step) < 0) {
269  mmode = removeNoise;
270  }
271  endPer = m_extrapolator->extrapolate(
272  ctx, *pntOnTrk, surfEnd, dir, true, prtType, mmode).release();
273  return endPer;
274  }
275  //
276  // Extrapolation for perigee and B-hit
277  //
278  if (Strategy == 2) {
279  double Border = 25.;
280  bool dirPositive = true;
281  if (pmom.dot(step) < 0.)
282  dirPositive = false;
283  if ((*endPoint).perp() > Border && iniPoint.perp() > Border) {
284  if (dirPositive) {
285  endPer = m_extrapolator->extrapolate(
286  ctx, *pntOnTrk, surfEnd, alongMomentum, true, prtType, addNoise).release();
287  } else {
288  endPer = m_extrapolator->extrapolate(ctx,
289  *pntOnTrk,
290  surfEnd,
292  true,
293  prtType,
294  removeNoise).release();
295  }
296  return endPer;
297  }
298  if ((*endPoint).perp() < Border && iniPoint.perp() < Border) {
299  if (dirPositive) {
300  endPer = m_extrapolator->extrapolate(ctx,
301  *pntOnTrk,
302  surfEnd,
304  true,
305  prtType,
306  removeNoise).release();
307  } else {
308  endPer = m_extrapolator->extrapolate(ctx,
309  *pntOnTrk,
310  surfEnd,
312  true,
313  prtType,
314  addNoise).release();
315  }
316  return endPer;
317  }
318  Amg::Transform3D trnsf;
319  trnsf.setIdentity();
320  CylinderSurface surfBorder(trnsf, Border, 3000.);
321  if (iniPoint.perp() < Border) {
322  tmpPer = m_extrapolator->extrapolate(ctx,
323  *pntOnTrk,
324  surfBorder,
326  true,
327  prtType,
328  removeNoise).release();
329  if (tmpPer == nullptr) {
330  return nullptr;
331  }
332  endPer = m_extrapolator->extrapolate(
333  ctx, *tmpPer, surfEnd, alongMomentum, true, prtType, addNoise).release();
334  } else {
335  endPer = m_extrapolator->extrapolate(
336  ctx, *pntOnTrk, surfEnd, oppositeMomentum, true, prtType, addNoise).release();
337  return endPer;
338  }
339  delete tmpPer;
340  return endPer;
341  }
342  return endPer;
343  }
344 
345 /*--------------------------------------------------------------------------------------*/
346 /* Extrapolation to line */
347 /* */
348 /* All points are in GLOBAL Atlas frame here! */
349 /* DOESN't WORK FOR COMBINED TRACKS! */
350 /*--------------------------------------------------------------------------------------*/
352  const TrackParameters *inpPer,
353  Amg::Vector3D * endPoint,
354  StraightLineSurface &lineTarget,
355  const IVKalState& istate) const
356  {
357  const TrkVKalVrtFitter::State& state = static_cast<const TrkVKalVrtFitter::State&> (istate);
358  const EventContext& ctx = (state.m_eventContext)
359  ? *(state.m_eventContext)
360  : Gaudi::Hive::currentContext();
361 
362 
363  const Trk::TrackParameters* endPer=nullptr;
364  ParticleHypothesis prtType = muon;
365 //Initial point
366  Amg::Vector3D iniPoint = inpPer->position();
367  Amg::Vector3D step = (*endPoint) - iniPoint;
368 //
369  int Strategy = 0; if(TrkID>=0) Strategy = state.m_trkControl[TrkID].extrapolationType;
370 //
371 // Extrapolation for new track - no material at all
372 //
373  const TrackParameters *pntOnTrk=nullptr;
374  if(TrkID<0){
375  return endPer;
376  }
377  pntOnTrk=dynamic_cast<const TrackParameters*> (state.m_trkControl.at(TrkID).TrkPnt);
378  if(!pntOnTrk) return endPer;
379  //double inpMass=state.m_trkControl[TrkID].prtMass;
380  //if( inpMass > 0. && inpMass< 20.) { prtType=electron; } //VK Disabled according to users request
381  //else if(inpMass > 20. && inpMass<120.) { prtType=muon; } // May be activated in future
382  //else if(inpMass >120. && inpMass<200.) { prtType=pion; }
383  //else if(inpMass >200. && inpMass<700.) { prtType=kaon; }
384  //else if(inpMass >700. && inpMass<999.) { prtType=proton; }
385  //else { prtType=undefined; }
386  prtType=muon; // Muon hypothesis is always used for extrapolation
387  iniPoint = pntOnTrk->position();
388  step = (*endPoint) - iniPoint;
389 
390  Amg::Vector3D pmom=pntOnTrk->momentum();
391 //
392 // Extrapolation for first measured point strategy. Start from it and always add material
393 //
394  if( Strategy == 0) {
396  endPer = m_extrapolator->extrapolate(ctx, *pntOnTrk, lineTarget, dir, true, prtType, addNoise).release();
397  if (!endPer)
398  endPer = m_extrapolator->extrapolateDirectly(ctx, *pntOnTrk, lineTarget, dir, true, prtType).release();
399  return endPer;
400  }
401 //
402 // Extrapolation for any measured point
403 //
404  if( Strategy == 1 || Strategy == 2) {
406  if(pmom.dot(step)<0){ dir=oppositeMomentum; mmode=removeNoise;}
407  endPer = m_extrapolator->extrapolate(ctx, *pntOnTrk, lineTarget, dir, true, prtType, mmode).release();
408  return endPer;
409  }
410  return endPer;
411  }
412 
413 /*--------------------------------------------------------------------------------------*/
414 /* Simple extrapolation of neutral tracks */
415 /*--------------------------------------------------------------------------------------*/
417  Amg::Vector3D * endPoint) const
418  {
419  const Trk::NeutralParameters* endPer=nullptr;
420 //End surface
421  PerigeeSurface surfEnd( *endPoint );
422  endPer = m_extrapolator->extrapolate( *inpPer, surfEnd, anyDirection, true).release();
423  return endPer;
424  }
425 
426 /*--------------------------------------------------------------------------------------*/
427 /* xAOD method to find FirstMeasuredPoint on TrackParticle. */
428 /* */
429 /* All points are in GLOBAL Atlas frame here! */
430 /* */
431 /* xAOD::TrackParticle Perigee is extrapolated to cylinder with radius of first hit */
432 /*--------------------------------------------------------------------------------------*/
434  {
435  static const SG::ConstAccessor<float> radiusOfFirstHitAcc ("radiusOfFirstHit");
436  if(!radiusOfFirstHitAcc.isAvailable (*xprt)) return nullptr; // No radiusOfFirstHit on track
437 
438  const EventContext& ctx = Gaudi::Hive::currentContext();
439  const Trk::Perigee* mPer = &(xprt->perigeeParameters());
440  Amg::Transform3D trnsf;
441  trnsf.setIdentity();
442  CylinderSurface surfacePntOnTrk( trnsf, xprt->radiusOfFirstHit(), 20000.);
443  ParticleHypothesis prtType = pion;
444 
445  const TrackParameters *hitOnTrk =
447  *mPer,
448  surfacePntOnTrk,
450  true, prtType, removeNoise).release();
451 //std::cout<<" Radius="<<xprt->radiusOfFirstHit()<<" extrap="<<hitOnTrk<<'\n';
452  if(hitOnTrk==nullptr)hitOnTrk=m_extrapolator->extrapolateDirectly(ctx,
453  *mPer,
454  surfacePntOnTrk,
456  true, prtType).release();
457  if(hitOnTrk==nullptr)return nullptr;
458 
459  //convert result to Perigee
460  PerigeeSurface surfacePerigee( hitOnTrk->position() );
461  const TrackParameters *hitOnTrkPerig = m_extrapolator->extrapolate(ctx,
462  *hitOnTrk,
463  surfacePerigee).release();
464  delete hitOnTrk; // Delete temporary results
465  if(hitOnTrkPerig==nullptr)return nullptr;
466 //std::cout<<" perig="<<(*hitOnTrkPerig)<<'\n';
467  return dynamic_cast<const Perigee* > (hitOnTrkPerig);
468  }
469 
470 //--------------------------------------------------------------------------
471 //&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&&
472 // Setting interface to setup VKalVrt InDet extrapolator
474  {
475  // Save external propagator in VKalExtPropagator object and send it to TrkVKalVrtCore
476 //
477  if(m_fitPropagator != nullptr) delete m_fitPropagator;
478  m_fitPropagator = new VKalExtPropagator( this );
480  m_InDetExtrapolator = Pnt; // Pointer to InDet extrapolator
481  }
482 
483 }
Trk::TrkVKalVrtFitter::State
Definition: TrkVKalVrtFitter.h:382
Trk::anyDirection
@ anyDirection
Definition: PropDirection.h:22
Trk::TrkVKalVrtFitter::State::m_refFrameZ
double m_refFrameZ
Definition: TrkVKalVrtFitter.h:390
StraightLineSurface.h
Trk::TrkVKalVrtFitter::State::m_refFrameX
double m_refFrameX
Definition: TrkVKalVrtFitter.h:388
Trk::TrkVKalVrtFitter::m_InDetExtrapolator
const IExtrapolator * m_InDetExtrapolator
Pointer to Extrapolator AlgTool.
Definition: TrkVKalVrtFitter.h:470
Trk::TrkVKalVrtFitter::CreatePerigee
virtual std::unique_ptr< Trk::Perigee > CreatePerigee(const std::vector< double > &VKPerigee, const std::vector< double > &VKCov, IVKalState &istate) const override final
Definition: CvtPerigee.cxx:159
Trk::TrkVKalVrtFitter::VKalTransform
void VKalTransform(double MAG, double A0V, double ZV, double PhiV, double ThetaV, double PInv, const double[15], long int &Charge, double[5], double[15]) const
Definition: VKalTransform.cxx:59
Trk::PerigeeSurface
Definition: PerigeeSurface.h:43
Trk::ParametersBase::position
const Amg::Vector3D & position() const
Access method for the position.
Trk::oppositeMomentum
@ oppositeMomentum
Definition: PropDirection.h:21
Trk::ParametersT
Dummy class used to allow special convertors to be called for surfaces owned by a detector element.
Definition: EMErrorDetail.h:25
Trig::MatchingStrategy::Strategy
Strategy
Definition: MatchingImplementation.h:26
Trk::VKalExtPropagator::VKalExtPropagator
VKalExtPropagator(TrkVKalVrtFitter *)
Definition: VKalExtPropagator.cxx:33
Trk::TrkVKalVrtFitter::m_fitPropagator
VKalExtPropagator * m_fitPropagator
Definition: TrkVKalVrtFitter.h:469
Trk::VKalAtlasMagFld::getMagFld
virtual void getMagFld(const double, const double, const double, double &, double &, double &) override
Definition: VKalAtlasMagFld.cxx:61
Trk::TrkVKalVrtFitter::VKalExtPropagator
friend class VKalExtPropagator
Definition: TrkVKalVrtFitter.h:68
Trk::MaterialUpdateMode
MaterialUpdateMode
This is a steering enum to force the material update it can be: (1) addNoise (-1) removeNoise Second ...
Definition: MaterialUpdateMode.h:18
Trk::alongMomentum
@ alongMomentum
Definition: PropDirection.h:20
Trk::VKalExtPropagator::myExtrapWithMatUpdate
const TrackParameters * myExtrapWithMatUpdate(long int TrkID, const TrackParameters *inpPer, Amg::Vector3D *endPoint, const IVKalState &istate) const
Definition: VKalExtPropagator.cxx:198
SG::ConstAccessor< float >
Trk::undefined
@ undefined
Definition: ParticleHypothesis.h:38
TrkVKalVrtFitter.h
Trk::IExtrapolator::extrapolateDirectly
virtual std::unique_ptr< TrackParameters > extrapolateDirectly(const EventContext &ctx, const TrackParameters &parm, const Surface &sf, PropDirection dir=anyDirection, const BoundaryCheck &bcheck=true, ParticleHypothesis particle=pion) const =0
Extrapolate directly: Forwards directly the call to the configured "Global" propagator.
Trk::VKalExtPropagator::myxAODFstPntOnTrk
const Perigee * myxAODFstPntOnTrk(const xAOD::TrackParticle *xprt) const
Definition: VKalExtPropagator.cxx:433
Trk::AmgSymMatrix
AmgSymMatrix(5) &GXFTrackState
Definition: GXFTrackState.h:156
Trk::ParticleHypothesis
ParticleHypothesis
Definition: ParticleHypothesis.h:25
xAOD::TrackParticle_v1::perigeeParameters
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
Definition: TrackParticle_v1.cxx:485
Trk::PropDirection
PropDirection
Definition: PropDirection.h:19
Trk::TrkVKalVrtFitter::m_IDsizeR
SimpleProperty< double > m_IDsizeR
Definition: TrkVKalVrtFitter.h:326
Trk::TrkVKalVrtFitter::setAthenaPropagator
void setAthenaPropagator(const Trk::IExtrapolator *)
Definition: VKalExtPropagator.cxx:473
Trk::CylinderSurface
Definition: CylinderSurface.h:55
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
Trk::TrkVKalVrtFitter::State::m_trkControl
std::vector< TrkMatControl > m_trkControl
Definition: TrkVKalVrtFitter.h:391
xAOD::TrackParticle_v1::radiusOfFirstHit
float radiusOfFirstHit() const
Returns the radius of the first hit.
Trk::pion
@ pion
Definition: ParticleHypothesis.h:29
Trk::IExtrapolator::extrapolate
virtual std::unique_ptr< NeutralParameters > extrapolate(const NeutralParameters &parameters, const Surface &sf, PropDirection dir=anyDirection, const BoundaryCheck &bcheck=true) const =0
Main extrapolation Interface starting from neutral parameters and aiming at surface.
Amg::Transform3D
Eigen::Affine3d Transform3D
Definition: GeoPrimitives.h:46
Trk::TrkVKalVrtFitter::m_IDsizeZ
SimpleProperty< double > m_IDsizeZ
Definition: TrkVKalVrtFitter.h:327
CylinderSurface.h
Trk::VKalExtPropagator::m_vkalFitSvc
TrkVKalVrtFitter * m_vkalFitSvc
Pointer to TrkVKalVrtFitter.
Definition: VKalExtPropagator.h:54
Trk::ParametersBase
Definition: ParametersBase.h:55
Scale
void Scale(TH1 *h, double d=1)
Definition: comparitor.cxx:76
Trk::muon
@ muon
Definition: ParticleHypothesis.h:28
Trk::VKalExtPropagator::Propagate
virtual void Propagate(long int trkID, long int Charge, double *ParOld, double *CovOld, double *RefStart, double *RefEnd, double *ParNew, double *CovNew, IVKalState &istate) const override
Definition: VKalExtPropagator.cxx:82
beamspotman.dir
string dir
Definition: beamspotman.py:623
Trk::VKalExtPropagator::checkTarget
virtual bool checkTarget(double *, const IVKalState &istate) const override
Definition: VKalExtPropagator.cxx:65
Trk
Ensure that the ATLAS eigen extensions are properly loaded.
Definition: FakeTrackBuilder.h:9
Trk::TrkVKalVrtFitter::State::m_fitField
VKalAtlasMagFld m_fitField
Definition: TrkVKalVrtFitter.h:401
Trk::VKalExtPropagator::setPropagator
void setPropagator(const IExtrapolator *)
Definition: VKalExtPropagator.cxx:41
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
Trk::IExtrapolator
Definition: IExtrapolator.h:62
Trk::VKalExtPropagator::myExtrapNeutral
const NeutralParameters * myExtrapNeutral(const NeutralParameters *inpPer, Amg::Vector3D *endPoint) const
Definition: VKalExtPropagator.cxx:416
Trk::TrkVKalVrtFitter::State::m_refFrameY
double m_refFrameY
Definition: TrkVKalVrtFitter.h:389
Trk::ParametersBase::momentum
const Amg::Vector3D & momentum() const
Access method for the momentum.
Trk::IVKalState
Definition: IVKalState.h:21
Trk::addNoise
@ addNoise
Definition: MaterialUpdateMode.h:19
Trk::VKalExtPropagator::Protection
double Protection(const double *, const IVKalState &istate) const
Definition: VKalExtPropagator.cxx:48
LArCellBinning.step
step
Definition: LArCellBinning.py:158
SG::ConstAccessor::isAvailable
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
Trk::VKalExtPropagator::~VKalExtPropagator
virtual ~VKalExtPropagator()
calibdata.copy
bool copy
Definition: calibdata.py:27
xAOD::TrackParticle_v1
Class describing a TrackParticle.
Definition: TrackParticle_v1.h:43
ConstAccessor.h
Helper class to provide constant type-safe access to aux data.
Trk::VKalExtPropagator::m_extrapolator
const IExtrapolator * m_extrapolator
Pointer to Extrapolator AlgTool.
Definition: VKalExtPropagator.h:53
Trk::SurfaceType::Line
@ Line
Trk::TrkVKalVrtFitter::State::m_eventContext
const EventContext * m_eventContext
Definition: TrkVKalVrtFitter.h:403
Trk::removeNoise
@ removeNoise
Definition: MaterialUpdateMode.h:20
Trk::VKalExtPropagator::myExtrapToLine
const TrackParameters * myExtrapToLine(long int TrkID, const TrackParameters *inpPer, Amg::Vector3D *endPoint, StraightLineSurface &lineTarget, const IVKalState &istate) const
Definition: VKalExtPropagator.cxx:351
VKalExtPropagator.h
Trk::StraightLineSurface
Definition: StraightLineSurface.h:51
Trk::TrkVKalVrtFitter
Definition: TrkVKalVrtFitter.h:67