ATLAS Offline Software
Loading...
Searching...
No Matches
RefitTracksAndVertex Class Reference

#include <RefitTracksAndVertex.h>

Inheritance diagram for RefitTracksAndVertex:
Collaboration diagram for RefitTracksAndVertex:

Public Member Functions

 RefitTracksAndVertex (const std::string &name, ISvcLocator *pSvcLocator)
 The RefitTracksAndVertex is an implementation to add the so-called TRT momentum constraint on a track.
 ~RefitTracksAndVertex ()
StatusCode initialize ()
 initialize method of this algorithm.
StatusCode execute (const EventContext &ctx)
 execute method of this algorithm that is called for each event
StatusCode finalize ()
 finalize method of this algorithm.
MsgStream & dump (MsgStream &outst) const
virtual StatusCode sysInitialize () override
 Override sysInitialize.
virtual const DataObjIDColl & extraOutputDeps () const override
 Return the list of extra output dependencies.
bool filterPassed (const EventContext &ctx) const
void setFilterPassed (bool state, const EventContext &ctx) const
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const
const EventContext & getContext () const
 Deprecated methods (use the ones with EventContext).
bool filterPassed () const
void setFilterPassed (bool state) const

Protected Member Functions

virtual bool isReEntrant () const override final
 Legacy algorithms are not thread-safe.
void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

bool accept (const Trk::Track &track)
 Verifies if the given track passes the track selection criteria specified via the jobOptions.
const Trk::MeasurementSet addPM (Trk::MeasurementSet &ms, const Trk::PseudoMeasurementOnTrack *pm)
 adds a PseudoMeasurement to a MeasurementSet
const Trk::PseudoMeasurementOnTrackcreatePMfromSi (const Trk::Perigee *mp)
 creates a PseudoMeasurement with (z0, theta) from extended track perigee parameters
Trk::TrackfitSCTOnlyTrack (const Trk::Track *track)
 Strips of all TRT hits from the track and replaces them with a TRT momentum constraint from a the TRT segment belonging to the (extended) track (returns NULL in case of failure).
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

const AtlasDetectorIDm_idHelper
 Detector ID helper.
std::string m_vertexListInput
 Name of the TrackCollection (input).
std::string m_trackListOutput
 Name of the TrackCollection (Output).
std::string m_outputVertexContainerName
 Name of vertex container.
ToolHandle< Trk::ITrackFitterm_trackFitter
 The TrackFitter to refit the tracks (segment, momentum constraint).
ToolHandle< Trk::IVertexFitterm_vertexFitter
 The TrackFitter to refit the tracks (segment, momentum constraint).
bool m_applyTrkSel
 apply a selection on tracks or not
bool m_refitTracks
 refitTracks
bool m_addPM
 apply a pseudo measurement based on the original track (theta,z0)
double m_selPtMin
 minimal pT cut value for the TrackSelection
double m_selEtaMin
 minimal eta cut value for the TrackSelection
double m_selEtaMax
 maximal eta cut value for the TrackSelection
int m_selNHitPIXMin
 minimal number of PIX hits cut value for the TrackSelection
int m_selNHitSCTMin
 minimal number of SCT hits cut value for the TrackSelection
size_t m_nTracksProcessed
 Counter for number of tracks processed.
size_t m_nTracksPresel
 Counter for number of tracks passing the preselection.
size_t m_nTracksAccepted
 Counter for number of tracks passing the preselection and with PM.
size_t m_nRejectPIX
 Counter for number of tracks failing the min number of PIX hits req.
size_t m_nRejectSCT
 Counter for number of tracks failing the min number of SCT hits req.
size_t m_nRejectTRT
 Counter for number of tracks failing the min number of TRT hits req.
size_t m_nRejectPM
 Counter for number of tracks failing the addition of a pseudo-measurement (PM).
DataObjIDColl m_extendedExtraObjects
StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default).
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default).
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 31 of file RefitTracksAndVertex.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ RefitTracksAndVertex()

RefitTracksAndVertex::RefitTracksAndVertex ( const std::string & name,
ISvcLocator * pSvcLocator )

The RefitTracksAndVertex is an implementation to add the so-called TRT momentum constraint on a track.

