CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Private Member Functions | Private Attributes
DeepCoreSeedGenerator Class Reference

#include <trackJet/DeepCoreSeedGenerator/plugins/DeepCoreSeedGenerator.cc>

Inheritance diagram for DeepCoreSeedGenerator:
edm::stream::EDProducer< edm::GlobalCache< DeepCoreCache > >

Public Member Functions

 DeepCoreSeedGenerator (const edm::ParameterSet &, const DeepCoreCache *)
 
 ~DeepCoreSeedGenerator () override
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DeepCoreCache > >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
static void globalEndJob (DeepCoreCache *)
 
static std::unique_ptr< DeepCoreCacheinitializeGlobalCache (const edm::ParameterSet &)
 

Public Attributes

double jetEta_
 
double jetPt_
 
double pitchX_ = 0.01
 
double pitchY_ = 0.015
 

Static Public Attributes

static constexpr int jetDimX = 30
 
static constexpr int jetDimY = 30
 
static constexpr int Nlayer = 4
 
static constexpr int Nover = 3
 
static constexpr int Npar = 5
 

Private Member Functions

void beginJob ()
 
const GeomDetDetectorSelector (int, const reco::Candidate &, GlobalVector, const reco::Vertex &, const TrackerTopology *const, const edmNew::DetSetVector< SiPixelCluster > &)
 
void endJob ()
 
void fillPixelMatrix (const SiPixelCluster &, int, Point3DBase< float, LocalTag >, const GeomDet *, tensorflow::NamedTensorList)
 
std::pair< bool, Basic3DVector< float > > findIntersection (const GlobalVector &, const reco::Candidate::Point &, const GeomDet *)
 
std::pair< int, int > local2Pixel (double, double, const GeomDet *)
 
LocalPoint pixel2Local (int, int, const GeomDet *)
 
int pixelFlipper (const GeomDet *)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
std::pair< double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> SeedEvaluation (tensorflow::NamedTensorList, std::vector< std::string >)
 
std::vector< GlobalVectorsplittedClusterDirections (const reco::Candidate &, const TrackerTopology *const, const PixelClusterParameterEstimator *, const reco::Vertex &, int, const edmNew::DetSetVector< SiPixelCluster > &)
 

Private Attributes

double centralMIPCharge_
 
double chargeFracMin_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
 
double deltaR_
 
edm::ESHandle< GlobalTrackingGeometrygeometry_
 
edm::Handle< edmNew::DetSetVector< SiPixelCluster > > inputPixelClusters_
 
std::vector< std::string > inputTensorName_
 
edm::ESHandle< MagneticFieldmagfield_
 
std::vector< std::string > outputTensorName_
 
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
 
std::string pixelCPE_
 
double probThr_
 
edm::ESHandle< Propagatorpropagator_
 
std::string propagatorName_
 
double ptMin_
 
tensorflow::Session * session_
 
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcdtopoToken_
 
edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_
 
std::string weightfilename_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DeepCoreCache > >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 79 of file DeepCoreSeedGenerator.cc.

Constructor & Destructor Documentation

◆ DeepCoreSeedGenerator()

DeepCoreSeedGenerator::DeepCoreSeedGenerator ( const edm::ParameterSet iConfig,
const DeepCoreCache cache 
)
explicit

Definition at line 164 of file DeepCoreSeedGenerator.cc.

165  : vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
167  consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelClusters"))),
168  cores_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("cores"))),
169  topoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>()),
170  ptMin_(iConfig.getParameter<double>("ptMin")),
171  deltaR_(iConfig.getParameter<double>("deltaR")),
172  chargeFracMin_(iConfig.getParameter<double>("chargeFractionMin")),
173  centralMIPCharge_(iConfig.getParameter<double>("centralMIPCharge")),
174  pixelCPE_(iConfig.getParameter<std::string>("pixelCPE")),
175  weightfilename_(iConfig.getParameter<edm::FileInPath>("weightFile").fullPath()),
176  inputTensorName_(iConfig.getParameter<std::vector<std::string>>("inputTensorName")),
177  outputTensorName_(iConfig.getParameter<std::vector<std::string>>("outputTensorName")),
178  probThr_(iConfig.getParameter<double>("probThr")),
180 
181 {
182  produces<TrajectorySeedCollection>();
183  produces<reco::TrackCollection>();
184 }

◆ ~DeepCoreSeedGenerator()

DeepCoreSeedGenerator::~DeepCoreSeedGenerator ( )
override

Definition at line 186 of file DeepCoreSeedGenerator.cc.

186 {}

Member Function Documentation

◆ beginJob()

void DeepCoreSeedGenerator::beginJob ( void  )
private

Definition at line 553 of file DeepCoreSeedGenerator.cc.

553 {}

◆ DetectorSelector()

const GeomDet * DeepCoreSeedGenerator::DetectorSelector ( int  llay,
const reco::Candidate jet,
GlobalVector  jetDir,
const reco::Vertex jetVertex,
const TrackerTopology * const  tTopo,
const edmNew::DetSetVector< SiPixelCluster > &  clusters 
)
private

Definition at line 485 of file DeepCoreSeedGenerator.cc.

490  {
491  double minDist = 0.0;
492  GeomDet* output = (GeomDet*)nullptr;
493  for (const auto& detset : clusters) {
494  auto aClusterID = detset.id();
495  if (DetId(aClusterID).subdetId() != 1)
496  continue;
497  const GeomDet* det = geometry_->idToDet(aClusterID);
498  int lay = tTopo->layer(det->geographicalId());
499  if (lay != llay)
500  continue;
501  std::pair<bool, Basic3DVector<float>> interPair =
502  findIntersection(jetDir, (reco::Candidate::Point)jetVertex.position(), det);
503  if (interPair.first == false)
504  continue;
505  Basic3DVector<float> inter = interPair.second;
506  auto localInter = det->specificSurface().toLocal((GlobalPoint)inter);
507  if ((minDist == 0.0 || std::abs(localInter.x()) < minDist) && std::abs(localInter.y()) < 3.35) {
508  minDist = std::abs(localInter.x());
509  output = (GeomDet*)det;
510  }
511  } //detset
512  return output;
513 }

