ATLAS Offline Software
Loading...
Searching...
No Matches
BeamspotVertexPreProcessor.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef TRKALIGNGENTOOLS_BEAMSPOTVERTEXPREPROCESSOR_H
6#define TRKALIGNGENTOOLS_BEAMSPOTVERTEXPREPROCESSOR_H
7
8#include "GaudiKernel/ToolHandle.h"
9#include "GaudiKernel/IAlgTool.h"
12
14
20
22
25
41
42
43namespace Trk {
44 class Track;
45 class AlignTrack;
46 class AlignVertex;
47 class VertexOnTrack;
48 class VxTrackAtVertex;
49
51 {
52
53 public:
54
55 BeamspotVertexPreProcessor(const std::string & type, const std::string & name, const IInterface * parent);
57
58 virtual StatusCode initialize() override;
59 virtual StatusCode finalize() override;
60
61 virtual DataVector<Track> * processTrackCollection(const DataVector<Track> * trks) override;
62
63 void accumulateVTX(AlignTrack* alignTrack) override;
64
65 void solveVTX() override;
66
68 virtual void printSummary() override;
69
70
71
72 private:
73
75
76 bool isAssociatedToPV(const Track * track, const xAOD::VertexContainer* vertices);
77 bool isAssociatedToVertex(const Track * track, const xAOD::Vertex * vertex);
78
79 bool selectVertices(const xAOD::Vertex * vtx) const;
80 bool selectUpdatedVertices(const xAOD::Vertex * updatedVtx) const;
81
82 const xAOD::Vertex* findVertexCandidate(const Track* track) const; // MD: changed to be now only a vertex
83 const VertexOnTrack* provideVotFromVertex(const Track* track, const xAOD::Vertex* &vtx) const;
84 const VertexOnTrack* provideVotFromBeamspot(const Track* track) const;
85 void provideVtxBeamspot(const AlignVertex* b, AmgSymMatrix(3)* q, Amg::Vector3D* v) const;
86
87
88 const Track* doConstraintRefit(ToolHandle<IGlobalTrackFitter>& fitter, const Track* track, const VertexOnTrack* vot, const ParticleHypothesis& particleHypothesis) const;
90 AlignTrack* doTrackRefit(const Track* track);
91
92
93 ToolHandle<IGlobalTrackFitter> m_trackFitter{
94 this, "TrackFitter", "Trk::GlobalChi2Fitter/InDetTrackFitter",
95 "normal track fitter"};
96 ToolHandle<IGlobalTrackFitter> m_SLTrackFitter{
97 this, "SLTrackFitter", "", "straight line track fitter"};
98 ToolHandle<IExtrapolator> m_extrapolator{
99 this, "Extrapolator", "Trk::Extrapolator/AtlasExtrapolator"};
100 ToolHandle<InDet::IInDetTrackSelectionTool> m_trkSelector{
101 this, "TrackSelector", "", "new track selector tool"};
102 ToolHandle<InDet::IInDetTrackSelectionTool> m_BSTrackSelector{
103 this, "BSConstraintTrackSelector", "",
104 "new track selector tool for tracks to be used with beam-spot constraint"};
105 ToolHandle<ITrackToVertexIPEstimator> m_trackToVertexIPEstimatorTool{
106 this, "TrackToVertexIPEstimatorTool", ""};
107
109 PublicToolHandle<IAlignModuleTool> m_alignModuleTool{
110 this, "AlignModuleTool", "InDet::InDetAlignModuleTool/InDetAlignModuleTool"};
111
113 this, "BeamSpotKey", "BeamSpotData", "SG key for beam spot" };
114
116 this, "PVContainerName", "PrimaryVertices"};
117
118 BooleanProperty m_runOutlierRemoval{this, "RunOutlierRemoval", false,
119 "switch whether to run outlier logics or not"};
120 IntegerProperty m_particleNumber{this, "ParticleNumber", 3,
121 "type of material interaction in extrapolation, 3=pion, 0=non-interacting"};
122 BooleanProperty m_doTrkSelection{this, "DoTrackSelection", true,
123 "to activate the preprocessor track selection"};
124 BooleanProperty m_doBSTrackSelection{this, "DoBSTrackSelection", false,
125 "the selection mechanism which is based on cutting the perigee parameters, pt, etc."};
127 this, "DoAssociatedToPVSelection", true,
128 "the selection mechanism that only use the tracks associated to PV"};
129
130 UnsignedIntegerProperty m_constraintMode{this, "ConstraintMode", 0};
131
132 std::vector< std::pair< const xAOD::Vertex*, std::vector<VxTrackAtVertex> > > m_allTracksVector;
133
134 BooleanProperty m_doBeamspotConstraint{this, "DoBSConstraint", true,
135 "Constrain tracks to the beamspot (x,y) position"};
136 BooleanProperty m_doPrimaryVertexConstraint{this, "DoPVConstraint", false,
137 "Constrain tracks to the associated primary vertex (x,y,z) position"};
138 BooleanProperty m_doFullVertexConstraint{this, "DoFullVertex", false,
139 "Full 3D vertex constraint. Note DoPVConstraint needs to be set to true to use this option. If DoBSConstraint vertex position will be constrained to the BS"};
140 BooleanProperty m_doNormalRefit{this, "doNormalRefit", true,
141 "provide tracks in the case failed BS, PV and FullVertex constraints."};
142
143 DoubleProperty m_maxPt{this, "maxPt", 0.,
144 "Max pT range for refitting tracks"};
145
146 BooleanProperty m_refitTracks{this, "RefitTracks", true,
147 "flag to refit tracks"};
148 BooleanProperty m_storeFitMatrices{this, "StoreFitMatrices", true,
149 "flag to store derivative and covariance matrices after refit"};
150 BooleanProperty m_useSingleFitter{this, "UseSingleFitter", false,
151 "only use 1 fitter for refitting track"};
152 DoubleProperty m_BSScalingFactor{this, "BeamspotScalingFactor", 1.,
153 "scaling factor on beasmpot width"};
154 DoubleProperty m_PVScalingFactor{this, "PrimaryVertexScalingFactor", 1.,
155 "scaling factor on primary vertex position error"};
156
157 IntegerProperty m_minTrksInVtx{this, "MinTrksInVtx", 3,
158 "requirement to the minimal number of tracks in the vertex"};
159
160 int m_nTracks = 0;
161 std::vector<int> m_trackTypeCounter{};
165
167
168
169 };
170
171
172
174
175 public:
176 CompareTwoTracks(const Track* track, const std::string& compareMethod)
177 :m_method(compareMethod)
178 ,m_track(track)
179 { //std::cout <<"compareMethod: "<< m_method <<std::endl;
180 }
181
182 bool operator()(VxTrackAtVertex vtxTrk);
183
184 private:
185 std::string m_method;
187 };
188
189
190}
191
192#endif // TRKALIGNGENTOOLS_BEAMSPOTVERTEXPREPROCESSOR_H
193
#define AmgSymMatrix(dim)
Property holding a SG store/key/clid from which a ReadHandle is made.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Derived DataVector<T>.
Definition DataVector.h:795
Property holding a SG store/key/clid from which a ReadHandle is made.
bool doBeamspotConstraintTrackSelection(const Track *track)
PublicToolHandle< IAlignModuleTool > m_alignModuleTool
Pointer to AlignModuleTool.
SG::ReadHandleKey< xAOD::VertexContainer > m_PVContainerName
void accumulateVTX(AlignTrack *alignTrack) override
methods added for the full VTX fit:
const VertexOnTrack * provideVotFromVertex(const Track *track, const xAOD::Vertex *&vtx) const
bool selectVertices(const xAOD::Vertex *vtx) const
bool isAssociatedToVertex(const Track *track, const xAOD::Vertex *vertex)
ToolHandle< InDet::IInDetTrackSelectionTool > m_trkSelector
bool isAssociatedToPV(const Track *track, const xAOD::VertexContainer *vertices)
ToolHandle< IGlobalTrackFitter > m_trackFitter
const xAOD::Vertex * findVertexCandidate(const Track *track) const
AlignTrack * doTrackRefit(const Track *track)
virtual DataVector< Track > * processTrackCollection(const DataVector< Track > *trks) override
Main processing of track collection.
void provideVtxBeamspot(const AlignVertex *b, AmgSymMatrix(3) *q, Amg::Vector3D *v) const
ToolHandle< IExtrapolator > m_extrapolator
DataVector< AlignVertex > m_AlignVertices
collection of AlignVertices used in FullVertex constraint option
std::vector< std::pair< const xAOD::Vertex *, std::vector< VxTrackAtVertex > > > m_allTracksVector
const Track * doConstraintRefit(ToolHandle< IGlobalTrackFitter > &fitter, const Track *track, const VertexOnTrack *vot, const ParticleHypothesis &particleHypothesis) const
ToolHandle< IGlobalTrackFitter > m_SLTrackFitter
virtual StatusCode initialize() override
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
ToolHandle< InDet::IInDetTrackSelectionTool > m_BSTrackSelector
ToolHandle< ITrackToVertexIPEstimator > m_trackToVertexIPEstimatorTool
virtual void printSummary() override
Print processing summary to logfile.
bool selectUpdatedVertices(const xAOD::Vertex *updatedVtx) const
BeamspotVertexPreProcessor(const std::string &type, const std::string &name, const IInterface *parent)
const VertexOnTrack * provideVotFromBeamspot(const Track *track) const
CompareTwoTracks(const Track *track, const std::string &compareMethod)
bool operator()(VxTrackAtVertex vtxTrk)
Class to handle Vertex On Tracks, it inherits from the common MeasurementBase.
The VxTrackAtVertex is a common class for all present TrkVertexFitters The VxTrackAtVertex is designe...
Eigen::Matrix< double, 3, 1 > Vector3D
Ensure that the ATLAS eigen extensions are properly loaded.
@ v
Definition ParamDefs.h:78
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.