CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Attributes
CombinedSVComputer Class Reference

#include <CombinedSVComputer.h>

Inheritance diagram for CombinedSVComputer:
CombinedSVSoftLeptonComputer

Classes

struct  IterationRange
 

Public Member Functions

 CombinedSVComputer (const edm::ParameterSet &params)
 
edm::ParameterSet dropDeltaR (const edm::ParameterSet &pset) const
 
template<class SVTI , class IPTI >
void fillCommonVariables (reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
 
IterationRange flipIterate (int size, bool vertex) const
 
double flipValue (double value, bool vertex) const
 
virtual reco::TaggingVariableList operator() (const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
 
virtual reco::TaggingVariableList operator() (const reco::CandIPTagInfo &ipInfo, const reco::CandSecondaryVertexTagInfo &svInfo) const
 
const reco::btag::TrackIPDatathreshTrack (const reco::CandIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
 
const reco::btag::TrackIPDatathreshTrack (const reco::TrackIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
 
virtual ~CombinedSVComputer ()=default
 

Private Attributes

double charmCut
 
double minTrackWeight
 
unsigned int pseudoMultiplicityMin
 
reco::V0Filter pseudoVertexV0Filter
 
reco::btag::SortCriteria sortCriterium
 
std::vector< reco::btau::TaggingVariableNametaggingVariables
 
bool trackFlip
 
unsigned int trackMultiplicityMin
 
reco::TrackSelector trackNoDeltaRSelector
 
reco::V0Filter trackPairV0Filter
 
reco::TrackSelector trackPseudoSelector
 
reco::TrackSelector trackSelector
 
bool useTrackWeights
 
bool vertexFlip
 
bool vertexMassCorrection
 

Detailed Description

Definition at line 42 of file CombinedSVComputer.h.

Constructor & Destructor Documentation

CombinedSVComputer::CombinedSVComputer ( const edm::ParameterSet params)
explicit

Definition at line 13 of file CombinedSVComputer.cc.

13  :
14  trackFlip(params.getParameter<bool>("trackFlip")),
15  vertexFlip(params.getParameter<bool>("vertexFlip")),
16  charmCut(params.getParameter<double>("charmCut")),
18  trackSelector(params.getParameter<edm::ParameterSet>("trackSelection")),
20  trackPseudoSelector(params.getParameter<edm::ParameterSet>("trackPseudoSelection")),
21  pseudoMultiplicityMin(params.getParameter<unsigned int>("pseudoMultiplicityMin")),
22  trackMultiplicityMin(params.getParameter<unsigned int>("trackMultiplicityMin")),
23  minTrackWeight(params.getParameter<double>("minimumTrackWeight")),
24  useTrackWeights(params.getParameter<bool>("useTrackWeights")),
25  vertexMassCorrection(params.getParameter<bool>("correctVertexMass")),
26  pseudoVertexV0Filter(params.getParameter<edm::ParameterSet>("pseudoVertexV0Filter")),
27  trackPairV0Filter(params.getParameter<edm::ParameterSet>("trackPairV0Filter"))
28 {
29 
30 }
reco::TrackSelector trackSelector
T getParameter(std::string const &) const
reco::btag::SortCriteria getCriterium(const std::string &name)
Definition: TrackSorting.cc:11
unsigned int pseudoMultiplicityMin
reco::TrackSelector trackNoDeltaRSelector
reco::V0Filter pseudoVertexV0Filter
reco::btag::SortCriteria sortCriterium
reco::TrackSelector trackPseudoSelector
reco::V0Filter trackPairV0Filter
edm::ParameterSet dropDeltaR(const edm::ParameterSet &pset) const
unsigned int trackMultiplicityMin
virtual CombinedSVComputer::~CombinedSVComputer ( )
virtualdefault

Member Function Documentation

edm::ParameterSet CombinedSVComputer::dropDeltaR ( const edm::ParameterSet pset) const
inline

Definition at line 6 of file CombinedSVComputer.cc.

References edm::ParameterSet::addParameter().

7 {
8  edm::ParameterSet psetCopy(pset);
9  psetCopy.addParameter<double>("jetDeltaRMax", 99999.0);
10  return psetCopy;
11 }
template<class SVTI , class IPTI >
void CombinedSVComputer::fillCommonVariables ( reco::TaggingVariableList vars,
reco::TrackKinematics vertexKinematics,
const IPTI &  ipInfo,
const SVTI &  svInfo,
double &  vtx_track_ptSum,
double &  vtx_track_ESum 
) const

Definition at line 95 of file CombinedSVComputer.h.

References reco::TrackKinematics::add(), reco::PFJet::chargedHadronEnergyFraction(), reco::btau::chargedHadronEnergyFraction, pat::Jet::chargedHadronEnergyFraction(), reco::btau::chargedHadronMultiplicity, reco::PFJet::chargedHadronMultiplicity(), pat::Jet::chargedHadronMultiplicity(), reco::btag::TrackIPData::closestToJetAxis, data, HiRegitMuonDetachedTripletStep_cff::DeltaR, boostedElectronIsolation_cff::deltaR, dir, reco::btag::TrackIPData::distanceToJetAxis, reco::btau::electronEnergyFraction, reco::PFJet::electronEnergyFraction(), pat::Jet::electronEnergyFraction(), reco::btau::electronMultiplicity, reco::PFJet::electronMultiplicity(), pat::Jet::electronMultiplicity(), reco::btau::etaRel(), reco::btau::flightDistance1dSig, reco::btau::flightDistance1dVal, reco::btau::flightDistance2dSig, reco::btau::flightDistance2dVal, reco::btau::flightDistance3dSig, reco::btau::flightDistance3dVal, flipIterate(), flipValue(), reco::btau::hadronMultiplicity, reco::btau::hadronPhotonMultiplicity, mps_fire::i, training_settings::idx, reco::TaggingVariableList::insert(), reco::btag::TrackIPData::ip2d, reco::btag::IP2DSig, reco::btag::TrackIPData::ip3d, electrons_cff::ip3d, reco::btag::IP3DSig, edm::RefToBase< T >::isNonnull(), pat::Jet::isPFJet(), metsig::jet, reco::btau::jetAbsEta, reco::btau::jetEta, reco::btau::jetNSecondaryVertices, reco::btau::jetNTracks, reco::btau::jetPt, PV3DBase< T, PVType, FrameType >::mag2(), reco::btau::massVertexEnergyFraction, reco::btau::muonEnergyFraction, reco::PFJet::muonEnergyFraction(), pat::Jet::muonEnergyFraction(), reco::btau::muonMultiplicity, reco::PFJet::muonMultiplicity(), pat::Jet::muonMultiplicity(), reco::PFJet::neutralHadronEnergyFraction(), reco::btau::neutralHadronEnergyFraction, pat::Jet::neutralHadronEnergyFraction(), reco::btau::neutralHadronMultiplicity, reco::PFJet::neutralHadronMultiplicity(), pat::Jet::neutralHadronMultiplicity(), reco::btag::Vertices::NoVertex, reco::TrackKinematics::numberOfTracks(), convertSQLiteXML::ok, reco::btau::photonEnergyFraction, reco::PFJet::photonEnergyFraction(), pat::Jet::photonEnergyFraction(), reco::btau::photonMultiplicity, reco::PFJet::photonMultiplicity(), pat::Jet::photonMultiplicity(), reco::ParticleMasses::piPlus, pseudoMultiplicityMin, reco::btag::Vertices::PseudoVertex, pseudoVertexV0Filter, MetAnalyzer::pv(), range_for, reco::btag::Vertices::RecoVertex, Measurement1D::significance(), sortCriterium, mathSSE::sqrt(), HGCalGeometryMode::Square, threshTrack(), reco::btau::totalMultiplicity, HiIsolationCommonParameters_cff::track, reco::btau::trackDecayLenVal, reco::btau::trackDeltaR, reco::btau::trackEta, reco::btau::trackEtaRel, reco::btau::trackJetDistVal, reco::btau::trackJetPt, reco::btau::trackMomentum, trackMultiplicityMin, trackPairV0Filter, reco::btau::trackPPar, reco::btau::trackPParRatio, trackPseudoSelector, reco::btau::trackPtRatio, reco::btau::trackPtRel, l1t::tracks, trackSelector, reco::btau::trackSip2dSig, reco::btau::trackSip2dSigAboveCharm, reco::btau::trackSip2dVal, reco::btau::trackSip2dValAboveCharm, reco::btau::trackSip3dSig, reco::btau::trackSip3dSigAboveCharm, reco::btau::trackSip3dVal, reco::btau::trackSip3dValAboveCharm, reco::btau::trackSumJetDeltaR, reco::btau::trackSumJetEtRatio, useTrackWeights, Measurement1D::value(), reco::TrackKinematics::vectorSum(), reco::btau::vertexBoostOverSqrtJetPt, reco::btau::vertexCategory, reco::btau::vertexEnergyRatio, vertexFlip, reco::btau::vertexJetDeltaR, reco::btau::vertexMass, vertexMassCorrection, reco::btau::vertexNTracks, extraflags_cff::vtx, reco::TrackKinematics::weightedVectorSum(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by operator()().

98 {
99  using namespace ROOT::Math;
100  using namespace reco;
101 
102  typedef typename IPTI::input_container Container;
103  typedef typename Container::value_type TrackRef;
104 
105  edm::RefToBase<Jet> jet = ipInfo.jet();
106  math::XYZVector jetDir = jet->momentum().Unit();
107  bool havePv = ipInfo.primaryVertex().isNonnull();
108  GlobalPoint pv;
109  if (havePv)
110  pv = GlobalPoint(ipInfo.primaryVertex()->x(),
111  ipInfo.primaryVertex()->y(),
112  ipInfo.primaryVertex()->z());
113 
115 
116 
117  vars.insert(btau::jetPt, jet->pt(), true);
118  vars.insert(btau::jetEta, jet->eta(), true);
119  vars.insert(btau::jetAbsEta, fabs(jet->eta()), true);
120 
121  if (ipInfo.selectedTracks().size() < trackMultiplicityMin)
122  return;
123 
124  vars.insert(btau::jetNTracks, ipInfo.selectedTracks().size(), true);
125 
126  TrackKinematics allKinematics;
127  TrackKinematics trackJetKinematics;
128 
129  double jet_track_ESum= 0.;
130 
131  int vtx = -1;
132 
133  IterationRange range = flipIterate(svInfo.nVertices(), true);
134  range_for(i , range) {
135  if (vtx < 0) vtx = i;
136  }
137 
138  if (vtx >= 0) {
139  vtxType = btag::Vertices::RecoVertex;
140  vars.insert(btau::flightDistance1dVal,flipValue(svInfo.flightDistance(vtx, 1).value(),true),true);
141  vars.insert(btau::flightDistance1dSig,flipValue(svInfo.flightDistance(vtx, 1).significance(),true),true);
142  vars.insert(btau::flightDistance2dVal,flipValue(svInfo.flightDistance(vtx, 2).value(),true),true);
143  vars.insert(btau::flightDistance2dSig,flipValue(svInfo.flightDistance(vtx, 2).significance(),true),true);
144  vars.insert(btau::flightDistance3dVal,flipValue(svInfo.flightDistance(vtx, 3).value(),true),true);
145  vars.insert(btau::flightDistance3dSig,flipValue(svInfo.flightDistance(vtx, 3).significance(),true),true);
146  vars.insert(btau::vertexJetDeltaR,Geom::deltaR(svInfo.flightDirection(vtx), vertexFlip ? -jetDir : jetDir),true);
147  vars.insert(btau::jetNSecondaryVertices, svInfo.nVertices(), true);
148  }
149 
150  std::vector<std::size_t> indices = ipInfo.sortedIndexes(sortCriterium);
151  const std::vector<reco::btag::TrackIPData> &ipData = ipInfo.impactParameterData();
152 
153  const Container &tracks = ipInfo.selectedTracks();
154  std::vector<TrackRef> pseudoVertexTracks;
155 
156  std::vector<TrackRef> trackPairV0Test(2);
157  range = flipIterate(indices.size(), false);
158  range_for(i, range) {
159  std::size_t idx = indices[i];
160  const reco::btag::TrackIPData &data = ipData[idx];
161  const TrackRef &track = tracks[idx];
162 
163  jet_track_ESum += std::sqrt(track->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
164 
165  // add track to kinematics for all tracks in jet
166  //allKinematics.add(track); // would make more sense for some variables, e.g. vertexEnergyRatio nicely between 0 and 1, but not necessarily the best option for the discriminating power...
167 
168  // filter track -> this track selection can be tighter than the vertex track selection (used to fill the track related variables...)
169  if (!trackSelector(track, data, *jet, pv))
170  continue;
171 
172  // add track to kinematics for all tracks in jet
173  allKinematics.add(track);
174 
175  // if no vertex was reconstructed, attempt pseudo vertex
176  if (vtxType == btag::Vertices::NoVertex && trackPseudoSelector(track, data, *jet, pv)) {
177  pseudoVertexTracks.push_back(track);
178  vertexKinematics.add(track);
179  }
180 
181  // check against all other tracks for V0 track pairs
182  trackPairV0Test[0] = track;
183  bool ok = true;
184  range_for(j, range) {
185  if (i == j)
186  continue;
187 
188  std::size_t pairIdx = indices[j];
189  const reco::btag::TrackIPData &pairTrackData = ipData[pairIdx];
190  const TrackRef &pairTrack = tracks[pairIdx];
191 
192  if (!trackSelector(pairTrack, pairTrackData, *jet, pv))
193  continue;
194 
195  trackPairV0Test[1] = pairTrack;
196  if (!trackPairV0Filter(trackPairV0Test)) {
197  ok = false;
198  break;
199  }
200  }
201  if (!ok)
202  continue;
203 
204  trackJetKinematics.add(track);
205 
206  // add track variables
207  math::XYZVector trackMom = track->momentum();
208  double trackMag = std::sqrt(trackMom.Mag2());
209 
210  vars.insert(btau::trackSip3dVal, flipValue(data.ip3d.value(), false), true);
211  vars.insert(btau::trackSip3dSig, flipValue(data.ip3d.significance(), false), true);
212  vars.insert(btau::trackSip2dVal, flipValue(data.ip2d.value(), false), true);
213  vars.insert(btau::trackSip2dSig, flipValue(data.ip2d.significance(), false), true);
215 // vars.insert(btau::trackJetDistSig, data.distanceToJetAxis.significance(), true);
216 // vars.insert(btau::trackFirstTrackDist, data.distanceToFirstTrack, true);
217 // vars.insert(btau::trackGhostTrackVal, data.distanceToGhostTrack.value(), true);
218 // vars.insert(btau::trackGhostTrackSig, data.distanceToGhostTrack.significance(), true);
219  vars.insert(btau::trackDecayLenVal, havePv ? (data.closestToJetAxis - pv).mag() : -1.0, true);
220 
221  vars.insert(btau::trackMomentum, trackMag, true);
222  vars.insert(btau::trackEta, trackMom.Eta(), true);
223 
224  // check for erroneous Perp^2 values before evaluating Perp (fix for DeepCSV NaN outputs)
225  double perp_trackMom_jetDir = VectorUtil::Perp2(trackMom, jetDir);
226  perp_trackMom_jetDir = (perp_trackMom_jetDir > 0. ? std::sqrt(perp_trackMom_jetDir ) : 0. );
227 
228  vars.insert(btau::trackPtRel, perp_trackMom_jetDir, true);
229  vars.insert(btau::trackPPar, jetDir.Dot(trackMom), true);
230  vars.insert(btau::trackDeltaR, VectorUtil::DeltaR(trackMom, jetDir), true);
231  vars.insert(btau::trackPtRatio, perp_trackMom_jetDir / trackMag, true);
232  vars.insert(btau::trackPParRatio, jetDir.Dot(trackMom) / trackMag, true);
233  }
234 
235  if (vtxType == btag::Vertices::NoVertex && vertexKinematics.numberOfTracks() >= pseudoMultiplicityMin && pseudoVertexV0Filter(pseudoVertexTracks))
236  {
238  for(typename std::vector<TrackRef>::const_iterator trkIt = pseudoVertexTracks.begin(); trkIt != pseudoVertexTracks.end(); ++trkIt)
239  {
240  vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir,(*trkIt)->momentum()), true);
241  vtx_track_ptSum += std::sqrt((*trkIt)->momentum().Perp2());
242  vtx_track_ESum += std::sqrt((*trkIt)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
243  }
244  }
245 
246  vars.insert(btau::vertexCategory, vtxType, true);
247 
248  vars.insert(btau::trackJetPt, trackJetKinematics.vectorSum().Pt(), true);
249  vars.insert(btau::trackSumJetDeltaR,VectorUtil::DeltaR(allKinematics.vectorSum(), jetDir), true);
250  vars.insert(btau::trackSumJetEtRatio,allKinematics.vectorSum().Et() / ipInfo.jet()->et(), true);
251 
252  vars.insert(btau::trackSip3dSigAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.significance(),false),true);
253  vars.insert(btau::trackSip3dValAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP3DSig, *jet, pv).ip3d.value(),false),true);
254  vars.insert(btau::trackSip2dSigAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.significance(),false),true);
255  vars.insert(btau::trackSip2dValAboveCharm, flipValue(threshTrack(ipInfo, reco::btag::IP2DSig, *jet, pv).ip2d.value(),false),true);
256 
257  if (vtxType != btag::Vertices::NoVertex) {
259  ? allKinematics.weightedVectorSum()
260  : allKinematics.vectorSum();
262  ? vertexKinematics.weightedVectorSum()
263  : vertexKinematics.vectorSum();
264 
265  if (vtxType != btag::Vertices::RecoVertex) {
266  vars.insert(btau::vertexNTracks,vertexKinematics.numberOfTracks(), true);
267  vars.insert(btau::vertexJetDeltaR,VectorUtil::DeltaR(vertexSum, jetDir), true);
268  }
269 
270  double vertexMass = vertexSum.M();
271  if (vtxType == btag::Vertices::RecoVertex &&
273  const GlobalVector& dir = svInfo.flightDirection(vtx);
274  double vertexPt2 = math::XYZVector(dir.x(), dir.y(), dir.z()).Cross(vertexSum).Mag2() / dir.mag2();
275  vertexMass = std::sqrt(vertexMass * vertexMass + vertexPt2) + std::sqrt(vertexPt2);
276  }
277  vars.insert(btau::vertexMass, vertexMass, true);
278 
279  double varPi = (vertexMass/5.2794) * (vtx_track_ESum /jet_track_ESum); // 5.2794 should be the average B meson mass of PDG! CHECK!!!
280  vars.insert(btau::massVertexEnergyFraction, varPi / (varPi + 0.04), true);
281  double varB = (std::sqrt(5.2794) * vtx_track_ptSum) / ( vertexMass * std::sqrt(jet->pt()));
282  vars.insert(btau::vertexBoostOverSqrtJetPt,varB*varB/(varB*varB + 10.), true);
283 
284  if (allKinematics.numberOfTracks()) {
285  vars.insert(btau::vertexEnergyRatio, vertexSum.E() / allSum.E(), true);
286  }
287  else {
288  vars.insert(btau::vertexEnergyRatio, 1, true);
289  }
290  }
291 
292  reco::PFJet const * pfJet = dynamic_cast<reco::PFJet const *>( &* jet ) ;
293  pat::Jet const * patJet = dynamic_cast<pat::Jet const *>( &* jet ) ;
294  if ( pfJet != nullptr )
295  {
305  vars.insert(btau::muonMultiplicity,pfJet->muonMultiplicity(), true);
309 
310  }
311  else if( patJet != nullptr && patJet->isPFJet() )
312  {
317  vars.insert(btau::muonEnergyFraction,patJet->muonEnergyFraction(), true);
320  vars.insert(btau::photonMultiplicity,patJet->photonMultiplicity(), true);
322  vars.insert(btau::muonMultiplicity,patJet->muonMultiplicity(), true);
326 
327  }
328  else
329  {
332  vars.insert(btau::photonEnergyFraction,0., true);
333  vars.insert(btau::electronEnergyFraction,0., true);
334  vars.insert(btau::muonEnergyFraction,0., true);
337  vars.insert(btau::photonMultiplicity,0, true);
338  vars.insert(btau::electronMultiplicity,0, true);
339  vars.insert(btau::muonMultiplicity,0, true);
340  vars.insert(btau::hadronMultiplicity,0, true);
342  vars.insert(btau::totalMultiplicity,0, true);
343  }
344 }
reco::TrackSelector trackSelector
const double piPlus
Definition: ParticleMasses.h:9
T mag2() const
Definition: PV3DBase.h:66
int photonMultiplicity() const
photonMultiplicity
Definition: PFJet.h:131
const math::XYZTLorentzVector & weightedVectorSum() const
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:377
float photonEnergyFraction() const
photonEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:387
int electronMultiplicity() const
electronMultiplicity
Definition: Jet.h:412
unsigned int pseudoMultiplicityMin
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction (relative to uncorrected jet energy)
Definition: Jet.h:375
void add(const reco::Track &track, double weight=1.0)
int photonMultiplicity() const
photonMultiplicity
Definition: Jet.h:410
float chargedHadronEnergyFraction() const
chargedHadronEnergyFraction
Definition: PFJet.h:100
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
Measurement1D ip2d
Definition: IPTagInfo.h:31
reco::V0Filter pseudoVertexV0Filter
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
int muonMultiplicity() const
muonMultiplicity
Definition: Jet.h:686
bool isNonnull() const
Checks for non-null.
Definition: RefToBase.h:340
int chargedHadronMultiplicity() const
chargedHadronMultiplicity
Definition: Jet.h:406
float photonEnergyFraction() const
photonEnergyFraction
Definition: PFJet.h:108
Jets made from PFObjects.
Definition: PFJet.h:21
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::btag::SortCriteria sortCriterium
reco::TrackSelector trackPseudoSelector
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
def pv(vc)
Definition: MetAnalyzer.py:7
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: PFJet.h:129
float neutralHadronEnergyFraction() const
neutralHadronEnergyFraction
Definition: PFJet.h:104
Measurement1D ip3d
Definition: IPTagInfo.h:32
IterationRange flipIterate(int size, bool vertex) const
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:258
reco::V0Filter trackPairV0Filter
double significance() const
Definition: Measurement1D.h:29
unsigned int numberOfTracks() const
float electronEnergyFraction() const
electronEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:391
#define range_for(i, x)
GlobalPoint closestToJetAxis
Definition: IPTagInfo.h:29
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
double flipValue(double value, bool vertex) const
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:21
float electronEnergyFraction() const
electronEnergyFraction
Definition: PFJet.h:112
Analysis-level calorimeter jet class.
Definition: Jet.h:80
double value() const
Definition: Measurement1D.h:25
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
Definition: PFJet.h:127
Measurement1D distanceToJetAxis
Definition: IPTagInfo.h:33
int muonMultiplicity() const
muonMultiplicity
Definition: PFJet.h:135
fixed size matrix
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
int neutralHadronMultiplicity() const
neutralHadronMultiplicity
Definition: Jet.h:408
const math::XYZTLorentzVector & vectorSum() const
float muonEnergyFraction() const
muonEnergyFraction (relative to corrected jet energy)
Definition: Jet.h:395
dbl *** dir
Definition: mlp_gen.cc:35
float muonEnergyFraction() const
muonEnergyFraction
Definition: PFJet.h:116
T x() const
Definition: PV3DBase.h:62
edm::AssociationVector< reco::JetRefBaseProd, Values > Container
int electronMultiplicity() const
electronMultiplicity
Definition: PFJet.h:133
void insert(const TaggingVariable &variable, bool delayed=false)
unsigned int trackMultiplicityMin
CombinedSVComputer::IterationRange CombinedSVComputer::flipIterate ( int  size,
bool  vertex 
) const
inline
double CombinedSVComputer::flipValue ( double  value,
bool  vertex 
) const
inline

Definition at line 32 of file CombinedSVComputer.cc.

References trackFlip, relativeConstraints::value, and vertexFlip.

Referenced by fillCommonVariables().

33 {
34  return (vertex ? vertexFlip : trackFlip) ? -value : value;
35 }
Definition: value.py:1
TaggingVariableList CombinedSVComputer::operator() ( const reco::TrackIPTagInfo ipInfo,
const reco::SecondaryVertexTagInfo svInfo 
) const
virtual

Definition at line 164 of file CombinedSVComputer.cc.

References reco::TrackKinematics::add(), reco::btau::etaRel(), fillCommonVariables(), reco::TaggingVariableList::finalize(), flipIterate(), reco::Vertex::hasRefittedTracks(), mps_fire::i, reco::TaggingVariableList::insert(), minTrackWeight, reco::TrackBase::momentum(), nTracks(), reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX >::nVertices(), reco::ParticleMasses::piPlus, range_for, reco::Vertex::refittedTrack(), reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX >::secondaryVertex(), mathSSE::sqrt(), HGCalGeometryMode::Square, HiIsolationCommonParameters_cff::track, reco::btau::trackEtaRel, reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::trackWeight(), reco::btau::vertexFitProb, reco::btau::vertexNTracks, extraflags_cff::vtx, and w.

Referenced by CombinedSVSoftLeptonComputer::operator()().

166 {
167  using namespace ROOT::Math;
168 
169  edm::RefToBase<Jet> jet = ipInfo.jet();
170  math::XYZVector jetDir = jet->momentum().Unit();
172 
173  TrackKinematics vertexKinematics;
174 
175  double vtx_track_ptSum = 0.;
176  double vtx_track_ESum = 0.;
177 
178  // the following is specific depending on the type of vertex
179  int vtx = -1;
180  unsigned int numberofvertextracks = 0;
181 
182  IterationRange range = flipIterate(svInfo.nVertices(), true);
183  range_for(i, range) {
184 
185  numberofvertextracks = numberofvertextracks + (svInfo.secondaryVertex(i)).nTracks();
186 
187  const Vertex &vertex = svInfo.secondaryVertex(i);
188  bool hasRefittedTracks = vertex.hasRefittedTracks();
189  for(reco::Vertex::trackRef_iterator track = vertex.tracks_begin(); track != vertex.tracks_end(); track++) {
190  double w = vertex.trackWeight(*track);
191  if (w < minTrackWeight)
192  continue;
193  if (hasRefittedTracks) {
194  const Track actualTrack = vertex.refittedTrack(*track);
195  vertexKinematics.add(actualTrack, w);
196  vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir,actualTrack.momentum()), true);
197  if(vtx < 0) // calculate this only for the first vertex
198  {
199  vtx_track_ptSum += std::sqrt(actualTrack.momentum().Perp2());
200  vtx_track_ESum += std::sqrt(actualTrack.momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
201  }
202  } else {
203  vertexKinematics.add(**track, w);
204  vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir,(*track)->momentum()), true);
205  if(vtx < 0) // calculate this only for the first vertex
206  {
207  vtx_track_ptSum += std::sqrt((*track)->momentum().Perp2());
208  vtx_track_ESum += std::sqrt((*track)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
209  }
210  }
211  }
212 
213  if (vtx < 0) vtx = i;
214  }
215  if(vtx>=0){
216  vars.insert(btau::vertexNTracks, numberofvertextracks, true);
217  vars.insert(btau::vertexFitProb,(svInfo.secondaryVertex(vtx)).normalizedChi2(), true);
218  }
219 
220  // after we collected vertex information we let the common code complete the job
221  fillCommonVariables(vars,vertexKinematics,ipInfo,svInfo,vtx_track_ptSum,vtx_track_ESum);
222 
223  vars.finalize();
224  return vars;
225 }
const double piPlus
Definition: ParticleMasses.h:9
const unsigned int nTracks(const reco::Vertex &sv)
const double w
Definition: UKUtility.cc:23
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
void add(const reco::Track &track, double weight=1.0)
const VTX & secondaryVertex(unsigned int index) const
Track refittedTrack(const TrackBaseRef &track) const
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
bool hasRefittedTracks() const
Checks whether refitted tracks are stored.
Definition: Vertex.h:149
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:714
T sqrt(T t)
Definition: SSEVec.h:18
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:81
IterationRange flipIterate(int size, bool vertex) const
#define range_for(i, x)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void fillCommonVariables(reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
vars
Definition: DeepTauId.cc:77
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
void insert(const TaggingVariable &variable, bool delayed=false)
TaggingVariableList CombinedSVComputer::operator() ( const reco::CandIPTagInfo ipInfo,
const reco::CandSecondaryVertexTagInfo svInfo 
) const
virtual