References funct::abs(), bsc_activity_cfg::clusters, findIntersection(), GeomDet::geographicalId(), geometry_, GlobalTrackingGeometry::idToDet(), inter, TrackerTopology::layer(), convertSQLitetoXML_cfg::output, reco::Vertex::position(), GeomDet::specificSurface(), and GloballyPositioned< T >::toLocal().

Referenced by produce().

◆ endJob()

void DeepCoreSeedGenerator::endJob ( void  )
private

Definition at line 556 of file DeepCoreSeedGenerator.cc.

556 {}

◆ fillDescriptions()

void DeepCoreSeedGenerator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 559 of file DeepCoreSeedGenerator.cc.

559  {
561  desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
562  desc.add<edm::InputTag>("pixelClusters", edm::InputTag("siPixelClustersPreSplitting"));
563  desc.add<edm::InputTag>("cores", edm::InputTag("jetsForCoreTracking"));
564  desc.add<double>("ptMin", 100);
565  desc.add<double>("deltaR", 0.25); // the current training makes use of 0.1
566  desc.add<double>("chargeFractionMin", 18000.0);
567  desc.add<double>("centralMIPCharge", 2);
568  desc.add<std::string>("pixelCPE", "PixelCPEGeneric");
569  desc.add<edm::FileInPath>(
570  "weightFile",
571  edm::FileInPath("RecoTracker/TkSeedGenerator/data/DeepCore/DeepCoreSeedGenerator_TrainedModel_barrel_2017.pb"));
572  desc.add<std::vector<std::string>>("inputTensorName", {"input_1", "input_2", "input_3"});
573  desc.add<std::vector<std::string>>("outputTensorName", {"output_node0", "output_node1"});
574  desc.add<double>("probThr", 0.85);
575  descriptions.add("deepCoreSeedGenerator", desc);
576 }

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, HLT_FULL_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

◆ fillPixelMatrix()

void DeepCoreSeedGenerator::fillPixelMatrix ( const SiPixelCluster cluster,
int  layer,
Point3DBase< float, LocalTag inter,
const GeomDet det,
tensorflow::NamedTensorList  input_tensors 
)
private

Definition at line 434 of file DeepCoreSeedGenerator.cc.

438  {
439  int flip = pixelFlipper(det); // 1=not flip, -1=flip
440 
441  for (int i = 0; i < cluster.size(); i++) {
442  SiPixelCluster::Pixel pix = cluster.pixel(i);
443  std::pair<int, int> pixInter = local2Pixel(inter.x(), inter.y(), det);
444  int nx = pix.x - pixInter.first;
445  int ny = pix.y - pixInter.second;
446  nx = flip * nx;
447 
448  if (abs(nx) < jetDimX / 2 && abs(ny) < jetDimY / 2) {
449  nx = nx + jetDimX / 2;
450  ny = ny + jetDimY / 2;
451  //14000 = normalization of ACD counts used in the NN
452  input_tensors[2].second.tensor<float, 4>()(0, nx, ny, layer - 1) += (pix.adc) / (14000.f);
453  }
454  }
455 }

References funct::abs(), SiPixelCluster::Pixel::adc, f, pfNegativeDeepFlavourTagInfos_cfi::flip, dqmMemoryStats::float, mps_fire::i, inter, jetDimX, jetDimY, phase1PixelTopology::layer, local2Pixel(), SiPixelCluster::pixel(), pixelFlipper(), SiPixelCluster::size(), SiPixelCluster::Pixel::x, and SiPixelCluster::Pixel::y.

Referenced by produce().

◆ findIntersection()

std::pair< bool, Basic3DVector< float > > DeepCoreSeedGenerator::findIntersection ( const GlobalVector dir,
const reco::Candidate::Point vertex,
const GeomDet det 
)
private

Definition at line 396 of file DeepCoreSeedGenerator.cc.

398  {
400  Basic3DVector<float>(dir.x(), dir.y(), dir.z()));
401 
402  std::pair<bool, Basic3DVector<float>> pos = vertexPlane.position(det->specificSurface());
403 
404  return pos;
405 }

References DeadROC_duringRun::dir, StraightLinePlaneCrossing::position(), GeomDet::specificSurface(), and bphysicsOniaDQM_cfi::vertex.

Referenced by DetectorSelector(), and produce().

◆ globalEndJob()

void DeepCoreSeedGenerator::globalEndJob ( DeepCoreCache cache)
static

Definition at line 550 of file DeepCoreSeedGenerator.cc.

550 { delete cache->graph_def; }

References utilities::cache().

◆ initializeGlobalCache()

std::unique_ptr< DeepCoreCache > DeepCoreSeedGenerator::initializeGlobalCache ( const edm::ParameterSet iConfig)
static

Definition at line 542 of file DeepCoreSeedGenerator.cc.

542  {
543  // this method is supposed to create, initialize and return a DeepCoreCache instance
544  std::unique_ptr<DeepCoreCache> cache = std::make_unique<DeepCoreCache>();
545  std::string graphPath = iConfig.getParameter<edm::FileInPath>("weightFile").fullPath();
546  cache->graph_def = tensorflow::loadGraphDef(graphPath);
547  return cache;
548 }

References utilities::cache(), contentValuesFiles::fullPath, edm::ParameterSet::getParameter(), tensorflow::loadGraphDef(), and AlCaHLTBitMon_QueryRunRegistry::string.

