CMS 3D CMS Logo

JetId.cc
Go to the documentation of this file.
3 #include <cmath>
4 
5 JetId::JetId(const std::string &iInput,
6  const std::string &iOutput,
7  const std::shared_ptr<hls4mlEmulator::Model> model,
8  int iNParticles)
9  : modelRef_(model) {
10  NNvectorVar_.clear();
11  fNParticles_ = iNParticles;
12 
13  fPt_ = std::make_unique<float[]>(fNParticles_);
14  fEta_ = std::make_unique<float[]>(fNParticles_);
15  fPhi_ = std::make_unique<float[]>(fNParticles_);
16  fId_ = std::make_unique<int[]>(fNParticles_);
17  fCharge_ = std::make_unique<int[]>(fNParticles_);
18  fDZ_ = std::make_unique<float[]>(fNParticles_);
19  fDX_ = std::make_unique<float[]>(fNParticles_);
20  fDY_ = std::make_unique<float[]>(fNParticles_);
21  fInput_ = iInput;
22  fOutput_ = iOutput;
23 }
24 
25 //--BJet algo specific constructor
26 JetId::JetId(const std::string &iInput, const std::string &iOutput, const BJetTFCache *cache, int iNParticles)
27  : sessionRef_(cache->session) {
28  NNvectorVar_.clear();
29  fNParticles_ = iNParticles;
30 
31  fPt_ = std::make_unique<float[]>(fNParticles_);
32  fEta_ = std::make_unique<float[]>(fNParticles_);
33  fPhi_ = std::make_unique<float[]>(fNParticles_);
34  fId_ = std::make_unique<int[]>(fNParticles_);
35  fCharge_ = std::make_unique<int[]>(fNParticles_);
36  fDZ_ = std::make_unique<float[]>(fNParticles_);
37  fDX_ = std::make_unique<float[]>(fNParticles_);
38  fDY_ = std::make_unique<float[]>(fNParticles_);
39  fInput_ = iInput;
40  fOutput_ = iOutput;
41 }
42 
44  NNvectorVar_.clear();
45  for (int i0 = 0; i0 < fNParticles_; i0++) {
46  if (fPt_.get()[i0] == 0) {
47  for (int i1 = 0; i1 < 13; i1++)
48  NNvectorVar_.push_back(0);
49  continue;
50  }
51  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Electron && fCharge_.get()[i0] < 0); // Electron
52  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Electron && fCharge_.get()[i0] > 0); // Positron
53  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Muon && fCharge_.get()[i0] < 0); // Muon
54  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Muon && fCharge_.get()[i0] > 0); // Anti-Muon
55  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Photon); // Photon
56  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::NeutralHadron); // Neutral Had
57  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::ChargedHadron && fCharge_.get()[i0] < 0); // Pion
58  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::ChargedHadron && fCharge_.get()[i0] > 0); // Anti-Pion
59  NNvectorVar_.push_back(fDZ_.get()[i0]); //dZ
60  NNvectorVar_.push_back(std::hypot(fDX_.get()[i0], fDY_.get()[i0])); //d0
61  NNvectorVar_.push_back(fPt_.get()[i0]); //pT as a fraction of jet pT
62  NNvectorVar_.push_back(fEta_.get()[i0]); //dEta from jet axis
63  NNvectorVar_.push_back(fPhi_.get()[i0]); //dPhi from jet axis
64  }
65 }
67  tensorflow::Tensor input(tensorflow::DT_FLOAT, {1, (unsigned int)NNvectorVar_.size(), 1});
68  for (unsigned int i = 0; i < NNvectorVar_.size(); i++) {
69  input.tensor<float, 3>()(0, i, 0) = float(NNvectorVar_[i]);
70  }
71  std::vector<tensorflow::Tensor> outputs;
73  return outputs[0].matrix<float>()(0, 0);
74 } //end EvaluateNN
75 
76 ap_fixed<16, 6> JetId::EvaluateNNFixed() {
77  ap_fixed<16, 6> modelInput[140] = {};
78  for (unsigned int i = 0; i < NNvectorVar_.size(); i++) {
79  modelInput[i] = NNvectorVar_[i];
80  }
81  ap_fixed<16, 6> modelResult[1] = {-1};
82 
83  modelRef_->prepare_input(modelInput);
84  modelRef_->predict();
85  modelRef_->read_result(modelResult);
86  ap_fixed<16, 6> modelResult_ = modelResult[0];
87  return modelResult_;
88 } //end EvaluateNNFixed
89 
90 float JetId::compute(const l1t::PFJet &iJet, float vz, bool useRawPt) {
91  for (int i0 = 0; i0 < fNParticles_; i0++) {
92  fPt_.get()[i0] = 0;
93  fEta_.get()[i0] = 0;
94  fPhi_.get()[i0] = 0;
95  fId_.get()[i0] = 0;
96  fCharge_.get()[i0] = 0;
97  fDZ_.get()[i0] = 0;
98  fDX_.get()[i0] = 0;
99  fDY_.get()[i0] = 0;
100  }
101  auto iParts = iJet.constituents();
102  std::sort(iParts.begin(), iParts.end(), [](edm::Ptr<l1t::PFCandidate> i, edm::Ptr<l1t::PFCandidate> j) {
103  return (i->pt() > j->pt());
104  });
105  float jetpt = useRawPt ? iJet.rawPt() : iJet.pt();
106  for (unsigned int i0 = 0; i0 < iParts.size(); i0++) {
107  if (i0 >= (unsigned int)fNParticles_)
108  break;
109  fPt_.get()[i0] = iParts[i0]->pt() / jetpt;
110  fEta_.get()[i0] = iParts[i0]->eta() - iJet.eta();
111  fPhi_.get()[i0] = deltaPhi(iParts[i0]->phi(), iJet.phi());
112  fId_.get()[i0] = iParts[i0]->id();
113  fCharge_.get()[i0] = iParts[i0]->charge();
114  if (iParts[i0]->pfTrack().isNonnull()) {
115  fDX_.get()[i0] = iParts[i0]->pfTrack()->vx();
116  fDY_.get()[i0] = iParts[i0]->pfTrack()->vy();
117  fDZ_.get()[i0] = iParts[i0]->pfTrack()->vz() - vz;
118  }
119  }
120  setNNVectorVar();
121  return EvaluateNN();
122 }
123 
124 ap_fixed<16, 6> JetId::computeFixed(const l1t::PFJet &iJet, float vz, bool useRawPt) {
125  for (int i0 = 0; i0 < fNParticles_; i0++) {
126  fPt_.get()[i0] = 0;
127  fEta_.get()[i0] = 0;
128  fPhi_.get()[i0] = 0;
129  fId_.get()[i0] = 0;
130  fCharge_.get()[i0] = 0;
131  fDZ_.get()[i0] = 0;
132  fDX_.get()[i0] = 0;
133  fDY_.get()[i0] = 0;
134  }
135  auto iParts = iJet.constituents();
136  std::sort(iParts.begin(), iParts.end(), [](edm::Ptr<l1t::PFCandidate> i, edm::Ptr<l1t::PFCandidate> j) {
137  return (i->pt() > j->pt());
138  });
139  float jetpt = useRawPt ? iJet.rawPt() : iJet.pt();
140  for (unsigned int i0 = 0; i0 < iParts.size(); i0++) {
141  if (i0 >= (unsigned int)fNParticles_)
142  break;
143  fPt_.get()[i0] = iParts[i0]->pt() / jetpt;
144  fEta_.get()[i0] = iParts[i0]->eta() - iJet.eta();
145  fPhi_.get()[i0] = deltaPhi(iParts[i0]->phi(), iJet.phi());
146  fId_.get()[i0] = iParts[i0]->id();
147  fCharge_.get()[i0] = iParts[i0]->charge();
148  if (iParts[i0]->pfTrack().isNonnull()) {
149  fDX_.get()[i0] = iParts[i0]->pfTrack()->vx();
150  fDY_.get()[i0] = iParts[i0]->pfTrack()->vy();
151  fDZ_.get()[i0] = iParts[i0]->pfTrack()->vz() - vz;
152  }
153  }
154  setNNVectorVar();
155  return EvaluateNNFixed();
156 }
ap_fixed< 16, 6 > EvaluateNNFixed()
Definition: JetId.cc:76
unique_ptr< float[]> fDZ_
Definition: JetId.h:48
MPlex< T, D1, D2, N > hypot(const MPlex< T, D1, D2, N > &a, const MPlex< T, D1, D2, N > &b)
Definition: Matriplex.h:436
double pt() const final
transverse momentum
unique_ptr< float[]> fDX_
Definition: JetId.h:49
ap_fixed< 16, 6 > computeFixed(const l1t::PFJet &iJet, float vz, bool useRawPt)
Definition: JetId.cc:124
JetId(const std::string &iInput, const std::string &iOutput, const std::shared_ptr< hls4mlEmulator::Model > model, int iNParticles)
Definition: JetId.cc:5
std::vector< float > NNvectorVar_
Definition: JetId.h:39
float compute(const l1t::PFJet &iJet, float vz, bool useRawPt)
Definition: JetId.cc:90
std::string fInput_
Definition: JetId.h:40
static std::string const input
Definition: EdmProvDump.cc:50
void setNNVectorVar()
Definition: JetId.cc:43
unique_ptr< int[]> fId_
Definition: JetId.h:46
float rawPt() const
Definition: PFJet.h:27
unique_ptr< float[]> fPhi_
Definition: JetId.h:45
unique_ptr< float[]> fEta_
Definition: JetId.h:44
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:271
const Constituents & constituents() const
constituent information. note that this is not going to be available in the hardware! ...
Definition: PFJet.h:30
std::string fOutput_
Definition: JetId.h:41
unique_ptr< int[]> fCharge_
Definition: JetId.h:47
float EvaluateNN()
Definition: JetId.cc:66
unique_ptr< float[]> fPt_
Definition: JetId.h:43
unique_ptr< float[]> fDY_
Definition: JetId.h:50
def cache(function)
Definition: utilities.py:3
std::shared_ptr< hls4mlEmulator::Model > modelRef_
Definition: JetId.h:52
double phi() const final
momentum azimuthal angle
tensorflow::Session * sessionRef_
Definition: JetId.h:51
int fNParticles_
Definition: JetId.h:42
double eta() const final
momentum pseudorapidity