ATLAS Offline Software
Loading...
Searching...
No Matches
VertexRetriever.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include <string>
7#include <vector>
8#include <map>
9#include <cmath>
10
13#include "VxVertex/RecVertex.h"
15
16#include "TrkTrack/TrackInfo.h"
19
20//#include "TrkVertexAnalysisUtils/V0Tools.h"
21
26
27
28namespace JiveXML {
29
36 VertexRetriever::VertexRetriever(const std::string& type,const std::string& name,const IInterface* parent):
37 AthAlgTool(type,name,parent),
38 m_typeName("RVx"){
39
40 //Declare the interface
41 declareInterface<IDataRetriever>(this);
42
43 //Declare AlgTool properties
44 declareProperty("PrimaryVertexCollection", m_primaryVertexKey = "VxPrimaryCandidate", "Vertices to use as primary vertex");
45 declareProperty("SecondaryVertexCollection", m_secondaryVertexKey = "SecVertices", "Vertices to use as secondary vertex");
46 declareProperty("ConversionVertexCollection", m_conversionVertexKey = "ConversionCandidate", "Vertices to use as conversion vertex");
47 declareProperty("DoWritePrimAndSecVertexOnly", m_doWritePrimAndSecVertexOnly = false,
48 "if true only write primary and secondary vertex, placeholder to be removed");
49 declareProperty("DoWriteHLT", m_doWriteHLT = false, "whether to write HLTAutoKey objects");
50 declareProperty("TrackCollection", m_trackCollection = "Tracks");
51 declareProperty("Chi2OverDOFCut", m_chi2Cut = 10.); // <10 is 'tight', <100 is 'loose'
52 }
53
54
60
62 const Trk::Perigee* per2, bool& matched){
63
64 double diff1 = ( per1->parameters()[Trk::d0] -
65 per2->parameters()[Trk::d0] );
66 double diff2 = ( per1->parameters()[Trk::z0] -
67 per2->parameters()[Trk::z0] );
68 double diff3 = ( per1->parameters()[Trk::phi0] -
69 per2->parameters()[Trk::phi0] );
70 double diff4 = ( per1->charge() -
71 per2->charge() );
72 double diff5 = ( per1->pT() -
73 per2->pT() );
74
75 matched = diff1+diff2+diff3+diff4+diff5 == 0.;
76 }
77
83
85
87// from:
88// http://atlas-sw.cern.ch/cgi-bin/viewcvs-all.cgi/users/jsvirzi/PhysicsAnalysis/Atlas_ESD_Analysis/src/Atlas_ESD_Analysis.cxx?root=atlas&view=co
90// InnerDetector/InDetMonitoring/InDetAlignmentMonitoring/TrackSplitterTool.cxx
92// const Trk::Perigee* trackPerigee = inputTrack->perigeeParameters();
93// const Trk::perigee* perigee = dynamic_cast<const Trk::perigee*>( trackPerigee );
94// double m_d0 = perigee->parameters()[Trk::d0];
95
96 std::string searchStr = "TrackParticle";
97 size_t found=m_trackCollection.find(searchStr);
98
99 m_perigeeVector.clear(); // need to clear, otherwise accumulates over events
100 if (found!=std::string::npos){ // User selected a Rec::TrackParticle Collection
102 if (!tracks.isValid()){
103 ATH_MSG_WARNING("Unable to retrieve track collection" << m_trackCollection << " for association ");
104 } else {
105 ATH_MSG_DEBUG("Retrieved " << m_trackCollection);
106
107 //xAOD::VertexContainer::const_iterator VertexItr = cont->begin();
108 //for ( ; VertexItr != cont->end(); ++VertexItr) {
109 for(const auto track : *tracks) {
110 const Trk::Perigee *perigee = track->perigee();
111 //if(perigee == 0) continue; // not skip ! need to keep order for index !
112 m_perigeeVector.push_back( perigee ); // this perigee match works !
113 }
114 }
115 }
116
117 else{ // User selected a Trk::Track collection
118 ATH_MSG_DEBUG(" User selected a Trk::Track collection ! ");
120 if (!tracks.isValid()){
121 ATH_MSG_WARNING("Unable to retrieve track collection" << m_trackCollection << " for association ");
122 } else {
123 ATH_MSG_DEBUG("Retrieved " << m_trackCollection);
124 for(const auto track : *tracks) {
125 const Trk::Perigee* trackPerigee = track->perigeeParameters();
126 //const Trk::Perigee *perigee = dynamic_cast<const Trk::Perigee*>( trackPerigee );
127 //if(perigee == 0) continue; // not skip ! need to keep order for index !
128 //m_perigeeVector.push_back( perigee );
129 m_perigeeVector.push_back( trackPerigee );
130 }
131 }
132 }
133 return StatusCode::SUCCESS;
134 }
135
143 StatusCode VertexRetriever::retrieve(ToolHandle<IFormatTool> &FormatTool){
144
145 ATH_MSG_DEBUG("In retrieve()");
146
147 //Get an iterator over all vertex collections,
148 //return if there are none
149 SG::ConstIterator<VxContainer> vtxCollectionItr, vtxCollectionsEnd;
150 if (evtStore()->retrieve(vtxCollectionItr,vtxCollectionsEnd).isFailure()) {
151 if (msgLvl(MSG::DEBUG )) msg(MSG::DEBUG ) << "No VxContainer containers found in this event" << endmsg;
152 return StatusCode::SUCCESS;
153 }
154
155 //See if we can find the requested secondary vertex collection
156 const VxContainer* secondaryVtxCollection;
157 if (evtStore()->retrieve(secondaryVtxCollection,m_secondaryVertexKey).isFailure()) {
158 if (msgLvl(MSG::DEBUG )) msg(MSG::DEBUG ) << "No Secondary vertex container found at SecVertices" << endmsg;
159 }else{
160 if (msgLvl(MSG::DEBUG )) msg(MSG::DEBUG ) << "Secondary vertex container size: " << secondaryVtxCollection->size() << endmsg;
161 }
162
163 //See if we can find the requested primary vertex collection
164 const VxContainer* primaryVtxCollection;
165 if ( evtStore()->retrieve(primaryVtxCollection,m_primaryVertexKey).isFailure()) {
166 if (msgLvl(MSG::WARNING)) msg(MSG::WARNING) << "Primary vertex container "
167 << m_primaryVertexKey << " not found" << endmsg;
168 }
169
170 //Declare all the data vectors we want to retrieve
171 DataVect x;
172 DataVect y;
173 DataVect z;
174 DataVect primVxCand;
176 DataVect sgkey;
177 DataVect covMatrix;
178 DataVect numTracks;
179 DataVect tracks;
180 DataVect vertexType;
181
182 //Loop over all vertex containers
183 for ( ; vtxCollectionItr != vtxCollectionsEnd; ++vtxCollectionItr ) {
184
185 //Check whether we should ignore HLTAutoKey collections
186 if ( (!m_doWriteHLT) && ( vtxCollectionItr.key().find("HLT") != std::string::npos)){
187 if (msgLvl(MSG::DEBUG)) msg(MSG::DEBUG) << "Ignoring HLT collection " << vtxCollectionItr.key() << endmsg;
188 continue;
189 }
190
191 //Get size of current container
192 VxContainer::size_type NVtx = vtxCollectionItr->size();
193
194 ATH_MSG_DEBUG("Reading vertex container " << vtxCollectionItr.key()
195 << " with " << NVtx << " entries");
196
197 //Declare all the data vectors we want to retrieve and reserve space
198 x.reserve(x.size()+NVtx);
199 y.reserve(y.size()+NVtx);
200 z.reserve(z.size()+NVtx);
201 primVxCand.reserve(primVxCand.size()+NVtx);
202 chi2.reserve(chi2.size()+NVtx);
203 sgkey.reserve(sgkey.size()+NVtx);
204 covMatrix.reserve(covMatrix.size()+NVtx);
205 numTracks.reserve(numTracks.size()+NVtx);
206 vertexType.reserve(numTracks.size()+NVtx);
207
208 StatusCode sc = fillPerigeeList();
209 if (!sc.isFailure()) {
210 ATH_MSG_DEBUG("Perigee list filled with " << m_perigeeVector.size() << " entries ");
211 }
212
213 //Loop over vertices
214 VxContainer::const_iterator vertexItr = vtxCollectionItr->begin();
215 for ( ; vertexItr != vtxCollectionItr->end(); ++vertexItr) {
217// if ( m_doWritePrimAndSecVertexOnly && vtxCollectionItr.key() == m_conversionVertexKey ){
218// if (msgLvl(MSG::DEBUG)){ msg(MSG::DEBUG) << " DoWritePrimAndSecVertexOnly switch true - skipping "
219// << vtxCollectionItr.key() << endmsg; }
220// continue; }
222
226 //Get fit quality object
227 Trk::FitQuality fitQuality = (*vertexItr)->recVertex().fitQuality();
228
229 //degrees of freedom might be zero - beware
230 float this_chi2 = -1;
231 if ( fitQuality.doubleNumberDoF() != 0 ){
232 this_chi2 = fitQuality.chiSquared()/fitQuality.doubleNumberDoF();
233 }
234 //Cut: on Chi^2 over NumberOfDegreesOfFreedom only for ConversionCandidate
235 if ( this_chi2 > m_chi2Cut && ( vtxCollectionItr.key() == m_conversionVertexKey )) continue;
236
237 float this_x = (*vertexItr)->recVertex().position().x()/10.; //Atlantis units are cm
238 float this_y = (*vertexItr)->recVertex().position().y()/10.;
239 float this_z = (*vertexItr)->recVertex().position().z()/10.;
240 float R = std::hypot (this_x, this_y); // distance from beamline
241
242 ATH_MSG_DEBUG(" Collection: " << vtxCollectionItr.key()
243 << ", this_chi2: " << this_chi2 << " - chi2: " << fitQuality.chiSquared()
244 //<< " ," << (*vertexItr)->recVertex().position().x()/10. << " ," << (*vertexItr)->recVertex().position().x()*CLHEP::cm
245 << ", R: " << R);
246
247 chi2.emplace_back( this_chi2 );
248 x.emplace_back( this_x );
249 y.emplace_back( this_y );
250 z.emplace_back( this_z );
251
252 // from: Tracking/TrkEvent/TrkEventPrimitives/VertexType.h
253 const Trk::VertexType vtx_type = (*vertexItr)->vertexType();
254 vertexType.emplace_back( vtx_type );
255 ATH_MSG_DEBUG(" collection " << vtxCollectionItr.key() << ": VertexType: " << vtx_type);
256
257 //Store primary vertex candidate flag
258 if ( &(*vtxCollectionItr) == primaryVtxCollection ){
259 if ( Trk::PriVtx == vtx_type ){ primVxCand.emplace_back( 1 ); // type 1 'real' primary vertex
260 }else{ primVxCand.emplace_back( 0 ); } // hack ! 'type 3 pileup' Should properly use 'vertexType'
261 }else if ( &(*vtxCollectionItr) == secondaryVtxCollection ){
262 primVxCand.emplace_back( 2 ); // normally those are 'type 9 Kshort'
263 }else{
264 primVxCand.emplace_back( 0 );
265 }
266
267// !!! This doesn't work after Migration. 'errorposition' not known anymore ? jpt Dec'13 !!!
269 //double scale = (CLHEP::cm*CLHEP::cm)/(CLHEP::mm*CLHEP::mm); // in order to convert the covMatrix elements from (mm)^2 to (cm)^2 divide them by 100
270 //double precision = 10000; //// otherwise you sometimes loose precision by putting only 6 decimal points into the xml
271 //for (int i = 0;i<3; i++)
272 //for (int j = 0; j<i+1; j++ )
273 //covMatrix.push_back (DataType((*vertexItr)->recVertex().errorPosition().covariance()[i][j]/scale*precision ));
274//-----------
279 covMatrix.emplace_back("2 -.1 .5 -.01 0.002 .01");
280
285
286//---- association from:
287//---- from http://alxr.usatlas.bnl.gov/lxr/source/atlas/Reconstruction/tauRec/src/PhotonConversionPID.cxx
288
289 const std::vector<Trk::VxTrackAtVertex*>* trklist = (*vertexItr)->vxTrackAtVertex();
290
291 ATH_MSG_DEBUG("Tracks at vertex: " << trklist->size());
292
293 numTracks.emplace_back( trklist->size() );
294 sgkey.emplace_back( m_trackCollection ); // sgkey in current scheme is _not_ a multiple !
295
296 tracks.emplace_back( -1 );
297 }
298 } // loop over vertex containers
299
300 //Finally add all retrieved data to the data map
301 DataMap dataMap;
302 dataMap["x"] = x;
303 dataMap["y"] = y;
304 dataMap["z"] = z;
305 dataMap["primVxCand"] = primVxCand;
306 dataMap["chi2"] = chi2;
307 dataMap["covMatrix multiple=\"6\""] = covMatrix;
308 dataMap["numTracks"] = numTracks;
309 dataMap["sgkey"] = sgkey;
310 dataMap["vertexType"] = vertexType;
311
312 //If there had been any tracks, add a tag
313 if (!numTracks.empty()){
314 //Calculate average number of tracks per vertex
315 double NTracksPerVertex = tracks.size()*1./numTracks.size();
316 std::string tag = "tracks multiple=\"" +DataType(NTracksPerVertex).toString()+"\"";
317 dataMap[tag] = tracks;
319// std::string tag2 = "sgkey multiple=\"" +DataType(NTracksPerVertex).toString()+"\"";
320// dataMap[tag2] = sgkey;
321 }
322
326 //return FormatTool->AddToEvent(dataTypeName(), "PrimSecConvRecVx", &dataMap);
327 return FormatTool->AddToEvent(dataTypeName(), "", &dataMap);
328
329 } // end retrieve
330
331} // namepspace JiveXML
#define endmsg
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
OFFLINE_FRAGMENTS_NAMESPACE::PointerType DataType
static Double_t sc
#define y
#define x
#define z
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ServiceHandle< StoreGateSvc > & evtStore()
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
typename DataVectorBase< Trk::VxCandidate >::Base::size_type size_type
Definition DataVector.h:814
size_type size() const noexcept
Returns the number of elements in the collection.
virtual std::string dataTypeName() const
Return the name of the data type.
std::string m_trackCollection
StoreGate key for track collection for association.
std::vector< const Trk::Perigee * > m_perigeeVector
std::string m_conversionVertexKey
StoreGate key for conversion candidate collection.
virtual StatusCode fillPerigeeList()
Retrieve measured perigee, automatically switch between Trk::Track and Rec::TrackParticle depending o...
std::string m_secondaryVertexKey
StoreGate key for secondary vertex candidate collection.
const std::string m_typeName
The data type that is generated by this retriever.
virtual StatusCode retrieve(ToolHandle< IFormatTool > &FormatTool)
Retrieve all the data.
VertexRetriever(const std::string &t, const std::string &n, const IInterface *p)
Standard constructor.
float m_chi2Cut
Chi^2 over NumberOfDegreesOfFreedom cut.
bool m_doWriteHLT
wether to write HLTAutoKey objects
std::string m_primaryVertexKey
StoreGate key for primary vertex candidate collection.
bool m_doWritePrimAndSecVertexOnly
write primary and secondary vertizes only - placeholder, to be removed
a const_iterator facade to DataHandle.
Definition SGIterator.h:164
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const std::string & key() const
Get the key string with which the current object was stored.
double chiSquared() const
returns the of the overall track fit
Definition FitQuality.h:56
double doubleNumberDoF() const
returns the number of degrees of freedom of the overall track or vertex fit as double
Definition FitQuality.h:68
Class to represent and store fit qualities from track reconstruction in terms of and number of degre...
Definition FitQuality.h:97
double charge() const
Returns the charge.
double pT() const
Access method for transverse momentum.
double chi2(TH1 *h0, TH1 *h1)
This header is shared inbetween the C-style server thread and the C++ Athena ServerSvc.
void manualPerigeeMatch(const Trk::Perigee *per1, const Trk::Perigee *per2, bool &matched)
Manually match measured perigee d0,z0,phi0 of Trk::Track and Rec::TrackParticle, some problem with di...
std::map< std::string, DataVect > DataMap
Definition DataType.h:59
std::vector< DataType > DataVect
Defines a map with a key and a vector of DataType objects e.g.
Definition DataType.h:58
VertexType
This file defines the enums in the Trk namespace for the different vertex types.
Definition VertexType.h:25
@ PriVtx
Primary Vertex.
Definition VertexType.h:27
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ phi0
Definition ParamDefs.h:65
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64