◆ local2Pixel()

std::pair< int, int > DeepCoreSeedGenerator::local2Pixel ( double  locX,
double  locY,
const GeomDet det 
)
private

Definition at line 407 of file DeepCoreSeedGenerator.cc.

407  {
408  LocalPoint locXY(locX, locY);
409  float pixX = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().pixel(locXY).first;
410  float pixY = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().pixel(locXY).second;
411  std::pair<int, int> out(pixX, pixY);
412  return out;
413 }

References MillePedeFileConverter_cfg::out.

Referenced by fillPixelMatrix(), and produce().

◆ pixel2Local()

LocalPoint DeepCoreSeedGenerator::pixel2Local ( int  pixX,
int  pixY,
const GeomDet det 
)
private

Definition at line 415 of file DeepCoreSeedGenerator.cc.

415  {
416  float locX = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().localX(pixX);
417  float locY = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().localY(pixY);
418  LocalPoint locXY(locX, locY);
419  return locXY;
420 }

Referenced by produce().

◆ pixelFlipper()

int DeepCoreSeedGenerator::pixelFlipper ( const GeomDet det)
private

Definition at line 422 of file DeepCoreSeedGenerator.cc.

422  {
423  int out = 1;
424  LocalVector locZdir(0, 0, 1);
425  GlobalVector globZdir = det->specificSurface().toGlobal(locZdir);
426  const GlobalPoint& globDetCenter = det->position();
427  float direction =
428  globZdir.x() * globDetCenter.x() + globZdir.y() * globDetCenter.y() + globZdir.z() * globDetCenter.z();
429  if (direction < 0)
430  out = -1;
431  return out;
432 }

