CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
CandidateBoostedDoubleSecondaryVertexComputer Class Reference

#include <CandidateBoostedDoubleSecondaryVertexComputer.h>

Inheritance diagram for CandidateBoostedDoubleSecondaryVertexComputer:
JetTagComputer

Public Member Functions

 CandidateBoostedDoubleSecondaryVertexComputer (const edm::ParameterSet &parameters)
 
float discriminator (const TagInfoHelper &tagInfos) const override
 
- Public Member Functions inherited from JetTagComputer
const std::vector< std::string > & getInputLabels () const
 
virtual void initialize (const JetTagComputerRecord &)
 
 JetTagComputer ()
 
 JetTagComputer (const edm::ParameterSet &configuration)
 
float operator() (const reco::BaseTagInfo &info) const
 
float operator() (const TagInfoHelper &helper) const
 
void setupDone ()
 
virtual ~JetTagComputer ()
 

Private Member Functions

void calcNsubjettiness (const reco::JetBaseRef &jet, float &tau1, float &tau2, std::vector< fastjet::PseudoJet > &currentAxes) const
 
void setTracksPV (const reco::CandidatePtr &trackRef, const reco::VertexRef &vertexRef, float &PVweight) const
 
void setTracksPVBase (const reco::TrackRef &trackRef, const reco::VertexRef &vertexRef, float &PVweight) const
 
void vertexKinematics (const reco::VertexCompositePtrCandidate &vertex, reco::TrackKinematics &vertexKinematics) const
 

Private Attributes

const double beta_
 
const double maxSVDeltaRToJet_
 
std::unique_ptr< TMVAEvaluatormvaID
 
fastjet::contrib::Njettiness njettiness_
 
const double R0_
 
const edm::FileInPath weightFile_
 

Additional Inherited Members

- Protected Member Functions inherited from JetTagComputer
virtual float discriminator (const reco::BaseTagInfo &) const
 
void uses (unsigned int id, const std::string &label)
 
void uses (const std::string &label)
 

Detailed Description

Definition at line 16 of file CandidateBoostedDoubleSecondaryVertexComputer.h.

Constructor & Destructor Documentation

CandidateBoostedDoubleSecondaryVertexComputer::CandidateBoostedDoubleSecondaryVertexComputer ( const edm::ParameterSet parameters)

Definition at line 11 of file CandidateBoostedDoubleSecondaryVertexComputer.cc.

References edm::FileInPath::fullPath(), mvaID, JetTagComputer::uses(), makeLayoutFileForGui::variables, and weightFile_.

11  :
12  beta_(parameters.getParameter<double>("beta")),
13  R0_(parameters.getParameter<double>("R0")),
14  njettiness_(fastjet::contrib::OnePass_KT_Axes(), fastjet::contrib::NormalizedMeasure(beta_,R0_)),
15  maxSVDeltaRToJet_(parameters.getParameter<double>("maxSVDeltaRToJet")),
16  weightFile_(parameters.getParameter<edm::FileInPath>("weightFile"))
17 {
18  uses(0, "ipTagInfos");
19  uses(1, "svTagInfos");
20  uses(2, "muonTagInfos");
21  uses(3, "elecTagInfos");
22 
23  mvaID.reset(new TMVAEvaluator());
24 
25  // variable order needs to be the same as in the training
26  std::vector<std::string> variables({"PFLepton_ptrel", "z_ratio1", "tau_dot", "SV_mass_0", "SV_vtx_EnergyRatio_0",
27  "SV_vtx_EnergyRatio_1","PFLepton_IP2D", "tau2/tau1", "nSL", "jetNTracksEtaRel"});
28  std::vector<std::string> spectators({"massGroomed", "flavour", "nbHadrons", "ptGroomed", "etaGroomed"});
29 
30  mvaID->initialize("Color:Silent:Error", "BDTG", weightFile_.fullPath(), variables, spectators,true,false);
31 }
T getParameter(std::string const &) const
void uses(unsigned int id, const std::string &label)
std::string fullPath() const
Definition: FileInPath.cc:165

Member Function Documentation

void CandidateBoostedDoubleSecondaryVertexComputer::calcNsubjettiness ( const reco::JetBaseRef jet,
float &  tau1,
float &  tau2,
std::vector< fastjet::PseudoJet > &  currentAxes 
) const
private

Definition at line 170 of file CandidateBoostedDoubleSecondaryVertexComputer.cc.

References reco::CompositePtrCandidate::daughterPtrVector(), and njettiness_.

Referenced by discriminator().