Definition at line 228 of file CombinedSVComputer.cc.

References reco::TrackKinematics::add(), reco::CompositePtrCandidate::daughterPtrVector(), reco::btau::etaRel(), fillCommonVariables(), reco::TaggingVariableList::finalize(), flipIterate(), mps_fire::i, reco::TaggingVariableList::insert(), reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX >::nVertices(), reco::ParticleMasses::piPlus, range_for, reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX >::secondaryVertex(), mathSSE::sqrt(), HGCalGeometryMode::Square, HiIsolationCommonParameters_cff::track, reco::btau::trackEtaRel, l1t::tracks, reco::btau::vertexFitProb, reco::btau::vertexNTracks, and extraflags_cff::vtx.

230 {
231  using namespace ROOT::Math;
232 
233  edm::RefToBase<Jet> jet = ipInfo.jet();
234  math::XYZVector jetDir = jet->momentum().Unit();
236 
237  TrackKinematics vertexKinematics;
238 
239  double vtx_track_ptSum = 0.;
240  double vtx_track_ESum = 0.;
241 
242  // the following is specific depending on the type of vertex
243  int vtx = -1;
244  unsigned int numberofvertextracks = 0;
245 
246  IterationRange range = flipIterate(svInfo.nVertices(), true);
247  range_for(i, range) {
248 
249  numberofvertextracks = numberofvertextracks + (svInfo.secondaryVertex(i)).numberOfSourceCandidatePtrs();
250 
251  const reco::VertexCompositePtrCandidate &vertex = svInfo.secondaryVertex(i);
252  const std::vector<CandidatePtr> & tracks = vertex.daughterPtrVector();
253  for(std::vector<CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track) {
254  vertexKinematics.add(*track);
255  vars.insert(btau::trackEtaRel, reco::btau::etaRel(jetDir,(*track)->momentum()), true);
256  if(vtx < 0) // calculate this only for the first vertex
257  {
258  vtx_track_ptSum += std::sqrt((*track)->momentum().Perp2());
259  vtx_track_ESum += std::sqrt((*track)->momentum().Mag2() + ROOT::Math::Square(ParticleMasses::piPlus));
260  }
261  }
262 
263  if (vtx < 0) vtx = i;
264  }
265  if(vtx>=0){
266  vars.insert(btau::vertexNTracks, numberofvertextracks, true);
267  vars.insert(btau::vertexFitProb,(svInfo.secondaryVertex(vtx)).vertexNormalizedChi2(), true);
268  }
269 
270  // after we collected vertex information we let the common code complete the job
271  fillCommonVariables(vars,vertexKinematics,ipInfo,svInfo,vtx_track_ptSum,vtx_track_ESum);
272 
273  vars.finalize();
274  return vars;
275 }
const double piPlus
Definition: ParticleMasses.h:9
void add(const reco::Track &track, double weight=1.0)
const VTX & secondaryVertex(unsigned int index) const
double etaRel(const math::XYZVector &dir, const math::XYZVector &track)
virtual const daughters & daughterPtrVector() const
references to daughtes
T sqrt(T t)
Definition: SSEVec.h:18
IterationRange flipIterate(int size, bool vertex) const
#define range_for(i, x)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void fillCommonVariables(reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
vars
Definition: DeepTauId.cc:77
void insert(const TaggingVariable &variable, bool delayed=false)
const reco::btag::TrackIPData& CombinedSVComputer::threshTrack ( const reco::CandIPTagInfo trackIPTagInfo,
const reco::btag::SortCriteria  sort,
const reco::Jet jet,
const GlobalPoint pv 
) const

