1 #ifndef RecoBTag_SecondaryVertex_CombinedSVComputer_h 2 #define RecoBTag_SecondaryVertex_CombinedSVComputer_h 10 #include <Math/VectorUtil.h> 38 #define range_for(i, x) \ 39 for(int i = (x).begin; i != (x).end; i += (x).increment) 71 template <
class SVTI,
class IPTI>
73 const IPTI & ipInfo,
const SVTI & svInfo,
74 double & vtx_track_ptSum,
double & vtx_track_ESum)
const;
94 template <
class SVTI,
class IPTI>
96 const IPTI & ipInfo,
const SVTI & svInfo,
97 double & vtx_track_ptSum,
double & vtx_track_ESum)
const 100 using namespace reco;
102 typedef typename IPTI::input_container
Container;
107 bool havePv = ipInfo.primaryVertex().
isNonnull();
111 ipInfo.primaryVertex()->y(),
112 ipInfo.primaryVertex()->z());
129 double jet_track_ESum= 0.;
135 if (vtx < 0) vtx =
i;
150 std::vector<std::size_t> indices = ipInfo.sortedIndexes(
sortCriterium);
151 const std::vector<reco::btag::TrackIPData> &ipData = ipInfo.impactParameterData();
153 const Container &
tracks = ipInfo.selectedTracks();
154 std::vector<TrackRef> pseudoVertexTracks;
156 std::vector<TrackRef> trackPairV0Test(2);
159 std::size_t
idx = indices[
i];
161 const TrackRef &
track = tracks[
idx];
173 allKinematics.
add(track);
177 pseudoVertexTracks.push_back(track);
178 vertexKinematics.
add(track);
182 trackPairV0Test[0] =
track;
188 std::size_t pairIdx = indices[j];
190 const TrackRef &pairTrack = tracks[pairIdx];
195 trackPairV0Test[1] = pairTrack;
204 trackJetKinematics.
add(track);
208 double trackMag =
std::sqrt(trackMom.Mag2());
225 double perp_trackMom_jetDir = VectorUtil::Perp2(trackMom, jetDir);
226 perp_trackMom_jetDir = (perp_trackMom_jetDir > 0. ?
std::sqrt(perp_trackMom_jetDir ) : 0. );
238 for(
typename std::vector<TrackRef>::const_iterator trkIt = pseudoVertexTracks.begin(); trkIt != pseudoVertexTracks.end(); ++trkIt)
241 vtx_track_ptSum +=
std::sqrt((*trkIt)->momentum().Perp2());
279 double varPi = (vertexMass/5.2794) * (vtx_track_ESum /jet_track_ESum);
281 double varB = (
std::sqrt(5.2794) * vtx_track_ptSum) / ( vertexMass *
std::sqrt(jet->pt()));
294 if ( pfJet !=
nullptr )
311 else if( patJet !=
nullptr && patJet->
isPFJet() )
347 #endif // RecoBTag_SecondaryVertex_CombinedSVComputer_h
reco::TrackSelector trackSelector
int photonMultiplicity() const
photonMultiplicity
const math::XYZTLorentzVector & weightedVectorSum() const
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction (relative to uncorrected jet energy)
float photonEnergyFraction() const
photonEnergyFraction (relative to corrected jet energy)
int electronMultiplicity() const
electronMultiplicity
virtual ~CombinedSVComputer()=default
unsigned int pseudoMultiplicityMin
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
void add(const reco::Track &track, double weight=1.0)
int photonMultiplicity() const
photonMultiplicity
bool vertexMassCorrection
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Base class for all types of Jets.
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
reco::TrackSelector trackNoDeltaRSelector
reco::V0Filter pseudoVertexV0Filter
Global3DPoint GlobalPoint
int muonMultiplicity() const
muonMultiplicity
bool isNonnull() const
Checks for non-null.
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
CombinedSVComputer(const edm::ParameterSet ¶ms)
float photonEnergyFraction() const
photonEnergyFraction
Jets made from PFObjects.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
reco::btag::SortCriteria sortCriterium
Container::value_type value_type
reco::TrackSelector trackPseudoSelector
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
IterationRange flipIterate(int size, bool vertex) const
bool isPFJet() const
check to see if the jet is a reco::PFJet
reco::V0Filter trackPairV0Filter
double deltaR(double eta1, double eta2, double phi1, double phi2)
std::vector< reco::btau::TaggingVariableName > taggingVariables
double significance() const
unsigned int numberOfTracks() const
float electronEnergyFraction() const
electronEnergyFraction (relative to corrected jet energy)
GlobalPoint closestToJetAxis
XYZVectorD XYZVector
spatial vector with cartesian internal representation
double flipValue(double value, bool vertex) const
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
float electronEnergyFraction() const
electronEnergyFraction
Analysis-level calorimeter jet class.
edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset) const
const reco::btag::TrackIPData & threshTrack(const reco::CandIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Measurement1D distanceToJetAxis
int muonMultiplicity() const
muonMultiplicity
void fillCommonVariables(reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
char data[epos_bytes_allocation]
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
const math::XYZTLorentzVector & vectorSum() const
float muonEnergyFraction() const
muonEnergyFraction (relative to corrected jet energy)
float muonEnergyFraction() const
muonEnergyFraction
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
int electronMultiplicity() const
electronMultiplicity
void insert(const TaggingVariable &variable, bool delayed=false)
unsigned int trackMultiplicityMin