CMS 3D CMS Logo

deep_helpers.cc
Go to the documentation of this file.
4 
5 namespace btagbtvdeep {
6 
7  constexpr static int qualityMap[8] = {1, 0, 1, 1, 4, 4, 5, 6};
8 
16  muonFlagsMask = 0x0600,
18  };
19 
20  // remove infs and NaNs with value (adapted from DeepNTuples)
21  const float catch_infs(const float in, const float replace_value) {
22  if (edm::isNotFinite(in))
23  return replace_value;
24  if (in < -1e32 || in > 1e32)
25  return replace_value;
26  return in;
27  }
28 
29  // remove infs/NaN and bound (adapted from DeepNTuples)
30  const float catch_infs_and_bound(const float in,
31  const float replace_value,
32  const float lowerbound,
33  const float upperbound,
34  const float offset,
35  const bool use_offsets) {
36  float withoutinfs = catch_infs(in, replace_value);
37  if (withoutinfs + offset < lowerbound)
38  return lowerbound;
39  if (withoutinfs + offset > upperbound)
40  return upperbound;
41  if (use_offsets)
42  withoutinfs += offset;
43  return withoutinfs;
44  }
45 
46  // 2D distance between SV and PV (adapted from DeepNTuples)
48  VertexDistanceXY dist;
50  svcand.fillVertexCovariance(csv);
51  reco::Vertex svtx(svcand.vertex(), csv);
52  return dist.distance(svtx, pv);
53  }
54 
55  //3D distance between SV and PV (adapted from DeepNTuples)
57  VertexDistance3D dist;
59  svcand.fillVertexCovariance(csv);
60  reco::Vertex svtx(svcand.vertex(), csv);
61  return dist.distance(svtx, pv);
62  }
63 
64  // dot product between SV and PV (adapted from DeepNTuples)
66  reco::Candidate::Vector p = sv.momentum();
67  reco::Candidate::Vector d(sv.vx() - pv.x(), sv.vy() - pv.y(), sv.vz() - pv.z());
68  return p.Unit().Dot(d.Unit());
69  }
70 
71  // compute minimum dr between SVs and a candidate (from DeepNTuples, now polymorphic)
72  float mindrsvpfcand(const std::vector<reco::VertexCompositePtrCandidate> &svs,
73  const reco::Candidate *cand,
74  float mindr) {
75  for (unsigned int i0 = 0; i0 < svs.size(); ++i0) {
76  float tempdr = reco::deltaR(svs[i0], *cand);
77  if (tempdr < mindr) {
78  mindr = tempdr;
79  }
80  }
81  return mindr;
82  }
83 
84  // instantiate template
85  template bool sv_vertex_comparator<reco::VertexCompositePtrCandidate, reco::Vertex>(
87 
88  float vtx_ass_from_pfcand(const reco::PFCandidate &pfcand, int pv_ass_quality, const reco::VertexRef &pv) {
89  float vtx_ass = pat::PackedCandidate::PVAssociationQuality(qualityMap[pv_ass_quality]);
90  if (pfcand.trackRef().isNonnull() && pv->trackWeight(pfcand.trackRef()) > 0.5 && pv_ass_quality == 7) {
92  }
93  return vtx_ass;
94  }
95 
97  const auto &pseudo_track = (pfcand.bestTrack()) ? *pfcand.bestTrack() : reco::Track();
98  // conditions from PackedCandidate producer
99  bool highPurity = pfcand.trackRef().isNonnull() && pseudo_track.quality(reco::Track::highPurity);
100  // do same bit operations than in PackedCandidate
101  uint16_t qualityFlags = 0;
102  qualityFlags = (qualityFlags & ~trackHighPurityMask) | ((highPurity << trackHighPurityShift) & trackHighPurityMask);
103  bool isHighPurity = (qualityFlags & trackHighPurityMask) >> trackHighPurityShift;
104  // to do as in TrackBase
105  uint8_t quality = (1 << reco::TrackBase::loose);
106  if (isHighPurity) {
108  }
109  return quality;
110  }
111 
113  const auto &pseudo_track = (pfcand.bestTrack()) ? *pfcand.bestTrack() : reco::Track();
114  // conditions from PackedCandidate producer
115  bool highPurity = pfcand.trackRef().isNonnull() && pseudo_track.quality(reco::Track::highPurity);
116  // do same bit operations than in PackedCandidate
117  uint16_t qualityFlags = 0;
118  qualityFlags = (qualityFlags & ~trackHighPurityMask) | ((highPurity << trackHighPurityShift) & trackHighPurityMask);
119  return int16_t((qualityFlags & lostInnerHitsMask) >> lostInnerHitsShift) - 1;
120  }
121 
122  std::pair<float, float> getDRSubjetFeatures(const reco::Jet &jet, const reco::Candidate *cand) {
123  const auto *patJet = dynamic_cast<const pat::Jet *>(&jet);
124  std::pair<float, float> features;
125  // Do Subjets
126  if (patJet) {
127  if (patJet->nSubjetCollections() > 0) {
128  auto subjets = patJet->subjets();
129  std::nth_element(
130  subjets.begin(),
131  subjets.begin() + 1,
132  subjets.end(),
133  [](const edm::Ptr<pat::Jet> &p1, const edm::Ptr<pat::Jet> &p2) { return p1->pt() > p2->pt(); });
134  features.first = !subjets.empty() ? reco::deltaR(*cand, *subjets[0]) : -1;
135  features.second = subjets.size() > 1 ? reco::deltaR(*cand, *subjets[1]) : -1;
136  } else {
137  features.first = -1;
138  features.second = -1;
139  }
140  } else {
141  features.first = -1;
142  features.second = -1;
143  }
144  return features;
145  }
146 
147  int center_norm_pad(const std::vector<float> &input,
148  float center,
149  float norm_factor,
150  unsigned min_length,
151  unsigned max_length,
152  std::vector<float> &datavec,
153  int startval,
154  float pad_value,
155  float replace_inf_value,
156  float min,
157  float max) {
158  // do variable shifting/scaling/padding/clipping in one go
159 
160  assert(min <= pad_value && pad_value <= max);
161  assert(min_length <= max_length);
162 
163  unsigned target_length = std::clamp((unsigned)input.size(), min_length, max_length);
164  for (unsigned i = 0; i < target_length; ++i) {
165  if (i < input.size()) {
166  datavec[i + startval] = std::clamp((catch_infs(input[i], replace_inf_value) - center) * norm_factor, min, max);
167  } else {
168  datavec[i + startval] = pad_value;
169  }
170  }
171  return target_length;
172  }
173 
174  int center_norm_pad_halfRagged(const std::vector<float> &input,
175  float center,
176  float norm_factor,
177  unsigned target_length,
178  std::vector<float> &datavec,
179  int startval,
180  float pad_value,
181  float replace_inf_value,
182  float min,
183  float max) {
184  // do variable shifting/scaling/padding/clipping in one go
185 
186  assert(min <= pad_value && pad_value <= max);
187 
188  for (unsigned i = 0; i < std::min(static_cast<unsigned int>(input.size()), target_length); ++i) {
189  datavec.push_back(std::clamp((catch_infs(input[i], replace_inf_value) - center) * norm_factor, min, max));
190  }
191  if (input.size() < target_length)
192  datavec.insert(datavec.end(), target_length - input.size(), pad_value);
193 
194  return target_length;
195  }
196 
198  bool doExtra,
199  std::vector<std::string> &input_names_,
200  std::unordered_map<std::string, PreprocessParams> &prep_info_map_,
201  std::vector<std::vector<int64_t>> &input_shapes_,
202  std::vector<unsigned> &input_sizes_,
203  cms::Ort::FloatArrays *data_) {
204  // load preprocessing info
205  auto json_path = Config_.getParameter<std::string>("preprocess_json");
206  if (!json_path.empty()) {
207  // use preprocessing json file if available
208  std::ifstream ifs(edm::FileInPath(json_path).fullPath());
210  js.at("input_names").get_to(input_names_);
211  for (const auto &group_name : input_names_) {
212  const auto &group_pset = js.at(group_name);
213  auto &prep_params = prep_info_map_[group_name];
214  group_pset.at("var_names").get_to(prep_params.var_names);
215  if (group_pset.contains("var_length")) {
216  prep_params.min_length = group_pset.at("var_length");
217  prep_params.max_length = prep_params.min_length;
218  } else {
219  prep_params.min_length = group_pset.at("min_length");
220  prep_params.max_length = group_pset.at("max_length");
221  input_shapes_.push_back({1, (int64_t)prep_params.var_names.size(), -1});
222  }
223  const auto &var_info_pset = group_pset.at("var_infos");
224  for (const auto &var_name : prep_params.var_names) {
225  const auto &var_pset = var_info_pset.at(var_name);
226  double median = var_pset.at("median");
227  double norm_factor = var_pset.at("norm_factor");
228  double replace_inf_value = var_pset.at("replace_inf_value");
229  double lower_bound = var_pset.at("lower_bound");
230  double upper_bound = var_pset.at("upper_bound");
231  double pad = var_pset.contains("pad") ? double(var_pset.at("pad")) : 0;
232  prep_params.var_info_map[var_name] =
234  }
235 
236  if (doExtra && data_ != nullptr) {
237  // create data storage with a fixed size vector initialized w/ 0
238  const auto &len = input_sizes_.emplace_back(prep_params.max_length * prep_params.var_names.size());
239  data_->emplace_back(len, 0);
240  }
241  }
242  } else {
243  // otherwise use the PSet in the python config file
244  const auto &prep_pset = Config_.getParameterSet("preprocessParams");
245  input_names_ = prep_pset.getParameter<std::vector<std::string>>("input_names");
246  for (const auto &group_name : input_names_) {
247  const edm::ParameterSet &group_pset = prep_pset.getParameterSet(group_name);
248  auto &prep_params = prep_info_map_[group_name];
249  prep_params.var_names = group_pset.getParameter<std::vector<std::string>>("var_names");
250  prep_params.min_length = group_pset.getParameter<unsigned>("var_length");
251  prep_params.max_length = prep_params.min_length;
252  const auto &var_info_pset = group_pset.getParameterSet("var_infos");
253  for (const auto &var_name : prep_params.var_names) {
254  const edm::ParameterSet &var_pset = var_info_pset.getParameterSet(var_name);
255  double median = var_pset.getParameter<double>("median");
256  double norm_factor = var_pset.getParameter<double>("norm_factor");
257  double replace_inf_value = var_pset.getParameter<double>("replace_inf_value");
258  double lower_bound = var_pset.getParameter<double>("lower_bound");
259  double upper_bound = var_pset.getParameter<double>("upper_bound");
260  prep_params.var_info_map[var_name] =
262  }
263 
264  if (doExtra && data_ != nullptr) {
265  // create data storage with a fixed size vector initialized w/ 0
266  const auto &len = input_sizes_.emplace_back(prep_params.max_length * prep_params.var_names.size());
267  data_->emplace_back(len, 0);
268  }
269  }
270  }
271  }
272 
273 } // namespace btagbtvdeep
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)
Definition: deep_helpers.cc:96
vector< string > parse(string line, const string &delimiter)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float vertexDdotP(const reco::VertexCompositePtrCandidate &sv, const reco::Vertex &pv)
Definition: deep_helpers.cc:65
Measurement1D vertexD3d(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv)
Definition: deep_helpers.cc:56
math::XYZVector Vector
point in the space
Definition: Candidate.h:42
const float catch_infs(const float in, const float replace_value=0.)
Definition: deep_helpers.cc:21
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
nlohmann::json json
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
void fillVertexCovariance(CovarianceMatrix &v) const override
fill SMatrix
Base class for all types of Jets.
Definition: Jet.h:20
ParameterSet const & getParameterSet(std::string const &) const
const Point & vertex() const override
vertex position (overwritten by PF...)
std::pair< float, float > getDRSubjetFeatures(const reco::Jet &jet, const reco::Candidate *cand)
std::vector< std::vector< float > > FloatArrays
Definition: ONNXRuntime.h:23
assert(be >=bs)
math::Error< dimension >::type CovarianceMatrix
covariance error matrix (3x3)
Definition: Vertex.h:46
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:30
string quality
qualityFlagsShiftsAndMasks
Definition: deep_helpers.cc:9
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:72
std::vector< float > features(const reco::PreId &ecal, const reco::PreId &hcal, double rho, const reco::BeamSpot &spot, noZS::EcalClusterLazyTools &ecalTools)
Measurement1D vertexDxy(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv)
Definition: deep_helpers.cc:47
float vtx_ass_from_pfcand(const reco::PFCandidate &pfcand, int pv_ass_quality, const reco::VertexRef &pv)
Definition: deep_helpers.cc:88
d
Definition: ztail.py:151
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
Analysis-level calorimeter jet class.
Definition: Jet.h:77
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
static constexpr int qualityMap[8]
Definition: deep_helpers.cc:7