It was developed in the context of alignment, and its main idea is that if the TRT is free from "weak modes", its curvature measurement can be added as a PseudoMeasurement (PM) to a Si-only track. This extra information, taken into account while refitting the Si-only part of the track or read in by the alignment algorithms directly, can serve to eliminate the alignment weak modes from the Si tracker of ATLAS. More details can be found at: https://twiki.cern.ch/twiki/bin/view/Atlas/TRTMomConstraint

Definition at line 27 of file RefitTracksAndVertex.cxx.

28 : AthAlgorithm (name, pSvcLocator)
29 , m_idHelper(nullptr)
30 , m_vertexListInput("PrimaryVertices")
31 , m_trackListOutput("SiTracksWithRefitTracksAndVertex")
32 , m_outputVertexContainerName("RefittedPrimaryVertex")
33 , m_trackFitter("Trk::GlobalChi2Fitter/InDetTrackFitter")
34 , m_vertexFitter("Trk::FullVertexFitter")
35 , m_applyTrkSel(false)
36 , m_refitTracks(true)
37 , m_addPM(false)
38 , m_selPtMin(0)
39 , m_selEtaMin(-2.5)
40 , m_selEtaMax(2.5)
46 , m_nRejectPIX(0)
47 , m_nRejectSCT(0)
48 , m_nRejectTRT(0)
49 , m_nRejectPM(0)
50{
51 // declare algorithm parameters
52 declareProperty("VertexListInput" , m_vertexListInput) ;
53 declareProperty("TrackListOutput" , m_trackListOutput) ;
54 declareProperty("OutputVertexContainer" , m_outputVertexContainerName) ;
55 declareProperty("TrackFitter" , m_trackFitter) ;
56 declareProperty("VertexFitterTool" , m_vertexFitter);
57 declareProperty("ApplyTrkSel" , m_applyTrkSel) ;
58 declareProperty("RefitTracks" , m_refitTracks) ;
59 declareProperty("AddPseudoMeasurement", m_addPM);
60 declareProperty("SelPtMin" , m_selPtMin) ;
61 declareProperty("SelEtaMin" , m_selEtaMin) ;
62 declareProperty("SelEtaMax" , m_selEtaMax) ;
63 declareProperty("SelNHitPIXMin" , m_selNHitPIXMin) ;
64 declareProperty("SelNHitSCTMin" , m_selNHitSCTMin) ;
65
66}
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
std::string m_vertexListInput
Name of the TrackCollection (input).
size_t m_nTracksPresel
Counter for number of tracks passing the preselection.
const AtlasDetectorID * m_idHelper
Detector ID helper.
double m_selEtaMax
maximal eta cut value for the TrackSelection
double m_selEtaMin
minimal eta cut value for the TrackSelection
int m_selNHitSCTMin
minimal number of SCT hits cut value for the TrackSelection
int m_selNHitPIXMin
minimal number of PIX hits cut value for the TrackSelection
size_t m_nRejectSCT
Counter for number of tracks failing the min number of SCT hits req.
ToolHandle< Trk::IVertexFitter > m_vertexFitter
The TrackFitter to refit the tracks (segment, momentum constraint).
size_t m_nTracksProcessed
Counter for number of tracks processed.
std::string m_outputVertexContainerName
Name of vertex container.
size_t m_nRejectPM
Counter for number of tracks failing the addition of a pseudo-measurement (PM).
size_t m_nRejectPIX
Counter for number of tracks failing the min number of PIX hits req.
bool m_refitTracks
refitTracks
double m_selPtMin
minimal pT cut value for the TrackSelection
bool m_applyTrkSel
apply a selection on tracks or not
size_t m_nRejectTRT
Counter for number of tracks failing the min number of TRT hits req.
ToolHandle< Trk::ITrackFitter > m_trackFitter
The TrackFitter to refit the tracks (segment, momentum constraint).
std::string m_trackListOutput
Name of the TrackCollection (Output).
bool m_addPM
apply a pseudo measurement based on the original track (theta,z0)
size_t m_nTracksAccepted
Counter for number of tracks passing the preselection and with PM.

◆ ~RefitTracksAndVertex()

RefitTracksAndVertex::~RefitTracksAndVertex ( )

Definition at line 68 of file RefitTracksAndVertex.cxx.

68{}

Member Function Documentation

◆ accept()

bool RefitTracksAndVertex::accept ( const Trk::Track & track)
private

Verifies if the given track passes the track selection criteria specified via the jobOptions.

◆ addPM()