171 {
172  std::vector<fastjet::PseudoJet> fjParticles;
173 
174  // loop over jet constituents and push them in the vector of FastJet constituents
175  for(const reco::CandidatePtr & daughter : jet->daughterPtrVector())
176  {
177  if ( daughter.isNonnull() && daughter.isAvailable() )
178  fjParticles.push_back( fastjet::PseudoJet( daughter->px(), daughter->py(), daughter->pz(), daughter->energy() ) );
179  else
180  edm::LogWarning("MissingJetConstituent") << "Jet constituent required for N-subjettiness computation is missing!";
181  }
182 
183  // calculate N-subjettiness
184  tau1 = njettiness_.getTau(1, fjParticles);
185  tau2 = njettiness_.getTau(2, fjParticles);
186  currentAxes = njettiness_.currentAxes();
187 }
const daughters & daughterPtrVector() const
references to daughtes
float CandidateBoostedDoubleSecondaryVertexComputer::discriminator ( const TagInfoHelper tagInfos) const
overridevirtual

Reimplemented from JetTagComputer.

Definition at line 34 of file CandidateBoostedDoubleSecondaryVertexComputer.cc.

References reco::TrackKinematics::add(), calcNsubjettiness(), generateEDF::cont, reco::deltaR(), reco::deltaR2(), JetTagComputer::TagInfoHelper::get(), SiPixelRawToDigiRegional_cfi::inputs, metsig::jet, PV3DBase< T, PVType, FrameType >::mag(), maxSVDeltaRToJet_, reco::LeafCandidate::momentum(), mvaID, reco::LeafCandidate::p4(), reco::IPTagInfo< Container, Base >::primaryVertex(), TrackCollections2monitor_cff::selectedTracks, reco::IPTagInfo< Container, Base >::selectedTracks(), setTracksPV(), mathSSE::sqrt(), reco::btag::toTrack(), relativeConstraints::value, vertexKinematics(), reco::btau::vertexNTracks, reco::TrackKinematics::weightedVectorSum(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

35 {
36  // get TagInfos
37  const reco::CandIPTagInfo & ipTagInfo = tagInfo.get<reco::CandIPTagInfo>(0);
41 
42  // default discriminator value
43  float value = -10.;
44 
45  // default variable values
46  float z_ratio = -1. , tau_dot = -1., SV_pt_0 = -1., SV_mass_0 = -1., SV_EnergyRatio_0 = -1., SV_EnergyRatio_1 = -1., tau21 = -1.;
47  int contSV = 0, vertexNTracks = 0;
48  int nSL = 0, nSM = 0, nSE = 0;
49 
50  // get the jet reference
51  const reco::JetBaseRef jet = svTagInfo.jet();
52 
53  std::vector<fastjet::PseudoJet> currentAxes;
54  float tau2, tau1;
55  // calculate N-subjettiness
56  calcNsubjettiness(jet, tau1, tau2, currentAxes);
57  if (tau1 != 0.) tau21 = tau2/tau1;
58 
59  const std::vector<reco::CandidatePtr> & selectedTracks( ipTagInfo.selectedTracks() );
60  size_t trackSize = selectedTracks.size();
61  const reco::VertexRef & vertexRef = ipTagInfo.primaryVertex();
62  reco::TrackKinematics allKinematics;
63 
64  for (size_t itt=0; itt < trackSize; ++itt)
65  {
66  const reco::Track & ptrack = *(reco::btag::toTrack(selectedTracks[itt]));
67  const reco::CandidatePtr ptrackRef = selectedTracks[itt];
68 
69  float track_PVweight = 0.;
70  setTracksPV(ptrackRef, vertexRef, track_PVweight);
71  if (track_PVweight>0.) { allKinematics.add(ptrack, track_PVweight); }
72  }
73 
74  math::XYZVector jetDir = jet->momentum().Unit();
75 
76  std::map<double, size_t> VTXmass;
77  for (size_t vtx = 0; vtx < svTagInfo.nVertices(); ++vtx)
78  {
79  vertexNTracks += (svTagInfo.secondaryVertex(vtx)).numberOfSourceCandidatePtrs();
80  GlobalVector flightDir = svTagInfo.flightDirection(vtx);
81  if (reco::deltaR2(flightDir, jetDir)<(maxSVDeltaRToJet_*maxSVDeltaRToJet_))
82  {
83  ++contSV;
84  VTXmass[svTagInfo.secondaryVertex(vtx).p4().mass()]=vtx;
85  }
86  }
87 
88  int cont=0;
89  GlobalVector flightDir_0, flightDir_1;
90  reco::Candidate::LorentzVector SV_p4_0 , SV_p4_1;
91  for ( std::map<double, size_t>::reverse_iterator iVtx=VTXmass.rbegin(); iVtx!=VTXmass.rend(); ++iVtx)
92  {
93  ++cont;
94  const reco::VertexCompositePtrCandidate &vertex = svTagInfo.secondaryVertex(iVtx->second);
95  reco::TrackKinematics vtxKinematics;
96  vertexKinematics(vertex, vtxKinematics);
97  math::XYZTLorentzVector allSum = allKinematics.weightedVectorSum();
98  math::XYZTLorentzVector vertexSum = vtxKinematics.weightedVectorSum();
99  if (cont==1)
100  {
101  SV_mass_0 = vertex.p4().mass() ;
102  SV_EnergyRatio_0 = vertexSum.E() / allSum.E();
103  SV_pt_0 = vertex.p4().pt();
104  flightDir_0 = svTagInfo.flightDirection(iVtx->second);
105  SV_p4_0 = vertex.p4();
106 
107  if (reco::deltaR2(flightDir_0,currentAxes[1])<reco::deltaR2(flightDir_0,currentAxes[0]))
108  tau_dot = (currentAxes[1].px()*flightDir_0.x()+currentAxes[1].py()*flightDir_0.y()+currentAxes[1].pz()*flightDir_0.z())/(sqrt(currentAxes[1].modp2())*flightDir_0.mag());
109  else
110  tau_dot = (currentAxes[0].px()*flightDir_0.x()+currentAxes[0].py()*flightDir_0.y()+currentAxes[0].pz()*flightDir_0.z())/(sqrt(currentAxes[0].modp2())*flightDir_0.mag());
111  }
112  if (cont==2)
113  {
114  SV_EnergyRatio_1= vertexSum.E() / allSum.E();
115  flightDir_1 = svTagInfo.flightDirection(iVtx->second);
116  SV_p4_1 = vertex.p4();
117  z_ratio = reco::deltaR(flightDir_0,flightDir_1)*SV_pt_0/(SV_p4_0+SV_p4_1).mass();
118  break;
119  }
120  }
121 
122  nSM = muonTagInfo.leptons();
123  nSE = elecTagInfo.leptons();
124  nSL = nSM + nSE;
125 
126  float PFLepton_ptrel = -1., PFLepton_IP2D = -1.;
127 
128  // PFMuon information
129  for (size_t leptIdx = 0; leptIdx < muonTagInfo.leptons() ; ++leptIdx)
130  {
131  float PFMuon_ptrel = (muonTagInfo.properties(leptIdx).ptRel);
132  if (PFMuon_ptrel > PFLepton_ptrel )
133  {
134  PFLepton_ptrel = PFMuon_ptrel;
135  PFLepton_IP2D = (muonTagInfo.properties(leptIdx).sip2d);
136  }
137  }
138 
139  // PFElectron information
140  for (size_t leptIdx = 0; leptIdx < elecTagInfo.leptons() ; ++leptIdx)
141  {
142  float PFElectron_ptrel = (elecTagInfo.properties(leptIdx).ptRel);
143  if (PFElectron_ptrel > PFLepton_ptrel )
144  {
145  PFLepton_ptrel = PFElectron_ptrel;
146  PFLepton_IP2D = (elecTagInfo.properties(leptIdx).sip2d);
147  }
148  }
149 
150  std::map<std::string,float> inputs;
151  inputs["z_ratio1"] = z_ratio;
152  inputs["tau_dot"] = tau_dot;
153  inputs["SV_mass_0"] = SV_mass_0;
154  inputs["SV_vtx_EnergyRatio_0"] = SV_EnergyRatio_0;
155  inputs["SV_vtx_EnergyRatio_1"] = SV_EnergyRatio_1;
156  inputs["jetNTracksEtaRel"] = vertexNTracks;
157  inputs["PFLepton_ptrel"] = PFLepton_ptrel;
158  inputs["PFLepton_IP2D"] = PFLepton_IP2D;
159  inputs["nSL"] = nSL;
160  inputs["tau2/tau1"] = tau21;
161 
162  // evaluate the MVA
163  value = mvaID->evaluate(inputs);
164 
165  // return the final discriminator value
166  return value;
167 }
const math::XYZTLorentzVector & weightedVectorSum() const
tuple cont
load Luminosity info ##
Definition: generateEDF.py:622
void add(const reco::Track &track, double weight=1.0)
void setTracksPV(const reco::CandidatePtr &trackRef, const reco::VertexRef &vertexRef, float &PVweight) const
virtual Vector momentum() const
spatial momentum vector
T y() const
Definition: PV3DBase.h:63
const Container & selectedTracks() const
Definition: IPTagInfo.h:101
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
const edm::Ref< VertexCollection > & primaryVertex() const
Definition: IPTagInfo.h:133
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T mag() const
Definition: PV3DBase.h:67
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
void vertexKinematics(const reco::VertexCompositePtrCandidate &vertex, reco::TrackKinematics &vertexKinematics) const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
void calcNsubjettiness(const reco::JetBaseRef &jet, float &tau1, float &tau2, std::vector< fastjet::PseudoJet > &currentAxes) const
T x() const
Definition: PV3DBase.h:62
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
void CandidateBoostedDoubleSecondaryVertexComputer::setTracksPV ( const reco::CandidatePtr trackRef,
const reco::VertexRef vertexRef,
float &  PVweight 
) const
private

Definition at line 213 of file CandidateBoostedDoubleSecondaryVertexComputer.cc.

References pat::PackedCandidate::fromPV(), edm::Ptr< T >::get(), pat::PackedCandidate::PVUsedInFit, setTracksPVBase(), and reco::PFCandidate::trackRef().

Referenced by discriminator().

214 {
215  PVweight = 0.;
216 
217  const pat::PackedCandidate * pcand = dynamic_cast<const pat::PackedCandidate *>(trackRef.get());
218 
219  if(pcand) // MiniAOD case
220  {
221  if( pcand->fromPV() == pat::PackedCandidate::PVUsedInFit )
222  {
223  PVweight = 1.;
224  }
225  }
226  else
227  {
228  const reco::PFCandidate * pfcand = dynamic_cast<const reco::PFCandidate *>(trackRef.get());
229 
230  setTracksPVBase(pfcand->trackRef(), vertexRef, PVweight);
231  }
232 }
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:160
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:433
const PVAssoc fromPV(size_t ipv=0) const
void setTracksPVBase(const reco::TrackRef &trackRef, const reco::VertexRef &vertexRef, float &PVweight) const
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
void CandidateBoostedDoubleSecondaryVertexComputer::setTracksPVBase ( const reco::TrackRef trackRef,
const reco::VertexRef vertexRef,
float &  PVweight 
) const
private

Definition at line 190 of file CandidateBoostedDoubleSecondaryVertexComputer.cc.

References reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), and reco::Vertex::trackWeight().

