ATLAS Offline Software
InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
5 // Header include
7 //-------------------------------------------------
9 #include "xAODTracking/Vertex.h"
10 #include <cmath>
11 #include <vector>
12 
13 namespace InDet{
14 
15 
16 
17  /* Technicalities */
18 
19 
20  double InDetSVWithMuonTool::VrtVrtDist(const xAOD::Vertex & PrimVrt, const Amg::Vector3D & SecVrt,
21  const std::vector<double>& SecVrtErr, double& Signif)
22 
23  {
24  double distx = PrimVrt.x()- SecVrt.x();
25  double disty = PrimVrt.y()- SecVrt.y();
26  double distz = PrimVrt.z()- SecVrt.z();
27 
28 
29  AmgSymMatrix(3) PrimCovMtx=PrimVrt.covariancePosition(); //Create
30  PrimCovMtx(0,0) += SecVrtErr[0];
31  PrimCovMtx(0,1) += SecVrtErr[1];
32  PrimCovMtx(1,0) += SecVrtErr[1];
33  PrimCovMtx(1,1) += SecVrtErr[2];
34  PrimCovMtx(0,2) += SecVrtErr[3];
35  PrimCovMtx(2,0) += SecVrtErr[3];
36  PrimCovMtx(1,2) += SecVrtErr[4];
37  PrimCovMtx(2,1) += SecVrtErr[4];
38  PrimCovMtx(2,2) += SecVrtErr[5];
39 
40  AmgSymMatrix(3) WgtMtx = PrimCovMtx.inverse();
41 
42  Signif = distx*WgtMtx(0,0)*distx
43  +disty*WgtMtx(1,1)*disty
44  +distz*WgtMtx(2,2)*distz
45  +2.*distx*WgtMtx(0,1)*disty
46  +2.*distx*WgtMtx(0,2)*distz
47  +2.*disty*WgtMtx(1,2)*distz;
48  Signif=(Signif<0) ? 0.:std::sqrt(Signif);
49  return std::sqrt(distx*distx+disty*disty+distz*distz);
50  }
51 
52 
53 
54  double InDetSVWithMuonTool::ConeDist(const AmgVector(5) & VectPerig, const TLorentzVector & Dir)
55 
56  {
57  double etaTr = -std::log(std::tan(VectPerig[3]*0.5));
58  double etaJet = Dir.PseudoRapidity();
59  double adphi = std::abs(Dir.Phi()-VectPerig[2]);
60  while(adphi> M_PI)adphi-=2.*M_PI;
61  return std::sqrt(adphi*adphi + (etaJet-etaTr)*(etaJet-etaTr));
62  }
63 
64 
65 
66 
67  int InDetSVWithMuonTool::FindMaxAfterFirst( std::vector<double>& Chi2PerTrk)
68 
69  {
70  double Chi2Ref=0.;
71  int Position=1;
72  if( Chi2PerTrk.size() < 2 ) return Position ;
73  for (int i=1; i<(int)Chi2PerTrk.size(); i++){
74  if( Chi2PerTrk[i] > Chi2Ref) { Chi2Ref=Chi2PerTrk[i]; Position=i;}
75  }
76  return Position;
77  }
78 
79 
80 
81 // Function returns a transverse momentum of track w/r some direction
82 //
83  double InDetSVWithMuonTool::pTvsDir(const Amg::Vector3D &Dir, const std::vector< double >& InpTrk)
84 
85  {
86  double Norm=std::sqrt(Dir.x()*Dir.x() + Dir.y()*Dir.y() + Dir.z()*Dir.z());
87  double sx=Dir.x()/Norm;
88  double sy=Dir.y()/Norm;
89  double sz=Dir.z()/Norm;
90 
91  double px=0.,py=0.,pz=0.;
92  double scale{};
93  const auto invDenom{1./std::abs(InpTrk[2])};
94  px = std::cos ( InpTrk[0]) * std::sin(InpTrk[1])*invDenom;
95  py = std::sin ( InpTrk[0]) * std::sin(InpTrk[1])*invDenom;
96  pz = std::cos(InpTrk[1])*invDenom;
97  scale = px*sx + py*sy + pz*sz;
98  px -= sx*scale;
99  py -= sy*scale;
100  pz -= sz*scale;
101  return std::sqrt( px*px +py*py + pz*pz );
102  }
103 
104 
105 
106  StatusCode InDetSVWithMuonTool::VKalVrtFitFastBase(const std::vector<const xAOD::TrackParticle*>& listTrk,
107  Amg::Vector3D & FitVertex,
108  Trk::IVKalState& istate)
109  const
110  { return m_fitSvc->VKalVrtFitFast(listTrk,FitVertex,istate); }
111 
112 
113 
114 
115  StatusCode InDetSVWithMuonTool::VKalVrtFitBase(const std::vector<const xAOD::TrackParticle*> & listPart,
117  TLorentzVector& Momentum,
118  long int& Charge,
119  std::vector<double>& ErrorMatrix,
120  std::vector<double>& Chi2PerTrk,
121  std::vector< std::vector<double> >& TrkAtVrt,
122  double& Chi2,
123  Trk::IVKalState& istate,
124  bool ifCovV0) const
125  {
126  std::vector<const xAOD::NeutralParticle*> netralPartDummy(0);
127  return m_fitSvc->VKalVrtFit( listPart, netralPartDummy,Vertex, Momentum, Charge,
128  ErrorMatrix, Chi2PerTrk, TrkAtVrt, Chi2,
129  istate, ifCovV0);
130 
131  }
132 
133  TLorentzVector InDetSVWithMuonTool::TotalMom(const std::vector<const xAOD::TrackParticle*>& InpTrk)
134 
135  {
136  TLorentzVector sum(0.,0.,0.,0.);
137  for (const auto *i : InpTrk) {
138  if( i == nullptr ) continue;
139  sum += i->p4();
140  }
141  return sum;
142  }
143 
144 
145 } //end namespace
Trk::TrkVKalVrtFitter::VKalVrtFit
virtual StatusCode VKalVrtFit(const std::vector< const xAOD::TrackParticle * > &, const std::vector< const xAOD::NeutralParticle * > &, Amg::Vector3D &Vertex, TLorentzVector &Momentum, long int &Charge, dvect &ErrorMatrix, dvect &Chi2PerTrk, std::vector< std::vector< double >> &TrkAtVrt, double &Chi2, IVKalState &istate, bool ifCovV0=false) const override final
xAOD::Vertex_v1::x
float x() const
Returns the x position.
fitman.sy
sy
Definition: fitman.py:524
fitman.sz
sz
Definition: fitman.py:527
test_pyathena.px
px
Definition: test_pyathena.py:18
InDetSVWithMuonTool.h
xAOD::Vertex
Vertex_v1 Vertex
Define the latest version of the vertex class.
Definition: Event/xAOD/xAODTracking/xAODTracking/Vertex.h:16
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
InDet
DUMMY Primary Vertex Finder.
Definition: VP1ErrorUtils.h:36
Trk::TrkVKalVrtFitter::VKalVrtFitFast
virtual StatusCode VKalVrtFitFast(std::span< const xAOD::TrackParticle *const >, Amg::Vector3D &Vertex, double &minDZ, IVKalState &istate) const
Definition: VKalVrtFitFastSvc.cxx:56
InDet::InDetSVWithMuonTool::m_fitSvc
Trk::TrkVKalVrtFitter * m_fitSvc
Definition: InDetSVWithMuonTool.h:139
python.compressB64.sx
string sx
Definition: compressB64.py:96
M_PI
#define M_PI
Definition: ActiveFraction.h:11
drawFromPickle.cos
cos
Definition: drawFromPickle.py:36
yodamerge_tmp.scale
scale
Definition: yodamerge_tmp.py:138
AmgSymMatrix
#define AmgSymMatrix(dim)
Definition: EventPrimitives.h:52
TrkVKalVrtFitter.h
convertTimingResiduals.sum
sum
Definition: convertTimingResiduals.py:55
lumiFormat.i
int i
Definition: lumiFormat.py:92
InDet::InDetSVWithMuonTool::VrtVrtDist
static double VrtVrtDist(const xAOD::Vertex &PrimVrt, const Amg::Vector3D &SecVrt, const std::vector< double > &VrtErr, double &Signif)
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:20
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
AmgVector
AmgVector(4) T2BSTrackFilterTool
Definition: T2BSTrackFilterTool.cxx:114
JetVoronoiDiagramHelpers::Norm
Point Norm(const Point &a)
Definition: JetVoronoiDiagramHelpers.cxx:79
Amg::pz
@ pz
Definition: GeoPrimitives.h:40
drawFromPickle.tan
tan
Definition: drawFromPickle.py:36
xAOD::Vertex_v1::z
float z() const
Returns the z position.
Vertex.h
Amg::py
@ py
Definition: GeoPrimitives.h:39
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
InDet::InDetSVWithMuonTool::ConeDist
static double ConeDist(const AmgVector(5) &, const TLorentzVector &)
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:54
InDet::InDetSVWithMuonTool::pTvsDir
static double pTvsDir(const Amg::Vector3D &Dir, const std::vector< double > &InpTrk)
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:83
xAOD::Vertex_v1
Class describing a Vertex.
Definition: Vertex_v1.h:42
Trk::IVKalState
Definition: IVKalState.h:21
InDet::InDetSVWithMuonTool::VKalVrtFitFastBase
StatusCode VKalVrtFitFastBase(const std::vector< const xAOD::TrackParticle * > &listPart, Amg::Vector3D &Vertex, Trk::IVKalState &istate) const
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:106
InDet::InDetSVWithMuonTool::VKalVrtFitBase
StatusCode VKalVrtFitBase(const std::vector< const xAOD::TrackParticle * > &listPart, Amg::Vector3D &Vertex, TLorentzVector &Momentum, long int &Charge, std::vector< double > &ErrorMatrix, std::vector< double > &Chi2PerTrk, std::vector< std::vector< double > > &TrkAtVrt, double &Chi2, Trk::IVKalState &istate, bool ifCovV0) const
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:115
xAOD::Vertex_v1::y
float y() const
Returns the y position.
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
drawFromPickle.sin
sin
Definition: drawFromPickle.py:36
InDet::InDetSVWithMuonTool::FindMaxAfterFirst
static int FindMaxAfterFirst(std::vector< double > &Chi2PerTrk)
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:67
python.exampleDriverScript.Dir
Dir
Definition: exampleDriverScript.py:20
InDet::InDetSVWithMuonTool::TotalMom
static TLorentzVector TotalMom(const std::vector< const xAOD::TrackParticle * > &InpTrk)
Definition: InnerDetector/InDetRecTools/InDetSVWithMuonTool/src/Utilities.cxx:133