const Trk::MeasurementSet RefitTracksAndVertex::addPM ( Trk::MeasurementSet & ms,
const Trk::PseudoMeasurementOnTrack * pm )
private

adds a PseudoMeasurement to a MeasurementSet

Definition at line 236 of file RefitTracksAndVertex.cxx.

236 {
237 Trk::MeasurementSet sortedMS ;
238 sortedMS.push_back( pm ) ;
239 for( int i=0, i_max=ms.size() ; i!=i_max ; ++i ) {
240 sortedMS.push_back( ms[i] ) ;
241 }
242 return sortedMS ;
243}
std::vector< const MeasurementBase * > MeasurementSet
vector of fittable measurements
Definition FitterTypes.h:30

◆ createPMfromSi()

const Trk::PseudoMeasurementOnTrack * RefitTracksAndVertex::createPMfromSi ( const Trk::Perigee * mp)
private

creates a PseudoMeasurement with (z0, theta) from extended track perigee parameters

Definition at line 212 of file RefitTracksAndVertex.cxx.

212 {
213 Trk::DefinedParameter z0 ( mp->parameters()[Trk::z0], Trk::z0 ) ;
214 Trk::DefinedParameter theta( mp->parameters()[Trk::theta], Trk::theta ) ;
215 std::vector<Trk::DefinedParameter> defPar ;
216 defPar.push_back( z0 ) ;
217 defPar.push_back( theta ) ;
218 if( !mp->covariance() ) return nullptr;
219 Trk::LocalParameters parFromSi( defPar ) ;
220 AmgSymMatrix(2) covFromSi;
221
222 covFromSi( 0, 0 ) = (*mp->covariance())( Trk::z0,Trk::z0 ) ;
223 covFromSi( 1, 1 ) = (*mp->covariance())( Trk::theta,Trk::theta ) ;
224 covFromSi( 1, 0 ) = (*mp->covariance())( Trk::z0, Trk::theta ) ;
225 covFromSi( 0, 1 ) = (*mp->covariance())( Trk::z0, Trk::theta ) ;
226
227 const Trk::Surface& mpSurf = mp->associatedSurface() ;
228
229 Trk::PseudoMeasurementOnTrack *pm = new Trk::PseudoMeasurementOnTrack( std::move(parFromSi)
230 , covFromSi
231 , mpSurf
232 ) ;
233 return pm ;
234}
Scalar theta() const
theta method
#define AmgSymMatrix(dim)
@ theta
Definition ParamDefs.h:66
@ z0
Definition ParamDefs.h:64
std::pair< double, ParamDefs > DefinedParameter
Typedef to of a std::pair<double, ParamDefs> to identify a passed-through double as a specific type o...

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ dump()

MsgStream & RefitTracksAndVertex::dump ( MsgStream & outst) const

Definition at line 344 of file RefitTracksAndVertex.cxx.

344 {
345 outst << std::endl ;
346 outst << "|-------------------------------------------------------------------";
347 outst << "-----------------------------|" << std::endl ;
348 outst << "| processed : "
349 << std::setw(7) << m_nTracksProcessed
350 << " tracks |"
351 << std::endl ;
352 outst << "| accepted by track presel. : "
353 << std::setw(7) << m_nTracksPresel
354 << " tracks |"
355 << std::endl ;
356 outst << "| accepted by track presel. + PM : "
357 << std::setw(7) << m_nTracksAccepted
358 << " tracks |"
359 << std::endl ;
360 outst << "| ------------------------------------------------------------------";
361 outst << "---------------------------- |" << std::endl ;
362 outst << "| reject by # PIX hits : "
363 << std::setw(7) << m_nRejectPIX
364 << " tracks |"
365 << std::endl ;
366 outst << "| reject by # SCT hits : "
367 << std::setw(7) << m_nRejectSCT
368 << " tracks |"
369 << std::endl ;
370 outst << "| reject by # TRT hits : "
371 << std::setw(7) << m_nRejectTRT
372 << " tracks |"
373 << std::endl ;
374 outst << "| ------------------------------------------------------------------";
375 outst << "---------------------------- |" << std::endl ;
376 outst << "| reject by exist. PM(TRT) : "
377 << std::setw(7) << m_nRejectPM
378 << " tracks |"
379 << std::endl ;
380 outst << "|-------------------------------------------------------------------";
381 outst << "-----------------------------|" << std::endl ;
382 return outst ;
383}

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