Referenced by fillCommonVariables(), and flipIterate().

const reco::btag::TrackIPData& CombinedSVComputer::threshTrack ( const reco::TrackIPTagInfo trackIPTagInfo,
const reco::btag::SortCriteria  sort,
const reco::Jet jet,
const GlobalPoint pv 
) const

Member Data Documentation

double CombinedSVComputer::charmCut
private

Definition at line 79 of file CombinedSVComputer.h.

Referenced by flipIterate().

double CombinedSVComputer::minTrackWeight
private

Definition at line 86 of file CombinedSVComputer.h.

Referenced by operator()().

unsigned int CombinedSVComputer::pseudoMultiplicityMin
private

Definition at line 84 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

reco::V0Filter CombinedSVComputer::pseudoVertexV0Filter
private

Definition at line 89 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

reco::btag::SortCriteria CombinedSVComputer::sortCriterium
private

Definition at line 80 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

std::vector<reco::btau::TaggingVariableName> CombinedSVComputer::taggingVariables
private

Definition at line 91 of file CombinedSVComputer.h.

bool CombinedSVComputer::trackFlip
private

Definition at line 77 of file CombinedSVComputer.h.

Referenced by flipIterate(), and flipValue().

unsigned int CombinedSVComputer::trackMultiplicityMin
private

Definition at line 85 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

reco::TrackSelector CombinedSVComputer::trackNoDeltaRSelector
private

Definition at line 82 of file CombinedSVComputer.h.

Referenced by flipIterate().

reco::V0Filter CombinedSVComputer::trackPairV0Filter
private

Definition at line 90 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

reco::TrackSelector CombinedSVComputer::trackPseudoSelector
private

Definition at line 83 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

reco::TrackSelector CombinedSVComputer::trackSelector
private

Definition at line 81 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

bool CombinedSVComputer::useTrackWeights
private

Definition at line 87 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().

bool CombinedSVComputer::vertexFlip
private

Definition at line 78 of file CombinedSVComputer.h.

Referenced by fillCommonVariables(), flipIterate(), and flipValue().

bool CombinedSVComputer::vertexMassCorrection
private

Definition at line 88 of file CombinedSVComputer.h.

Referenced by fillCommonVariables().