CMS 3D CMS Logo

TauNNId.cc
Go to the documentation of this file.
4 #include <cmath>
5 
6 static constexpr unsigned int n_particles_max = 10;
7 
9  const tensorflow::Session *session,
10  const std::string &iWeightFile,
11  int iNParticles)
12  : session_(session) {
13  NNvectorVar_.clear();
14  edm::FileInPath fp(iWeightFile);
15  fNParticles_ = iNParticles;
16 
17  fPt_ = std::make_unique<float[]>(fNParticles_);
18  fEta_ = std::make_unique<float[]>(fNParticles_);
19  fPhi_ = std::make_unique<float[]>(fNParticles_);
20  fId_ = std::make_unique<float[]>(fNParticles_);
21  fInput_ = iInput;
22 }
23 
25  NNvectorVar_.clear();
26  for (int i0 = 0; i0 < fNParticles_; i0++) {
27  NNvectorVar_.push_back(fPt_.get()[i0]); //pT
28  NNvectorVar_.push_back(fEta_.get()[i0]); //dEta from jet axis
29  NNvectorVar_.push_back(fPhi_.get()[i0]); //dPhi from jet axis
30  if (fPt_.get()[i0] == 0) {
31  for (int i1 = 0; i1 < 5; i1++)
32  NNvectorVar_.push_back(0);
33  continue;
34  }
35  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Photon); // Photon
36  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Electron); // Electron
37  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::Muon); // Muon
38  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::NeutralHadron); // Neutral Had
39  NNvectorVar_.push_back(fId_.get()[i0] == l1t::PFCandidate::ChargedHadron); // Charged Had
40  }
41 }
43  tensorflow::Tensor input(tensorflow::DT_FLOAT,
44  {1, (unsigned int)NNvectorVar_.size()}); //was {1,35} but get size mismatch, CHECK
45  for (unsigned int i = 0; i < NNvectorVar_.size(); i++) {
46  input.matrix<float>()(0, i) = float(NNvectorVar_[i]);
47  }
48  std::vector<tensorflow::Tensor> outputs;
49  tensorflow::run(session_, {{fInput_, input}}, {"dense_4/Sigmoid:0"}, &outputs);
50  return outputs[0].matrix<float>()(0, 0);
51 } //end EvaluateNN
52 
54  for (int i0 = 0; i0 < fNParticles_; i0++) {
55  fPt_.get()[i0] = 0;
56  fEta_.get()[i0] = 0;
57  fPhi_.get()[i0] = 0;
58  fId_.get()[i0] = 0;
59  }
60  std::sort(iParts.begin(), iParts.end(), [](l1t::PFCandidate i, l1t::PFCandidate j) { return (i.pt() > j.pt()); });
61  for (unsigned int i0 = 0; i0 < iParts.size(); i0++) {
62  if (i0 > n_particles_max || i0 >= (unsigned int)fNParticles_)
63  break;
64  fPt_.get()[i0] = iParts[i0].pt();
65  fEta_.get()[i0] = iSeed.eta() - iParts[i0].eta();
66  fPhi_.get()[i0] = deltaPhi(iSeed.phi(), iParts[i0].phi());
67  fId_.get()[i0] = iParts[i0].id();
68  }
70  return EvaluateNN();
71 }
const tensorflow::Session * session_
Definition: TauNNId.h:22
unique_ptr< float[]> fPhi_
Definition: TauNNId.h:28
static constexpr unsigned int n_particles_max
Definition: TauNNId.cc:6
std::vector< l1t::PFCandidate > PFCandidateCollection
Definition: PFCandidate.h:86
TauNNId(const std::string &iInput, const tensorflow::Session *session, const std::string &iWeightFile, int iNParticles)
Definition: TauNNId.cc:8
int fNParticles_
Definition: TauNNId.h:25
unique_ptr< float[]> fId_
Definition: TauNNId.h:29
static std::string const input
Definition: EdmProvDump.cc:50
std::vector< float > NNvectorVar_
Definition: TauNNId.h:23
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:281
void setNNVectorVar()
Definition: TauNNId.cc:24
unique_ptr< float[]> fPt_
Definition: TauNNId.h:26
std::string fInput_
Definition: TauNNId.h:24
float EvaluateNN()
Definition: TauNNId.cc:42
double phi() const final
momentum azimuthal angle
float compute(const l1t::PFCandidate &iSeed, l1t::PFCandidateCollection &iParts)
Definition: TauNNId.cc:53
unique_ptr< float[]> fEta_
Definition: TauNNId.h:27
double eta() const final
momentum pseudorapidity