StatusCode RefitTracksAndVertex::execute ( const EventContext & ctx)
virtual

execute method of this algorithm that is called for each event

Implements AthAlgorithm.

Definition at line 123 of file RefitTracksAndVertex.cxx.

123 {
124
125 TrackCollection* outputtracks = new TrackCollection() ;
126
127 //---------------------------
128 // Primary Vertex
129 //---------------------------
130 const xAOD::VertexContainer* vertices = nullptr;
131 if ( !evtStore()->retrieve( vertices, m_vertexListInput ).isSuccess() ){ // retrieve arguments: container type, container key
132 ATH_MSG_WARNING("execute() Failed to retrieve Reconstructed vertex container. " << m_vertexListInput );
133 //ATH_MSG_ERROR(evtStore()->dump());
134 delete outputtracks;
135 return StatusCode::SUCCESS; //?? really
136 }
137
138 if( evtStore()->record( outputtracks, m_trackListOutput ).isFailure() ) {
139 ATH_MSG_ERROR( "Failed to record trackcollection with name " << m_trackListOutput );
140 }
141
142 xAOD::VertexContainer* theVertexContainer = new xAOD::VertexContainer;
143 xAOD::VertexAuxContainer* theVertexAuxContainer = new xAOD::VertexAuxContainer;
144 theVertexContainer->setStore( theVertexAuxContainer );
145
146 CHECK( evtStore()->record(theVertexContainer, m_outputVertexContainerName) );
147 CHECK( evtStore()->record(theVertexAuxContainer, m_outputVertexContainerName + "Aux.") );
148
149 const xAOD::Vertex* primaryVertex = nullptr;
150 for( const xAOD::Vertex* vertex: *vertices ) {
151 if( vertex->vertexType() == xAOD::VxType::PriVtx ) {
152 primaryVertex = vertex;
153 break;
154 }
155 }
156
157 std::vector< const Trk::TrackParameters* > trackParametersToFit;
158
159 // Check we have found a PV
160 if( primaryVertex ){
161 if(primaryVertex->nTrackParticles() < 2)
162 return StatusCode::SUCCESS;
163
164 //Loop over all tracks in the PV
165 for(unsigned int i =0; i< primaryVertex->nTrackParticles(); ++i ){
166 auto trackParticle = primaryVertex->trackParticle( i );
167 auto track = trackParticle->track();
168 // Check that there is Trk::Track
169 if(track){
170 //Refit Track using on the SCT hits
171 if(m_refitTracks){
172 auto newTrack = fitSCTOnlyTrack( track ) ;
173 if(newTrack){
174 outputtracks->push_back( newTrack ) ;
175 if(newTrack->perigeeParameters())
176 trackParametersToFit.push_back( newTrack->perigeeParameters() );
177 }
178 } else {
179 if(track->perigeeParameters())
180 trackParametersToFit.push_back( track->perigeeParameters() );
181 }
182 }
183 }
184 m_nTracksProcessed += primaryVertex->nTrackParticles() ;
185 m_nTracksAccepted += outputtracks->size() ;
186 // Make sure there are at least 2 tracks to fit in the vertex fit
187
188 ATH_MSG_DEBUG(" From: " << primaryVertex->nTrackParticles() << " tracks " << trackParametersToFit.size() << " will be used" );
189
190 if(trackParametersToFit.size() > 1){
191
192 Amg::Vector3D position = primaryVertex->position();
193 // Fit and store vertex
194 auto vertex = m_vertexFitter->fit(ctx, trackParametersToFit, position );
195 if(vertex){
196 ATH_MSG_DEBUG("Old Vtx " << primaryVertex->position());
197 ATH_MSG_DEBUG("New Vtx " << vertex->position());
198 vertex->setVertexType(xAOD::VxType::PriVtx);
199 theVertexContainer->push_back (std::move(vertex));
200 }
201 }
202 }
203
204
205
206
207
208 return StatusCode::SUCCESS;
209} // execute(.)
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define CHECK(...)
Evaluate an expression and check for errors.
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
Trk::Track * fitSCTOnlyTrack(const Trk::Track *track)
Strips of all TRT hits from the track and replaces them with a TRT momentum constraint from a the TRT...
size_t nTrackParticles() const
Get the number of tracks associated with this vertex.
const TrackParticle * trackParticle(size_t i) const
Get the pointer to a given track that was used in vertex reco.
const Amg::Vector3D & position() const
Returns the 3-pos.
Eigen::Matrix< double, 3, 1 > Vector3D
const Amg::Vector3D & position() const
Method to retrieve the position of the Intersection.
@ PriVtx
Primary vertex.
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ extraOutputDeps()

