ATLAS Offline Software
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 #include "Utils.h"
6 #include <numeric>
7 
8 
10 
11 
12 // Detector region: Separation into barrel and endcaps
13 bool inBarrel(double eta){
14  return std::abs(eta) <= fidVol_barrel_etaCut;
15 }
16 
17 
18 bool inBarrel(const Amg::Vector3D &vtx){
19  return inBarrel(vtx.eta());
20 }
21 
22 
23 bool inEndcaps(double eta){
24  return (std::abs(eta) > fidVol_endcaps_etaCut_low && std::abs(eta) < fidVol_endcaps_etaCut_up);
25 }
26 
27 
28 bool inEndcaps(const Amg::Vector3D &vtx){
29  return inEndcaps(vtx.eta());
30 }
31 
32 
34  return inBarrel(vtx) || inEndcaps(vtx);
35 }
36 
37 
38 int getNvtxBarrel(const std::vector<Amg::Vector3D> &vertices){
39  // find the number of vertices in the barrel region
40  return std::accumulate(vertices.begin(), vertices.end(), 0, [](unsigned int N, const Amg::Vector3D& vtx){return N+inBarrel(vtx);});
41 }
42 
43 
44 int getNvtxEndcaps(const std::vector<Amg::Vector3D> &vertices){
45  // find the number of vertices in the endcap region
46  return std::accumulate(vertices.begin(), vertices.end(), 0, [](unsigned int N, const Amg::Vector3D& vtx){return N+inEndcaps(vtx);});
47 }
48 
49 
50 int getNvtxDetectorRegion(const std::vector<Amg::Vector3D> &vertices){
51  // find the number of vertices in the either the barrel or endcaps
52  return getNvtxBarrel(vertices) + getNvtxEndcaps(vertices);
53 }
54 // --- //
55 
56 
57 // Muon spectrometer displaced vertex search fiducial region //
59  bool eta = inBarrel(vtx);
60  bool Lxy = vtx.perp() > fidVol_Lxy_low && vtx.perp() < fidVol_Lxy_up;
61  return eta && Lxy;
62 }
63 
64 
66  bool eta = inEndcaps(vtx);
67  bool z = std::abs(vtx.z()) > fidVol_z_low && std::abs(vtx.z()) < fidVol_z_up;
68  return eta && z;
69 }
70 
71 
72 bool inFiducialVol(const Amg::Vector3D &vtx){
73  return inFiducialVolBarrel(vtx) || inFiducialVolEndcaps(vtx);
74 }
75 
76 
77 int NvtxFiducialVol(const std::vector<Amg::Vector3D> &vertices){
78  // find the number ofvertices in the fiducial volume
79  unsigned int N = 0;
80  for (const Amg::Vector3D &vtx : vertices) if (inFiducialVol(vtx)) ++N;
81  return N;
82 }
83 
84 
85 bool isGoodVtx(const Amg::Vector3D &vtx){
86  // add further conditions for a vertex to be considered in efficiency calculations here
87  bool is_good = true;
88  if (!inFiducialVol(vtx)) is_good = false;
89  return is_good;
90 }
91 // --- //
92 
93 
94 // matching of vertices
95 double getMatchMetric(const Amg::Vector3D &vtx1, const Amg::Vector3D &vtx2){
96  // Computes the metric on which matches are made between two vertices
97  return xAOD::P4Helpers::deltaR(vtx1.eta(), vtx1.phi(), vtx2.eta(), vtx2.phi());
98 }
99 
100 
101 Amg::Vector3D findBestMatch(const Amg::Vector3D &vtx, const std::vector<Amg::Vector3D> &candidates){
102  // finds the best match for vtx among candidates
103  // The match metric is assumed to decreases for better matches. If no match could be found, a zero vector is returned
104  Amg::Vector3D best{Amg::Vector3D::Zero()}; // initialised with zeros
105  double aux(0), min(999);
106 
107  for (unsigned int i=0; i<candidates.size(); ++i){
108  aux = getMatchMetric(vtx, candidates[i]);
109  if (aux > match_max) continue;
110  if (aux < min){
111  min = aux;
112  best = candidates[i];
113  }
114  }
115 
116  return best;
117 }
118 
119 
120 bool isValidMatch(const Amg::Vector3D &match_candidate){
121  return match_candidate.mag() > 0.001;
122 }
123 
124 
125 bool hasMatch(const Amg::Vector3D &vtx1, const std::vector<Amg::Vector3D> &vtx2_vec){
126  // Checks if vtx1 has a match in vtx2_vec
127  Amg::Vector3D candidate = findBestMatch(vtx1, vtx2_vec);
128  return isValidMatch(candidate);
129 }
130 
131 
132 std::vector<Amg::Vector3D> getVertexPos(const std::vector<double> &vtx_x, const std::vector<double> &vtx_y, const std::vector<double> &vtx_z){
133  std::vector<Amg::Vector3D> vertices;
134  for (unsigned int i=0; i<vtx_x.size(); ++i){
135  Amg::Vector3D vtx_pos(vtx_x[i], vtx_y[i], vtx_z[i]);
136  vertices.push_back(vtx_pos);
137  }
138 
139  return vertices;
140 }
141 
142 
143 std::vector<std::vector<Amg::Vector3D>> getConstituentPos(int Nvtx, const std::vector<int> &obj_vtx_link,
144  const std::vector<double> &obj_x, const std::vector<double> &obj_y, const std::vector<double> &obj_z){
145  // returns the constituent position vectors for each reconstructed vertex
146  std::vector<std::vector<Amg::Vector3D>> consti_vec;
147  for (int j=0; j<Nvtx; ++j){
148  std::vector<Amg::Vector3D> consti;
149  for (unsigned int k=0; k<obj_x.size(); ++k){
150  if (obj_vtx_link[k] == j){consti.push_back(Amg::Vector3D(obj_x[k], obj_y[k], obj_z[k]));}
151  }
152  consti_vec.push_back(consti);
153  }
154 
155  return consti_vec;
156 }
157 
158 } // namespace MuonVertexValidationMacroUtils
drawFromPickle.candidates
candidates
Definition: drawFromPickle.py:271
MuonVertexValidationMacroUtils::getMatchMetric
double getMatchMetric(const Amg::Vector3D &vtx1, const Amg::Vector3D &vtx2)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:95
min
constexpr double min()
Definition: ap_fixedTest.cxx:26
MuonVertexValidationMacroUtils::getNvtxDetectorRegion
int getNvtxDetectorRegion(const std::vector< Amg::Vector3D > &vertices)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:50
MuonVertexValidationMacroUtils::fidVol_endcaps_etaCut_low
constexpr double fidVol_endcaps_etaCut_low
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:22
accumulate
bool accumulate(AccumulateMap &map, std::vector< module_t > const &modules, FPGATrackSimMatrixAccumulator const &acc)
Accumulates an accumulator (e.g.
Definition: FPGATrackSimMatrixAccumulator.cxx:22
MuonVertexValidationMacroUtils::getVertexPos
std::vector< Amg::Vector3D > getVertexPos(const std::vector< double > &vtx_x, const std::vector< double > &vtx_y, const std::vector< double > &vtx_z)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:132
MuonVertexValidationMacroUtils::fidVol_z_low
constexpr double fidVol_z_low
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:24
MuonVertexValidationMacroUtils::isGoodVtx
bool isGoodVtx(const Amg::Vector3D &vtx)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:85
MuonVertexValidationMacroUtils::getConstituentPos
std::vector< std::vector< Amg::Vector3D > > getConstituentPos(int Nvtx, const std::vector< int > &obj_vtx_link, const std::vector< double > &obj_x, const std::vector< double > &obj_y, const std::vector< double > &obj_z)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:143
JetTiledMap::N
@ N
Definition: TiledEtaPhiMap.h:44
MuonVertexValidationMacroUtils::fidVol_Lxy_low
constexpr double fidVol_Lxy_low
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:19
MuonVertexValidationMacroUtils
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:9
MuonVertexValidationMacroUtils::fidVol_endcaps_etaCut_up
constexpr double fidVol_endcaps_etaCut_up
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:23
MuonVertexValidationMacroUtils::match_max
constexpr double match_max
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:27
lumiFormat.i
int i
Definition: lumiFormat.py:85
z
#define z
xAOD::P4Helpers::deltaR
double deltaR(double rapidity1, double phi1, double rapidity2, double phi2)
from bare bare rapidity,phi
Definition: xAODP4Helpers.h:150
Utils.h
python.utils.best
def best(iterable, priorities=[3, 2, 1, -1, 0])
Definition: DataQuality/DQUtils/python/utils.py:50
MuonVertexValidationMacroUtils::inDetectorRegion
bool inDetectorRegion(const Amg::Vector3D &vtx)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:33
MuonVertexValidationMacroUtils::getNvtxBarrel
int getNvtxBarrel(const std::vector< Amg::Vector3D > &vertices)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:38
MuonVertexValidationMacroUtils::isValidMatch
bool isValidMatch(const Amg::Vector3D &match_candidate)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:120
MuonVertexValidationMacroUtils::inFiducialVol
bool inFiducialVol(const Amg::Vector3D &vtx)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:72
MuonVertexValidationMacroUtils::inEndcaps
bool inEndcaps(double eta)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:23
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
MuonVertexValidationMacroUtils::inFiducialVolEndcaps
bool inFiducialVolEndcaps(const Amg::Vector3D &vtx)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:65
MuonVertexValidationMacroUtils::fidVol_z_up
constexpr double fidVol_z_up
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:25
MuonVertexValidationMacroUtils::inFiducialVolBarrel
bool inFiducialVolBarrel(const Amg::Vector3D &vtx)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:58
MuonVertexValidationMacroUtils::inBarrel
bool inBarrel(double eta)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:13
MuonVertexValidationMacroUtils::NvtxFiducialVol
int NvtxFiducialVol(const std::vector< Amg::Vector3D > &vertices)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:77
MuonVertexValidationMacroUtils::getNvtxEndcaps
int getNvtxEndcaps(const std::vector< Amg::Vector3D > &vertices)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:44
MuonVertexValidationMacroUtils::fidVol_Lxy_up
constexpr double fidVol_Lxy_up
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:20
MuonVertexValidationMacroUtils::findBestMatch
Amg::Vector3D findBestMatch(const Amg::Vector3D &vtx, const std::vector< Amg::Vector3D > &candidates)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:101
MuonVertexValidationMacroUtils::hasMatch
bool hasMatch(const Amg::Vector3D &vtx1, const std::vector< Amg::Vector3D > &vtx2_vec)
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.cxx:125
MuonVertexValidationMacroUtils::fidVol_barrel_etaCut
constexpr double fidVol_barrel_etaCut
Definition: MuonSpectrometer/MuonValidation/MuonVertexValidation/util/Utils.h:18
fitman.k
k
Definition: fitman.py:528
generate::Zero
void Zero(TH1D *hin)
Definition: generate.cxx:32