References MillePedeFileConverter_cfg::out, GeomDet::position(), GeomDet::specificSurface(), Surface::toGlobal(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by fillPixelMatrix(), and produce().

◆ produce()

void DeepCoreSeedGenerator::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 188 of file DeepCoreSeedGenerator.cc.

188  {
189  auto result = std::make_unique<TrajectorySeedCollection>();
190  auto resultTracks = std::make_unique<reco::TrackCollection>();
191 
192  const tensorflow::TensorShape input_size_eta({1, 1});
193  const tensorflow::TensorShape input_size_pt({1, 1});
194  const tensorflow::TensorShape input_size_cluster({1, jetDimX, jetDimY, Nlayer});
195  std::vector<std::string> output_names;
196  output_names.push_back(outputTensorName_[0]);
197  output_names.push_back(outputTensorName_[1]);
198 
199  using namespace edm;
200  using namespace reco;
201 
204  iSetup.get<TrackingComponentsRecord>().get("AnalyticalPropagator", propagator_);
205 
206  const auto& inputPixelClusters_ = iEvent.get(pixelClusters_);
207  const auto& vertices = iEvent.get(vertices_);
208  const auto& cores = iEvent.get(cores_);
209 
212  iSetup.get<TkPixelCPERecord>().get(pixelCPE_, pixelCPEhandle);
213  pixelCPE = pixelCPEhandle.product();
214 
215  const TrackerTopology* const tTopo = &iSetup.getData(topoToken_);
216  auto output = std::make_unique<edmNew::DetSetVector<SiPixelCluster>>();
217 
218  for (const auto& jet : cores) { // jet loop
219 
220  if (jet.pt() > ptMin_) {
221  std::set<long long int> ids;
222  const reco::Vertex& jetVertex = vertices[0];
223 
224  std::vector<GlobalVector> splitClustDirSet =
226  bool l2off = (splitClustDirSet.empty());
227  if (splitClustDirSet.empty()) { //if layer 1 is broken find direcitons on layer 2
228  splitClustDirSet = splittedClusterDirections(jet, tTopo, pixelCPE, jetVertex, 2, inputPixelClusters_);
229  }
230  splitClustDirSet.emplace_back(GlobalVector(jet.px(), jet.py(), jet.pz()));
231  for (int cc = 0; cc < (int)splitClustDirSet.size(); cc++) {
232  tensorflow::NamedTensorList input_tensors; //TensorFlow tensors init
233  input_tensors.resize(3);
234  input_tensors[0] =
235  tensorflow::NamedTensor(inputTensorName_[0], tensorflow::Tensor(tensorflow::DT_FLOAT, input_size_eta));
236  input_tensors[1] =
237  tensorflow::NamedTensor(inputTensorName_[1], tensorflow::Tensor(tensorflow::DT_FLOAT, input_size_pt));
238  input_tensors[2] = tensorflow::NamedTensor(inputTensorName_[2],
239  tensorflow::Tensor(tensorflow::DT_FLOAT, {input_size_cluster}));
240 
241  //put all the input tensor to 0
242  input_tensors[0].second.matrix<float>()(0, 0) = 0.0;
243  input_tensors[1].second.matrix<float>()(0, 0) = 0.0;
244  for (int x = 0; x < jetDimX; x++) {
245  for (int y = 0; y < jetDimY; y++) {
246  for (int l = 0; l < Nlayer; l++) {
247  input_tensors[2].second.tensor<float, 4>()(0, x, y, l) = 0.0;
248  }
249  }
250  } //end of TensorFlow tensors init
251 
252  GlobalVector bigClustDir = splitClustDirSet[cc];
253 
254  jetEta_ = jet.eta();
255  jetPt_ = jet.pt();
256  input_tensors[0].second.matrix<float>()(0, 0) = jet.eta();
257  input_tensors[1].second.matrix<float>()(0, 0) = jet.pt();
258 
259  const GeomDet* globDet = DetectorSelector(
260  2, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_); //det. aligned to bigClustDir on L2
261 
262  if (globDet == nullptr) //no intersection between bigClustDir and pixel detector modules found
263  continue;
264 
265  const GeomDet* goodDet1 = DetectorSelector(1, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_);
266  const GeomDet* goodDet3 = DetectorSelector(3, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_);
267  const GeomDet* goodDet4 = DetectorSelector(4, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_);
268 
269  for (const auto& detset : inputPixelClusters_) {
270  const GeomDet* det = geometry_->idToDet(detset.id());
271 
272  for (const auto& aCluster : detset) {
273  det_id_type aClusterID = detset.id();
274  if (DetId(aClusterID).subdetId() != PixelSubdetector::PixelBarrel)
275  continue;
276 
277  int lay = tTopo->layer(det->geographicalId());
278 
279  std::pair<bool, Basic3DVector<float>> interPair =
280  findIntersection(bigClustDir, (reco::Candidate::Point)jetVertex.position(), det);
281  if (interPair.first == false)
282  continue;
283  Basic3DVector<float> inter = interPair.second;
284  auto localInter = det->specificSurface().toLocal((GlobalPoint)inter);
285 
286  GlobalPoint pointVertex(jetVertex.position().x(), jetVertex.position().y(), jetVertex.position().z());
287 
288  LocalPoint clustPos_local =
289  pixelCPE->localParametersV(aCluster, (*geometry_->idToDetUnit(detset.id())))[0].first;
290 
291  if (std::abs(clustPos_local.x() - localInter.x()) / pitchX_ <= jetDimX / 2 &&
292  std::abs(clustPos_local.y() - localInter.y()) / pitchY_ <=
293  jetDimY / 2) { //used the baricenter, better description maybe useful
294 
295  if (det == goodDet1 || det == goodDet3 || det == goodDet4 || det == globDet) {
296  fillPixelMatrix(aCluster, lay, localInter, det, input_tensors);
297  }
298  } //cluster in ROI
299  } //cluster
300  } //detset
301 
302  //here the NN produce the seed from the filled input
303  std::pair<double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> seedParamNN =
304  DeepCoreSeedGenerator::SeedEvaluation(input_tensors, output_names);
305 
306  for (int i = 0; i < jetDimX; i++) {
307  for (int j = 0; j < jetDimY; j++) {
308  for (int o = 0; o < Nover; o++) {
309  if (seedParamNN.second[i][j][o] > (probThr_ - o * 0.1 - (l2off ? 0.35 : 0))) {
310  std::pair<bool, Basic3DVector<float>> interPair =
311  findIntersection(bigClustDir, (reco::Candidate::Point)jetVertex.position(), globDet);
312  auto localInter = globDet->specificSurface().toLocal((GlobalPoint)interPair.second);
313 
314  int flip = pixelFlipper(globDet); // 1=not flip, -1=flip
315  int nx = i - jetDimX / 2;
316  int ny = j - jetDimY / 2;
317  nx = flip * nx;
318  std::pair<int, int> pixInter = local2Pixel(localInter.x(), localInter.y(), globDet);
319  nx = nx + pixInter.first;
320  ny = ny + pixInter.second;
321  LocalPoint xyLocal = pixel2Local(nx, ny, globDet);
322 
323  double xx = xyLocal.x() + seedParamNN.first[i][j][o][0] * 0.01; //0.01=internal normalization of NN
324  double yy = xyLocal.y() + seedParamNN.first[i][j][o][1] * 0.01; //0.01=internal normalization of NN
325  LocalPoint localSeedPoint = LocalPoint(xx, yy, 0);
326 
327  double track_eta =
328  seedParamNN.first[i][j][o][2] * 0.01 + bigClustDir.eta(); //pay attention to this 0.01
329  double track_theta = 2 * std::atan(std::exp(-track_eta));
330  double track_phi =
331  seedParamNN.first[i][j][o][3] * 0.01 + bigClustDir.phi(); //pay attention to this 0.01
332 
333  double pt = 1. / seedParamNN.first[i][j][o][4];
334  double normdirR = pt / sin(track_theta);
335 
336  const GlobalVector globSeedDir(
337  GlobalVector::Polar(Geom::Theta<double>(track_theta), Geom::Phi<double>(track_phi), normdirR));
338  LocalVector localSeedDir = globDet->surface().toLocal(globSeedDir);
339  int64_t seedid = (int64_t(xx * 200.) << 0) + (int64_t(yy * 200.) << 16) +
340  (int64_t(track_eta * 400.) << 32) + (int64_t(track_phi * 400.) << 48);
341  if (ids.count(seedid) != 0) {
342  continue;
343  }
344  //to disable jetcore iteration comment from here to probability block end--->
345  //seeing iteration skipped, useful to total eff and FakeRate comparison
346  ids.insert(seedid);
347 
348  //Covariance matrix, the hadrcoded variances = NN residuals width
349  //(see https://twiki.cern.ch/twiki/bin/view/CMSPublic/NNJetCoreAtCtD2019)
350  float em[15] = {0,
351  0,
352  0,
353  0,
354  0,
355  0,
356  0,
357  0,
358  0,
359  0,
360  0,
361  0,
362  0,
363  0,
364  0}; // (see LocalTrajectoryError for details), order as follow:
365  em[0] = 0.15 * 0.15; // q/pt
366  em[2] = 0.5e-5; // dxdz
367  em[5] = 0.5e-5; // dydz
368  em[9] = 2e-5; // x
369  em[14] = 2e-5; // y
370  long int detId = globDet->geographicalId();
371  LocalTrajectoryParameters localParam(localSeedPoint, localSeedDir, TrackCharge(1));
372  result->emplace_back(TrajectorySeed(PTrajectoryStateOnDet(localParam, pt, em, detId, /*surfaceSide*/ 0),
375 
376  GlobalPoint globalSeedPoint = globDet->surface().toGlobal(localSeedPoint);
378  resultTracks->emplace_back(
379  reco::Track(1,
380  1,
381  reco::Track::Point(globalSeedPoint.x(), globalSeedPoint.y(), globalSeedPoint.z()),
382  reco::Track::Vector(globSeedDir.x(), globSeedDir.y(), globSeedDir.z()),
383  1,
384  mm));
385  } //probability
386  }
387  }
388  }
389  } //bigcluster
390  } //jet > pt
391  } //jet
392  iEvent.put(std::move(result));
393  iEvent.put(std::move(resultTracks));
394 }

References funct::abs(), alongMomentum, HLT_FULL_cff::cores, cores_, DetectorSelector(), MillePedeFileConverter_cfg::e, PV3DBase< T, PVType, FrameType >::eta(), JetChargeProducer_cfi::exp, fillPixelMatrix(), findIntersection(), first, pfNegativeDeepFlavourTagInfos_cfi::flip, dqmMemoryStats::float, GeomDet::geographicalId(), geometry_, edm::EventSetup::get(), get, edm::EventSetup::getData(), mps_fire::i, GlobalTrackingGeometry::idToDet(), GlobalTrackingGeometry::idToDetUnit(), iEvent, inputPixelClusters_, inputTensorName_, createfilelist::int, inter, dqmiolumiharvest::j, metsig::jet, jetDimX, jetDimY, jetEta_, jetPt_, cmsLHEtoEOSManager::l, TrackerTopology::layer(), local2Pixel(), magfield_, eostools::move(), Nlayer, Nover, EcalTangentSkim_cfg::o, convertSQLitetoXML_cfg::output, outputTensorName_, PV3DBase< T, PVType, FrameType >::phi(), pitchX_, pitchY_, pixel2Local(), PixelSubdetector::PixelBarrel, pixelClusters_, HLT_FULL_cff::pixelCPE, pixelCPE_, pixelFlipper(), reco::Vertex::position(), probThr_, edm::ESHandle< T >::product(), propagator_, DiDispStaMuonMonitor_cfi::pt, ptMin_, mps_fire::result, SeedEvaluation(), funct::sin(), GeomDet::specificSurface(), splittedClusterDirections(), GeomDet::surface(), Surface::toGlobal(), GloballyPositioned< T >::toLocal(), topoToken_, AlignmentTracksFromVertexSelector_cfi::vertices, vertices_, x, PV3DBase< T, PVType, FrameType >::x(), geometryCSVtoXML::xx, y, PV3DBase< T, PVType, FrameType >::y(), geometryCSVtoXML::yy, and PV3DBase< T, PVType, FrameType >::z().

◆ SeedEvaluation()

std::pair< double[DeepCoreSeedGenerator::jetDimX][DeepCoreSeedGenerator::jetDimY][DeepCoreSeedGenerator::Nover][DeepCoreSeedGenerator::Npar], double[DeepCoreSeedGenerator::jetDimX][DeepCoreSeedGenerator::jetDimY][DeepCoreSeedGenerator::Nover]> DeepCoreSeedGenerator::SeedEvaluation ( tensorflow::NamedTensorList  input_tensors,
std::vector< std::string >  output_names 
)
private

Definition at line 460 of file DeepCoreSeedGenerator.cc.

461  {
462  std::vector<tensorflow::Tensor> outputs;
463  tensorflow::run(session_, input_tensors, output_names, &outputs);
464  auto matrix_output_par = outputs.at(0).tensor<float, 5>();
465  auto matrix_output_prob = outputs.at(1).tensor<float, 5>();
466 
467  std::pair<double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> output_combined;
468 
469  for (int x = 0; x < jetDimX; x++) {
470  for (int y = 0; y < jetDimY; y++) {
471  for (int trk = 0; trk < Nover; trk++) {
472  output_combined.second[x][y][trk] =
473  matrix_output_prob(0, x, y, trk, 0); //outputs.at(1).matrix<double>()(0,x,y,trk);
474 
475  for (int p = 0; p < Npar; p++) {
476  output_combined.first[x][y][trk][p] =
477  matrix_output_par(0, x, y, trk, p); //outputs.at(0).matrix<double>()(0,x,y,trk,p);
478  }
479  }
480  }
481  }
482  return output_combined;
483 }

References dqmMemoryStats::float, jetDimX, jetDimY, Nover, Npar, PatBasicFWLiteJetAnalyzer_Selector_cfg::outputs, AlCaHLTBitMon_ParallelJobs::p, tensorflow::run(), session_, x, and y.

Referenced by produce().

◆ splittedClusterDirections()

std::vector< GlobalVector > DeepCoreSeedGenerator::splittedClusterDirections ( const reco::Candidate jet,
const TrackerTopology * const  tTopo,
const PixelClusterParameterEstimator pixelCPE,
const reco::Vertex jetVertex,
int  layer,
const edmNew::DetSetVector< SiPixelCluster > &  clusters 
)
private

Definition at line 515 of file DeepCoreSeedGenerator.cc.

521  {
522  std::vector<GlobalVector> clustDirs;
523  for (const auto& detset_int : clusters) {
524  const GeomDet* det_int = geometry_->idToDet(detset_int.id());
525  int lay = tTopo->layer(det_int->geographicalId());
526  if (lay != layer)
527  continue; //NB: saved bigClusters on all the layers!!
528  auto detUnit = geometry_->idToDetUnit(detset_int.id());
529  for (const auto& aCluster : detset_int) {
530  GlobalPoint clustPos = det_int->surface().toGlobal(pixelCPE->localParametersV(aCluster, (*detUnit))[0].first);
531  GlobalPoint vertexPos(jetVertex.position().x(), jetVertex.position().y(), jetVertex.position().z());
532  GlobalVector clusterDir = clustPos - vertexPos;
533  GlobalVector jetDir(jet.px(), jet.py(), jet.pz());
534  if (Geom::deltaR(jetDir, clusterDir) < deltaR_) {
535  clustDirs.emplace_back(clusterDir);
536  }
537  }
538  }
539  return clustDirs;
540 }

References bsc_activity_cfg::clusters, PbPb_ZMuSkimMuonDPG_cff::deltaR, deltaR_, GeomDet::geographicalId(), geometry_, GlobalTrackingGeometry::idToDet(), GlobalTrackingGeometry::idToDetUnit(), metsig::jet, phase1PixelTopology::layer, TrackerTopology::layer(), HLT_FULL_cff::pixelCPE, reco::Vertex::position(), GeomDet::surface(), and Surface::toGlobal().

Referenced by produce().

Member Data Documentation

◆ centralMIPCharge_

double DeepCoreSeedGenerator::centralMIPCharge_
private

Definition at line 121 of file DeepCoreSeedGenerator.cc.

◆ chargeFracMin_

double DeepCoreSeedGenerator::chargeFracMin_
private

Definition at line 120 of file DeepCoreSeedGenerator.cc.

◆ cores_

edm::EDGetTokenT<edm::View<reco::Candidate> > DeepCoreSeedGenerator::cores_
private

Definition at line 115 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ deltaR_

double DeepCoreSeedGenerator::deltaR_
private

Definition at line 119 of file DeepCoreSeedGenerator.cc.

Referenced by splittedClusterDirections().

◆ geometry_

edm::ESHandle<GlobalTrackingGeometry> DeepCoreSeedGenerator::geometry_
private

Definition at line 109 of file DeepCoreSeedGenerator.cc.

Referenced by DetectorSelector(), produce(), and splittedClusterDirections().

◆ inputPixelClusters_

edm::Handle<edmNew::DetSetVector<SiPixelCluster> > DeepCoreSeedGenerator::inputPixelClusters_
private

Definition at line 114 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ inputTensorName_

std::vector<std::string> DeepCoreSeedGenerator::inputTensorName_
private

Definition at line 124 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ jetDimX

constexpr int DeepCoreSeedGenerator::jetDimX = 30
staticconstexpr

Definition at line 95 of file DeepCoreSeedGenerator.cc.

Referenced by fillPixelMatrix(), produce(), and SeedEvaluation().

◆ jetDimY

constexpr int DeepCoreSeedGenerator::jetDimY = 30
staticconstexpr

Definition at line 96 of file DeepCoreSeedGenerator.cc.

Referenced by fillPixelMatrix(), produce(), and SeedEvaluation().

◆ jetEta_

double DeepCoreSeedGenerator::jetEta_

Definition at line 92 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ jetPt_

double DeepCoreSeedGenerator::jetPt_

Definition at line 91 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ magfield_

edm::ESHandle<MagneticField> DeepCoreSeedGenerator::magfield_
private

Definition at line 108 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ Nlayer

constexpr int DeepCoreSeedGenerator::Nlayer = 4
staticconstexpr

Definition at line 97 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ Nover

constexpr int DeepCoreSeedGenerator::Nover = 3
staticconstexpr

Definition at line 98 of file DeepCoreSeedGenerator.cc.

Referenced by produce(), and SeedEvaluation().

◆ Npar

constexpr int DeepCoreSeedGenerator::Npar = 5
staticconstexpr

Definition at line 99 of file DeepCoreSeedGenerator.cc.

Referenced by SeedEvaluation().

◆ outputTensorName_

std::vector<std::string> DeepCoreSeedGenerator::outputTensorName_
private

Definition at line 125 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pitchX_

double DeepCoreSeedGenerator::pitchX_ = 0.01

Definition at line 93 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pitchY_

double DeepCoreSeedGenerator::pitchY_ = 0.015

Definition at line 94 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pixelClusters_

edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > DeepCoreSeedGenerator::pixelClusters_
private

Definition at line 113 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pixelCPE_

std::string DeepCoreSeedGenerator::pixelCPE_
private

Definition at line 122 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ probThr_

double DeepCoreSeedGenerator::probThr_
private

Definition at line 126 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ propagator_

edm::ESHandle<Propagator> DeepCoreSeedGenerator::propagator_
private

Definition at line 110 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ propagatorName_

std::string DeepCoreSeedGenerator::propagatorName_
private

Definition at line 107 of file DeepCoreSeedGenerator.cc.

◆ ptMin_

double DeepCoreSeedGenerator::ptMin_
private

Definition at line 118 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ session_

tensorflow::Session* DeepCoreSeedGenerator::session_
private

Definition at line 127 of file DeepCoreSeedGenerator.cc.

Referenced by SeedEvaluation().

◆ topoToken_

const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> DeepCoreSeedGenerator::topoToken_
private

Definition at line 116 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ vertices_

edm::EDGetTokenT<std::vector<reco::Vertex> > DeepCoreSeedGenerator::vertices_
private

Definition at line 112 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ weightfilename_

std::string DeepCoreSeedGenerator::weightfilename_
private

Definition at line 123 of file DeepCoreSeedGenerator.cc.

Vector3DBase< float, LocalTag >
GeomDet::position
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
tensorflow::createSession
Session * createSession(SessionOptions &sessionOptions)
Definition: TensorFlow.cc:85
PixelClusterParameterEstimator
Definition: PixelClusterParameterEstimator.h:15
DDAxes::y
DeepCoreSeedGenerator::local2Pixel
std::pair< int, int > local2Pixel(double, double, const GeomDet *)
Definition: DeepCoreSeedGenerator.cc:407
mps_fire.i
i
Definition: mps_fire.py:428
DeepCoreSeedGenerator::outputTensorName_
std::vector< std::string > outputTensorName_
Definition: DeepCoreSeedGenerator.cc:125
DeepCoreSeedGenerator::weightfilename_
std::string weightfilename_
Definition: DeepCoreSeedGenerator.cc:123
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
DeepCoreSeedGenerator::Npar
static constexpr int Npar
Definition: DeepCoreSeedGenerator.cc:99
TrackCharge
int TrackCharge
Definition: TrackCharge.h:4
GeomDet
Definition: GeomDet.h:27
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
DeepCoreSeedGenerator::geometry_
edm::ESHandle< GlobalTrackingGeometry > geometry_
Definition: DeepCoreSeedGenerator.cc:109
f
double f[11][100]
Definition: MuScleFitUtils.cc:78
DeepCoreSeedGenerator::pixel2Local
LocalPoint pixel2Local(int, int, const GeomDet *)
Definition: DeepCoreSeedGenerator.cc:415
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
contentValuesFiles.fullPath
fullPath
Definition: contentValuesFiles.py:64
edm
HLT enums.
Definition: AlignableModifier.h:19
TrackerTopology
Definition: TrackerTopology.h:16
tensorflow::NamedTensor
std::pair< std::string, Tensor > NamedTensor
Definition: TensorFlow.h:29
Geom::Theta
Definition: Theta.h:12
pos
Definition: PixelAliasList.h:18
PatBasicFWLiteJetAnalyzer_Selector_cfg.outputs
outputs
Definition: PatBasicFWLiteJetAnalyzer_Selector_cfg.py:48
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
TkPixelCPERecord
Definition: TkPixelCPERecord.h:18
DeepCoreSeedGenerator::centralMIPCharge_
double centralMIPCharge_
Definition: DeepCoreSeedGenerator.cc:121
Geom::Spherical2Cartesian
Definition: CoordinateSets.h:58
GlobalTrackingGeometryRecord
Definition: GlobalTrackingGeometryRecord.h:17
TrackerTopology::layer
unsigned int layer(const DetId &id) const
Definition: TrackerTopology.cc:47
reco::Vertex::position
const Point & position() const
position
Definition: Vertex.h:127
DDAxes::x
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
align::LocalPoint
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
GlobalVector
Global3DVector GlobalVector
Definition: GlobalVector.h:10
DeepCoreSeedGenerator::vertices_
edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_
Definition: DeepCoreSeedGenerator.cc:112
DeepCoreSeedGenerator::splittedClusterDirections
std::vector< GlobalVector > splittedClusterDirections(const reco::Candidate &, const TrackerTopology *const, const PixelClusterParameterEstimator *, const reco::Vertex &, int, const edmNew::DetSetVector< SiPixelCluster > &)
Definition: DeepCoreSeedGenerator.cc:515
DeepCoreSeedGenerator::fillPixelMatrix
void fillPixelMatrix(const SiPixelCluster &, int, Point3DBase< float, LocalTag >, const GeomDet *, tensorflow::NamedTensorList)
Definition: DeepCoreSeedGenerator.cc:434
SiPixelCluster::Pixel
Definition: SiPixelCluster.h:31
DeepCoreSeedGenerator::pixelCPE_
std::string pixelCPE_
Definition: DeepCoreSeedGenerator.cc:122
AlignmentTracksFromVertexSelector_cfi.vertices
vertices
Definition: AlignmentTracksFromVertexSelector_cfi.py:5
EcalTangentSkim_cfg.o
o
Definition: EcalTangentSkim_cfg.py:42
cc
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
LocalTrajectoryParameters
Definition: LocalTrajectoryParameters.h:25
DeepCoreSeedGenerator::pitchX_
double pitchX_
Definition: DeepCoreSeedGenerator.cc:93
IdealMagneticFieldRecord
Definition: IdealMagneticFieldRecord.h:11
DetId
Definition: DetId.h:17
edm::FileInPath
Definition: FileInPath.h:61
GeomDet::surface
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
SiPixelCluster::Pixel::adc
uint16_t adc
Definition: SiPixelCluster.h:37
HLT_FULL_cff.pixelCPE
pixelCPE
Definition: HLT_FULL_cff.py:50744
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
pfNegativeDeepFlavourTagInfos_cfi.flip
flip
Definition: pfNegativeDeepFlavourTagInfos_cfi.py:8
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
DeepCoreSeedGenerator::Nlayer
static constexpr int Nlayer
Definition: DeepCoreSeedGenerator.cc:97
DeepCoreSeedGenerator::deltaR_
double deltaR_
Definition: DeepCoreSeedGenerator.cc:119
Surface::toGlobal
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
reco::Track
Definition: Track.h:27
edm::ESHandle
Definition: DTSurvey.h:22
DeepCoreSeedGenerator::ptMin_
double ptMin_
Definition: DeepCoreSeedGenerator.cc:118
Point3DBase< float, GlobalTag >
first
auto first
Definition: CAHitNtupletGeneratorKernelsImpl.h:125
utilities.cache
def cache(function)
Definition: utilities.py:3
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
DeepCoreSeedGenerator::findIntersection
std::pair< bool, Basic3DVector< float > > findIntersection(const GlobalVector &, const reco::Candidate::Point &, const GeomDet *)
Definition: DeepCoreSeedGenerator.cc:396
edm::View
Definition: CaloClusterFwd.h:14
GeomDet::geographicalId
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
bphysicsOniaDQM_cfi.vertex
vertex
Definition: bphysicsOniaDQM_cfi.py:7
DeepCoreSeedGenerator::pixelClusters_
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
Definition: DeepCoreSeedGenerator.cc:113
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
geometryCSVtoXML.yy
yy
Definition: geometryCSVtoXML.py:19
DeepCoreSeedGenerator::cores_
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
Definition: DeepCoreSeedGenerator.cc:115
DeepCoreSeedGenerator::session_
tensorflow::Session * session_
Definition: DeepCoreSeedGenerator.cc:127
DeepCoreSeedGenerator::inputTensorName_
std::vector< std::string > inputTensorName_
Definition: DeepCoreSeedGenerator.cc:124
tensorflow::NamedTensorList
std::vector< NamedTensor > NamedTensorList
Definition: TensorFlow.h:30
PV3DBase::eta
T eta() const
Definition: PV3DBase.h:73
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
DeepCoreSeedGenerator::magfield_
edm::ESHandle< MagneticField > magfield_
Definition: DeepCoreSeedGenerator.cc:108
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
Geom::Phi
Definition: Phi.h:52
SiPixelCluster::Pixel::y
uint16_t y
Definition: SiPixelCluster.h:36
DeepCoreSeedGenerator::jetEta_
double jetEta_
Definition: DeepCoreSeedGenerator.cc:92
SiPixelCluster::size
int size() const
Definition: SiPixelCluster.h:134
SiPixelCluster::pixel
Pixel pixel(int i) const
Definition: SiPixelCluster.h:170
GeomDet::specificSurface
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
DeepCoreSeedGenerator::inputPixelClusters_
edm::Handle< edmNew::DetSetVector< SiPixelCluster > > inputPixelClusters_
Definition: DeepCoreSeedGenerator.cc:114
get
#define get
DeepCoreSeedGenerator::propagator_
edm::ESHandle< Propagator > propagator_
Definition: DeepCoreSeedGenerator.cc:110
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
cmsLHEtoEOSManager.l
l
Definition: cmsLHEtoEOSManager.py:204
DeepCoreSeedGenerator::Nover
static constexpr int Nover
Definition: DeepCoreSeedGenerator.cc:98
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
reco::TrackBase::Point
math::XYZPoint Point
point in the space
Definition: TrackBase.h:80
tensorflow::loadGraphDef
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:68
DeepCoreSeedGenerator::pitchY_
double pitchY_
Definition: DeepCoreSeedGenerator.cc:94
edmNew::DetSetVector
Definition: DetSetNew.h:13
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
DeepCoreSeedGenerator::DetectorSelector
const GeomDet * DetectorSelector(int, const reco::Candidate &, GlobalVector, const reco::Vertex &, const TrackerTopology *const, const edmNew::DetSetVector< SiPixelCluster > &)
Definition: DeepCoreSeedGenerator.cc:485
DeepCoreSeedGenerator::jetPt_
double jetPt_
Definition: DeepCoreSeedGenerator.cc:91
GlobalTrackingGeometry::idToDetUnit
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: GlobalTrackingGeometry.cc:34
DeepCoreSeedGenerator::probThr_
double probThr_
Definition: DeepCoreSeedGenerator.cc:126
DeepCoreSeedGenerator::jetDimX
static constexpr int jetDimX
Definition: DeepCoreSeedGenerator.cc:95
metsig::jet
Definition: SignAlgoResolutions.h:47
GloballyPositioned::toLocal
LocalPoint toLocal(const GlobalPoint &gp) const
Definition: GloballyPositioned.h:98
GlobalTrackingGeometry::idToDet
const GeomDet * idToDet(DetId) const override
Definition: GlobalTrackingGeometry.cc:44
tensorflow::run
void run(Session *session, const NamedTensorList &inputs, const std::vector< std::string > &outputNames, std::vector< Tensor > *outputs, const thread::ThreadPoolOptions &threadPoolOptions)
Definition: TensorFlow.cc:213
DeepCoreSeedGenerator::pixelFlipper
int pixelFlipper(const GeomDet *)
Definition: DeepCoreSeedGenerator.cc:422
HLT_FULL_cff.cores
cores
Definition: HLT_FULL_cff.py:104308
TrajectorySeed
Definition: TrajectorySeed.h:18
DeepCoreSeedGenerator::chargeFracMin_
double chargeFracMin_
Definition: DeepCoreSeedGenerator.cc:120
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MillePedeFileConverter_cfg.out
out
Definition: MillePedeFileConverter_cfg.py:31
inter
int inter
Definition: CascadeWrapper.h:27
mps_fire.result
result
Definition: mps_fire.py:311
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::det_id_type
uint32_t det_id_type
Definition: DetSet.h:20
PTrajectoryStateOnDet
Definition: PTrajectoryStateOnDet.h:10
SiPixelCluster::Pixel::x
uint16_t x
Definition: SiPixelCluster.h:35
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
reco::Candidate::Point
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6
DeepCoreSeedGenerator::jetDimY
static constexpr int jetDimY
Definition: DeepCoreSeedGenerator.cc:96
reco::TrackBase::CovarianceMatrix
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:74
StraightLinePlaneCrossing
Definition: StraightLinePlaneCrossing.h:14
edm::InputTag
Definition: InputTag.h:15
alongMomentum
Definition: PropagationDirection.h:4
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Basic3DVector< float >
geometryCSVtoXML.xx
xx
Definition: geometryCSVtoXML.py:19
reco::Vertex
Definition: Vertex.h:35
edm::FileInPath::fullPath
std::string fullPath() const
Definition: FileInPath.cc:161
edm::OwnVector< TrackingRecHit >
reco::TrackBase::Vector
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:77
DeepCoreSeedGenerator::topoToken_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
Definition: DeepCoreSeedGenerator.cc:116
DeadROC_duringRun.dir
dir
Definition: DeadROC_duringRun.py:23
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
DeepCoreSeedGenerator::SeedEvaluation
std::pair< double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> SeedEvaluation(tensorflow::NamedTensorList, std::vector< std::string >)
Definition: DeepCoreSeedGenerator.cc:460
TrackingComponentsRecord
Definition: TrackingComponentsRecord.h:12