|  | ATLAS Offline Software
    | 
 
 
 
Go to the documentation of this file.
   26 #ifndef _VKalVrt_NewVrtSecInclusiveTool_H 
   27 #define _VKalVrt_NewVrtSecInclusiveTool_H 
   32 #include "GaudiKernel/ToolHandle.h" 
   33 #include "GaudiKernel/ServiceHandle.h" 
   36 #define BOOST_ALLOW_DEPRECATED_HEADERS 
   37 #include "boost/graph/adjacency_list.hpp" 
   49 #include "Math/LorentzVector.h" 
   59   class TrkVKalVrtFitter;
 
   70   struct workVectorArrxAOD{
 
   73         std::vector<const xAOD::TrackParticle*> 
inpTrk;         
 
  104       void lockDecorations (
const std::vector<const xAOD::TrackParticle*> & inpTrk) 
const;
 
  144       std::unique_ptr<Hists> 
m_h;
 
  145       Gaudi::Property<bool> 
m_fillHist{
this, 
"FillHist",   
false, 
"Fill debugging and development histograms+ntuple"  };
 
  148       Gaudi::Property<int> 
m_cutSctHits{
this,   
"CutSctHits",   4, 
"Remove track if it has less SCT hits" };
 
  149       Gaudi::Property<int> 
m_cutPixelHits{
this, 
"CutPixelHits", 2, 
"Remove track if it has less Pixel hits"};
 
  150       Gaudi::Property<int> 
m_cutTRTHits{
this,   
"CutTRTHits",  10, 
"Remove track if it has less TRT hits"};
 
  151       Gaudi::Property<int> 
m_cutSiHits{
this,    
"CutSiHits",    8, 
"Remove track if it has less Pixel+SCT hits"  };
 
  152       Gaudi::Property<int> 
m_cutBLayHits{
this,  
"CutBLayHits",  0, 
"Remove track if it has less B-layer hits"   };
 
  153       Gaudi::Property<int> 
m_cutSharedHits{
this,
"CutSharedHits",1, 
"Reject final 2tr vertices if tracks have shared hits" };
 
  154       Gaudi::Property<float> 
m_cutPt{
this,      
"CutPt",     500.,  
"Track Pt selection cut"  };
 
  155       Gaudi::Property<float> 
m_cutD0Min{
this,   
"CutD0Min",    0.,  
"Track minimal D0 selection cut"  };
 
  156       Gaudi::Property<float> 
m_cutD0Max{
this,   
"CutD0Max",   10.,  
"Track maximal D0 selection cut"  };
 
  157       Gaudi::Property<float> 
m_maxZVrt{
this,    
"MaxZVrt",    15.,  
"Track Z impact selection max"};
 
  158       Gaudi::Property<float> 
m_minZVrt{
this,    
"MinZVrt",     0.,  
"Track Z impact selection min"};
 
  159       Gaudi::Property<float> 
m_cutChi2{
this,    
"CutChi2",     5.,  
"Track Chi2 selection cut" };
 
  160       Gaudi::Property<float> 
m_antiPileupSigRCut{
this, 
"AntiPileupSigRCut", 2.0,  
"Upper cut on significance of 2D distance between beam and perigee"  };
 
  163       Gaudi::Property<float> 
m_trkSigCut{
this,    
"TrkSigCut",     2.0, 
"Track 3D impact significance w/r primary vertex. Should be >=AntiPileupSigRCut" };
 
  164       Gaudi::Property<float> 
m_dRdZRatioCut{
this, 
"dRdZRatioCut", 0.25, 
"Cut on dR/dZ ratio to remove pileup tracks"  };
 
  167       Gaudi::Property<float> 
m_fastZSVCut{
this, 
"FastZSVCut",        10.,  
"Cut to remove SV candidates based on fast SV estimation. To save full fit CPU."  };
 
  170       Gaudi::Property<float> 
m_removeTrkMatSignif{
this, 
"removeTrkMatSignif", 0., 
"Significance of Vertex-TrackingMaterial distance for removal. No removal if <=0." };
 
  171       Gaudi::Property<float> 
m_beampipeR{
this, 
"BeampipeR",   24.3, 
"Radius of the beampipe material for aggressive material rejection" };
 
  174       Gaudi::Property<float> 
m_vrtMassLimit{
this,    
"VrtMassLimit",  5500., 
"Maximal allowed mass for found vertices" };
 
  175       Gaudi::Property<float> 
m_globVrtProbCut{
this, 
"GlobVrtProbCut", 0.005, 
"Cut on probability of any vertex for final selection"  };
 
  176       Gaudi::Property<float> 
m_maxSVRadiusCut{
this, 
"MaxSVRadiusCut",  140., 
"Cut on maximal radius of SV (def = Pixel detector size)"  };
 
  177       Gaudi::Property<float> 
m_selVrtSigCut{
this,   
"SelVrtSigCut",     3.0, 
"Cut on significance of 3D distance between vertex and PV"  };
 
  178       Gaudi::Property<float> 
m_vertexMergeCut{
this,  
"VertexMergeCut",   4., 
"To allow vertex merging for MultiVertex Finder" };
 
  179       Gaudi::Property<bool> 
m_multiWithOneTrkVrt{
this, 
"MultiWithOneTrkVrt", 
true, 
"Allow one-track-vertex addition to already found secondary vertices"};
 
  184       ToolHandle<Trk::IExtrapolator>  
m_extrapolator{
this,
"ExtrapolatorName",
"Trk::Extrapolator/Extrapolator", 
"Name of the extrapolator tool"};
 
  185       ToolHandle<Trk::TrkVKalVrtFitter>  
m_fitSvc{
this, 
"VertexFitterTool", 
"Trk::TrkVKalVrtFitter/VertexFitterTool", 
"Name of the Vertex Fitter tool"};
 
  186       ToolHandle<Reco::ITrackToVertex>  
m_trackToVertexTool{
this, 
"TrackToVertexTool", 
"Reco::TrackToVertex/TrackToVertex", 
"Name of the TrackToVertex tool"};
 
  188       ToolHandle<Rec::ITwoTrackVertexSelector>  
m_ini_v2trselector{
this, 
"TwoTrkVtxSelectorIni", 
"Rec::TwoTrackVrtBDTSelector/V2TrBDTSelectorIni",
 
  189                                                                          "Name of the initial 2-track vertex selector"};
 
  190       ToolHandle<Rec::ITwoTrackVertexSelector>  
m_fin_v2trselector{
this, 
"TwoTrkVtxSelectorFinal", 
"Rec::TwoTrackVrtBDTSelector/V2TrBDTSelectorFin",
 
  191                                                                          "Name of the final 2-track vertex selector"};
 
  193       Gaudi::Property<std::string> 
m_augString {
this, 
"AugmentingVersionString", 
"_NVSI", 
"Augmentation version string"};
 
  336       void printWrkSet(
const std::vector<WrkVrt> * WrkSet, 
const std::string &
name ) 
const;
 
  340       static double massV0(
const std::vector< std::vector<double> >& TrkAtVrt, 
double massP, 
double massPi ) ;
 
  343       ROOT::Math::PxPyPzEVector 
momAtVrt(
const std::vector<double>& inpTrk) 
const;
 
  347       static int   nTrkCommon( std::vector<WrkVrt> *WrkVrtSet, 
int indexV1, 
int indexV2) ;
 
  348       double minVrtVrtDist( std::vector<WrkVrt> *WrkVrtSet, 
int & indexV1, 
int & indexV2, std::vector<double> & 
check) 
const;
 
  349       static bool isPart( 
const std::deque<long int>& 
test, std::deque<long int> 
base) ;
 
  350       static std::vector<double> 
estimVrtPos( 
int nTrk, std::deque<long int> &selTrk, std::map<
long int,std::vector<double>> & vrt) ;
 
  353                                   const std::vector<double>& vrtErr,
double& signif ) ;
 
  355                                   const std::vector<double>& vrtErr,
