ATLAS Offline Software
Loading...
Searching...
No Matches
JetFitterV0FinderTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef JETFITTER_V0FINDER_TOOL_H
6#define JETFITTER_V0FINDER_TOOL_H 1
7
9#include "GaudiKernel/ToolHandle.h"
10
16
17#include <TH2D.h>
18
20
21namespace InDetDD {
23}
24
25namespace InDet {
26
27 static const InterfaceID IID_JetFitterV0FinderTool("JetFitterV0FinderTool", 1, 0);
28
30 public:
31
32 static const InterfaceID& interfaceID() {
34 }
35
36 StatusCode initialize();
37 StatusCode finalize();
38
39 JetFitterV0FinderTool(const std::string &t, const std::string &n, const IInterface *p);
41
43 const TLorentzVector&,
44 std::vector< const Trk::ITrackLink* >&,
45 const std::vector< const xAOD::Vertex* >&,
46 std::vector< const Trk::ITrackLink* >&,
47 std::vector< const Trk::ITrackLink* >&,
48 Amg::Vector3D& ) const;
49
50 private:
51 std::vector< const xAOD::Vertex* > findV0candidates( const xAOD::Vertex&,
52 const TLorentzVector&,
53 std::vector< const Trk::ITrackLink* >&,
54 const std::vector< const xAOD::Vertex* >& ) const;
55
57 const TLorentzVector&,
58 const xAOD::Vertex& ) const;
59
61 const Trk::ITrackLink* ) const;
62
64 const TLorentzVector&,
65 const std::vector< Trk::PositionAndWeight >& ) const;
66
67 private:
68 ToolHandle< InDet::InDetJetFitterUtils > m_jetFitterUtils {this,"InDetJetFitterUtils","InDet::InDetJetFitterUtils/InDetJetFitterUtils",""};
69 ToolHandle< Trk::IMode3dFinder > m_mode3dfinder {this,"Mode3dFinder","Trk::Mode3dTo1dFinder/Mode3dTo1dFinder",""};
70
71 private:
74
75 private:
76 // @brief
77 Gaudi::Property< bool > m_revertFromPositiveToNegativeTags {this,"revertFromPositiveToNegativeTags",false,""};
78
79 Gaudi::Property< double > m_cutTwoTrkVtxVertexProbForBFirstSelectionFirstCriterium {this,"cutTwoTrkVtxVtxProbForBFirstSelectCriteriumA",0.05,""};
80 Gaudi::Property< double > m_cutTwoTrkVtxVertexProbForBFirstSelectionSecondCriterium {this,"cutTwoTrkVtxVtxProbForBFirstSelectCriteriumB",0.034,""};
81
82 Gaudi::Property< double > m_cutCompatibilityPrimaryVertexSingleTrackForBFirstSelection {this,"cutCompatibilityPrimaryVertexSingleTrackForBFirstSelection",1e-1,""};
83 Gaudi::Property< double > m_cutCompatibilityPrimaryVertexBothTracksForBFirstSelection {this,"cutCompatibilityPrimaryVertexBothTracksForBFirstSelection",1e-2,""};
84 Gaudi::Property< double > m_cutIPD0BothTracksForBFirstSelection {this,"cutIPD0BothTracksForBFirstSelection",3.5,""};
85 Gaudi::Property< double > m_cutIPZ0BothTracksForBFirstSelection {this,"cutIPZ0BothTracksForBFirstSelection",5.,""};
86 Gaudi::Property< double > m_cutPtBothTracksForBFirstSelection {this,"cutPtBothTracksForBFirstSelection",500.,""};
87
88 Gaudi::Property< double > m_cutTwoTrkVtxLifetimeSignificanceForBFirstSelectionFirstCriterium {this,"cutTwoTrkVtxLifeSignForBFirstSelectCriteriumA",1.,""};
89 Gaudi::Property< double > m_cutTwoTrkVtxLifetimeSignificanceForBFirstSelectionSecondCriterium {this,"cutTwoTrkVtxLifeSignForBFirstSelectCriteriumB",1.5,""};
90 Gaudi::Property< double > m_cutCompatibilityToPrimarySingleTrackForMatInteractions {this,"cutCompToPrimarySingleTrackForMatInterac",1e-4,""};
91 Gaudi::Property< double > m_cutCompatibilityToPrimaryBothTracksForMatInteractions {this,"cutCompToPrimaryBothTracksForMatInterac",1e-6,""};
92
93
94 Gaudi::Property< double > m_firstLayer_min {this,"firstLayer_min",34.0-2.5,""};
95 Gaudi::Property< double > m_firstLayer_max {this,"firstLayer_max",34.0+2.5,""};
96 Gaudi::Property< double > m_secondLayer_min {this,"secondLayer_min",51.5-3,""};
97 Gaudi::Property< double > m_secondLayer_max {this,"secondLayer_max",51.5+3,""};
98
99 Gaudi::Property< double > m_cutCompatibilityPrimaryVertexSinglePositiveLifetimeTrackForBSecondSelection {this,"cutCompPVSinglePosLifeTrackForBSecondSelect",5e-2,""};
100 Gaudi::Property< double > m_cutCompatibilityPrimaryVertexSingleNegativeLifetimeTrackForBSecondSelection {this,"cutCompPVSingleNegLifeTrackForBSecondSelect",1e-2,""};
101 Gaudi::Property< double > m_cutIPD0SigBoxSingleTrackForBSecondSelection {this,"cutIPD0SigBoxSingleTrackForBSecondSelection",2.,""};
102 Gaudi::Property< double > m_cutIPZ0SigBoxSingleTrackForBSecondSelection {this,"cutIPZ0SigBoxSingleTrackForBSecondSelection",5.,""};
103
104 Gaudi::Property< double > m_cutIPD0SingleTrackForBSecondSelection {this,"cutIPD0SingleTrackForBSecondSelection",1.5,""};
105 Gaudi::Property< double > m_cutIPZ0SingleTrackForBSecondSelection {this,"cutIPZ0SingleTrackForBSecondSelection",3.,""};
106 Gaudi::Property< double > m_cutPtSingleTrackForBSecondSelection {this,"cutPtSingleTrackForBSecondSelection",750,""};
107
108 Gaudi::Property< bool > m_useITkMaterialRejection {this,"useITkMaterialRejection",false,"Reject vertices from hadronic interactions in detector material using ITk layout"};
111 std::unique_ptr<TH2D> m_ITkPixMaterialMap;
112
113 };
114
115}
116
117#endif
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Dedicated detector manager extending the functionality of the SiDetectorManager with dedicated pixel ...
Gaudi::Property< double > m_cutIPZ0SigBoxSingleTrackForBSecondSelection
static const InterfaceID & interfaceID()
Gaudi::Property< double > m_cutCompatibilityToPrimarySingleTrackForMatInteractions
Gaudi::Property< double > m_cutTwoTrkVtxLifetimeSignificanceForBFirstSelectionSecondCriterium
ToolHandle< Trk::IMode3dFinder > m_mode3dfinder
Gaudi::Property< double > m_cutCompatibilityPrimaryVertexSingleNegativeLifetimeTrackForBSecondSelection
Gaudi::Property< double > m_firstLayer_min
Gaudi::Property< double > m_cutIPZ0BothTracksForBFirstSelection
SG::AuxElement::Accessor< float > m_compatibilityAccessor
SG::AuxElement::Accessor< std::vector< const Trk::ITrackLink * > > m_tracksAccessor
Gaudi::Property< double > m_cutPtBothTracksForBFirstSelection
const InDetDD::PixelDetectorManager * m_pixelManager
Gaudi::Property< bool > m_revertFromPositiveToNegativeTags
std::unique_ptr< TH2D > m_ITkPixMaterialMap
Gaudi::Property< double > m_cutTwoTrkVtxLifetimeSignificanceForBFirstSelectionFirstCriterium
Gaudi::Property< double > m_cutIPD0SingleTrackForBSecondSelection
bool checkCriteriaFirstFit(const xAOD::Vertex &, const TLorentzVector &, const xAOD::Vertex &) const
Gaudi::Property< bool > m_useITkMaterialRejection
Gaudi::Property< double > m_cutCompatibilityPrimaryVertexBothTracksForBFirstSelection
Gaudi::Property< double > m_secondLayer_max
JetFitterV0FinderTool(const std::string &t, const std::string &n, const IInterface *p)
Gaudi::Property< double > m_cutIPZ0SingleTrackForBSecondSelection
Gaudi::Property< double > m_cutPtSingleTrackForBSecondSelection
const Trk::TwoTrackVerticesInJet * doV0Finding(const xAOD::Vertex &, const TLorentzVector &, std::vector< const Trk::ITrackLink * > &, const std::vector< const xAOD::Vertex * > &, std::vector< const Trk::ITrackLink * > &, std::vector< const Trk::ITrackLink * > &, Amg::Vector3D &) const
Gaudi::Property< double > m_cutIPD0SigBoxSingleTrackForBSecondSelection
Gaudi::Property< double > m_secondLayer_min
Gaudi::Property< double > m_cutTwoTrkVtxVertexProbForBFirstSelectionSecondCriterium
std::vector< const xAOD::Vertex * > findV0candidates(const xAOD::Vertex &, const TLorentzVector &, std::vector< const Trk::ITrackLink * > &, const std::vector< const xAOD::Vertex * > &) const
bool checkCriteriaSecondFit(const xAOD::Vertex &, const Trk::ITrackLink *) const
Gaudi::Property< double > m_cutCompatibilityPrimaryVertexSinglePositiveLifetimeTrackForBSecondSelection
Gaudi::Property< double > m_cutCompatibilityPrimaryVertexSingleTrackForBFirstSelection
Gaudi::Property< double > m_cutIPD0BothTracksForBFirstSelection
Gaudi::Property< double > m_cutCompatibilityToPrimaryBothTracksForMatInteractions
Gaudi::Property< double > m_cutTwoTrkVtxVertexProbForBFirstSelectionFirstCriterium
const BeamPipeDetectorManager * m_beamPipeMgr
Amg::Vector3D computeSeedDirection(const xAOD::Vertex &, const TLorentzVector &, const std::vector< Trk::PositionAndWeight > &) const
Gaudi::Property< double > m_firstLayer_max
ToolHandle< InDet::InDetJetFitterUtils > m_jetFitterUtils
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
Eigen::Matrix< double, 3, 1 > Vector3D
Message Stream Member.
Primary Vertex Finder.
static const InterfaceID IID_JetFitterV0FinderTool("JetFitterV0FinderTool", 1, 0)
Vertex_v1 Vertex
Define the latest version of the vertex class.