CMS 3D CMS Logo

ChargedCandidateConverter.h
Go to the documentation of this file.
1 #ifndef RecoBTag_FeatureTools_ChargedCandidateConverter_h
2 #define RecoBTag_FeatureTools_ChargedCandidateConverter_h
3 
7 
11 
12 namespace btagbtvdeep {
13 
14  template <typename CandidateType>
15  void commonCandidateToFeatures(const CandidateType* c_pf,
16  const reco::Jet& jet,
17  const TrackInfoBuilder& track_info,
18  const float& drminpfcandsv,
19  const float& jetR,
20  ChargedCandidateFeatures& c_pf_features,
21  const bool flip = false) {
22  float trackSip2dVal = track_info.getTrackSip2dVal();
23  float trackSip2dSig = track_info.getTrackSip2dSig();
24  float trackSip3dVal = track_info.getTrackSip3dVal();
25  float trackSip3dSig = track_info.getTrackSip3dSig();
26  if (flip == true) {
31  }
32 
33  c_pf_features.deltaR = reco::deltaR(*c_pf, jet);
34  c_pf_features.ptrel = catch_infs_and_bound(c_pf->pt() / jet.pt(), 0, -1, 0, -1);
35  c_pf_features.ptrel_noclip = c_pf->pt() / jet.pt();
36  c_pf_features.erel = c_pf->energy() / jet.energy();
37  const float etasign = jet.eta() > 0 ? 1 : -1;
38  c_pf_features.etarel = etasign * (c_pf->eta() - jet.eta());
39 
40  c_pf_features.btagPf_trackEtaRel = catch_infs_and_bound(track_info.getTrackEtaRel(), 0, -5, 15);
41  c_pf_features.btagPf_trackPtRel = catch_infs_and_bound(track_info.getTrackPtRel(), 0, -1, 4);
42  c_pf_features.btagPf_trackPPar = catch_infs_and_bound(track_info.getTrackPPar(), 0, -1e5, 1e5);
43  c_pf_features.btagPf_trackDeltaR = catch_infs_and_bound(track_info.getTrackDeltaR(), 0, -5, 5);
44  c_pf_features.btagPf_trackPtRatio = catch_infs_and_bound(track_info.getTrackPtRatio(), 0, -1, 10);
45  c_pf_features.btagPf_trackPParRatio = catch_infs_and_bound(track_info.getTrackPParRatio(), 0, -10, 100);
46  c_pf_features.btagPf_trackSip3dVal = catch_infs_and_bound(trackSip3dVal, 0, -1, 1e5);
48  c_pf_features.btagPf_trackSip2dVal = catch_infs_and_bound(trackSip2dVal, 0, -1, 70);
50  c_pf_features.btagPf_trackJetDistVal = catch_infs_and_bound(track_info.getTrackJetDistVal(), 0, -20, 1);
51 
52  c_pf_features.drminsv = catch_infs_and_bound(drminpfcandsv, 0, -1. * jetR, 0, -1. * jetR);
53 
54  // subjet related
55  const auto* patJet = dynamic_cast<const pat::Jet*>(&jet);
56  if (!patJet) {
57  throw edm::Exception(edm::errors::InvalidReference) << "Input is not a pat::Jet.";
58  }
59 
60  if (patJet->nSubjetCollections() > 0) {
61  auto subjets = patJet->subjets();
62  // sort by pt
63  std::sort(subjets.begin(), subjets.end(), [](const edm::Ptr<pat::Jet>& p1, const edm::Ptr<pat::Jet>& p2) {
64  return p1->pt() > p2->pt();
65  });
66  c_pf_features.drsubjet1 = !subjets.empty() ? reco::deltaR(*c_pf, *subjets.at(0)) : -1;
67  c_pf_features.drsubjet2 = subjets.size() > 1 ? reco::deltaR(*c_pf, *subjets.at(1)) : -1;
68  } else {
69  // AK4 jets don't have subjets
70  c_pf_features.drsubjet1 = -1;
71  c_pf_features.drsubjet2 = -1;
72  }
73  }
74 
76  const pat::Jet& jet,
77  const TrackInfoBuilder& track_info,
78  const float drminpfcandsv,
79  const float jetR,
80  ChargedCandidateFeatures& c_pf_features,
81  const bool flip = false);
82 
84  const reco::Jet& jet,
85  const TrackInfoBuilder& track_info,
86  const float drminpfcandsv,
87  const float jetR,
88  const float puppiw,
89  const int pv_ass_quality,
90  const reco::VertexRef& pv,
91  ChargedCandidateFeatures& c_pf_features,
92  const bool flip = false);
93 
94 } // namespace btagbtvdeep
95 
96 #endif //RecoBTag_FeatureTools_ChargedCandidateConverter_h
btagbtvdeep::TrackInfoBuilder::getTrackSip3dVal
const float getTrackSip3dVal() const
Definition: TrackInfoBuilder.h:33
btagbtvdeep
Definition: BoostedDoubleSVTagInfoFeatures.h:4
btagbtvdeep::recoCandidateToFeatures
void recoCandidateToFeatures(const reco::PFCandidate *c_pf, const reco::Jet &jet, const TrackInfoBuilder &track_info, const float drminpfcandsv, const float jetR, const float puppiw, const int pv_ass_quality, const reco::VertexRef &pv, ChargedCandidateFeatures &c_pf_features, const bool flip=false)
Definition: ChargedCandidateConverter.cc:37
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
btagbtvdeep::TrackInfoBuilder::getTrackSip2dSig
const float getTrackSip2dSig() const
Definition: TrackInfoBuilder.h:30
btagbtvdeep::commonCandidateToFeatures
void commonCandidateToFeatures(const CandidateType *c_pf, const reco::Jet &jet, const TrackInfoBuilder &track_info, const float &drminpfcandsv, const float &jetR, ChargedCandidateFeatures &c_pf_features, const bool flip=false)
Definition: ChargedCandidateConverter.h:15
edm::errors::InvalidReference
Definition: EDMException.h:39
PFCandidate.h
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackPParRatio
float btagPf_trackPParRatio
Definition: ChargedCandidateFeatures.h:20
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackSip2dVal
float btagPf_trackSip2dVal
Definition: ChargedCandidateFeatures.h:23
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackEtaRel
float btagPf_trackEtaRel
Definition: ChargedCandidateFeatures.h:15
btagbtvdeep::ChargedCandidateFeatures::etarel
float etarel
Definition: ChargedCandidateFeatures.h:11
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackPtRel
float btagPf_trackPtRel
Definition: ChargedCandidateFeatures.h:16
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackJetDistVal
float btagPf_trackJetDistVal
Definition: ChargedCandidateFeatures.h:26
LHEJetFilter_cfi.jetR
jetR
Definition: LHEJetFilter_cfi.py:5
btagbtvdeep::TrackInfoBuilder
Definition: TrackInfoBuilder.h:12
btagbtvdeep::packedCandidateToFeatures
void packedCandidateToFeatures(const pat::PackedCandidate *c_pf, const pat::Jet &jet, const TrackInfoBuilder &track_info, const float drminpfcandsv, const float jetR, ChargedCandidateFeatures &c_pf_features, const bool flip=false)
Definition: ChargedCandidateConverter.cc:5
reco::btau::trackSip2dVal
Definition: TaggingVariable.h:54
edm::Ref< VertexCollection >
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackSip2dSig
float btagPf_trackSip2dSig
Definition: ChargedCandidateFeatures.h:24
btagbtvdeep::ChargedCandidateFeatures::ptrel
float ptrel
Definition: ChargedCandidateFeatures.h:8
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackPtRatio
float btagPf_trackPtRatio
Definition: ChargedCandidateFeatures.h:17
pat::Jet
Analysis-level calorimeter jet class.
Definition: Jet.h:77
pfNegativeDeepFlavourTagInfos_cfi.flip
flip
Definition: pfNegativeDeepFlavourTagInfos_cfi.py:8
btagbtvdeep::TrackInfoBuilder::getTrackPParRatio
const float getTrackPParRatio() const
Definition: TrackInfoBuilder.h:27
p2
double p2[4]
Definition: TauolaWrapper.h:90
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackPPar
float btagPf_trackPPar
Definition: ChargedCandidateFeatures.h:18
btagbtvdeep::catch_infs_and_bound
const float catch_infs_and_bound(const float in, const float replace_value, const float lowerbound, const float upperbound, const float offset=0., const bool use_offsets=true)
Definition: deep_helpers.cc:32
btagbtvdeep::TrackInfoBuilder::getTrackEtaRel
const float getTrackEtaRel() const
Definition: TrackInfoBuilder.h:22
btagbtvdeep::TrackInfoBuilder::getTrackSip2dVal
const float getTrackSip2dVal() const
Definition: TrackInfoBuilder.h:31
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackDeltaR
float btagPf_trackDeltaR
Definition: ChargedCandidateFeatures.h:19
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackSip3dSig
float btagPf_trackSip3dSig
Definition: ChargedCandidateFeatures.h:22
btagbtvdeep::ChargedCandidateFeatures::ptrel_noclip
float ptrel_noclip
Definition: ChargedCandidateFeatures.h:9
ChargedCandidateFeatures.h
btagbtvdeep::ChargedCandidateFeatures::drsubjet2
float drsubjet2
Definition: ChargedCandidateFeatures.h:29
btagbtvdeep::ChargedCandidateFeatures::drminsv
float drminsv
Definition: ChargedCandidateFeatures.h:39
btagbtvdeep::TrackInfoBuilder::getTrackPtRel
const float getTrackPtRel() const
Definition: TrackInfoBuilder.h:29
reco::btau::trackSip2dSig
Definition: TaggingVariable.h:55
reco::btau::trackSip3dSig
Definition: TaggingVariable.h:57
pat::PackedCandidate
Definition: PackedCandidate.h:22
PackedCandidate.h
TrackInfoBuilder.h
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
p1
double p1[4]
Definition: TauolaWrapper.h:89
btagbtvdeep::ChargedCandidateFeatures::deltaR
float deltaR
Definition: ChargedCandidateFeatures.h:34
btagbtvdeep::TrackInfoBuilder::getTrackPPar
const float getTrackPPar() const
Definition: TrackInfoBuilder.h:26
Jet.h
reco::btau::trackSip3dVal
Definition: TaggingVariable.h:56
edm::Ptr< pat::Jet >
btagbtvdeep::ChargedCandidateFeatures::drsubjet1
float drsubjet1
Definition: ChargedCandidateFeatures.h:28
btagbtvdeep::ChargedCandidateFeatures
Definition: ChargedCandidateFeatures.h:6
metsig::jet
Definition: SignAlgoResolutions.h:47
btagbtvdeep::TrackInfoBuilder::getTrackJetDistVal
const float getTrackJetDistVal() const
Definition: TrackInfoBuilder.h:24
Exception
Definition: hltDiff.cc:246
reco::deltaR
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
btagbtvdeep::ChargedCandidateFeatures::erel
float erel
Definition: ChargedCandidateFeatures.h:10
deep_helpers.h
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
btagbtvdeep::ChargedCandidateFeatures::btagPf_trackSip3dVal
float btagPf_trackSip3dVal
Definition: ChargedCandidateFeatures.h:21
btagbtvdeep::TrackInfoBuilder::getTrackDeltaR
const float getTrackDeltaR() const
Definition: TrackInfoBuilder.h:20
btagbtvdeep::TrackInfoBuilder::getTrackSip3dSig
const float getTrackSip3dSig() const
Definition: TrackInfoBuilder.h:32
vertexPlots.e4
e4
Definition: vertexPlots.py:64
btagbtvdeep::TrackInfoBuilder::getTrackPtRatio
const float getTrackPtRatio() const
Definition: TrackInfoBuilder.h:28