const DataObjIDColl & AthAlgorithm::extraOutputDeps ( ) const
overridevirtualinherited

Return the list of extra output dependencies.

This list is extended to include symlinks implied by inheritance relations.

Definition at line 50 of file AthAlgorithm.cxx.

51{
52 // If we didn't find any symlinks to add, just return the collection
53 // from the base class. Otherwise, return the extended collection.
54 if (!m_extendedExtraObjects.empty()) {
56 }
57 return Algorithm::extraOutputDeps();
58}
DataObjIDColl m_extendedExtraObjects

◆ filterPassed() [1/2]

bool AthAlgorithm::filterPassed ( ) const
inherited

Definition at line 94 of file AthAlgorithm.cxx.

94 {
95 return filterPassed( Gaudi::Hive::currentContext() );
96}
bool filterPassed() const

◆ filterPassed() [2/2]

bool AthAlgorithm::filterPassed ( const EventContext & ctx) const
inherited

Definition at line 98 of file AthAlgorithm.cxx.

98 {
99 return execState( ctx ).filterPassed();
100}

◆ finalize()

StatusCode RefitTracksAndVertex::finalize ( )

finalize method of this algorithm.

Prints out a summary of all events

Definition at line 118 of file RefitTracksAndVertex.cxx.

118 {
119 ATH_MSG_DEBUG( (*this) ) ;
120 return StatusCode::SUCCESS;
121}

◆ fitSCTOnlyTrack()

Trk::Track * RefitTracksAndVertex::fitSCTOnlyTrack ( const Trk::Track * track)
private

Strips of all TRT hits from the track and replaces them with a TRT momentum constraint from a the TRT segment belonging to the (extended) track (returns NULL in case of failure).

Definition at line 245 of file RefitTracksAndVertex.cxx.

