ATLAS Offline Software
Loading...
Searching...
No Matches
VertexParametersPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
9
12#include "../TrackParametersHelper.h" // also includes VertexParametersHelper.h
13
14
19 PlotMgr* pParent,
20 const std::string& dirName,
21 const std::string& anaTag,
22 const std::string& vertexType,
23 bool doTrackPlots,
24 bool doGlobalPlots,
25 bool doTruthMuPlots ) :
26 PlotMgr( dirName, anaTag, pParent ),
28 m_doTrackPlots( doTrackPlots ),
29 m_doGlobalPlots( doGlobalPlots ),
30 m_doTruthMuPlots( doTruthMuPlots ) { }
31
32
37{
38 StatusCode sc = bookPlots();
39 if( sc.isFailure() ) {
40 ATH_MSG_ERROR( "Failed to book vertex parameters plots" );
41 }
42}
43
44
46{
47 ATH_MSG_DEBUG( "Booking vertex parameters plots in " << getDirectory() );
48
49 if( not m_doGlobalPlots ) {
55 if( m_vertexType != "truth" ) {
63 }
64 if( m_doTrackPlots ) {
76 }
77 } else {
79 ATH_CHECK( retrieveAndBook( m_nVtx_vs_actualMu_2D, m_vertexType+"_nVtx_vs_actualMu_2D" ) );
81 if( m_doTruthMuPlots ) {
84 }
85 }
86
87 return StatusCode::SUCCESS;
88}
89
90
95template< typename VERTEX, typename PARTICLE >
97 const VERTEX& vertex,
98 const std::vector< const PARTICLE* >& associatedTracks,
99 const std::vector< float >& associatedTrackWeights,
100 float weight )
101{
102 if( m_doGlobalPlots ) return StatusCode::SUCCESS;
103
105 float vx = posX( vertex );
106 float vy = posY( vertex );
107 float vz = posZ( vertex );
108 float vtime = time( vertex );
109 bool vhasValidTime = bool( hasValidTime( vertex ) );
110
111 ATH_CHECK( fill( m_vtx_x, vx, weight ) );
112 ATH_CHECK( fill( m_vtx_y, vy, weight ) );
113 ATH_CHECK( fill( m_vtx_z, vz, weight ) );
114 if( vhasValidTime ) ATH_CHECK( fill( m_vtx_time, vtime, weight ) );
115
116 if( m_vertexType != "truth" ) {
118 float vxErr = error( vertex, VposX );
119 float vyErr = error( vertex, VposY );
120 float vzErr = error( vertex, VposZ );
121 float vtimeErr = timeErr( vertex );
122 float vchi2 = chiSquared( vertex );
123 float vndof = ndof( vertex );
124 float vchi2OverNdof = ( vndof > 0 ) ? vchi2 / vndof : -9999.;
125 float vtype = vertexType( vertex );
126
127 ATH_CHECK( fill( m_vtx_x_err, vxErr, weight ) );
128 ATH_CHECK( fill( m_vtx_y_err, vyErr, weight ) );
129 ATH_CHECK( fill( m_vtx_z_err, vzErr, weight ) );
130 if( vhasValidTime ) ATH_CHECK( fill( m_vtx_time_err, vtimeErr, weight ) );
131 ATH_CHECK( fill( m_vtx_chi2OverNdof, vchi2OverNdof, weight ) );
132 ATH_CHECK( fill( m_vtx_type, vtype, weight ) );
133 } // close if m_vertexType != "truth"
134
136 if( m_doTrackPlots ) {
137 size_t nvTracks = associatedTracks.size();
138 if( nvTracks != associatedTrackWeights.size() ) {
139 ATH_MSG_ERROR( "Vertex-associated track number mismatch" );
140 return StatusCode::FAILURE;
141 }
142
143 ATH_CHECK( fill( m_vtx_nTracks, nvTracks, weight ) );
144
145 for( size_t it=0 ; it<nvTracks ; it++ ) {
146 ATH_CHECK( fill( m_vtx_track_weight, associatedTrackWeights[it], weight ) );
147
148 float tpt = pT( *associatedTracks[it] ) / Gaudi::Units::GeV;
149 float teta = eta( *associatedTracks[it] );
150 float tnSiHits = nSiHits( *associatedTracks[it] );
151 float tnSiHoles = nSiHoles( *associatedTracks[it] );
152 float td0 = d0( *associatedTracks[it] );
153 float tz0 = z0( *associatedTracks[it] ) - vz;
154 float td0Err = error( *associatedTracks[it], Trk::d0 );
155 // TODO: should this also account for the vertex z error?
156 float tz0Err = error( *associatedTracks[it], Trk::z0 );
157
158 ATH_CHECK( fill( m_vtx_track_pt, tpt, weight ) );
159 ATH_CHECK( fill( m_vtx_track_eta, teta, weight ) );
160 ATH_CHECK( fill( m_vtx_track_nSiHits, tnSiHits, weight ) );
161 ATH_CHECK( fill( m_vtx_track_nSiHoles, tnSiHoles, weight ) );
162 ATH_CHECK( fill( m_vtx_track_d0, td0, weight ) );
163 ATH_CHECK( fill( m_vtx_track_z0, tz0, weight ) );
164 ATH_CHECK( fill( m_vtx_track_d0_err, td0Err, weight ) );
165 ATH_CHECK( fill( m_vtx_track_z0_err, tz0Err, weight ) );
166 } // close loop over tracks
167 } // close if m_doTrackPlots
168
169 return StatusCode::SUCCESS;
170}
171
173 const xAOD::Vertex& vertex,
174 const std::vector< const xAOD::TrackParticle* >& associatedTracks,
175 const std::vector< float >& associatedTrackWeights,
176 float weight );
177
179 const xAOD::TruthVertex& vertex,
180 const std::vector< const xAOD::TruthParticle* >& associatedTracks,
181 const std::vector< float >& associatedTrackWeights,
182 float weight );
183
184
187 int nVertices,
188 float truthMu,
189 float actualMu,
190 float weight )
191{
192 if( not m_doGlobalPlots ) return StatusCode::SUCCESS;
193
195 ATH_CHECK( fill( m_nVtx_vs_actualMu_2D, actualMu, nVertices, weight ) );
196 ATH_CHECK( fill( m_nVtx_vs_actualMu, actualMu, nVertices, weight ) );
197 if( m_doTruthMuPlots ) {
198 ATH_CHECK( fill( m_nVtx_vs_truthMu_2D, truthMu, nVertices, weight ) );
199 ATH_CHECK( fill( m_nVtx_vs_truthMu, truthMu, nVertices, weight ) );
200 }
201
202 return StatusCode::SUCCESS;
203}
204
205
210{
211 ATH_MSG_DEBUG( "Finalising vertex parameters plots" );
213}
Scalar eta() const
pseudorapidity method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
static Double_t sc
Utility methods to access track/truth particles parmeters in a consitent way in this package.
StatusCode retrieveAndBook(P *&pHisto, const std::string &identifier, const std::string &folderOverride="", const std::string &nameOverride="")
Definition PlotMgr.h:64
PlotMgr(const std::string &dirName, const std::string &anaTag, PlotMgr *pParent=nullptr)
Constructor taking parent node and directory name for plots pParent = nullptr by default to book plot...
Definition PlotMgr.cxx:25
void finalizePlots()
Print out final stats on histograms.
void initializePlots()
Book the histograms.
TH1 * m_vtx_x
vertex parameters plots
TH2 * m_nVtx_vs_truthMu_2D
vertex multiplicity vs pileup plots
VertexParametersPlots(PlotMgr *pParent, const std::string &dirName, const std::string &anaTag, const std::string &vertexType, bool doTrackPlots=false, bool doGlobalPlots=false, bool doTruthMuPlots=false)
Constructor.
TH1 * m_vtx_nTracks
vertex-associated tracks plots
StatusCode fillPlots(const VERTEX &vertex, const std::vector< const PARTICLE * > &associatedTracks, const std::vector< float > &associatedTrackWeights, float weight)
Dedicated fill method (for reco and/or truth vertices)
const std::string & getDirectory()
Definition PlotBase.h:88
float nSiHoles(const U &p)
Accessor utility function for getting the value of nSiHoles.
float pT(const U &p)
Accessor utility function for getting the value of pT.
float posX(const V &v)
Accessor utility function for getting the value of vertex position x.
float posZ(const V &v)
Accessor utility function for getting the value of vertex position z.
float chiSquared(const U &p)
float posY(const V &v)
Accessor utility function for getting the value of vertex position y.
float z0(const U &p)
float d0(const U &p)
float timeErr(const V &v)
float ndof(const U &p)
float time(const U &p)
float nSiHits(const U &p)
uint8_t hasValidTime(const U &p)
int vertexType(const V &v)
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
Vertex_v1 Vertex
Define the latest version of the vertex class.
void fill(H5::Group &out_file, size_t iterations)