CMS 3D CMS Logo

DisplacedRegionSeedingVertexProducer.cc
Go to the documentation of this file.
1 #include <list>
2 #include <vector>
3 #include <limits>
4 #include <string>
5 #include <atomic>
6 
17 
24 
26 
28 
29 using namespace std;
32 
34 public:
37  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
39 
40 private:
41  // clustering parameters
42  const double rParam_;
43 
44  // selection parameters
45  const double minRadius_;
46  const double nearThreshold_;
47  const double farThreshold_;
48  const double discriminatorCut_;
49  const vector<string> input_names_;
50  const vector<string> output_names_;
51 
54 
55  tensorflow::Session *session_;
56 
57  double getDiscriminatorValue(const DisplacedVertexCluster &, const reco::BeamSpot &) const;
58 };
59 
61  : rParam_(cfg.getParameter<double>("rParam")),
62  minRadius_(cfg.getParameter<double>("minRadius")),
63  nearThreshold_(cfg.getParameter<double>("nearThreshold")),
64  farThreshold_(cfg.getParameter<double>("farThreshold")),
65  discriminatorCut_(cfg.getParameter<double>("discriminatorCut")),
66  input_names_(cfg.getParameter<vector<string> >("input_names")),
67  output_names_(cfg.getParameter<vector<string> >("output_names")),
68  beamSpotToken_(consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamSpot"))),
69  trackClustersToken_(
70  consumes<edm::View<reco::VertexCompositeCandidate> >(cfg.getParameter<edm::InputTag>("trackClusters"))) {
71  unsigned nThreads = cfg.getUntrackedParameter<unsigned>("nThreads");
73  options.setThreading(nThreads);
74  string pbFile = cfg.getParameter<edm::FileInPath>("graph_path").fullPath();
75  auto graphDef = tensorflow::loadGraphDef(pbFile);
77 
78  produces<vector<reco::Vertex> >("nearRegionsOfInterest");
79  produces<vector<reco::Vertex> >("farRegionsOfInterest");
80 }
81 
83  if (session_ != nullptr)
85 }
86 
89  const edm::EventSetup &setup) const {
90  const auto &beamSpot = event.get(beamSpotToken_);
91  const math::XYZVector bs(beamSpot.position());
92 
93  const auto &trackClusters = event.get(trackClustersToken_);
94 
95  // Initialize distances.
96  list<DisplacedVertexCluster> pseudoROIs;
97  list<Distance> distances;
98  const double minTrackClusterRadius = minRadius_ - rParam_;
99  for (unsigned i = 0; i < trackClusters.size(); i++) {
100  const reco::VertexCompositeCandidate &trackCluster = trackClusters[i];
101  const math::XYZVector x(trackCluster.vertex());
102  if (minRadius_ < 0.0 || minTrackClusterRadius < 0.0 || (x - bs).rho() > minTrackClusterRadius)
103  pseudoROIs.emplace_back(&trackClusters.at(i), rParam_);
104  }
105  if (pseudoROIs.size() > 1) {
106  DisplacedVertexClusterItr secondToLast = pseudoROIs.end();
107  secondToLast--;
108  for (DisplacedVertexClusterItr i = pseudoROIs.begin(); i != secondToLast; i++) {
110  j++;
111  for (; j != pseudoROIs.end(); j++) {
112  distances.emplace_back(i, j);
113 
114  // Track clusters farther apart than 4 times rParam_ (i.e., 16 times
115  // rParam_^2) cannot wind up in the same ROI, so remove these pairs.
116  if (distances.back().distance2() > 16.0 * rParam_ * rParam_)
117  distances.pop_back();
118  }
119  }
120  }
121 
122  // Do clustering.
123  while (!distances.empty()) {
124  const auto comp = [](const Distance &a, const Distance &b) { return a.distance2() <= b.distance2(); };
125  distances.sort(comp);
126  DistanceItr dBest = distances.begin();
127  if (dBest->distance2() > rParam_ * rParam_)
128  break;
129 
130  dBest->entities().first->merge(*dBest->entities().second);
131  dBest->entities().second->setInvalid();
132 
133  const auto distancePred = [](const Distance &a) {
134  return (!a.entities().first->valid() || !a.entities().second->valid());
135  };
136  const auto pseudoROIPred = [](const DisplacedVertexCluster &a) { return !a.valid(); };
137  distances.remove_if(distancePred);
138  pseudoROIs.remove_if(pseudoROIPred);
139  }
140 
141  // Remove invalid ROIs.
142  const auto roiPred = [&](const DisplacedVertexCluster &roi) {
143  if (!roi.valid())
144  return true;
145  const auto &x(roi.centerOfMass());
146  if ((x - bs).rho() < minRadius_)
147  return true;
148  const double discriminatorValue = ((discriminatorCut_ > 0.0) ? getDiscriminatorValue(roi, beamSpot) : 1.0);
149  if (discriminatorValue < discriminatorCut_)
150  return true;
151  return false;
152  };
153  pseudoROIs.remove_if(roiPred);
154 
155  auto nearRegionsOfInterest = make_unique<vector<reco::Vertex> >();
156  auto farRegionsOfInterest = make_unique<vector<reco::Vertex> >();
157 
158  constexpr std::array<double, 6> errorA{{1.0, 0.0, 1.0, 0.0, 0.0, 1.0}};
159  static const reco::Vertex::Error errorRegion(errorA.begin(), errorA.end(), true, true);
160 
161  for (const auto &roi : pseudoROIs) {
162  const auto &x(roi.centerOfMass());
163  if ((x - bs).rho() < nearThreshold_)
164  nearRegionsOfInterest->emplace_back(reco::Vertex::Point(roi.centerOfMass()), errorRegion);
165  if ((x - bs).rho() > farThreshold_)
166  farRegionsOfInterest->emplace_back(reco::Vertex::Point(roi.centerOfMass()), errorRegion);
167  }
168 
169  event.put(std::move(nearRegionsOfInterest), "nearRegionsOfInterest");
170  event.put(std::move(farRegionsOfInterest), "farRegionsOfInterest");
171 }
172 
175 
176  desc.add<double>("rParam", 1.0);
177  desc.add<double>("minRadius", -1.0);
178  desc.add<double>("nearThreshold", 9999.0);
179  desc.add<double>("farThreshold", -1.0);
180  desc.add<double>("discriminatorCut", -1.0);
181  desc.add<vector<string> >("input_names", {"phi_0", "phi_1"});
182  desc.add<vector<string> >("output_names", {"model_5/activation_10/Softmax"});
183  desc.addUntracked<unsigned>("nThreads", 1);
184  desc.add<edm::FileInPath>(
185  "graph_path",
187  "RecoTracker/DisplacedRegionalTracking/data/FullData_Phi-64-128-256_16-32-64_F-128-64-32_Model.pb"));
188  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
189  desc.add<edm::InputTag>("trackClusters", edm::InputTag("generalV0Candidates", "Kshort"));
190 
191  descriptions.add("displacedRegionProducer", desc);
192 }
193 
195  const reco::BeamSpot &bs) const {
196  // The network takes in two maps of data features, one with information
197  // related to the pairwise track vertices and one with information related to
198  // the tracks in an isolation annulus.
199 
200  constexpr int maxNVertices = 40; // maximum number of pairwise track vertices per ROI
201  constexpr int nVertexFeatures = 23; // number of features per pairwise track vertex
202  constexpr int nDimVertices = 3; // number of tensor dimensions for the map of pairwise track vertices
203 
204  constexpr int maxNAnnulusTracks = 10; // maximum number of annulus tracks per ROI
205  constexpr int nAnnulusTrackFeatures = 8; // number of features per annulus track
206  constexpr int nDimAnnulusTracks = 3; // number of tensor dimensions for the map of annulus tracks
207 
208  tensorflow::Tensor vertexTensor(tensorflow::DT_FLOAT, tensorflow::TensorShape({1, maxNVertices, nVertexFeatures}));
209  auto vertex_map = vertexTensor.tensor<float, nDimVertices>();
210  tensorflow::Tensor annulusTensor(tensorflow::DT_FLOAT,
211  tensorflow::TensorShape({1, maxNAnnulusTracks, nAnnulusTrackFeatures}));
212  auto annulus_map = annulusTensor.tensor<float, nDimAnnulusTracks>();
213 
214  for (int i = 0, map_i = 0; map_i < maxNVertices; i++, map_i++) {
215  if (i >= static_cast<int>(roi.nConstituents()))
216  for (unsigned j = 0; j < nVertexFeatures; j++)
217  vertex_map(0, map_i, j) = 0.0;
218  else {
219  const auto &trackCluster = *roi.constituent(i);
220  const auto &track0 = *trackCluster.daughter(0)->bestTrack();
221  const auto &track1 = *trackCluster.daughter(1)->bestTrack();
222 
223  vertex_map(0, map_i, 0) = trackCluster.vx() - bs.x0();
224  vertex_map(0, map_i, 1) = trackCluster.vy() - bs.y0();
225  vertex_map(0, map_i, 2) = trackCluster.vz() - bs.z0();
226 
227  vertex_map(0, map_i, 3) = trackCluster.vertexCovariance()(0, 0);
228  vertex_map(0, map_i, 4) = trackCluster.vertexCovariance()(0, 1);
229  vertex_map(0, map_i, 5) = trackCluster.vertexCovariance()(0, 2);
230  vertex_map(0, map_i, 6) = trackCluster.vertexCovariance()(1, 1);
231  vertex_map(0, map_i, 7) = trackCluster.vertexCovariance()(1, 2);
232  vertex_map(0, map_i, 8) = trackCluster.vertexCovariance()(2, 2);
233 
234  vertex_map(0, map_i, 9) = track0.charge() * track0.pt();
235  vertex_map(0, map_i, 10) = track0.eta();
236  vertex_map(0, map_i, 11) = track0.phi();
237  vertex_map(0, map_i, 12) = track0.dxy(bs);
238  vertex_map(0, map_i, 13) = track0.dz(bs.position());
239  vertex_map(0, map_i, 14) = track0.normalizedChi2();
240  vertex_map(0, map_i, 15) = track0.quality(reco::Track::highPurity) ? 1 : 0;
241 
242  vertex_map(0, map_i, 16) = track1.charge() * track1.pt();
243  vertex_map(0, map_i, 17) = track1.eta();
244  vertex_map(0, map_i, 18) = track1.phi();
245  vertex_map(0, map_i, 19) = track1.dxy(bs);
246  vertex_map(0, map_i, 20) = track1.dz(bs.position());
247  vertex_map(0, map_i, 21) = track1.normalizedChi2();
248  vertex_map(0, map_i, 22) = track1.quality(reco::Track::highPurity) ? 1 : 0;
249  }
250  }
251 
252  for (int i = 0; i < maxNAnnulusTracks; i++)
253  for (unsigned j = 0; j < nAnnulusTrackFeatures; j++)
254  annulus_map(0, i, j) = 0.0;
255 
256  tensorflow::NamedTensorList input_tensors;
257  input_tensors.resize(2);
258  input_tensors[0] = tensorflow::NamedTensor(input_names_.at(0), vertexTensor);
259  input_tensors[1] = tensorflow::NamedTensor(input_names_.at(1), annulusTensor);
260  vector<tensorflow::Tensor> outputs;
261  tensorflow::run(session_, input_tensors, output_names_, &outputs);
262 
263  return (outputs.at(0).flat<float>()(1));
264 }
265 
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
std::vector< NamedTensor > NamedTensorList
Definition: TensorFlow.h:31
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:655
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:129
std::list< Distance >::iterator DistanceItr
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:45
const Point & vertex() const override
vertex position (overwritten by PF...)
std::pair< std::string, Tensor > NamedTensor
Definition: TensorFlow.h:30
const reco::VertexCompositeCandidate * constituent(const unsigned i) const
DisplacedVertexCluster::DistanceItr DistanceItr
DisplacedVertexCluster::Distance Distance
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:268
bool closeSession(Session *&session)
Definition: TensorFlow.cc:243
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Session * createSession()
Definition: TensorFlow.cc:146
static void fillDescriptions(edm::ConfigurationDescriptions &)
const edm::EDGetTokenT< edm::View< reco::VertexCompositeCandidate > > trackClustersToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double getDiscriminatorValue(const DisplacedVertexCluster &, const reco::BeamSpot &) const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
std::list< DisplacedVertexCluster >::iterator DisplacedVertexClusterItr
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
virtual const Track * bestTrack() const
Definition: Candidate.h:268
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1