245 {
246 Trk::MeasurementSet setSCT ;
247 Trk::MeasurementSet setPix ;
248 Trk::MeasurementSet setTRT ;
249
250 const Trk::Perigee* perTrk = track->perigeeParameters();
251 if( !perTrk ) {
252 ATH_MSG_WARNING("RefitTracksAndVertex() : No Perigee parameter on track!");
253 return nullptr ;
254 }
255
256 if( perTrk->eta() < m_selEtaMin || perTrk->eta() > m_selEtaMax)
257 return nullptr ;
258
259 if( perTrk->pT() < m_selPtMin)
260 return nullptr ;
261
262
263 //store all silicon measurements into the measurementset
264 DataVector<const Trk::MeasurementBase>::const_iterator it = track->measurementsOnTrack()->begin();
265 DataVector<const Trk::MeasurementBase>::const_iterator itEnd = track->measurementsOnTrack()->end();
266 for ( ; it!=itEnd; ++it){
267 if( !(*it) ) {
268 // log (MSG::WARNING) << "The MeasurementBase set has a void"
269 // << " member! Skip it.." << endmsg;
270 } else {
271 const Trk::RIO_OnTrack* rio = dynamic_cast <const Trk::RIO_OnTrack*>(*it);
272 if (rio != nullptr) {
273 const Identifier& surfaceID = (rio->identify()) ;
274 if( m_idHelper->is_sct(surfaceID) ) {
275 setSCT.push_back ( *it ) ;
276 } else if( m_idHelper->is_trt(surfaceID) ) {
277 setTRT.push_back ( *it ) ;
278 } else if( m_idHelper->is_pixel(surfaceID) ){
279 setPix.push_back ( *it ) ;
280 }
281 }
282 }
283 }
284
285 if( (int)setSCT.size() < m_selNHitSCTMin )
286 return nullptr;
287
288 ATH_MSG_DEBUG("RefitTracksAndVertex() : Found " << setSCT.size() << " SCT measurm's!" ) ;
289 ATH_MSG_DEBUG("RefitTracksAndVertex() : Found " << setPix.size() << " Pix measurm's!" ) ;
290 ATH_MSG_DEBUG("RefitTracksAndVertex() : Found " << setTRT.size() << " TRT measurm's!") ;
291
292
293 ATH_MSG_DEBUG( std::setiosflags( std::ios::scientific )) ;
294 ATH_MSG_DEBUG (std::setprecision( 7 )) ;
295 ATH_MSG_VERBOSE( std::setiosflags( std::ios::scientific )) ;
296 ATH_MSG_VERBOSE( std::setprecision( 7 )) ;
297
298 // now add z_0 and theta from original track as PseudoMeas to the TRT MeasurementSet
299 if(m_addPM){
300 const Trk::PseudoMeasurementOnTrack *pmFromSi = createPMfromSi( perTrk ) ;
301 if( !pmFromSi ) {
302 ATH_MSG_ERROR("RefitTracksAndVertex() : PseudoMeasurementOnTrack creation failed! " );
303 return nullptr ;
304 }
305 ATH_MSG_DEBUG( "RefitTracksAndVertex() : pmFromSi " << *pmFromSi) ;
306 Trk::MeasurementSet setSCT = addPM( setSCT, pmFromSi ) ;
307 }
308
309 // Add TRT hits as they do no harm
310 for( int i=0, i_max=setTRT.size() ; i!=i_max ; ++i ) {
311 setSCT.push_back( setTRT[i] ) ;
312 }
313
314 ATH_MSG_VERBOSE ( "RefitTracksAndVertex() : Si+PM MeasurementSet : " );
315 for( int i=0, i_max=setSCT.size() ; i!=i_max ; ++i ) {
316 ATH_MSG_VERBOSE ("============== i=" << i << " =============");
317 ATH_MSG_VERBOSE ( *(setSCT[i]));
318 }
319 ATH_MSG_VERBOSE ("==========================================");
320
321 // fit TRT part of the track with PseudoMeas on z_0, theta
322 Trk::Track* trkSCT=(m_trackFitter->fit(Gaudi::Hive::currentContext()
323 ,setSCT
324 , *perTrk
325 , true
326 , Trk::pion
327 //, Trk::muon
328 //, Trk::nonInteracting
329 )).release() ;
330 if( !trkSCT ) {
331 ATH_MSG_DEBUG( "RefitTracksAndVertex() : Fit of SCT part of the track failed! " ) ;
332 return nullptr ;
333 }
334 const Trk::Perigee* perSCT = trkSCT->perigeeParameters();
335 ATH_MSG_VERBOSE( "RefitTracksAndVertex() : perSCT " << *perSCT) ;
336
337 return trkSCT;
338}
#define ATH_MSG_VERBOSE(x)
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
const Trk::MeasurementSet addPM(Trk::MeasurementSet &ms, const Trk::PseudoMeasurementOnTrack *pm)
adds a PseudoMeasurement to a MeasurementSet
const Trk::PseudoMeasurementOnTrack * createPMfromSi(const Trk::Perigee *mp)
creates a PseudoMeasurement with (z0, theta) from extended track perigee parameters
double eta() const
Access method for pseudorapidity - from momentum.
double pT() const
Access method for transverse momentum.
Identifier identify() const
return the identifier -extends MeasurementBase
const Perigee * perigeeParameters() const
return Perigee.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee

◆ getContext()

const EventContext & AthAlgorithm::getContext ( ) const
inherited

Deprecated methods (use the ones with EventContext).

Definition at line 90 of file AthAlgorithm.cxx.

90 {
91 return Gaudi::Hive::currentContext();
92}

◆ initialize()

StatusCode RefitTracksAndVertex::initialize ( )

initialize method of this algorithm.

Definition at line 70 of file RefitTracksAndVertex.cxx.

