ATLAS Offline Software
Loading...
Searching...
No Matches
xAODVertexRetriever.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
11
12#include "AthenaKernel/Units.h"
13using Athena::Units::cm;
14
15namespace JiveXML {
16
24 xAODVertexRetriever::xAODVertexRetriever(const std::string& type,const std::string& name,const IInterface* parent):
25 AthAlgTool(type,name,parent){}
26
27
29 ATH_CHECK(m_keys.initialize());
30 return StatusCode::SUCCESS;
31 }
32
33
38 StatusCode xAODVertexRetriever::retrieve(ToolHandle<IFormatTool> &FormatTool) {
39
40 ATH_MSG_DEBUG("In retrieve()");
41
42 // Loop through the keys and retrieve the corresponding data
43 for (const auto& key : m_keys) {
45 if (cont.isValid()) {
46 DataMap data = getData(&(*cont));
47 if (FormatTool->AddToEvent(dataTypeName(), key.key() + "_xAOD", &data).isFailure()) {
48 ATH_MSG_WARNING("Failed to add collection " << key.key());
49 } else {
50 ATH_MSG_DEBUG(" (" << key.key() << ") retrieved");
51 }
52 } else {
53 ATH_MSG_WARNING("Collection " << key.key() << " not found in SG");
54 }
55 }
56
57 return StatusCode::SUCCESS;
58 }
59
65
66 ATH_MSG_DEBUG("in getData()");
67
68 float chi2val = 0.;
69
71
72 DataVect x;
73 DataVect y;
74 DataVect z;
76 DataVect vertexType;
77 DataVect primVxCand;
78 DataVect covMatrix;
79 DataVect numTracks;
80 DataVect tracks;
81 DataVect sgkey;
82
83 //Get size of current container
85
86 x.reserve(x.size()+NVtx);
87 y.reserve(y.size()+NVtx);
88 z.reserve(z.size()+NVtx);
89 chi2.reserve(chi2.size()+NVtx);
90 vertexType.reserve(vertexType.size()+NVtx);
91 primVxCand.reserve(primVxCand.size()+NVtx);
92 covMatrix.reserve(covMatrix.size()+NVtx);
93 numTracks.reserve(numTracks.size()+NVtx);
94 tracks.reserve(tracks.size()+NVtx);
95 sgkey.reserve(sgkey.size()+NVtx);
96
97 int counter = 0;
98
99 //Loop over vertices
101 for ( ; VertexItr != cont->end(); ++VertexItr) {
102
103 ATH_MSG_DEBUG(" Vertex #" << counter++ << " : x = " << (*VertexItr)->x()/cm << ", y = "
104 << (*VertexItr)->y()/cm << ", z[GeV] = " << (*VertexItr)->z()/cm
105 << ", vertexType = " << (*VertexItr)->vertexType()
106 << ", chiSquared = " << (*VertexItr)->chiSquared()
107 << ", numberDoF = " << (*VertexItr)->numberDoF());
108
109 x.emplace_back(DataType((*VertexItr)->x()/cm));
110 y.emplace_back(DataType((*VertexItr)->y()/cm));
111 z.emplace_back(DataType((*VertexItr)->z()/cm));
112
113 vertexType.emplace_back( DataType((*VertexItr)->vertexType()));
114
115 if ((*VertexItr)->vertexType() == 1 ){
116 primVxCand.emplace_back( 1 );
117 }else{
118 primVxCand.emplace_back( 0 );
119 }
120 sgkey.emplace_back (m_tracksName.value());
121
123 covMatrix.emplace_back(DataType("2 -.1 .5 -.01 0.002 .01"));
124
125 //degrees of freedom might be zero - beware
126 if ( (*VertexItr)->numberDoF() != 0 ){
127 chi2val = (*VertexItr)->chiSquared()/(*VertexItr)->numberDoF() ;
128 }else{
129 chi2val = -1.;
130 }
131 chi2.emplace_back(DataType( chi2val ));
132
133 // track-vertex association code in xAOD from Nick Styles, Apr14:
134 // InnerDetector/InDetRecAlgs/InDetPriVxFinder/InDetVxLinksToTrackParticles
135
136 int trkCnt = 0;
137 const std::vector< ElementLink< xAOD::TrackParticleContainer > > tpLinks = (*VertexItr)->trackParticleLinks();
138
139 //iterating over the links
140 unsigned int tp_size = tpLinks.size();
141 numTracks.emplace_back(DataType( tp_size ));
142 if(tp_size){ // links exist
143 for(unsigned int tp = 0; tp<tp_size; ++tp)
144 {
146
147 ATH_MSG_DEBUG(" Vertex #" << counter << " track association index: " << tpl.index()
148 << ", collection : " << tpl.key()
149 << ", Tracks : " << tp << " out of " << tp_size << ", own count: " << trkCnt++);
150
151 if ( tpl.index() < 1000 ){ // sanity check, this can be huge number
152 tracks.emplace_back(DataType( tpl.index() ));
153 }else{
154 tracks.emplace_back(DataType( 0 ));
155 }
156 } //links exist
157 }//end of track particle collection size check
158
159 } // end VertexIterator
160
161 // four-vectors
162 DataMap["x"] = x;
163 DataMap["y"] = y;
164 DataMap["z"] = z;
165 DataMap["chi2"] = chi2;
166 DataMap["vertexType"] = vertexType;
167 DataMap["primVxCand"] = primVxCand;
168 DataMap["covMatrix multiple=\"6\""] = covMatrix;
169 DataMap["numTracks"] = numTracks;
170 DataMap["sgkey"] = sgkey;
171
172 //This is needed once we know numTracks and associations:
173 //If there had been any tracks, add a tag
174 if ((numTracks.size()) != 0){
175 //Calculate average number of tracks per vertex
176 double NTracksPerVertex = tracks.size()*1./numTracks.size();
177 std::string tag = "tracks multiple=\"" +DataType(NTracksPerVertex).toString()+"\"";
178 DataMap[tag] = tracks;
179 }
180
181 ATH_MSG_DEBUG(dataTypeName() << " retrieved with " << x.size() << " entries");
182
183 return DataMap;
184
185 }
186
187
188} // JiveXML namespace
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
OFFLINE_FRAGMENTS_NAMESPACE::PointerType DataType
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
Wrapper to avoid constant divisions when using units.
#define y
#define x
#define z
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
DataModel_detail::const_iterator< DataVector > const_iterator
Definition DataVector.h:838
typename DataVectorBase< xAOD::Vertex_v1 >::Base::size_type size_type
Definition DataVector.h:814
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
size_type size() const noexcept
Returns the number of elements in the collection.
xAODVertexRetriever(const std::string &type, const std::string &name, const IInterface *parent)
Standard Constructor.
const DataMap getData(const xAOD::VertexContainer *)
Puts the variables into a DataMap.
virtual StatusCode retrieve(ToolHandle< IFormatTool > &FormatTool)
For each Vertex collections retrieve basic parameters.
virtual std::string dataTypeName() const
Return the name of the data type that is generated by this retriever.
Gaudi::Property< std::string > m_tracksName
SG::ReadHandleKeyArray< xAOD::VertexContainer > m_keys
virtual bool isValid() override final
Can the handle be successfully dereferenced?
double chi2(TH1 *h0, TH1 *h1)
This header is shared inbetween the C-style server thread and the C++ Athena ServerSvc.
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
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".