FCCAnalyses
Loading...
Searching...
No Matches
VertexFinderLCFIPlus.h
Go to the documentation of this file.
1#ifndef VERTEXFINDERLCFIPLUS_ANALYZERS_H
2#define VERTEXFINDERLCFIPLUS_ANALYZERS_H
3
4#include <cmath>
5#include <vector>
6
7#include "ROOT/RVec.hxx"
8#include "edm4hep/ReconstructedParticleData.h"
9#include "edm4hep/TrackState.h"
10
14
15#include "fastjet/JetDefinition.hh"
16
17namespace FCCAnalyses{
26namespace VertexFinderLCFIPlus{
27
32 ROOT::VecOps::RVec<ROOT::VecOps::RVec<VertexingUtils::FCCAnalysesVertex>> get_SV_jets( ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> recoparticles,
33 ROOT::VecOps::RVec<edm4hep::TrackState> thetracks,
35 ROOT::VecOps::RVec<bool> isInPrimary,
36 ROOT::VecOps::RVec<fastjet::PseudoJet> jets,
37 std::vector<std::vector<int>> jet_consti,
38 bool V0_rej=true,
39 double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5. ) ;
40
45 ROOT::VecOps::RVec<VertexingUtils::FCCAnalysesVertex> get_SV_event( ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> recoparticles,
46 ROOT::VecOps::RVec<edm4hep::TrackState> thetracks,
48 ROOT::VecOps::RVec<bool> isInPrimary,
49 bool V0_rej=true,
50 double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5. ) ;
51
55 ROOT::VecOps::RVec<VertexingUtils::FCCAnalysesVertex> get_SV_event( ROOT::VecOps::RVec<edm4hep::TrackState> np_tracks,
56 ROOT::VecOps::RVec<edm4hep::TrackState> thetracks,
58 bool V0_rej=true,
59 double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5. ) ;
60
64 ROOT::VecOps::RVec<int> VertexSeed_best( ROOT::VecOps::RVec<edm4hep::TrackState> tracks,
66 double chi2_cut=9., double invM_cut=10.) ;
67
72 ROOT::VecOps::RVec<int> addTrack_best( ROOT::VecOps::RVec<edm4hep::TrackState> tracks,
73 ROOT::VecOps::RVec<int> vtx_tr,
75 double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.) ;
76
81 ROOT::VecOps::RVec<edm4hep::TrackState> V0rejection_tight( ROOT::VecOps::RVec<edm4hep::TrackState> tracks,
83 bool V0_rej=true ) ;
84
88 ROOT::VecOps::RVec<VertexingUtils::FCCAnalysesVertex> findSVfromTracks( ROOT::VecOps::RVec<edm4hep::TrackState> tracks_fin,
89 const ROOT::VecOps::RVec<edm4hep::TrackState>& alltracks,
91 double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.) ;
92
99 ROOT::VecOps::RVec<edm4hep::TrackState> tracks,
101 bool seed=true,
102 double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.) ;
103
110 ROOT::VecOps::RVec<bool> isV0( ROOT::VecOps::RVec<edm4hep::TrackState> np_tracks,
112 bool tight = false ) ;
113
115
119 VertexingUtils::FCCAnalysesV0 get_V0s( ROOT::VecOps::RVec<edm4hep::TrackState> np_tracks,
121 bool tight,
122 double chi2_cut=9. ) ;
123
127 VertexingUtils::FCCAnalysesV0 get_V0s( ROOT::VecOps::RVec<edm4hep::TrackState> np_tracks,
129 double Ks_invM_low=0.493, double Ks_invM_high=0.503, double Ks_dis=0.5, double Ks_cosAng=0.999,
130 double Lambda_invM_low=1.111, double Lambda_invM_high=1.121, double Lambda_dis=0.5, double Lambda_cosAng=0.99995,
131 double Gamma_invM_low=0., double Gamma_invM_high=0.005, double Gamma_dis=9, double Gamma_cosAng=0.99995,
132 double chi2_cut=9. ) ;
133
137 VertexingUtils::FCCAnalysesV0 get_V0s_jet( ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData> recoparticles,
138 ROOT::VecOps::RVec<edm4hep::TrackState> thetracks,
139 ROOT::VecOps::RVec<bool> isInPrimary,
140 ROOT::VecOps::RVec<fastjet::PseudoJet> jets,
141 std::vector<std::vector<int>> jet_consti,
143 bool tight = true,
144 double chi2_cut=9. );
145
156 ROOT::VecOps::RVec<double> get_V0candidate( VertexingUtils::FCCAnalysesVertex &V0_vtx,
157 ROOT::VecOps::RVec<edm4hep::TrackState> tr_pair,
159 bool chi2,
160 double chi2_cut=9. );
161
172 ROOT::VecOps::RVec<double> constraints_Ks(bool tight) ;
173 ROOT::VecOps::RVec<double> constraints_Lambda0(bool tight) ;
174 ROOT::VecOps::RVec<double> constraints_Gamma(bool tight) ;
175 //
176 ROOT::VecOps::RVec<double> constraints_Ks(double invM_low, double invM_high, double dis, double cosAng) ;
177 ROOT::VecOps::RVec<double> constraints_Lambda0(double invM_low, double invM_high, double dis, double cosAng) ;
178 ROOT::VecOps::RVec<double> constraints_Gamma(double invM_low, double invM_high, double dis, double cosAng) ;
182 //ROOT::VecOps::RVec<ROOT::VecOps::RVec<int>> VertexSeed_all( ROOT::VecOps::RVec<edm4hep::TrackState> tracks,
183 // VertexingUtils::FCCAnalysesVertex PV,
184 // double chi2_cut=9., double invM_cut=10.) ;
185
190 //ROOT::VecOps::RVec<int> addTrack_multi( ROOT::VecOps::RVec<edm4hep::TrackState> tracks,
191 // ROOT::VecOps::RVec<int> vtx_tr,
192 // VertexingUtils::FCCAnalysesVertex PV,
193 // double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.) ;
194
195
196}//end NS VertexFinderLCFIPlus
197
198}//end NS FCCAnalyses
199#endif
ROOT::VecOps::RVec< edm4hep::TrackState > V0rejection_tight(ROOT::VecOps::RVec< edm4hep::TrackState > tracks, VertexingUtils::FCCAnalysesVertex PV, bool V0_rej=true)
V0 rejection (tight) takes all (non-primary tracks) & removes tracks coming from V0s if user chooses ...
Definition VertexFinderLCFIPlus.cc:257
VertexingUtils::FCCAnalysesV0 get_V0s_jet(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > recoparticles, ROOT::VecOps::RVec< edm4hep::TrackState > thetracks, ROOT::VecOps::RVec< bool > isInPrimary, ROOT::VecOps::RVec< fastjet::PseudoJet > jets, std::vector< std::vector< int > > jet_consti, VertexingUtils::FCCAnalysesVertex PV, bool tight=true, double chi2_cut=9.)
returns V0s reconstructed in each jet of the event (as an FCCAnalysesV0 object) need to perform jet c...
Definition VertexFinderLCFIPlus.cc:638
VertexingUtils::FCCAnalysesV0 get_V0s(ROOT::VecOps::RVec< edm4hep::TrackState > np_tracks, VertexingUtils::FCCAnalysesVertex PV, bool tight, double chi2_cut=9.)
returns V0s reconstructed from a set of tracks (as an FCCAnalysesV0 object) constraint thresholds can...
Definition VertexFinderLCFIPlus.cc:443
ROOT::VecOps::RVec< ROOT::VecOps::RVec< VertexingUtils::FCCAnalysesVertex > > get_SV_jets(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > recoparticles, ROOT::VecOps::RVec< edm4hep::TrackState > thetracks, VertexingUtils::FCCAnalysesVertex PV, ROOT::VecOps::RVec< bool > isInPrimary, ROOT::VecOps::RVec< fastjet::PseudoJet > jets, std::vector< std::vector< int > > jet_consti, bool V0_rej=true, double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.)
returns SVs reconstructed from non-primary tracks of jets non-primary separated from all tracks using...
Definition VertexFinderLCFIPlus.cc:18
ROOT::VecOps::RVec< double > constraints_Lambda0(bool tight)
Definition VertexFinderLCFIPlus.cc:857
ROOT::VecOps::RVec< int > VertexSeed_best(ROOT::VecOps::RVec< edm4hep::TrackState > tracks, VertexingUtils::FCCAnalysesVertex PV, double chi2_cut=9., double invM_cut=10.)
returns indices of the best pair of tracks from a vector of (non-primary) tracks default chi2 thresho...
Definition VertexFinderLCFIPlus.cc:154
ROOT::VecOps::RVec< VertexingUtils::FCCAnalysesVertex > findSVfromTracks(ROOT::VecOps::RVec< edm4hep::TrackState > tracks_fin, const ROOT::VecOps::RVec< edm4hep::TrackState > &alltracks, VertexingUtils::FCCAnalysesVertex PV, double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.)
find SVs from a set of tracks default values of thresholds for the constraints are set
Definition VertexFinderLCFIPlus.cc:274
ROOT::VecOps::RVec< bool > isV0(ROOT::VecOps::RVec< edm4hep::TrackState > np_tracks, VertexingUtils::FCCAnalysesVertex PV, bool tight=false)
V0 rejection/identification takes all (non-primary) tracks & assigns "true" to pairs that form a V0 i...
Definition VertexFinderLCFIPlus.cc:368
ROOT::VecOps::RVec< double > get_V0candidate(VertexingUtils::FCCAnalysesVertex &V0_vtx, ROOT::VecOps::RVec< edm4hep::TrackState > tr_pair, VertexingUtils::FCCAnalysesVertex PV, bool chi2, double chi2_cut=9.)
returns invariant mass, distance from PV, and colliniarity variables for all V0 candidates [0] -> inv...
Definition VertexFinderLCFIPlus.cc:782
ROOT::VecOps::RVec< double > constraints_Ks(bool tight)
functions to fill constraint thresholds tight -> tight constraints !tight -> loose constraints also a...
Definition VertexFinderLCFIPlus.cc:836
ROOT::VecOps::RVec< double > constraints_Gamma(bool tight)
Definition VertexFinderLCFIPlus.cc:878
ROOT::VecOps::RVec< int > addTrack_best(ROOT::VecOps::RVec< edm4hep::TrackState > tracks, ROOT::VecOps::RVec< int > vtx_tr, VertexingUtils::FCCAnalysesVertex PV, double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.)
adds index of the best track (from the remaining tracks) to the (seed) vtx default chi2 threshold is ...
Definition VertexFinderLCFIPlus.cc:205
bool check_constraints(VertexingUtils::FCCAnalysesVertex vtx, ROOT::VecOps::RVec< edm4hep::TrackState > tracks, VertexingUtils::FCCAnalysesVertex PV, bool seed=true, double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.)
check constraints of vertex candidates default values of thresholds for the constraints are set defau...
Definition VertexFinderLCFIPlus.cc:327
ROOT::VecOps::RVec< VertexingUtils::FCCAnalysesVertex > get_SV_event(ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > recoparticles, ROOT::VecOps::RVec< edm4hep::TrackState > thetracks, VertexingUtils::FCCAnalysesVertex PV, ROOT::VecOps::RVec< bool > isInPrimary, bool V0_rej=true, double chi2_cut=9., double invM_cut=10., double chi2Tr_cut=5.)
returns SVs reconstructed from non-primary tracks of the event SV finding done before jet clustering ...
Definition VertexFinderLCFIPlus.cc:80
FCC analyzers collection.
Definition Algorithms.h:15
Structure to keep useful information that is related to the V0.
Definition VertexingUtils.h:50
Structure to keep useful track information that is related to the vertex.
Definition VertexingUtils.h:38