ATLAS Offline Software
Loading...
Searching...
No Matches
TrkVKalVrtFitter.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Header include
7#include "GaudiKernel/EventContext.h"
8#include "GaudiKernel/IChronoStatSvc.h"
13#include "TrkTrack/Track.h"
19
20
21//-------------------------------------------------
22// Other stuff
23#include<iostream>
24
25namespace Trk{
26
27
28//
29//Constructor--------------------------------------------------------------
31 const std::string& name,
32 const IInterface* parent):
33 base_class(type,name,parent)
34 {
35 declareInterface<IVertexFitter>(this);
36 declareInterface<ITrkVKalVrtFitter>(this);
37 declareInterface<IVertexCascadeFitter>(this);
38
39/*--------------------------------------------------------------------------*/
40/* New propagator object is created. It's provided to VKalVrtCore. */
41/* VKalVrtFitter must set up Core BEFORE any call required propagation!!! */
42/* This object is created ONLY if IExtrapolator pointer is provideded. */
43/* see VKalExtPropagator.cxx for details */
44/*--------------------------------------------------------------------------*/
45 m_fitPropagator = nullptr; //Pointer to VKalVrtFitter propagator object to supply to VKalVrtCore (specific interface)
46 m_InDetExtrapolator = nullptr; //Direct pointer to Athena propagator
47}
48
49
50//Destructor---------------------------------------------------------------
52 //log << MSG::DEBUG << "TrkVKalVrtFitter destructor called" << endmsg;
53 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<<"TrkVKalVrtFitter destructor called" << endmsg;
55}
56
57std::unique_ptr<IVKalState>
58TrkVKalVrtFitter::makeState(const EventContext& ctx) const
59{
60 auto state = std::make_unique<State>();
61 initState(ctx, *state);
62 return state;
63}
64
66{
67 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG) <<"TrkVKalVrtFitter finalize() successful" << endmsg;
68 return StatusCode::SUCCESS;
69}
70
71
73{
74
75// Checking ROBUST algoritms
77
78
79 if(!m_useFixedField){
80 // Read handle for AtlasFieldCacheCondObj
81 if (!m_fieldCacheCondObjInputKey.key().empty()){
82 if( (m_fieldCacheCondObjInputKey.initialize()).isSuccess() ){
83 m_isAtlasField = true;
84 ATH_MSG_DEBUG( "Found AtlasFieldCacheCondObj with key ="<< m_fieldCacheCondObjInputKey.key());
85 }else{
86 ATH_MSG_INFO( "No AtlasFieldCacheCondObj with key ="<< m_fieldCacheCondObjInputKey.key());
87 ATH_MSG_INFO( "Use fixed magnetic field instead");
88 }
89 }
90 }
91//
92// Only here the VKalVrtFitter propagator object is created if ATHENA propagator is provided (see setAthenaPropagator)
93// In this case the ATHENA propagator can be used via pointers:
94// m_InDetExtrapolator - direct access
95// m_fitPropagator - via VKalVrtFitter object VKalExtPropagator
96// If ATHENA propagator is not provided, only defined object is
97// myPropagator - extern propagator from TrkVKalVrtCore
98//
99 if (m_extPropagator.empty()){
100 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<< "External propagator is not supplied - use internal one"<<endmsg;
101 m_extPropagator.disable();
102 }else{
103 if (m_extPropagator.retrieve().isFailure()) {
104 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<< "Could not find external propagator=" <<m_extPropagator<<endmsg;
105 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<< "TrkVKalVrtFitter will uses internal propagator" << endmsg;
106 m_extPropagator.disable();
107 }else{
108 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<< "External propagator="<<m_extPropagator<<" retrieved" << endmsg;
110 }
111 }
112
113//
114//
115//
116 if(msgLvl(MSG::DEBUG))msg(MSG::DEBUG)<< "TrkVKalVrtFitter initialize() successful" << endmsg;
117 if(msgLvl(MSG::DEBUG)){
118 msg(MSG::DEBUG)<< "TrkVKalVrtFitter configuration:" << endmsg;
119 msg(MSG::DEBUG)<< " Frozen version for BTagging: "<< m_frozenVersionForBTagging <<endmsg;
120 msg(MSG::DEBUG)<< " Allow ultra displaced vertices: "<< m_allowUltraDisplaced <<endmsg;
121 msg(MSG::DEBUG)<< " A priori vertex constraint: "<< m_useAprioriVertex <<endmsg;
122 msg(MSG::DEBUG)<< " Angle dTheta=0 constraint: "<< m_useThetaCnst <<endmsg;
123 msg(MSG::DEBUG)<< " Angle dPhi=0 constraint: "<< m_usePhiCnst <<endmsg;
124 msg(MSG::DEBUG)<< " Pointing to other vertex constraint: "<< m_usePointingCnst <<endmsg;
125 msg(MSG::DEBUG)<< " ZPointing to other vertex constraint: "<< m_useZPointingCnst <<endmsg;
126 msg(MSG::DEBUG)<< " Comb. particle pass near other vertex:"<< m_usePassNear <<endmsg;
127 msg(MSG::DEBUG)<< " Pass near with comb.particle errors: "<< m_usePassWithTrkErr <<endmsg;
129 msg(MSG::DEBUG)<< " Mass constraint M="<< m_massForConstraint <<endmsg;
130 msg(MSG::DEBUG)<< " with particles M=";
131 for(int i=0; i<(int)m_c_MassInputParticles.size(); i++) msg(MSG::DEBUG)<<m_c_MassInputParticles[i]<<", ";
132 msg(MSG::DEBUG)<<endmsg; ;
133 }
134 if(m_IterationNumber==0){
135 msg(MSG::DEBUG)<< " Default iteration number limit 50 is used " <<endmsg;
136 } else {
137 msg(MSG::DEBUG)<< " Iteration number limit: "<< m_IterationNumber <<endmsg;
138 }
139
140 if(m_isAtlasField){ msg(MSG::DEBUG)<< " ATLAS magnetic field is used!"<<endmsg; }
141 else { msg(MSG::DEBUG)<< " Constant magnetic field is used! B="<<m_BMAG<<endmsg; }
142
143 if(m_InDetExtrapolator){ msg(MSG::DEBUG)<< " InDet extrapolator is used!"<<endmsg; }
144 else { msg(MSG::DEBUG)<< " Internal VKalVrt extrapolator is used!"<<endmsg;}
145
146 if(m_Robustness) { msg(MSG::DEBUG)<< " VKalVrt uses robust algorithm! Type="<<m_Robustness<<" with Scale="<<m_RobustScale<<endmsg; }
147
148 if(m_firstMeasuredPoint){ msg(MSG::DEBUG)<< " VKalVrt will use FirstMeasuredPoint strategy in fits with InDetExtrapolator"<<endmsg; }
149 else { msg(MSG::DEBUG)<< " VKalVrt will use Perigee strategy in fits with InDetExtrapolator"<<endmsg; }
150 if(m_firstMeasuredPointLimit){ msg(MSG::DEBUG)<< " VKalVrt will use FirstMeasuredPointLimit strategy "<<endmsg; }
151 }
152
153
154 return StatusCode::SUCCESS;
155}
156
157
159{
160 initState(Gaudi::Hive::currentContext(), state);
161}
162
163
164void TrkVKalVrtFitter::initState (const EventContext& ctx, State& state) const
165
166{
167 //----------------------------------------------------------------------
168 // New magnetic field object is created. It's provided to VKalVrtCore.
169 // VKalVrtFitter must set up Core BEFORE any call required propagation!!!
170 if (m_isAtlasField) {
171 // For the moment, use Gaudi Hive for the event context - would need to be passed in from clients
173 const AtlasFieldCacheCondObj* fieldCondObj{*readHandle};
174 if (fieldCondObj == nullptr) {
175 ATH_MSG_ERROR("Failed to retrieve AtlasFieldCacheCondObj with key " << m_fieldCacheCondObjInputKey.key());
176 return;
177 }
178 fieldCondObj->getInitializedCache (state.m_fieldCache);
180 } else {
182 }
183 state.m_eventContext = &ctx;
200}
201
203xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const TrackParameters*> & perigeeListC,
204 const Amg::Vector3D & startingPoint) const
205{
206 State state;
207 initState (state);
208 setApproximateVertex(startingPoint.x(),
209 startingPoint.y(),
210 startingPoint.z(),
211 state);
212 std::vector<const NeutralParameters*> perigeeListN(0);
214 TLorentzVector Momentum;
215 long int Charge;
216 std::vector<double> ErrorMatrix;
217 std::vector<double> Chi2PerTrk;
218 std::vector< std::vector<double> > TrkAtVrt;
219 double Chi2;
220 StatusCode sc=VKalVrtFit( perigeeListC, perigeeListN,
221 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
222
223 xAOD::Vertex * tmpVertex = nullptr;
224 if(sc.isSuccess()) {
225 tmpVertex = makeXAODVertex( 0, Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
226 }
227 return tmpVertex;
228}
229
230
231xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const TrackParameters*> & perigeeListC,
232 const std::vector<const NeutralParameters*> & perigeeListN,
233 const Amg::Vector3D & startingPoint) const
234{
235 State state;
236 initState (state);
237 setApproximateVertex(startingPoint.x(),
238 startingPoint.y(),
239 startingPoint.z(),
240 state);
242 TLorentzVector Momentum;
243 long int Charge;
244 std::vector<double> ErrorMatrix;
245 std::vector<double> Chi2PerTrk;
246 std::vector< std::vector<double> > TrkAtVrt;
247 double Chi2;
248 StatusCode sc=VKalVrtFit( perigeeListC,perigeeListN,
249 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
250
251 xAOD::Vertex * tmpVertex = nullptr;
252 if(sc.isSuccess()) {
253 tmpVertex = makeXAODVertex( (int)perigeeListN.size(), Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
254 }
255 return tmpVertex;
256}
257
258
259
260
261
262
265xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const TrackParameters*> & perigeeListC,
266 const xAOD::Vertex & constraint) const
267{
268 State state;
269 initState (state);
270 if(msgLvl(MSG::DEBUG)) msg(MSG::DEBUG)<< "A priori vertex constraint is activated in VKalVrt fitter!" << endmsg;
271 Amg::Vector3D VertexIni(0.,0.,0.);
272 StatusCode sc=VKalVrtFitFast(perigeeListC, VertexIni, state);
273 if( sc.isSuccess()){
274 setApproximateVertex(VertexIni.x(),VertexIni.y(),VertexIni.z(),state);
275 }else{
276 setApproximateVertex(constraint.position().x(),
277 constraint.position().y(),
278 constraint.position().z(),
279 state);
280 }
281 setVertexForConstraint(constraint.position().x(),
282 constraint.position().y(),
283 constraint.position().z(),
284 state);
285 setCovVrtForConstraint(constraint.covariancePosition()(Trk::x,Trk::x),
286 constraint.covariancePosition()(Trk::y,Trk::x),
287 constraint.covariancePosition()(Trk::y,Trk::y),
288 constraint.covariancePosition()(Trk::z,Trk::x),
289 constraint.covariancePosition()(Trk::z,Trk::y),
290 constraint.covariancePosition()(Trk::z,Trk::z),
291 state);
292 state.m_useAprioriVertex=true;
293 std::vector<const NeutralParameters*> perigeeListN(0);
295 TLorentzVector Momentum;
296 long int Charge;
297 std::vector<double> ErrorMatrix;
298 std::vector<double> Chi2PerTrk;
299 std::vector< std::vector<double> > TrkAtVrt;
300 double Chi2;
301 sc=VKalVrtFit( perigeeListC, perigeeListN,
302 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
303
304
305 xAOD::Vertex * tmpVertex = nullptr;
306 if(sc.isSuccess()) {
307 tmpVertex = makeXAODVertex( 0, Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
308 }
309 return tmpVertex;
310}
311
312
313xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const TrackParameters*> & perigeeListC,
314 const std::vector<const NeutralParameters*> & perigeeListN,
315 const xAOD::Vertex & constraint) const
316{
317 State state;
318 initState (state);
319
320 if(msgLvl(MSG::DEBUG)) msg(MSG::DEBUG)<< "A priori vertex constraint is activated in VKalVrt fitter!" << endmsg;
321 Amg::Vector3D VertexIni(0.,0.,0.);
322 StatusCode sc=VKalVrtFitFast(perigeeListC, VertexIni, state);
323 if( sc.isSuccess()){
324 setApproximateVertex(VertexIni.x(),VertexIni.y(),VertexIni.z(),state);
325 }else{
326 setApproximateVertex(constraint.position().x(),
327 constraint.position().y(),
328 constraint.position().z(),
329 state);
330 }
331 setVertexForConstraint(constraint.position().x(),
332 constraint.position().y(),
333 constraint.position().z(),
334 state);
335 setCovVrtForConstraint(constraint.covariancePosition()(Trk::x,Trk::x),
336 constraint.covariancePosition()(Trk::y,Trk::x),
337 constraint.covariancePosition()(Trk::y,Trk::y),
338 constraint.covariancePosition()(Trk::z,Trk::x),
339 constraint.covariancePosition()(Trk::z,Trk::y),
340 constraint.covariancePosition()(Trk::z,Trk::z),
341 state);
342 state.m_useAprioriVertex=true;
344 TLorentzVector Momentum;
345 long int Charge;
346 std::vector<double> ErrorMatrix;
347 std::vector<double> Chi2PerTrk;
348 std::vector< std::vector<double> > TrkAtVrt;
349 double Chi2;
350 sc=VKalVrtFit( perigeeListC, perigeeListN,
351 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
352
353
354 xAOD::Vertex * tmpVertex = nullptr;
355 if(sc.isSuccess()) {
356 tmpVertex = makeXAODVertex( (int)perigeeListN.size(), Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
357 }
358 return tmpVertex;
359}
360
364std::unique_ptr<xAOD::Vertex>
365TrkVKalVrtFitter::fit(const EventContext& ctx,
366 const std::vector<const xAOD::TrackParticle*>& xtpListC,
367 const Amg::Vector3D& startingPoint) const
368{
369 State state;
370 initState(ctx, state);
371 return std::unique_ptr<xAOD::Vertex>(fit(xtpListC, startingPoint, state));
372}
373
374xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const xAOD::TrackParticle*> & xtpListC,
375 const Amg::Vector3D & startingPoint,
376 IVKalState& istate) const
377{
378 assert(dynamic_cast<State*> (&istate)!=nullptr);
379 State& state = static_cast<State&> (istate);
380
381 xAOD::Vertex * tmpVertex = nullptr;
382 setApproximateVertex(startingPoint.x(),
383 startingPoint.y(),
384 startingPoint.z(),
385 state);
386 std::vector<const xAOD::NeutralParticle*> xtpListN(0);
388 TLorentzVector Momentum;
389 long int Charge;
390 std::vector<double> ErrorMatrix;
391 std::vector<double> Chi2PerTrk;
392 std::vector< std::vector<double> > TrkAtVrt;
393 double Chi2;
394 StatusCode sc=VKalVrtFit( xtpListC, xtpListN,
395 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
396 if(sc.isSuccess()) {
397 tmpVertex = makeXAODVertex( 0, Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
398 dvect fittrkwgt;
399 sc=VKalGetTrkWeights(fittrkwgt, state); if(sc.isFailure())fittrkwgt.clear();
400 for(int ii=0; ii<state.m_FitStatus; ii++) {
402 if(!fittrkwgt.empty()) tmpVertex->addTrackAtVertex(TEL,fittrkwgt[ii]);
403 else tmpVertex->addTrackAtVertex(TEL,1.);
404 }
405 }
406
407 return tmpVertex;
408}
409
410xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const xAOD::TrackParticle*> & xtpListC,
411 const std::vector<const xAOD::NeutralParticle*> & xtpListN,
412 const Amg::Vector3D & startingPoint) const
413{
414 State state;
415 initState (state);
416 xAOD::Vertex * tmpVertex = nullptr;
417 setApproximateVertex(startingPoint.x(),
418 startingPoint.y(),
419 startingPoint.z(),
420 state);
422 TLorentzVector Momentum;
423 long int Charge;
424 std::vector<double> ErrorMatrix;
425 std::vector<double> Chi2PerTrk;
426 std::vector< std::vector<double> > TrkAtVrt;
427 double Chi2;
428 StatusCode sc=VKalVrtFit( xtpListC, xtpListN,
429 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
430 if(sc.isSuccess()) {
431 tmpVertex = makeXAODVertex( (int)xtpListN.size(), Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
432 dvect fittrkwgt;
433 sc=VKalGetTrkWeights(fittrkwgt, state); if(sc.isFailure())fittrkwgt.clear();
434 for(int ii=0; ii<state.m_FitStatus; ii++) {
435 if(ii<(int)xtpListC.size()) {
437 if(!fittrkwgt.empty()) tmpVertex->addTrackAtVertex(TEL,fittrkwgt[ii]);
438 else tmpVertex->addTrackAtVertex(TEL,1.);
439 }else{
441 if(!fittrkwgt.empty()) tmpVertex->addNeutralAtVertex(TEL,fittrkwgt[ii]);
442 else tmpVertex->addNeutralAtVertex(TEL,1.);
443 }
444 }
445 }
446
447 return tmpVertex;
448}
449
452xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const xAOD::TrackParticle*> & xtpListC,
453 const xAOD::Vertex & constraint) const
454{
455 State state;
456 initState (state);
457 return fit (xtpListC, constraint, state);
458}
459xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const xAOD::TrackParticle*> & xtpListC,
460 const xAOD::Vertex & constraint,
461 IVKalState& istate) const
462{
463 assert(dynamic_cast<State*> (&istate)!=nullptr);
464 State& state = static_cast<State&> (istate);
465
466 if(msgLvl(MSG::DEBUG)) msg(MSG::DEBUG)<< "A priori vertex constraint is activated in VKalVrt fitter!" << endmsg;
467 xAOD::Vertex * tmpVertex = nullptr;
468 setApproximateVertex(constraint.position().x(), constraint.position().y(),constraint.position().z(),state);
469 setVertexForConstraint(constraint.position().x(),
470 constraint.position().y(),
471 constraint.position().z(),
472 state);
473 setCovVrtForConstraint(constraint.covariancePosition()(Trk::x,Trk::x),
474 constraint.covariancePosition()(Trk::y,Trk::x),
475 constraint.covariancePosition()(Trk::y,Trk::y),
476 constraint.covariancePosition()(Trk::z,Trk::x),
477 constraint.covariancePosition()(Trk::z,Trk::y),
478 constraint.covariancePosition()(Trk::z,Trk::z),
479 state);
480 state.m_useAprioriVertex=true;
481 std::vector<const xAOD::NeutralParticle*> xtpListN(0);
483 TLorentzVector Momentum;
484 long int Charge;
485 std::vector<double> ErrorMatrix;
486 std::vector<double> Chi2PerTrk;
487 std::vector< std::vector<double> > TrkAtVrt;
488 double Chi2;
489 StatusCode sc=VKalVrtFit( xtpListC, xtpListN,
490 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
491 if(sc.isSuccess()) {
492 tmpVertex = makeXAODVertex( 0, Vertex, ErrorMatrix,Chi2PerTrk, TrkAtVrt, Chi2, state );
493 dvect fittrkwgt;
494 sc=VKalGetTrkWeights(fittrkwgt, state); if(sc.isFailure())fittrkwgt.clear();
495 for(int ii=0; ii<state.m_FitStatus; ii++) {
497 if(!fittrkwgt.empty()) tmpVertex->addTrackAtVertex(TEL,fittrkwgt[ii]);
498 else tmpVertex->addTrackAtVertex(TEL,1.);
499 }
500 }
501
502 return tmpVertex;
503}
504
505xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const xAOD::TrackParticle*> & xtpListC,
506 const std::vector<const xAOD::NeutralParticle*> & xtpListN,
507 const xAOD::Vertex & constraint) const
508{
509 State state;
510 initState (state);
511
512 if(msgLvl(MSG::DEBUG)) msg(MSG::DEBUG)<< "A priori vertex constraint is activated in VKalVrt fitter!" << endmsg;
513 xAOD::Vertex * tmpVertex = nullptr;
514 setApproximateVertex(constraint.position().x(), constraint.position().y(),constraint.position().z(),state);
515 setVertexForConstraint(constraint.position().x(),
516 constraint.position().y(),
517 constraint.position().z(),
518 state);
519 setCovVrtForConstraint(constraint.covariancePosition()(Trk::x,Trk::x),
520 constraint.covariancePosition()(Trk::y,Trk::x),
521 constraint.covariancePosition()(Trk::y,Trk::y),
522 constraint.covariancePosition()(Trk::z,Trk::x),
523 constraint.covariancePosition()(Trk::z,Trk::y),
524 constraint.covariancePosition()(Trk::z,Trk::z),
525 state);
526 state.m_useAprioriVertex=true;
528 TLorentzVector Momentum;
529 long int Charge;
530 std::vector<double> ErrorMatrix;
531 std::vector<double> Chi2PerTrk;
532 std::vector< std::vector<double> > TrkAtVrt;
533 double Chi2;
534 StatusCode sc=VKalVrtFit( xtpListC, xtpListN,
535 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
536 if(sc.isSuccess()){
537 tmpVertex = makeXAODVertex( (int)xtpListN.size(), Vertex, ErrorMatrix,Chi2PerTrk, TrkAtVrt, Chi2, state );
538 dvect fittrkwgt;
539 sc=VKalGetTrkWeights(fittrkwgt, state); if(sc.isFailure())fittrkwgt.clear();
540 for(int ii=0; ii<state.m_FitStatus; ii++) {
541 if(ii<(int)xtpListC.size()) {
543 if(!fittrkwgt.empty()) tmpVertex->addTrackAtVertex(TEL,fittrkwgt[ii]);
544 else tmpVertex->addTrackAtVertex(TEL,1.);
545 }else{
547 if(!fittrkwgt.empty()) tmpVertex->addNeutralAtVertex(TEL,fittrkwgt[ii]);
548 else tmpVertex->addNeutralAtVertex(TEL,1.);
549 }
550 }
551 }
552
553 return tmpVertex;
554}
555
556
557xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const TrackParameters*> & perigeeListC) const
558{
559 State state;
560 initState (state);
561 Amg::Vector3D VertexIni(0.,0.,0.);
562 StatusCode sc=VKalVrtFitFast(perigeeListC, VertexIni, state);
563 if( sc.isSuccess()) setApproximateVertex(VertexIni.x(),VertexIni.y(),VertexIni.z(),state);
564 std::vector<const NeutralParameters*> perigeeListN(0);
566 TLorentzVector Momentum;
567 long int Charge;
568 std::vector<double> ErrorMatrix;
569 std::vector<double> Chi2PerTrk;
570 std::vector< std::vector<double> > TrkAtVrt;
571 double Chi2;
572 sc=VKalVrtFit( perigeeListC, perigeeListN,
573 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
574
575 xAOD::Vertex * tmpVertex = nullptr;
576 if(sc.isSuccess()) {
577 tmpVertex = makeXAODVertex( 0, Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
578 }
579 return tmpVertex;
580}
581
582xAOD::Vertex * TrkVKalVrtFitter::fit(const std::vector<const TrackParameters*> & perigeeListC,
583 const std::vector<const NeutralParameters*> & perigeeListN) const
584{
585 State state;
586 initState (state);
587 Amg::Vector3D VertexIni(0.,0.,0.);
588 StatusCode sc=VKalVrtFitFast(perigeeListC, VertexIni, state);
589 if( sc.isSuccess()) setApproximateVertex(VertexIni.x(),VertexIni.y(),VertexIni.z(),state);
591 TLorentzVector Momentum;
592 long int Charge;
593 std::vector<double> ErrorMatrix;
594 std::vector<double> Chi2PerTrk;
595 std::vector< std::vector<double> > TrkAtVrt;
596 double Chi2;
597 sc=VKalVrtFit( perigeeListC, perigeeListN,
598 Vertex, Momentum, Charge, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state, true );
599
600 xAOD::Vertex * tmpVertex = nullptr;
601 if(sc.isSuccess()) {
602 tmpVertex = makeXAODVertex( (int)perigeeListN.size(), Vertex, ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2, state );
603 }
604 return tmpVertex;
605}
606
607
608
609/* Filling of 3x3 HepSymMatrix with content of symmetric matrix
610 in packed form vector<double> (6x6 - 21 elem)
611 (VxVyVzPxPyPz) */
612
613// Fills 5x5 matrix. Input Matrix is track covariance only.
614void TrkVKalVrtFitter::FillMatrixP(AmgSymMatrix(5)& CovMtx, std::vector<double> & Matrix)
615{
616 CovMtx.setIdentity();
617 if( Matrix.size() < 21) return;
618 CovMtx(0,0) = 0;
619 CovMtx(1,1) = 0;
620 CovMtx(2,2)= Matrix[ 9];
621 CovMtx.fillSymmetric(2,3,Matrix[13]);
622 CovMtx(3,3)= Matrix[14];
623 CovMtx.fillSymmetric(2,4,Matrix[18]);
624 CovMtx.fillSymmetric(3,4,Matrix[19]);
625 CovMtx(4,4)= Matrix[20];
626}
627
628// Fills 5x5 matrix. Input Matrix is a full covariance
629void TrkVKalVrtFitter::FillMatrixP(int iTrk, AmgSymMatrix(5)& CovMtx, std::vector<double> & Matrix)
630{
631 int iTmp=(iTrk+1)*3;
632 int NContent = Matrix.size();
633 CovMtx.setIdentity(); //Clean matrix for the beginning, then fill needed elements
634 CovMtx(0,0) = 0;
635 CovMtx(1,1) = 0;
636 int pnt = (iTmp+1)*iTmp/2 + iTmp; if( pnt > NContent ) return;
637 CovMtx(2,2) = Matrix[pnt];
638 pnt = (iTmp+1+1)*(iTmp+1)/2 + iTmp; if( pnt+1 > NContent ){ CovMtx.setIdentity(); return; }
639 CovMtx.fillSymmetric(2,3,Matrix[pnt]);
640 CovMtx(3,3) = Matrix[pnt+1];
641 pnt = (iTmp+2+1)*(iTmp+2)/2 + iTmp; if( pnt+2 > NContent ){ CovMtx.setIdentity(); return; }
642 CovMtx.fillSymmetric(2,4,Matrix[pnt]);
643 CovMtx.fillSymmetric(3,4,Matrix[pnt+1]);
644 CovMtx(4,4) = Matrix[pnt+2];
645}
646
647
648
649Amg::MatrixX * TrkVKalVrtFitter::GiveFullMatrix(int NTrk, std::vector<double> & Matrix)
650{
651 Amg::MatrixX * mtx = new Amg::MatrixX(3+3*NTrk,3+3*NTrk);
652 long int ij=0;
653 for(int i=1; i<=(3+3*NTrk); i++){
654 for(int j=1; j<=i; j++){
655 if(i==j){ (*mtx)(i-1,j-1)=Matrix[ij];}
656 else { (*mtx).fillSymmetric(i-1,j-1,Matrix[ij]);}
657 ij++;
658 }
659 }
660 return mtx;
661}
662
663
664
666 const Amg::Vector3D& Vertex, const std::vector<double> & fitErrorMatrix,
667 const std::vector<double> & Chi2PerTrk, const std::vector< std::vector<double> >& TrkAtVrt,
668 double Chi2,
669 State& state) const
670{
671 long int NTrk = state.m_FitStatus;
672 long int Ndf = VKalGetNDOF(state)+state.m_planeCnstNDOF;
673
674 xAOD::Vertex * tmpVertex=new xAOD::Vertex();
675 tmpVertex->makePrivateStore();
676 tmpVertex->setPosition(Vertex);
677 tmpVertex->setFitQuality(Chi2, (float)Ndf);
678
679 std::vector<VxTrackAtVertex> & tmpVTAV=tmpVertex->vxTrackAtVertex();
680 tmpVTAV.clear();
681 std::vector <double> CovFull;
682 StatusCode sc = VKalGetFullCov( NTrk, CovFull, state);
683 int covarExist=0; if( sc.isSuccess() ) covarExist=1;
684
685 std::vector<float> floatErrMtx;
686 if( m_makeExtendedVertex && covarExist ) {
687 floatErrMtx.resize(CovFull.size());
688 for(int i=0; i<(int)CovFull.size(); i++) {
689 if( CovFull[i] < std::numeric_limits<float>::max() &&
690 CovFull[i] > std::numeric_limits<float>::lowest() ){
691 floatErrMtx[i]=static_cast<float>(CovFull[i]);
692 } else {
693 floatErrMtx[i]=std::numeric_limits<float>::max();
694 }
695 }
696 }else{
697 floatErrMtx.resize(fitErrorMatrix.size());
698 for(int i=0; i<(int)fitErrorMatrix.size(); i++) {
699 if( fitErrorMatrix[i] < std::numeric_limits<float>::max() &&
700 fitErrorMatrix[i] > std::numeric_limits<float>::lowest() ){
701 floatErrMtx[i]=static_cast<float>(fitErrorMatrix[i]);
702 } else {
703 floatErrMtx[i]=std::numeric_limits<float>::max();
704 }
705 }
706 }
707 tmpVertex->setCovariance(floatErrMtx);
708
709 for(int ii=0; ii<NTrk ; ii++) {
710 AmgSymMatrix(5) CovMtxP;
711 if(covarExist){ FillMatrixP( ii, CovMtxP, CovFull );}
712 else { CovMtxP.setIdentity();}
713 Perigee * tmpChargPer=nullptr;
714 NeutralPerigee * tmpNeutrPer=nullptr;
715 if(ii<NTrk-Neutrals){
716 tmpChargPer = new Perigee( 0.,0., TrkAtVrt[ii][0],
717 TrkAtVrt[ii][1],
718 TrkAtVrt[ii][2],
719 PerigeeSurface(Vertex), std::move(CovMtxP) );
720 }else{
721 tmpNeutrPer = new NeutralPerigee( 0.,0., TrkAtVrt[ii][0],
722 TrkAtVrt[ii][1],
723 TrkAtVrt[ii][2],
725 std::move(CovMtxP) );
726 }
727 tmpVTAV.emplace_back(Chi2PerTrk[ii], tmpChargPer, tmpNeutrPer );
728 }
729
730 return tmpVertex;
731}
732
733
734} // End Of Namespace
#define endmsg
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
static Double_t sc
void getInitializedCache(MagField::AtlasFieldCache &cache) const
get B field cache for evaluation as a function of 2-d or 3-d position.
void makePrivateStore()
Create a new (empty) private store for this object.
Class describing the Line to which the Perigee refers to.
std::vector< double > m_MassInputParticles
MagField::AtlasFieldCache m_fieldCache
std::vector< double > m_VertexForConstraint
std::vector< double > m_CovVrtForConstraint
const EventContext * m_eventContext
virtual std::unique_ptr< IVKalState > makeState() const
TrkVKalVrtFitter(const std::string &t, const std::string &name, const IInterface *parent)
const IExtrapolator * m_InDetExtrapolator
Pointer to Extrapolator AlgTool.
Gaudi::Property< bool > m_makeExtendedVertex
SG::ReadCondHandleKey< AtlasFieldCacheCondObj > m_fieldCacheCondObjInputKey
static Amg::MatrixX * GiveFullMatrix(int NTrk, std::vector< double > &)
Gaudi::Property< bool > m_usePhiCnst
virtual void setCovVrtForConstraint(double XX, double XY, double YY, double XZ, double YZ, double ZZ, IVKalState &istate) const override final
virtual StatusCode VKalVrtFitFast(std::span< const xAOD::TrackParticle *const >, Amg::Vector3D &Vertex, double &minDZ, IVKalState &istate) const
virtual void setVertexForConstraint(const xAOD::Vertex &, IVKalState &istate) const override final
Gaudi::Property< bool > m_usePointingCnst
ToolHandle< IExtrapolator > m_extPropagator
Gaudi::Property< int > m_IterationNumber
Gaudi::Property< bool > m_allowUltraDisplaced
virtual xAOD::Vertex * fit(const std::vector< const TrackParameters * > &perigeeList, const Amg::Vector3D &startingPoint) const override final
Interface for MeasuredPerigee with starting point.
Gaudi::Property< double > m_RobustScale
void initState(const EventContext &ctx, State &state) const
virtual StatusCode initialize() override final
Gaudi::Property< std::vector< double > > m_c_CovVrtForConstraint
xAOD::Vertex * makeXAODVertex(int, const Amg::Vector3D &, const dvect &, const dvect &, const std::vector< dvect > &, double, State &state) const
virtual StatusCode VKalGetFullCov(long int, dvect &CovMtx, IVKalState &istate, bool=false) const override final
Gaudi::Property< double > m_massForConstraint
virtual StatusCode VKalVrtFit(const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::NeutralParticle * > &, Amg::Vector3D &Vertex, TLorentzVector &Momentum, long int &Charge, dvect &ErrorMatrix, dvect &Chi2PerTrk, std::vector< std::vector< double > > &TrkAtVrt, double &Chi2, IVKalState &istate, bool ifCovV0=false) const override final
Gaudi::Property< bool > m_frozenVersionForBTagging
Gaudi::Property< bool > m_firstMeasuredPoint
Gaudi::Property< int > m_Robustness
virtual StatusCode VKalGetTrkWeights(dvect &Weights, const IVKalState &istate) const override final
static int VKalGetNDOF(const State &state)
Gaudi::Property< bool > m_usePassWithTrkErr
static void FillMatrixP(AmgSymMatrix(5)&, std::vector< double > &)
Gaudi::Property< bool > m_useZPointingCnst
Gaudi::Property< bool > m_useAprioriVertex
Gaudi::Property< bool > m_useFixedField
virtual StatusCode finalize() override final
Gaudi::Property< bool > m_usePassNear
VKalExtPropagator * m_fitPropagator
virtual void setApproximateVertex(double X, double Y, double Z, IVKalState &istate) const override final
Gaudi::Property< std::vector< double > > m_c_MassInputParticles
Gaudi::Property< std::vector< double > > m_c_VertexForConstraint
Gaudi::Property< bool > m_firstMeasuredPointLimit
void setAthenaPropagator(const Trk::IExtrapolator *)
Gaudi::Property< bool > m_useThetaCnst
void setAtlasField(MagField::AtlasFieldCache *)
const basePropagator * vk_objProp
This class is a simplest representation of a vertex candidate.
void addTrackAtVertex(const ElementLink< TrackParticleContainer > &tr, float weight=1.0)
Add a new track to the vertex.
void addNeutralAtVertex(const ElementLink< NeutralParticleContainer > &tr, float weight=1.0)
Add a new neutral to the vertex.
void setCovariance(const std::vector< float > &value)
Sets the covariance matrix as a simple vector of values.
void setPosition(const Amg::Vector3D &position)
Sets the 3-position.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
void setFitQuality(float chiSquared, float numberDoF)
Set the 'Fit Quality' information.
const Amg::Vector3D & position() const
Returns the 3-pos.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the ATLAS eigen extensions are properly loaded.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
ParametersT< NeutralParametersDim, Neutral, PerigeeSurface > NeutralPerigee
@ x
Definition ParamDefs.h:55
@ z
global position (cartesian)
Definition ParamDefs.h:57
@ y
Definition ParamDefs.h:56
std::vector< double > dvect
Vertex_v1 Vertex
Define the latest version of the vertex class.
MsgStream & msg
Definition testRead.cxx:32