70 {
71 //retrieve the DetectorStore service
72 StatusCode status=detStore().retrieve() ;
73 if( status.isFailure() ) {
74 ATH_MSG_FATAL( "DetectorStore service not found !" );
75 return status;
76 }else{
77 ATH_MSG_DEBUG( "DetectorStore retrieved!" );
78 }
79
80 if (detStore()->retrieve(m_idHelper, "AtlasID").isFailure()) {
81 ATH_MSG_FATAL( "Could not get AtlasDetectorID helper" );
82 return StatusCode::FAILURE;
83 }
84
85 // Set up fitter
86 status = m_trackFitter.retrieve();
87 if (status.isFailure()) {
88 ATH_MSG_FATAL( "Could not find tool " << m_trackFitter << ". Exiting.");
89 return status ;
90 } else {
91 ATH_MSG_DEBUG( " Got " << m_trackFitter << " as TrackFitter. " ) ;
92 }
93
94 // Set up fitter
95 status = m_vertexFitter.retrieve();
96 if (status.isFailure()) {
97 ATH_MSG_FATAL( "Could not find tool " << m_vertexFitter << ". Exiting.");
98 return status ;
99 } else {
100 ATH_MSG_DEBUG( " Got " << m_vertexFitter << " as VertexFitter. " ) ;
101 }
102
103 // Print input properties
104 if( m_applyTrkSel ) {
106 "Track selection will be applied:"
107 << "\n " << m_selEtaMin << " < eta < " << m_selEtaMax
108 << "\n pT > " << m_selPtMin/Gaudi::Units::GeV << " GeV"
109 << "\n number_of_PIXEL_hits >= " << m_selNHitPIXMin
110 << "\n number_of_SCT_hits >= " << m_selNHitSCTMin
111 ) ;
112 } else {
113 ATH_MSG_DEBUG( "NO selection will be applied to tracks" ) ;
114 }
115 return StatusCode::SUCCESS;
116} // initialize(.)
#define ATH_MSG_FATAL(x)
const ServiceHandle< StoreGateSvc > & detStore() const
::StatusCode StatusCode
StatusCode definition for legacy code.
status
Definition merge.py:16

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ isReEntrant()

virtual bool AthAlgorithm::isReEntrant ( ) const
inlinefinaloverrideprotectedvirtualinherited

Legacy algorithms are not thread-safe.

Definition at line 111 of file AthAlgorithm.h.

111{ return false; }

◆ msg()

MsgStream & AthCommonMsg< Gaudi::Algorithm >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msgLvl()

bool AthCommonMsg< Gaudi::Algorithm >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ setFilterPassed() [1/2]

void AthAlgorithm::setFilterPassed ( bool state) const
inherited

Definition at line 102 of file AthAlgorithm.cxx.

102 {
103 setFilterPassed( state, Gaudi::Hive::currentContext() );
104}
void setFilterPassed(bool state) const

◆ setFilterPassed() [2/2]

void AthAlgorithm::setFilterPassed ( bool state,
const EventContext & ctx ) const
inherited

Definition at line 106 of file AthAlgorithm.cxx.

106 {
107 execState( ctx ).setFilterPassed(state);
108}

◆ sysInitialize()

StatusCode AthAlgorithm::sysInitialize ( )
overridevirtualinherited

Override sysInitialize.

Override sysInitialize from the base class.

Loop through all output handles, and if they're WriteCondHandles, automatically register them and this Algorithm with the CondSvc.

Scan through all outputHandles, and if they're WriteCondHandles, register them with the CondSvc

Reimplemented from AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >.

Reimplemented in AthAnalysisAlgorithm, AthFilterAlgorithm, AthHistogramAlgorithm, and PyAthena::Alg.

Definition at line 66 of file AthAlgorithm.cxx.

66 {
68
69 if (sc.isFailure()) {
70 return sc;
71 }
72 ServiceHandle<ICondSvc> cs("CondSvc",name());
73 for (auto h : outputHandles()) {
74 if (h->isCondition() && h->mode() == Gaudi::DataHandle::Writer) {
75 // do this inside the loop so we don't create the CondSvc until needed
76 if ( cs.retrieve().isFailure() ) {
77 ATH_MSG_WARNING("no CondSvc found: won't autoreg WriteCondHandles");
78 return StatusCode::SUCCESS;
79 }
80 if (cs->regHandle(this,*h).isFailure()) {
81 sc = StatusCode::FAILURE;
82 ATH_MSG_ERROR("unable to register WriteCondHandle " << h->fullKey()
83 << " with CondSvc");
84 }
85 }
86 }
87 return sc;
88}
static Double_t sc
virtual StatusCode sysInitialize() override
Override sysInitialize.
AthCommonDataStore(const std::string &name, T... args)
virtual std::vector< Gaudi::DataHandle * > outputHandles() const override

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }

Member Data Documentation

◆ m_addPM

bool RefitTracksAndVertex::m_addPM
private

apply a pseudo measurement based on the original track (theta,z0)

Definition at line 73 of file RefitTracksAndVertex.h.

◆ m_applyTrkSel

bool RefitTracksAndVertex::m_applyTrkSel
private

apply a selection on tracks or not

