ATLAS Offline Software
Loading...
Searching...
No Matches
InDetGlobalPrimaryVertexMonAlg.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5
16
17//main header
19
23
24
25//Standard c++
26#include <vector>
27
28
29InDetGlobalPrimaryVertexMonAlg::InDetGlobalPrimaryVertexMonAlg( const std::string& name, ISvcLocator* pSvcLocator ) :
30 AthMonitorAlgorithm(name, pSvcLocator),
35{
36 declareProperty("splitVertexTrkInvFraction", m_splitVertexTrkInvFraction, "inverse fraction to split tracks (1:N)");
37 declareProperty("distanceSplitVertexMatch", m_distanceSplitVxMatch, "Distance for matching split-original Vertex in selection efficiency");
38 declareProperty("splitMatchingMetric", m_splitMatchingMetric, "Determines which function to use to calculate matching between split vertices and original input vertex -- used in selection efficiency");
39 declareProperty("doEnhancedMonitoring" , m_doEnhancedMonitoring, "turn on the enhanced vertex monitoring, it is triggered by the same InDetFlag that also triggers the creation of no beam constraint and split vertices");
40}
41
42
44
45
47
48
49 ATH_CHECK( m_vxContainerName.initialize() );
50 //ATH_CHECK( m_vxContainerNameSplit.initialize() );
51
52
54}
55
56
57StatusCode InDetGlobalPrimaryVertexMonAlg::fillHistograms( const EventContext& ctx ) const {
58 using namespace Monitored;
59
60 //*******************************************************************************
61 //************************** Begin of filling Track Histograms ******************
62 //*******************************************************************************
63
64 ATH_MSG_DEBUG("Filling InDetGlobalPrimaryVertexMonAlg");
65
66 // For histogram naming
67 auto pvGroup = getGroup("PrimaryVertex");
68
69 // retrieving vertices
70 auto handle_vxContainer = SG::makeHandle(m_vxContainerName, ctx); // another way to access ??
71
72 if (!handle_vxContainer.isPresent()) {
73 ATH_MSG_DEBUG ("InDetGlobalPrimaryVertexMonAlg: StoreGate doesn't contain primary vertex container with key "+m_vxContainerName.key());
74 return StatusCode::SUCCESS;
75 }
76 if (!handle_vxContainer.isValid()) {
77 ATH_MSG_ERROR ("InDetGlobalPrimaryVertexMonAlg: Could not retrieve primary vertex container with key "+m_vxContainerName.key());
78 return StatusCode::RECOVERABLE;
79 }
80
81 auto vertexContainer = handle_vxContainer.cptr();
82
83
84 // Total number of vertices (primary and pile up)
85 int pvN = vertexContainer->size()-1; // exclude dummy vertex
86 auto pvN_m = Monitored::Scalar<int>( "m_PvN", pvN);
87 fill(pvGroup, pvN_m);
88
89 int nPriVtx = 0;
90 int nPileupVtx = 0;
91
92 for(const auto vtx : *vertexContainer) {
93
94 if ( !vtx ) continue;
95
96 // Count different types of vertices
97 if (vtx->vertexType() == xAOD::VxType::PriVtx) nPriVtx++;
98 if (vtx->vertexType() == xAOD::VxType::PileUp) nPileupVtx++;
99
100
101 // Select primary vertex
102 if (vtx->vertexType() != xAOD::VxType::PriVtx) continue;
103 if (vtx->numberDoF() <= 0) continue;
104
105 float pvX = vtx->position().x();
106 auto pvX_m = Monitored::Scalar<float>( "m_PvX", pvX);
107 fill(pvGroup, pvX_m);
108
109 float pvY = vtx->position().y();
110 auto pvY_m = Monitored::Scalar<float>( "m_PvY", pvY);
111 fill(pvGroup, pvY_m);
112
113 float pvZ = vtx->position().z();
114 auto pvZ_m = Monitored::Scalar<float>( "m_PvZ", pvZ);
115 fill(pvGroup, pvZ_m);
116
117 float pvErrX = Amg::error( vtx->covariancePosition(), Trk::x);
118 auto pvErrX_m = Monitored::Scalar<float>( "m_PvErrX", pvErrX);
119 fill(pvGroup, pvErrX_m);
120
121 float pvErrY = Amg::error( vtx->covariancePosition(), Trk::y);
122 auto pvErrY_m = Monitored::Scalar<float>( "m_PvErrY", pvErrY);
123 fill(pvGroup, pvErrY_m);
124
125 float pvErrZ = Amg::error( vtx->covariancePosition(), Trk::z);
126 auto pvErrZ_m = Monitored::Scalar<float>( "m_PvErrZ", pvErrZ);
127 fill(pvGroup, pvErrZ_m);
128
129 float pvChiSqDoF = vtx->chiSquared() / vtx->numberDoF() ;
130 auto pvChiSqDoF_m = Monitored::Scalar<float>( "m_PvChiSqDoF", pvChiSqDoF);
131 fill(pvGroup, pvChiSqDoF_m);
132
133
134 auto & trackparticles = vtx->trackParticleLinks();
135
136 int pvNTracks = trackparticles.size() ;
137 auto pvNTracks_m = Monitored::Scalar<int>( "m_PvNTracks", pvNTracks);
138 fill(pvGroup, pvNTracks_m);
139
140
141 // original tracks used for primary vertex
142 for (const auto & trackparticle : trackparticles)
143 {
144 const Trk::Perigee & measuredPerigee = (*trackparticle)->perigeeParameters();
145
146 float pvTrackEta = measuredPerigee.eta() ;
147 auto pvTrackEta_m = Monitored::Scalar<float>( "m_PvTrackEta", pvTrackEta);
148 fill(pvGroup, pvTrackEta_m);
149
150 float pvTrackPt = measuredPerigee.pT()/1000. ; // Histo is in GeV
151 auto pvTrackPt_m = Monitored::Scalar<float>( "m_PvTrackPt", pvTrackPt);
152 fill(pvGroup, pvTrackPt_m);
153
154 }
155
156 } // vxContainer
157
158
159 auto nPriVtx_m = Monitored::Scalar<int>( "m_nPriVtx", nPriVtx);
160 fill(pvGroup, nPriVtx_m);
161
162 auto nPileupVtx_m = Monitored::Scalar<int>( "m_nPileupVtx", nPileupVtx);
163 fill(pvGroup, nPileupVtx_m);
164
165
166 // EnhancedMonitoring is OFF
167
168
169 //*******************************************************************************
170 //**************************** End of filling Track Histograms ******************
171 //*******************************************************************************
172
173 return StatusCode::SUCCESS;
174}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
Implementation of inner detector global Primary Vertex monitoring tool.
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const ToolHandle< GenericMonitoringTool > & getGroup(const std::string &name) const
Get a specific monitoring tool from the tool handle array.
virtual StatusCode initialize() override
initialize
AthMonitorAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
virtual StatusCode initialize() override
initialize
int m_splitMatchingMetric
store metric to be used for split vertex matching in selection efficiency Values currently implemente...
InDetGlobalPrimaryVertexMonAlg(const std::string &name, ISvcLocator *pSvcLocator)
SG::ReadHandleKey< xAOD::VertexContainer > m_vxContainerName
float m_distanceSplitVxMatch
store maximum distance for matching split vertices to original non-BC vertex
virtual StatusCode fillHistograms(const EventContext &ctx) const override
adds event to the monitoring histograms
int m_splitVertexTrkInvFraction
store inverse of the fraction of input tracks used for probe vertex (1:N)
Declare a monitored scalar variable.
double eta() const
Access method for pseudorapidity - from momentum.
double pT() const
Access method for transverse momentum.
void fill(const ToolHandle< GenericMonitoringTool > &groupHandle, std::vector< std::reference_wrapper< Monitored::IMonitoredVariable > > &&variables) const
Fills a vector of variables to a group by reference.
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Generic monitoring tool for athena components.
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ x
Definition ParamDefs.h:55
@ z
global position (cartesian)
Definition ParamDefs.h:57
@ y
Definition ParamDefs.h:56
@ PileUp
Pile-up vertex.
@ PriVtx
Primary vertex.