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