ATLAS Offline Software
Loading...
Searching...
No Matches
SecVtxValidationPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include "TVector3.h"
10
11using Athena::Units::GeV;
12
14 const std::string& sDir)
15 : PlotBase(pParent, sDir)
16{
17 // position
18 m_vertex_x = Book1D("x", "Vertex_x;x [mm];Entries;", 300, -300, 300, false);
19 m_vertex_y = Book1D("y", "Vertex_y;y [mm];Entries;", 300, -300, 300, false);
20 m_vertex_z = Book1D("z", "Vertex_z;z [mm];Entries;", 300, -300, 300, false);
21 m_vertex_r = Book1D("r", "Vertex_r;r [mm];Entries;", 300, 0, 300, false);
22
23 // four vector
24 m_vertex_pt = Book1D("pt", "Vertex_pt;p_{T} [GeV];Entries;", 100, 0, 100, false);
25 m_vertex_eta = Book1D("eta", "Vertex_eta;#eta;Entries;", 50, -2.5, 2.5, false);
26 m_vertex_phi = Book1D("phi", "Vertex_phi;#phi;Entries;", 64, -3.2, 3.2, false);
27 m_vertex_m = Book1D("m", "Vertex_m;m [GeV];Entries;", 100, 0, 100, false);
28
29 // misc
30 m_vertex_ntrk = Book1D("ntrk", "Vertex_ntrk;n_{trk};Entries;", 20, 0, 20, false);
31 m_vertex_chi2 = Book1D("chi2", "Vertex_chi2;#chi^{2}/n_{DoF};Entries;", 50, 0, 10, false);
32 m_vertex_charge = Book1D("charge", "Vertex_charge;charge;Entries;", 20, -10, 10, false);
33 m_vertex_mind0 = Book1D("mind0", "Vertex_mind0;d_{0,min} [mm];Entries;", 100, 0, 10, false);
34 m_vertex_maxd0 = Book1D("maxd0", "Vertex_maxd0;d_{0,max} [mm];Entries;", 300, 0, 300, false);
35}
36
38
39 TVector3 vtx_pos(secVtx->x(), secVtx->y(), secVtx->z());
40 float vtx_r = vtx_pos.Perp();
41
42 xAOD::Vertex::ConstAccessor<float> Chi2("chiSquared");
43 xAOD::Vertex::ConstAccessor<float> nDoF("numberDoF");
45
46 float minD0 = 1.0 * 1.e4;
47 float maxD0 = 0.0;
48 int charge = 0;
49
50 const xAOD::Vertex::TrackParticleLinks_t & trkParts = trkAcc( *secVtx );
51
52 size_t ntrk = trkParts.size();
53 TLorentzVector sumP4(0,0,0,0);
54
55 for (const auto &trklink : trkParts) {
56 if (!trklink.isValid()) {
57 continue;
58 }
59 const xAOD::TrackParticle & trk = **trklink;
60
61 double trk_d0 = std::abs(trk.definingParameters()[0]);
62 if(trk_d0 < minD0){ minD0 = trk_d0; }
63 if(trk_d0 > maxD0){ maxD0 = trk_d0; }
64
68
69 TLorentzVector trk_P4;
70 trk_P4.SetPtEtaPhiM(ptAcc(trk), etaAcc(trk), phiAcc(trk), trk.p4().M());
71
72 sumP4 += trk_P4;
73 charge += trk.charge();
74 }
75
76 m_vertex_x->Fill(secVtx->x());
77 m_vertex_y->Fill(secVtx->y());
78 m_vertex_z->Fill(secVtx->z());
79 m_vertex_r->Fill(vtx_r);
80 m_vertex_pt->Fill(sumP4.Pt() / GeV); // GeV
81 m_vertex_eta->Fill(sumP4.Eta());
82 m_vertex_phi->Fill(sumP4.Phi());
83 m_vertex_m->Fill(sumP4.M() / GeV); // GeV
84 m_vertex_ntrk->Fill(ntrk);
85 m_vertex_chi2->Fill(Chi2(*secVtx)/nDoF(*secVtx));
87 m_vertex_mind0->Fill(minD0);
88 m_vertex_maxd0->Fill(maxD0);
89
90}
double charge(const T &p)
Definition AtlasPID.h:997
Wrapper to avoid constant divisions when using units.
TH1D * Book1D(const std::string &name, const std::string &labels, int nBins, float start, float end, bool prependDir=true)
Book a TH1D histogram.
Definition PlotBase.cxx:94
PlotBase(PlotBase *parent, const std::string &sDir)
Definition PlotBase.cxx:29
SG::ConstAccessor< T, ALLOC > ConstAccessor
Helper class to provide type-safe access to aux data.
Definition AuxElement.h:131
void fill(const xAOD::Vertex *secVtx)
fill the histograms
SecVtxValidationPlots(PlotBase *pParent, const std::string &sDir)
Standard Constructor.
DefiningParameters_t definingParameters() const
Returns a SVector of the Perigee track parameters.
virtual FourMom_t p4() const override final
The full 4-momentum of the particle.
float charge() const
Returns the charge.
float z() const
Returns the z position.
std::vector< ElementLink< xAOD::TrackParticleContainer > > TrackParticleLinks_t
Type for the associated track particles.
Definition Vertex_v1.h:128
float y() const
Returns the y position.
float x() const
Returns the x position.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.