CMS 3D CMS Logo

MLPFModel.h
Go to the documentation of this file.
1 #ifndef RecoParticleFlow_PFProducer_interface_MLPFModel
2 #define RecoParticleFlow_PFProducer_interface_MLPFModel
3 
7 
8 namespace reco::mlpf {
9  //The model takes the following number of features for each input PFElement
10  static constexpr unsigned int NUM_ELEMENT_FEATURES = 25;
11  static constexpr unsigned int NUM_OUTPUT_FEATURES = 14;
12 
13  //these are defined at model creation time and set the random LSH codebook size
14  static constexpr int LSH_BIN_SIZE = 64;
16 
17  //In CPU mode, we want to evaluate each event separately
18  static constexpr int BATCH_SIZE = 1;
19 
20  //The model has 14 outputs for each particle:
21  // out[0-7]: particle classification logits for each pdgId
22  // out[8]: regressed charge
23  // out[9]: regressed pt
24  // out[10]: regressed eta
25  // out[11]: regressed sin phi
26  // out[12]: regressed cos phi
27  // out[13]: regressed energy
28  static constexpr unsigned int IDX_CLASS = 7;
29 
30  static constexpr unsigned int IDX_CHARGE = 8;
31 
32  static constexpr unsigned int IDX_PT = 9;
33  static constexpr unsigned int IDX_ETA = 10;
34  static constexpr unsigned int IDX_SIN_PHI = 11;
35  static constexpr unsigned int IDX_COS_PHI = 12;
36  static constexpr unsigned int IDX_ENERGY = 13;
37 
38  //for consistency with the baseline PFAlgo
39  static constexpr float PI_MASS = 0.13957;
40 
41  //index [0, N_pdgids) -> PDGID
42  //this maps the absolute values of the predicted PDGIDs to an array of ascending indices
43  static const std::vector<int> pdgid_encoding = {0, 211, 130, 1, 2, 22, 11, 13};
44 
45  //PFElement::type -> index [0, N_types)
46  //this maps the type of the PFElement to an ascending index that is used by the model to distinguish between different elements
47  static const std::map<int, int> elem_type_encoding = {
48  {0, 0},
49  {1, 1},
50  {2, 2},
51  {3, 3},
52  {4, 4},
53  {5, 5},
54  {6, 6},
55  {7, 7},
56  {8, 8},
57  {9, 9},
58  {10, 10},
59  {11, 11},
60  };
61 
62  std::array<float, NUM_ELEMENT_FEATURES> getElementProperties(const reco::PFBlockElement& orig);
63  float normalize(float in);
64 
65  int argMax(std::vector<float> const& vec);
66 
67  reco::PFCandidate makeCandidate(int pred_pid,
68  int pred_charge,
69  float pred_pt,
70  float pred_eta,
71  float pred_sin_phi,
72  float pred_cos_phi,
73  float pred_e);
74 
75  const std::vector<const reco::PFBlockElement*> getPFElements(const reco::PFBlockCollection& blocks);
76 
78  const std::vector<const reco::PFBlockElement*> elems,
79  size_t ielem_originator);
80 }; // namespace reco::mlpf
81 
82 #endif
Abstract base class for a PFBlock element (track, cluster...)
static constexpr unsigned int IDX_ETA
Definition: MLPFModel.h:33
std::array< float, NUM_ELEMENT_FEATURES > getElementProperties(const reco::PFBlockElement &orig)
Definition: MLPFModel.cc:17
int argMax(std::vector< float > const &vec)
Definition: MLPFModel.cc:216
static constexpr unsigned int NUM_ELEMENT_FEATURES
Definition: MLPFModel.h:10
static constexpr unsigned int IDX_PT
Definition: MLPFModel.h:32
static constexpr unsigned int IDX_ENERGY
Definition: MLPFModel.h:36
static constexpr int BATCH_SIZE
Definition: MLPFModel.h:18
const std::vector< const reco::PFBlockElement * > getPFElements(const reco::PFBlockCollection &blocks)
Definition: MLPFModel.cc:263
static constexpr unsigned int NUM_OUTPUT_FEATURES
Definition: MLPFModel.h:11
static constexpr unsigned int IDX_CHARGE
Definition: MLPFModel.h:30
std::vector< PFBlock > PFBlockCollection
collection of PFBlock objects
Definition: PFBlockFwd.h:10
reco::PFCandidate makeCandidate(int pred_pid, int pred_charge, float pred_pt, float pred_eta, float pred_sin_phi, float pred_cos_phi, float pred_e)
Definition: MLPFModel.cc:220
static const std::map< int, int > elem_type_encoding
Definition: MLPFModel.h:47
static constexpr unsigned int IDX_CLASS
Definition: MLPFModel.h:28
static constexpr unsigned int IDX_SIN_PHI
Definition: MLPFModel.h:34
float normalize(float in)
Definition: MLPFModel.cc:207
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
static constexpr unsigned int IDX_COS_PHI
Definition: MLPFModel.h:35
static constexpr float PI_MASS
Definition: MLPFModel.h:39
static const std::vector< int > pdgid_encoding
Definition: MLPFModel.h:43
void setCandidateRefs(reco::PFCandidate &cand, const std::vector< const reco::PFBlockElement *> elems, size_t ielem_originator)
Definition: MLPFModel.cc:284
static constexpr int LSH_BIN_SIZE
Definition: MLPFModel.h:14
static constexpr int NUM_MAX_ELEMENTS_BATCH
Definition: MLPFModel.h:15