CMS 3D CMS Logo

deep_helpers.h
Go to the documentation of this file.
1 #ifndef RecoBTag_FeatureTools_deep_helpers_h
2 #define RecoBTag_FeatureTools_deep_helpers_h
3 
6 
9 
11 
14 //#include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
15 
18 
21 
23 
24 #include <iostream>
25 #include <fstream>
26 #include <algorithm>
27 #include <numeric>
28 #include <nlohmann/json.hpp>
29 
30 namespace btagbtvdeep {
31 
32  // remove infs and NaNs with value (adapted from DeepNTuples)
33  const float catch_infs(const float in, const float replace_value = 0.);
34 
35  // remove infs/NaN and bound (adapted from DeepNTuples)
36  const float catch_infs_and_bound(const float in,
37  const float replace_value,
38  const float lowerbound,
39  const float upperbound,
40  const float offset = 0.,
41  const bool use_offsets = true);
42 
43  // 2D distance between SV and PV (adapted from DeepNTuples)
45 
46  //3D distance between SV and PV (adapted from DeepNTuples)
48 
49  // dot product between SV and PV (adapted from DeepNTuples)
51 
52  // helper to order vertices by significance (adapted from DeepNTuples)
53  template <typename SVType, typename PVType>
54  bool sv_vertex_comparator(const SVType &sva, const SVType &svb, const PVType &pv) {
55  auto adxy = vertexDxy(sva, pv);
56  auto bdxy = vertexDxy(svb, pv);
57  float aval = adxy.value();
58  float bval = bdxy.value();
59  float aerr = adxy.error();
60  float berr = bdxy.error();
61 
62  float asig = catch_infs(aval / aerr, 0.);
63  float bsig = catch_infs(bval / berr, 0.);
64  return bsig < asig;
65  }
66 
67  // write tagging variables to vector (adapted from DeepNTuples)
68  template <typename T>
70  std::vector<T> vals = from.getList(name, false);
71  size_t size = std::min(vals.size(), max);
72  if (size > 0) {
73  for (size_t i = 0; i < vals.size(); i++) {
74  to[i] = catch_infs(vals[i], -0.1);
75  }
76  }
77  return size;
78  }
79 
80  // compute minimum dr between SVs and a candidate (from DeepNTuples, now polymorphic)
81  float mindrsvpfcand(const std::vector<reco::VertexCompositePtrCandidate> &svs,
82  const reco::Candidate *cand,
83  float mindr = 0.4);
84 
85  // compute minimum distance between SVs and a candidate (from DeepNTuples, now polymorphic)
86  float mindistsvpfcand(const std::vector<reco::VertexCompositePtrCandidate> &svs, const reco::TransientTrack track);
87 
88  // mimic the calculation in PackedCandidate
89  float vtx_ass_from_pfcand(const reco::PFCandidate &pfcand, int pv_ass_quality, const reco::VertexRef &pv);
92 
93  std::pair<float, float> getDRSubjetFeatures(const reco::Jet &jet, const reco::Candidate *cand);
94 
95  // struct to hold preprocessing parameters
97  struct VarInfo {
98  VarInfo() {}
99  VarInfo(float median, float norm_factor, float replace_inf_value, float lower_bound, float upper_bound, float pad)
100  : center(median),
105  pad(pad) {}
106  float center = 0;
107  float norm_factor = 1;
108  float replace_inf_value = 0;
109  float lower_bound = -5;
110  float upper_bound = 5;
111  float pad = 0;
112  };
113 
114  unsigned min_length = 0;
115  unsigned max_length = 0;
116  std::vector<std::string> var_names;
117  std::unordered_map<std::string, VarInfo> var_info_map;
118 
119  VarInfo info(const std::string &name) const { return var_info_map.at(name); }
120  };
121 
122  int center_norm_pad(const std::vector<float> &input,
123  float center,
124  float scale,
125  unsigned min_length,
126  unsigned max_length,
127  std::vector<float> &datavec,
128  int startval,
129  float pad_value = 0,
130  float replace_inf_value = 0,
131  float min = 0,
132  float max = -1);
133 
134  int center_norm_pad_halfRagged(const std::vector<float> &input,
135  float center,
136  float scale,
137  unsigned target_length,
138  std::vector<float> &datavec,
139  int startval,
140  float pad_value = 0,
141  float replace_inf_value = 0,
142  float min = 0,
143  float max = -1);
144 
145  void ParticleNetConstructor(const edm::ParameterSet &Config_,
146  bool doExtra,
147  std::vector<std::string> &input_names_,
148  std::unordered_map<std::string, PreprocessParams> &prep_info_map_,
149  std::vector<std::vector<int64_t>> &input_shapes_,
150  std::vector<unsigned> &input_sizes_,
151  cms::Ort::FloatArrays *data_);
152 
153 } // namespace btagbtvdeep
154 #endif //RecoBTag_FeatureTools_deep_helpers_h
int center_norm_pad_halfRagged(const std::vector< float > &input, float center, float scale, unsigned target_length, std::vector< float > &datavec, int startval, float pad_value=0, float replace_inf_value=0, float min=0, float max=-1)
float quality_from_pfcand(const reco::PFCandidate &pfcand)
float vertexDdotP(const reco::VertexCompositePtrCandidate &sv, const reco::Vertex &pv)
Definition: deep_helpers.cc:78
std::vector< std::string > var_names
Definition: deep_helpers.h:116
Measurement1D vertexD3d(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv)
Definition: deep_helpers.cc:69
VarInfo info(const std::string &name) const
Definition: deep_helpers.h:119
std::unordered_map< std::string, VarInfo > var_info_map
Definition: deep_helpers.h:117
const float catch_infs(const float in, const float replace_value=0.)
Definition: deep_helpers.cc:34
Base class for all types of Jets.
Definition: Jet.h:20
int dump_vector(reco::TaggingVariableList &from, T *to, reco::btau::TaggingVariableName name, const size_t max)
Definition: deep_helpers.h:69
std::pair< float, float > getDRSubjetFeatures(const reco::Jet &jet, const reco::Candidate *cand)
std::vector< std::vector< float > > FloatArrays
Definition: ONNXRuntime.h:23
static std::string const input
Definition: EdmProvDump.cc:50
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:43
void ParticleNetConstructor(const edm::ParameterSet &Config_, bool doExtra, std::vector< std::string > &input_names_, std::unordered_map< std::string, PreprocessParams > &prep_info_map_, std::vector< std::vector< int64_t >> &input_shapes_, std::vector< unsigned > &input_sizes_, cms::Ort::FloatArrays *data_)
float mindrsvpfcand(const std::vector< reco::VertexCompositePtrCandidate > &svs, const reco::Candidate *cand, float mindr=0.4)
Definition: deep_helpers.cc:85
std::vector< TaggingValue > getList(TaggingVariableName tag, bool throwOnEmptyList=true) const
Measurement1D vertexDxy(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv)
Definition: deep_helpers.cc:60
float vtx_ass_from_pfcand(const reco::PFCandidate &pfcand, int pv_ass_quality, const reco::VertexRef &pv)
float mindistsvpfcand(const std::vector< reco::VertexCompositePtrCandidate > &svs, const reco::TransientTrack track)
Definition: deep_helpers.cc:98
bool sv_vertex_comparator(const SVType &sva, const SVType &svb, const PVType &pv)
Definition: deep_helpers.h:54
int center_norm_pad(const std::vector< float > &input, float center, float scale, unsigned min_length, unsigned max_length, std::vector< float > &datavec, int startval, float pad_value=0, float replace_inf_value=0, float min=0, float max=-1)
float lost_inner_hits_from_pfcand(const reco::PFCandidate &pfcand)
T median(std::vector< T > values)
Definition: median.h:16
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
VarInfo(float median, float norm_factor, float replace_inf_value, float lower_bound, float upper_bound, float pad)
Definition: deep_helpers.h:99
long double T