CMS 3D CMS Logo

NeutralCandidateConverter.h
Go to the documentation of this file.
1 #ifndef RecoBTag_FeatureTools_NeutralCandidateConverter_h
2 #define RecoBTag_FeatureTools_NeutralCandidateConverter_h
3 
6 
11 
12 
13 namespace btagbtvdeep {
14 
15 
16 
18  const pat::Jet & jet,
19  const float drminpfcandsv, const float jetR,
20  NeutralCandidateFeatures & n_pf_features) ;
21 
22 
24  const reco::Jet & jet,
25  const float drminpfcandsv, const float jetR, const float puppiw,
26  NeutralCandidateFeatures & n_pf_features) ;
27 
28 
29  template <typename CandidateType>
30  static void commonCandidateToFeatures(const CandidateType* n_pf,
31  const reco::Jet& jet,
32  const float& drminpfcandsv,
33  const float& jetR,
34  NeutralCandidateFeatures& n_pf_features) {
35  const auto* patJet = dynamic_cast<const pat::Jet*>(&jet);
36 
37  if (!patJet) {
38  throw edm::Exception(edm::errors::InvalidReference) << "Input is not a pat::Jet.";
39  }
40  // Do Subjets
41  if (patJet->nSubjetCollections() > 0) {
42  auto subjets = patJet->subjets();
43  // sort by pt
44  std::sort(subjets.begin(), subjets.end(), [](const edm::Ptr<pat::Jet>& p1, const edm::Ptr<pat::Jet>& p2) {
45  return p1->pt() > p2->pt();
46  });
47  n_pf_features.drsubjet1 = !subjets.empty() ? reco::deltaR(*n_pf, *subjets.at(0)) : -1;
48  n_pf_features.drsubjet2 = subjets.size() > 1 ? reco::deltaR(*n_pf, *subjets.at(1)) : -1;
49  } else {
50  n_pf_features.drsubjet1 = -1;
51  n_pf_features.drsubjet2 = -1;
52  }
53 
54  // Jet relative vars
55  n_pf_features.ptrel = catch_infs_and_bound(n_pf->pt() / jet.pt(), 0, -1, 0, -1);
56  n_pf_features.ptrel_noclip = n_pf->pt() / jet.pt();
57  n_pf_features.deltaR = catch_infs_and_bound(reco::deltaR(*n_pf, jet), 0, -0.6, 0, -0.6);
58  n_pf_features.deltaR_noclip = reco::deltaR(*n_pf, jet);
59  n_pf_features.erel = n_pf->energy() / jet.energy();
60  n_pf_features.isGamma = 0;
61  if(std::abs(n_pf->pdgId())==22) n_pf_features.isGamma = 1;
62 
63 
64  n_pf_features.drminsv = catch_infs_and_bound(drminpfcandsv,
65  0,-1.*jetR,0,-1.*jetR);
66 
67  }
68 
69 
70 
71 }
72 
73 #endif //RecoBTag_FeatureTools_NeutralCandidateConverter_h
Base class for all types of Jets.
Definition: Jet.h:20
double pt() const final
transverse momentum
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:34
double energy() const final
energy
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
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)
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)
Analysis-level calorimeter jet class.
Definition: Jet.h:80
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
double p1[4]
Definition: TauolaWrapper.h:89
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)