Definition at line 71 of file RefitTracksAndVertex.h.

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default).

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default).

Definition at line 390 of file AthCommonDataStore.h.

◆ m_extendedExtraObjects

DataObjIDColl AthAlgorithm::m_extendedExtraObjects
privateinherited

Definition at line 114 of file AthAlgorithm.h.

◆ m_idHelper

const AtlasDetectorID* RefitTracksAndVertex::m_idHelper
private

Detector ID helper.

Definition at line 60 of file RefitTracksAndVertex.h.

◆ m_nRejectPIX

size_t RefitTracksAndVertex::m_nRejectPIX
private

Counter for number of tracks failing the min number of PIX hits req.

Definition at line 82 of file RefitTracksAndVertex.h.

◆ m_nRejectPM

size_t RefitTracksAndVertex::m_nRejectPM
private

Counter for number of tracks failing the addition of a pseudo-measurement (PM).

Definition at line 85 of file RefitTracksAndVertex.h.

◆ m_nRejectSCT

size_t RefitTracksAndVertex::m_nRejectSCT
private

Counter for number of tracks failing the min number of SCT hits req.

Definition at line 83 of file RefitTracksAndVertex.h.

◆ m_nRejectTRT

size_t RefitTracksAndVertex::m_nRejectTRT
private

Counter for number of tracks failing the min number of TRT hits req.

Definition at line 84 of file RefitTracksAndVertex.h.

◆ m_nTracksAccepted

size_t RefitTracksAndVertex::m_nTracksAccepted
private

Counter for number of tracks passing the preselection and with PM.

Definition at line 81 of file RefitTracksAndVertex.h.

◆ m_nTracksPresel

size_t RefitTracksAndVertex::m_nTracksPresel
private

Counter for number of tracks passing the preselection.

Definition at line 80 of file RefitTracksAndVertex.h.

◆ m_nTracksProcessed

size_t RefitTracksAndVertex::m_nTracksProcessed
private

Counter for number of tracks processed.

Definition at line 79 of file RefitTracksAndVertex.h.

◆ m_outputVertexContainerName

std::string RefitTracksAndVertex::m_outputVertexContainerName
private

Name of vertex container.

Definition at line 66 of file RefitTracksAndVertex.h.

◆ m_refitTracks

bool RefitTracksAndVertex::m_refitTracks
private

refitTracks

Definition at line 72 of file RefitTracksAndVertex.h.

◆ m_selEtaMax

double RefitTracksAndVertex::m_selEtaMax
private

maximal eta cut value for the TrackSelection

Definition at line 76 of file RefitTracksAndVertex.h.

◆ m_selEtaMin

double RefitTracksAndVertex::m_selEtaMin
private

minimal eta cut value for the TrackSelection

Definition at line 75 of file RefitTracksAndVertex.h.

◆ m_selNHitPIXMin

int RefitTracksAndVertex::m_selNHitPIXMin
private

minimal number of PIX hits cut value for the TrackSelection

Definition at line 77 of file RefitTracksAndVertex.h.

◆ m_selNHitSCTMin

int RefitTracksAndVertex::m_selNHitSCTMin
private

minimal number of SCT hits cut value for the TrackSelection

Definition at line 78 of file RefitTracksAndVertex.h.

◆ m_selPtMin

double RefitTracksAndVertex::m_selPtMin
private

minimal pT cut value for the TrackSelection

Definition at line 74 of file RefitTracksAndVertex.h.

◆ m_trackFitter

ToolHandle<Trk::ITrackFitter> RefitTracksAndVertex::m_trackFitter
private

The TrackFitter to refit the tracks (segment, momentum constraint).

Definition at line 67 of file RefitTracksAndVertex.h.

◆ m_trackListOutput

std::string RefitTracksAndVertex::m_trackListOutput
private

Name of the TrackCollection (Output).

Definition at line 65 of file RefitTracksAndVertex.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexFitter

ToolHandle<Trk::IVertexFitter> RefitTracksAndVertex::m_vertexFitter
private

The TrackFitter to refit the tracks (segment, momentum constraint).

Definition at line 68 of file RefitTracksAndVertex.h.

◆ m_vertexListInput

std::string RefitTracksAndVertex::m_vertexListInput
private

Name of the TrackCollection (input).

Definition at line 61 of file RefitTracksAndVertex.h.

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< Gaudi::Algorithm > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


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