Referenced by setTracksPV().

191 {
192  PVweight = 0.;
193 
194  const reco::TrackBaseRef trackBaseRef( trackRef );
195 
197 
198  const reco::Vertex & vtx = *(vertexRef);
199  // loop over tracks in vertices
200  for(IT it=vtx.tracks_begin(); it!=vtx.tracks_end(); ++it)
201  {
202  const reco::TrackBaseRef & baseRef = *it;
203  // one of the tracks in the vertex is the same as the track considered in the function
204  if( baseRef == trackBaseRef )
205  {
206  PVweight = vtx.trackWeight(baseRef);
207  break;
208  }
209  }
210 }
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:44
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
std::vector< LinkConnSpec >::const_iterator IT
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:37
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:39
void CandidateBoostedDoubleSecondaryVertexComputer::vertexKinematics ( const reco::VertexCompositePtrCandidate vertex,
reco::TrackKinematics vertexKinematics 
) const
private

Definition at line 235 of file CandidateBoostedDoubleSecondaryVertexComputer.cc.

References reco::TrackKinematics::add(), reco::CompositePtrCandidate::daughterPtrVector(), and testEve_cfg::tracks.

Referenced by discriminator().

236 {
237  const std::vector<reco::CandidatePtr> & tracks = vertex.daughterPtrVector();
238 
239  for(std::vector<reco::CandidatePtr>::const_iterator track = tracks.begin(); track != tracks.end(); ++track) {
240  const reco::Track& mytrack = *(*track)->bestTrack();
241  vtxKinematics.add(mytrack, 1.0);
242  }
243 }
tuple tracks
Definition: testEve_cfg.py:39
const daughters & daughterPtrVector() const
references to daughtes

Member Data Documentation

const double CandidateBoostedDoubleSecondaryVertexComputer::beta_
private
const double CandidateBoostedDoubleSecondaryVertexComputer::maxSVDeltaRToJet_
private

Definition at line 34 of file CandidateBoostedDoubleSecondaryVertexComputer.h.

Referenced by discriminator().

std::unique_ptr<TMVAEvaluator> CandidateBoostedDoubleSecondaryVertexComputer::mvaID
private
fastjet::contrib::Njettiness CandidateBoostedDoubleSecondaryVertexComputer::njettiness_
private

Definition at line 32 of file CandidateBoostedDoubleSecondaryVertexComputer.h.

Referenced by calcNsubjettiness().

const double CandidateBoostedDoubleSecondaryVertexComputer::R0_
private
const edm::FileInPath CandidateBoostedDoubleSecondaryVertexComputer::weightFile_
private