Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:25

0001 #ifndef RecoBTag_FeatureTools_deep_helpers_h
0002 #define RecoBTag_FeatureTools_deep_helpers_h
0003 
0004 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
0005 #include "DataFormats/BTauReco/interface/TaggingVariable.h"
0006 
0007 #include "TrackingTools/IPTools/interface/IPTools.h"
0008 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0009 
0010 #include "DataFormats/BTauReco/interface/CandIPTagInfo.h"
0011 
0012 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0013 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0014 //#include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0015 
0016 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0017 #include "DataFormats/PatCandidates/interface/Jet.h"
0018 
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 
0022 #include "PhysicsTools/ONNXRuntime/interface/ONNXRuntime.h"
0023 
0024 #include <iostream>
0025 #include <fstream>
0026 #include <algorithm>
0027 #include <numeric>
0028 #include <nlohmann/json.hpp>
0029 
0030 namespace btagbtvdeep {
0031 
0032   // remove infs and NaNs with value  (adapted from DeepNTuples)
0033   const float catch_infs(const float in, const float replace_value = 0.);
0034 
0035   // remove infs/NaN and bound (adapted from DeepNTuples)
0036   const float catch_infs_and_bound(const float in,
0037                                    const float replace_value,
0038                                    const float lowerbound,
0039                                    const float upperbound,
0040                                    const float offset = 0.,
0041                                    const bool use_offsets = true);
0042 
0043   // 2D distance between SV and PV (adapted from DeepNTuples)
0044   Measurement1D vertexDxy(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv);
0045 
0046   //3D distance between SV and PV (adapted from DeepNTuples)
0047   Measurement1D vertexD3d(const reco::VertexCompositePtrCandidate &svcand, const reco::Vertex &pv);
0048 
0049   // dot product between SV and PV (adapted from DeepNTuples)
0050   float vertexDdotP(const reco::VertexCompositePtrCandidate &sv, const reco::Vertex &pv);
0051 
0052   // helper to order vertices by significance (adapted from DeepNTuples)
0053   template <typename SVType, typename PVType>
0054   bool sv_vertex_comparator(const SVType &sva, const SVType &svb, const PVType &pv) {
0055     auto adxy = vertexDxy(sva, pv);
0056     auto bdxy = vertexDxy(svb, pv);
0057     float aval = adxy.value();
0058     float bval = bdxy.value();
0059     float aerr = adxy.error();
0060     float berr = bdxy.error();
0061 
0062     float asig = catch_infs(aval / aerr, 0.);
0063     float bsig = catch_infs(bval / berr, 0.);
0064     return bsig < asig;
0065   }
0066 
0067   // write tagging variables to vector (adapted from DeepNTuples)
0068   template <typename T>
0069   int dump_vector(reco::TaggingVariableList &from, T *to, reco::btau::TaggingVariableName name, const size_t max) {
0070     std::vector<T> vals = from.getList(name, false);
0071     size_t size = std::min(vals.size(), max);
0072     if (size > 0) {
0073       for (size_t i = 0; i < vals.size(); i++) {
0074         to[i] = catch_infs(vals[i], -0.1);
0075       }
0076     }
0077     return size;
0078   }
0079 
0080   // compute minimum dr between SVs and a candidate (from DeepNTuples, now polymorphic)
0081   float mindrsvpfcand(const std::vector<reco::VertexCompositePtrCandidate> &svs,
0082                       const reco::Candidate *cand,
0083                       float mindr = 0.4);
0084 
0085   // compute minimum distance between SVs and a candidate (from DeepNTuples, now polymorphic)
0086   float mindistsvpfcand(const std::vector<reco::VertexCompositePtrCandidate> &svs, const reco::TransientTrack track);
0087 
0088   // mimic the calculation in PackedCandidate
0089   float vtx_ass_from_pfcand(const reco::PFCandidate &pfcand, int pv_ass_quality, const reco::VertexRef &pv);
0090   float quality_from_pfcand(const reco::PFCandidate &pfcand);
0091   float lost_inner_hits_from_pfcand(const reco::PFCandidate &pfcand);
0092 
0093   std::pair<float, float> getDRSubjetFeatures(const reco::Jet &jet, const reco::Candidate *cand);
0094 
0095   // struct to hold preprocessing parameters
0096   struct PreprocessParams {
0097     struct VarInfo {
0098       VarInfo() {}
0099       VarInfo(float median, float norm_factor, float replace_inf_value, float lower_bound, float upper_bound, float pad)
0100           : center(median),
0101             norm_factor(norm_factor),
0102             replace_inf_value(replace_inf_value),
0103             lower_bound(lower_bound),
0104             upper_bound(upper_bound),
0105             pad(pad) {}
0106       float center = 0;
0107       float norm_factor = 1;
0108       float replace_inf_value = 0;
0109       float lower_bound = -5;
0110       float upper_bound = 5;
0111       float pad = 0;
0112     };
0113 
0114     unsigned min_length = 0;
0115     unsigned max_length = 0;
0116     std::vector<std::string> var_names;
0117     std::unordered_map<std::string, VarInfo> var_info_map;
0118 
0119     VarInfo info(const std::string &name) const { return var_info_map.at(name); }
0120   };
0121 
0122   int center_norm_pad(const std::vector<float> &input,
0123                       float center,
0124                       float scale,
0125                       unsigned min_length,
0126                       unsigned max_length,
0127                       std::vector<float> &datavec,
0128                       int startval,
0129                       float pad_value = 0,
0130                       float replace_inf_value = 0,
0131                       float min = 0,
0132                       float max = -1);
0133 
0134   int center_norm_pad_halfRagged(const std::vector<float> &input,
0135                                  float center,
0136                                  float scale,
0137                                  unsigned target_length,
0138                                  std::vector<float> &datavec,
0139                                  int startval,
0140                                  float pad_value = 0,
0141                                  float replace_inf_value = 0,
0142                                  float min = 0,
0143                                  float max = -1);
0144 
0145   void ParticleNetConstructor(const edm::ParameterSet &Config_,
0146                               bool doExtra,
0147                               std::vector<std::string> &input_names_,
0148                               std::unordered_map<std::string, PreprocessParams> &prep_info_map_,
0149                               std::vector<std::vector<int64_t>> &input_shapes_,
0150                               std::vector<unsigned> &input_sizes_,
0151                               cms::Ort::FloatArrays *data_);
0152 
0153 }  // namespace btagbtvdeep
0154 #endif  //RecoBTag_FeatureTools_deep_helpers_h