double& signif ) ;
 
  357                         const Amg::Vector3D & vrt2, 
const std::vector<double>& vrtErr2) ;
 
  366       double refitVertex( 
WrkVrt &Vrt,std::vector<const xAOD::TrackParticle*> & SelectedTracks,
 
  370       static int mostHeavyTrk(
WrkVrt V, std::vector<const xAOD::TrackParticle*> AllTracks) ;
 
  374                                     std::vector<const xAOD::TrackParticle*> & AllTrackList,
 
  387                         std::map<
long int,std::vector<double>> & vrt,
 
  403     template <
typename Clique, 
typename Graph>
 
  406       std::vector<int> new_clique(0);
 
  407       for(
auto i = clq.begin(); 
i != clq.end(); ++
i) new_clique.push_back(*
i);
 
  
std::vector< const xAOD::TrackParticle * > listSelTracks
std::vector< const xAOD::TrackParticle * > tmpListTracks
void clique(const Clique &clq, Graph &)
Property holding a SG store/key/clid/attr name from which a WriteDecorHandle is made.
Helper class to provide type-safe access to aux data.
def TProfile(*args, **kwargs)
::StatusCode StatusCode
StatusCode definition for legacy code.
Class describing a truth particle in the MC record.
std::vector< std::vector< int > > & m_allCliques
Ensure that the ATLAS eigen extensions are properly loaded.
std::vector< const xAOD::TrackParticle * > inpTrk
Eigen::Matrix< double, 3, 1 > Vector3D
Class describing a Vertex.
clique_visitor(std::vector< std::vector< int > > &input)
def TH1F(name, title, nxbins, bins_par2, bins_par3=None, path='', **kwargs)
Class describing a TrackParticle.