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< tensorflow::SessionCache > >

Public Member Functions

 DeepCoreSeedGenerator (const edm::ParameterSet &, const tensorflow::SessionCache *)
 
 ~DeepCoreSeedGenerator () override
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< tensorflow::SessionCache > >
 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 (tensorflow::SessionCache *)
 
static std::unique_ptr< tensorflow::SessionCacheinitializeGlobalCache (const edm::ParameterSet &)
 

Public Attributes

double dth [3] = {0.0, 0.15, 0.3}
 
double dthl [3] = {0.1, 0.05, 0}
 
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

const GeomDetDetectorSelector (int, const reco::Candidate &, GlobalVector, const reco::Vertex &, const TrackerTopology *const, const edmNew::DetSetVector< SiPixelCluster > &)
 
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_
 
const GlobalTrackingGeometrygeometry_ = nullptr
 
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecordgeometryToken_
 
edm::Handle< edmNew::DetSetVector< SiPixelCluster > > inputPixelClusters_
 
std::vector< std::string > inputTensorName_
 
std::vector< std::string > outputTensorName_
 
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
 
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecordpixelCPEToken_
 
double probThr_
 
std::string propagatorName_
 
double ptMin_
 
const 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< tensorflow::SessionCache > >
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 76 of file DeepCoreSeedGenerator.cc.

Constructor & Destructor Documentation

◆ DeepCoreSeedGenerator()

DeepCoreSeedGenerator::DeepCoreSeedGenerator ( const edm::ParameterSet iConfig,
const tensorflow::SessionCache cache 
)
explicit

Definition at line 162 of file DeepCoreSeedGenerator.cc.

163  : vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
165  consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelClusters"))),
166  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  weightfilename_(iConfig.getParameter<edm::FileInPath>("weightFile").fullPath()),
175  inputTensorName_(iConfig.getParameter<std::vector<std::string>>("inputTensorName")),
176  outputTensorName_(iConfig.getParameter<std::vector<std::string>>("outputTensorName")),
177  probThr_(iConfig.getParameter<double>("probThr")),
178  session_(cache->getSession())
179 
180 {
181  produces<TrajectorySeedCollection>();
182  produces<reco::TrackCollection>();
183 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > geometryToken_
std::string fullPath() const
Definition: FileInPath.cc:161
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > pixelCPEToken_
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
std::vector< std::string > inputTensorName_
std::vector< std::string > outputTensorName_
edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_
def cache(function)
Definition: utilities.py:3
const tensorflow::Session * session_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_

◆ ~DeepCoreSeedGenerator()

DeepCoreSeedGenerator::~DeepCoreSeedGenerator ( )
override

Definition at line 185 of file DeepCoreSeedGenerator.cc.

185 {}

Member Function Documentation

◆ 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 486 of file DeepCoreSeedGenerator.cc.

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

Referenced by produce().

491  {
492  double minDist = 0.0;
493  GeomDet* output = (GeomDet*)nullptr;
494  for (const auto& detset : clusters) {
495  auto aClusterID = detset.id();
496  if (DetId(aClusterID).subdetId() != 1)
497  continue;
498  const GeomDet* det = geometry_->idToDet(aClusterID);
499  int lay = tTopo->layer(det->geographicalId());
500  if (lay != llay)
501  continue;
502  std::pair<bool, Basic3DVector<float>> interPair =
503  findIntersection(jetDir, (reco::Candidate::Point)jetVertex.position(), det);
504  if (interPair.first == false)
505  continue;
506  Basic3DVector<float> inter = interPair.second;
507  auto localInter = det->specificSurface().toLocal((GlobalPoint)inter);
508  if ((minDist == 0.0 || std::abs(localInter.x()) < minDist) && std::abs(localInter.y()) < 3.35) {
509  minDist = std::abs(localInter.x());
510  output = (GeomDet*)det;
511  }
512  } //detset
513  return output;
514 }
const Point & position() const
position
Definition: Vertex.h:128
LocalPoint toLocal(const GlobalPoint &gp) const
const GeomDet * idToDet(DetId) const override
unsigned int layer(const DetId &id) const
std::pair< bool, Basic3DVector< float > > findIntersection(const GlobalVector &, const reco::Candidate::Point &, const GeomDet *)
const GlobalTrackingGeometry * geometry_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
Definition: DetId.h:17
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
Definition: output.py:1
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40

◆ fillDescriptions()

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

Definition at line 554 of file DeepCoreSeedGenerator.cc.

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

554  {
556  desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
557  desc.add<edm::InputTag>("pixelClusters", edm::InputTag("siPixelClustersPreSplitting"));
558  desc.add<edm::InputTag>("cores", edm::InputTag("jetsForCoreTracking"));
559  desc.add<double>("ptMin", 100); // the current training uses 500 GeV
560  desc.add<double>("deltaR", 0.25); // the current training uses 0.1
561  desc.add<double>("chargeFractionMin", 18000.0);
562  desc.add<double>("centralMIPCharge", 2);
563  desc.add<std::string>("pixelCPE", "PixelCPEGeneric");
564  desc.add<edm::FileInPath>(
565  "weightFile",
566  edm::FileInPath("RecoTracker/TkSeedGenerator/data/DeepCore/DeepCoreSeedGenerator_TrainedModel_barrel_2023.pb"));
567  desc.add<std::vector<std::string>>("inputTensorName", {"input_1", "input_2", "input_3"});
568  desc.add<std::vector<std::string>>("outputTensorName", {"output_node0", "output_node1"});
569  desc.add<double>("probThr", 0.7); // DeepCore 2.2.1 baseline threshold
570  descriptions.add("deepCoreSeedGenerator", desc);
571 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ fillPixelMatrix()

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

Definition at line 435 of file DeepCoreSeedGenerator.cc.

References funct::abs(), SiPixelCluster::Pixel::adc, f, HLT_2024v14_cff::flip, nano_mu_digi_cff::float, mps_fire::i, jetDimX, jetDimY, nano_mu_digi_cff::layer, local2Pixel(), SiPixelCluster::pixel(), pixelFlipper(), SiPixelCluster::size(), SiPixelCluster::Pixel::x, PV3DBase< T, PVType, FrameType >::x(), SiPixelCluster::Pixel::y, and PV3DBase< T, PVType, FrameType >::y().

Referenced by produce().

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

◆ findIntersection()

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

Definition at line 397 of file DeepCoreSeedGenerator.cc.

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

Referenced by DetectorSelector(), and produce().

399  {
401  Basic3DVector<float>(dir.x(), dir.y(), dir.z()));
402 
403  std::pair<bool, Basic3DVector<float>> pos = vertexPlane.position(det->specificSurface());
404 
405  return pos;
406 }
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40

◆ globalEndJob()

void DeepCoreSeedGenerator::globalEndJob ( tensorflow::SessionCache cache)
static

Definition at line 551 of file DeepCoreSeedGenerator.cc.

551 {}

◆ initializeGlobalCache()

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

Definition at line 543 of file DeepCoreSeedGenerator.cc.

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

544  {
545  // this method is supposed to create, initialize and return a Tensorflow:SessionCache instance
546  std::string graphPath = iConfig.getParameter<edm::FileInPath>("weightFile").fullPath();
547  std::unique_ptr<tensorflow::SessionCache> cache = std::make_unique<tensorflow::SessionCache>(graphPath);
548  return cache;
549 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
def cache(function)
Definition: utilities.py:3

◆ local2Pixel()

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

Definition at line 408 of file DeepCoreSeedGenerator.cc.

References MillePedeFileConverter_cfg::out.

Referenced by fillPixelMatrix(), and produce().

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

◆ pixel2Local()

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

Definition at line 416 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

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

◆ pixelFlipper()

int DeepCoreSeedGenerator::pixelFlipper ( const GeomDet det)
private

Definition at line 423 of file DeepCoreSeedGenerator.cc.

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().

423  {
424  int out = 1;
425  LocalVector locZdir(0, 0, 1);
426  GlobalVector globZdir = det->specificSurface().toGlobal(locZdir);
427  const GlobalPoint& globDetCenter = det->position();
428  float direction =
429  globZdir.x() * globDetCenter.x() + globZdir.y() * globDetCenter.y() + globZdir.z() * globDetCenter.z();
430  if (direction < 0)
431  out = -1;
432  return out;
433 }
T z() const
Definition: PV3DBase.h:61
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:43
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40

◆ produce()

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

Definition at line 187 of file DeepCoreSeedGenerator.cc.

References funct::abs(), alongMomentum, gpuPixelDoublets::cc, HLT_FULL_cff::cores, cores_, DetectorSelector(), hcalRecHitTable_cff::detId, dth, dthl, MillePedeFileConverter_cfg::e, PV3DBase< T, PVType, FrameType >::eta(), JetChargeProducer_cfi::exp, fillPixelMatrix(), findIntersection(), dqmdumpme::first, HLT_2024v14_cff::flip, nano_mu_digi_cff::float, GeomDet::geographicalId(), geometry_, geometryToken_, edm::EventSetup::getData(), mps_fire::i, GlobalTrackingGeometry::idToDet(), GlobalTrackingGeometry::idToDetUnit(), iEvent, inputPixelClusters_, inputTensorName_, createfilelist::int, dqmiolumiharvest::j, metsig::jet, jetDimX, jetDimY, jetEta_, jetPt_, MainPageGenerator::l, TrackerTopology::layer(), local2Pixel(), eostools::move(), Nlayer, Nover, EcalTangentSkim_cfg::o, HLT_2024v14_cff::output_names, outputTensorName_, PV3DBase< T, PVType, FrameType >::phi(), pitchX_, pitchY_, pixel2Local(), PixelSubdetector::PixelBarrel, pixelClusters_, HLT_2024v14_cff::pixelCPE, pixelCPEToken_, pixelFlipper(), reco::Vertex::position(), probThr_, DiDispStaMuonMonitor_cfi::pt, ptMin_, mps_fire::result, SeedEvaluation(), funct::sin(), GeomDet::specificSurface(), splittedClusterDirections(), GeomDet::surface(), Surface::toGlobal(), GloballyPositioned< T >::toLocal(), topoToken_, leptonTimeLifeInfo_common_cff::track_phi, 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().

187  {
188  auto result = std::make_unique<TrajectorySeedCollection>();
189  auto resultTracks = std::make_unique<reco::TrackCollection>();
190 
191  const tensorflow::TensorShape input_size_eta({1, 1});
192  const tensorflow::TensorShape input_size_pt({1, 1});
193  const tensorflow::TensorShape input_size_cluster({1, jetDimX, jetDimY, Nlayer});
194  std::vector<std::string> output_names;
195  output_names.push_back(outputTensorName_[0]);
196  output_names.push_back(outputTensorName_[1]);
197 
198  using namespace edm;
199  using namespace reco;
200 
201  geometry_ = &iSetup.getData(geometryToken_);
202 
203  const auto& inputPixelClusters_ = iEvent.get(pixelClusters_);
204  const auto& vertices = iEvent.get(vertices_);
205  const auto& cores = iEvent.get(cores_);
206 
208 
209  const TrackerTopology* const tTopo = &iSetup.getData(topoToken_);
210  auto output = std::make_unique<edmNew::DetSetVector<SiPixelCluster>>();
211 
212  for (const auto& jet : cores) { // jet loop
213 
214  if (jet.pt() > ptMin_) {
215  std::set<unsigned long long int> ids;
216  const reco::Vertex& jetVertex = vertices[0];
217 
218  std::vector<GlobalVector> splitClustDirSet =
220  bool l2off = (splitClustDirSet.empty()); // no adc BPIX2
221  splitClustDirSet = splittedClusterDirections(jet, tTopo, pixelCPE, jetVertex, 3, inputPixelClusters_);
222  bool l134off = (splitClustDirSet.empty()); // no adc BPIX1 or BPIX3 or BPIX4
223  splitClustDirSet = splittedClusterDirections(jet, tTopo, pixelCPE, jetVertex, 4, inputPixelClusters_);
224  l134off = (splitClustDirSet.empty() && l134off);
225  splitClustDirSet = splittedClusterDirections(jet, tTopo, pixelCPE, jetVertex, 1, inputPixelClusters_);
226  l134off = (splitClustDirSet.empty() && l134off);
227 
228  if (splitClustDirSet.empty()) { //if layer 1 is broken find direcitons on layer 2
229  splitClustDirSet = splittedClusterDirections(jet, tTopo, pixelCPE, jetVertex, 2, inputPixelClusters_);
230  }
231  splitClustDirSet.emplace_back(GlobalVector(jet.px(), jet.py(), jet.pz()));
232  for (int cc = 0; cc < (int)splitClustDirSet.size(); cc++) {
233  tensorflow::NamedTensorList input_tensors; //TensorFlow tensors init
234  input_tensors.resize(3);
235  input_tensors[0] =
236  tensorflow::NamedTensor(inputTensorName_[0], tensorflow::Tensor(tensorflow::DT_FLOAT, input_size_eta));
237  input_tensors[1] =
238  tensorflow::NamedTensor(inputTensorName_[1], tensorflow::Tensor(tensorflow::DT_FLOAT, input_size_pt));
239  input_tensors[2] = tensorflow::NamedTensor(inputTensorName_[2],
240  tensorflow::Tensor(tensorflow::DT_FLOAT, {input_size_cluster}));
241 
242  //put all the input tensor to 0
243  input_tensors[0].second.matrix<float>()(0, 0) = 0.0;
244  input_tensors[1].second.matrix<float>()(0, 0) = 0.0;
245  for (int x = 0; x < jetDimX; x++) {
246  for (int y = 0; y < jetDimY; y++) {
247  for (int l = 0; l < Nlayer; l++) {
248  input_tensors[2].second.tensor<float, 4>()(0, x, y, l) = 0.0;
249  }
250  }
251  } //end of TensorFlow tensors init
252 
253  GlobalVector bigClustDir = splitClustDirSet[cc];
254 
255  jetEta_ = jet.eta();
256  jetPt_ = jet.pt();
257  input_tensors[0].second.matrix<float>()(0, 0) = jet.eta();
258  input_tensors[1].second.matrix<float>()(0, 0) = jet.pt();
259 
260  const GeomDet* globDet = DetectorSelector(
261  2, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_); //det. aligned to bigClustDir on L2
262 
263  if (globDet == nullptr) //no intersection between bigClustDir and pixel detector modules found
264  continue;
265 
266  const GeomDet* goodDet1 = DetectorSelector(1, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_);
267  const GeomDet* goodDet3 = DetectorSelector(3, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_);
268  const GeomDet* goodDet4 = DetectorSelector(4, jet, bigClustDir, jetVertex, tTopo, inputPixelClusters_);
269 
270  for (const auto& detset : inputPixelClusters_) {
271  const GeomDet* det = geometry_->idToDet(detset.id());
272 
273  for (const auto& aCluster : detset) {
274  det_id_type aClusterID = detset.id();
275  if (DetId(aClusterID).subdetId() != PixelSubdetector::PixelBarrel)
276  continue;
277 
278  int lay = tTopo->layer(det->geographicalId());
279 
280  std::pair<bool, Basic3DVector<float>> interPair =
281  findIntersection(bigClustDir, (reco::Candidate::Point)jetVertex.position(), det);
282  if (interPair.first == false)
283  continue;
284  Basic3DVector<float> inter = interPair.second;
285  auto localInter = det->specificSurface().toLocal((GlobalPoint)inter);
286 
287  GlobalPoint pointVertex(jetVertex.position().x(), jetVertex.position().y(), jetVertex.position().z());
288 
289  LocalPoint clustPos_local =
290  pixelCPE->localParametersV(aCluster, (*geometry_->idToDetUnit(detset.id())))[0].first;
291 
292  if (std::abs(clustPos_local.x() - localInter.x()) / pitchX_ <= jetDimX / 2 &&
293  std::abs(clustPos_local.y() - localInter.y()) / pitchY_ <=
294  jetDimY / 2) { //used the baricenter, better description maybe useful
295 
296  if (det == goodDet1 || det == goodDet3 || det == goodDet4 || det == globDet) {
297  fillPixelMatrix(aCluster, lay, localInter, det, input_tensors);
298  }
299  } //cluster in ROI
300  } //cluster
301  } //detset
302 
303  //here the NN produce the seed from the filled input
304  std::pair<double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> seedParamNN =
306 
307  for (int i = 0; i < jetDimX; i++) {
308  for (int j = 0; j < jetDimY; j++) {
309  for (int o = 0; o < Nover; o++) {
310  if (seedParamNN.second[i][j][o] >
311  (probThr_ + dth[o] + (l2off ? 0.3 : 0) + (l134off ? dthl[o] : 0))) { // DeepCore 2.2.1 Threshold
312  std::pair<bool, Basic3DVector<float>> interPair =
313  findIntersection(bigClustDir, (reco::Candidate::Point)jetVertex.position(), globDet);
314  auto localInter = globDet->specificSurface().toLocal((GlobalPoint)interPair.second);
315 
316  int flip = pixelFlipper(globDet); // 1=not flip, -1=flip
317  int nx = i - jetDimX / 2;
318  int ny = j - jetDimY / 2;
319  nx = flip * nx;
320  std::pair<int, int> pixInter = local2Pixel(localInter.x(), localInter.y(), globDet);
321  nx = nx + pixInter.first;
322  ny = ny + pixInter.second;
323  LocalPoint xyLocal = pixel2Local(nx, ny, globDet);
324 
325  double xx = xyLocal.x() + seedParamNN.first[i][j][o][0] * 0.01; //0.01=internal normalization of NN
326  double yy = xyLocal.y() + seedParamNN.first[i][j][o][1] * 0.01; //0.01=internal normalization of NN
327  LocalPoint localSeedPoint = LocalPoint(xx, yy, 0);
328 
329  double track_eta =
330  seedParamNN.first[i][j][o][2] * 0.01 + bigClustDir.eta(); //pay attention to this 0.01
331  double track_theta = 2 * std::atan(std::exp(-track_eta));
332  double track_phi =
333  seedParamNN.first[i][j][o][3] * 0.01 + bigClustDir.phi(); //pay attention to this 0.01
334  double pt = seedParamNN.first[i][j][o][4];
335  double normdirR = pt / sin(track_theta);
336 
337  const GlobalVector globSeedDir(
339  LocalVector localSeedDir = globDet->surface().toLocal(globSeedDir);
340  uint64_t seedid = (uint64_t(xx * 200.) << 0) + (uint64_t(yy * 200.) << 16) +
341  (uint64_t(track_eta * 400.) << 32) + (uint64_t(track_phi * 400.) << 48);
342  if (ids.count(seedid) != 0) {
343  continue;
344  }
345  //to disable jetcore iteration comment from here to probability block end--->
346  //seeing iteration skipped, useful to total eff and FakeRate comparison
347  ids.insert(seedid);
348 
349  //Covariance matrix, the hadrcoded variances = NN residuals width
350  //(see https://twiki.cern.ch/twiki/bin/view/CMSPublic/NNJetCoreAtCtD2019)
351  float em[15] = {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,
365  0}; // (see LocalTrajectoryError for details), order as follow:
366  em[0] = 0.15 * 0.15; // q/pt
367  em[2] = 0.5e-5; // dxdz
368  em[5] = 0.5e-5; // dydz
369  em[9] = 2e-5; // x
370  em[14] = 2e-5; // y
371  long int detId = globDet->geographicalId();
372  LocalTrajectoryParameters localParam(localSeedPoint, localSeedDir, TrackCharge(1));
373  result->emplace_back(TrajectorySeed(PTrajectoryStateOnDet(localParam, pt, em, detId, /*surfaceSide*/ 0),
376 
377  GlobalPoint globalSeedPoint = globDet->surface().toGlobal(localSeedPoint);
379  resultTracks->emplace_back(
380  reco::Track(1,
381  1,
382  reco::Track::Point(globalSeedPoint.x(), globalSeedPoint.y(), globalSeedPoint.z()),
383  reco::Track::Vector(globSeedDir.x(), globSeedDir.y(), globSeedDir.z()),
384  1,
385  mm));
386  } //probability
387  }
388  }
389  }
390  } //bigcluster
391  } //jet > pt
392  } //jet
393  iEvent.put(std::move(result));
394  iEvent.put(std::move(resultTracks));
395 }
std::pair< double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> SeedEvaluation(tensorflow::NamedTensorList, std::vector< std::string >)
std::vector< NamedTensor > NamedTensorList
Definition: TensorFlow.h:31
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > geometryToken_
std::vector< GlobalVector > splittedClusterDirections(const reco::Candidate &, const TrackerTopology *const, const PixelClusterParameterEstimator *, const reco::Vertex &, int, const edmNew::DetSetVector< SiPixelCluster > &)
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > pixelCPEToken_
const GeomDet * DetectorSelector(int, const reco::Candidate &, GlobalVector, const reco::Vertex &, const TrackerTopology *const, const edmNew::DetSetVector< SiPixelCluster > &)
T z() const
Definition: PV3DBase.h:61
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClusters_
edm::EDGetTokenT< edm::View< reco::Candidate > > cores_
const Point & position() const
position
Definition: Vertex.h:128
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
Definition: PV3DBase.h:73
int pixelFlipper(const GeomDet *)
edm::Handle< edmNew::DetSetVector< SiPixelCluster > > inputPixelClusters_
LocalPoint toLocal(const GlobalPoint &gp) const
const GeomDet * idToDet(DetId) const override
void fillPixelMatrix(const SiPixelCluster &, int, Point3DBase< float, LocalTag >, const GeomDet *, tensorflow::NamedTensorList)
unsigned int layer(const DetId &id) const
std::pair< bool, Basic3DVector< float > > findIntersection(const GlobalVector &, const reco::Candidate::Point &, const GeomDet *)
std::pair< std::string, Tensor > NamedTensor
Definition: TensorFlow.h:30
static constexpr int Nlayer
int TrackCharge
Definition: TrackCharge.h:4
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
std::vector< std::string > inputTensorName_
std::pair< int, int > local2Pixel(double, double, const GeomDet *)
const GlobalTrackingGeometry * geometry_
LocalPoint pixel2Local(int, int, const GeomDet *)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static constexpr int jetDimY
math::XYZPoint Point
point in the space
Definition: TrackBase.h:80
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
Definition: DetId.h:17
std::vector< std::string > outputTensorName_
uint32_t det_id_type
Definition: DetSet.h:20
edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
unsigned long long uint64_t
Definition: Time.h:13
static constexpr int Nover
static constexpr int jetDimX
fixed size matrix
HLT enums.
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
Definition: output.py:1
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:77
Definition: Phi.h:52
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
def move(src, dest)
Definition: eostools.py:511
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:74
Global3DVector GlobalVector
Definition: GlobalVector.h:10
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.

◆ 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 461 of file DeepCoreSeedGenerator.cc.

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

Referenced by produce().

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

◆ 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 516 of file DeepCoreSeedGenerator.cc.

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

Referenced by produce().

522  {
523  std::vector<GlobalVector> clustDirs;
524  for (const auto& detset_int : clusters) {
525  const GeomDet* det_int = geometry_->idToDet(detset_int.id());
526  int lay = tTopo->layer(det_int->geographicalId());
527  if (lay != layer)
528  continue; //NB: saved bigClusters on all the layers!!
529  auto detUnit = geometry_->idToDetUnit(detset_int.id());
530  for (const auto& aCluster : detset_int) {
531  GlobalPoint clustPos = det_int->surface().toGlobal(pixelCPE->localParametersV(aCluster, (*detUnit))[0].first);
532  GlobalPoint vertexPos(jetVertex.position().x(), jetVertex.position().y(), jetVertex.position().z());
533  GlobalVector clusterDir = clustPos - vertexPos;
534  GlobalVector jetDir(jet.px(), jet.py(), jet.pz());
535  if (Geom::deltaR(jetDir, clusterDir) < deltaR_) {
536  clustDirs.emplace_back(clusterDir);
537  }
538  }
539  }
540  return clustDirs;
541 }
const Point & position() const
position
Definition: Vertex.h:128
const GeomDet * idToDet(DetId) const override
unsigned int layer(const DetId &id) const
const GlobalTrackingGeometry * geometry_
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.

Member Data Documentation

◆ centralMIPCharge_

double DeepCoreSeedGenerator::centralMIPCharge_
private

Definition at line 120 of file DeepCoreSeedGenerator.cc.

◆ chargeFracMin_

double DeepCoreSeedGenerator::chargeFracMin_
private

Definition at line 119 of file DeepCoreSeedGenerator.cc.

◆ cores_

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

Definition at line 112 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ deltaR_

double DeepCoreSeedGenerator::deltaR_
private

Definition at line 118 of file DeepCoreSeedGenerator.cc.

Referenced by splittedClusterDirections().

◆ dth

double DeepCoreSeedGenerator::dth[3] = {0.0, 0.15, 0.3}

Definition at line 99 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ dthl

double DeepCoreSeedGenerator::dthl[3] = {0.1, 0.05, 0}

Definition at line 100 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ geometry_

const GlobalTrackingGeometry* DeepCoreSeedGenerator::geometry_ = nullptr
private

Definition at line 107 of file DeepCoreSeedGenerator.cc.

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

◆ geometryToken_

const edm::ESGetToken<GlobalTrackingGeometry, GlobalTrackingGeometryRecord> DeepCoreSeedGenerator::geometryToken_
private

Definition at line 113 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ inputPixelClusters_

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

Definition at line 111 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ inputTensorName_

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

Definition at line 122 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ jetDimX

constexpr int DeepCoreSeedGenerator::jetDimX = 30
static

Definition at line 92 of file DeepCoreSeedGenerator.cc.

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

◆ jetDimY

constexpr int DeepCoreSeedGenerator::jetDimY = 30
static

Definition at line 93 of file DeepCoreSeedGenerator.cc.

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

◆ jetEta_

double DeepCoreSeedGenerator::jetEta_

Definition at line 89 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ jetPt_

double DeepCoreSeedGenerator::jetPt_

Definition at line 88 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ Nlayer

constexpr int DeepCoreSeedGenerator::Nlayer = 4
static

Definition at line 94 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ Nover

constexpr int DeepCoreSeedGenerator::Nover = 3
static

Definition at line 95 of file DeepCoreSeedGenerator.cc.

Referenced by produce(), and SeedEvaluation().

◆ Npar

constexpr int DeepCoreSeedGenerator::Npar = 5
static

Definition at line 96 of file DeepCoreSeedGenerator.cc.

Referenced by SeedEvaluation().

◆ outputTensorName_

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

Definition at line 123 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pitchX_

double DeepCoreSeedGenerator::pitchX_ = 0.01

Definition at line 90 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pitchY_

double DeepCoreSeedGenerator::pitchY_ = 0.015

Definition at line 91 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pixelClusters_

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

Definition at line 110 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pixelCPEToken_

const edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> DeepCoreSeedGenerator::pixelCPEToken_
private

Definition at line 114 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ probThr_

double DeepCoreSeedGenerator::probThr_
private

Definition at line 124 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ propagatorName_

std::string DeepCoreSeedGenerator::propagatorName_
private

Definition at line 106 of file DeepCoreSeedGenerator.cc.

◆ ptMin_

double DeepCoreSeedGenerator::ptMin_
private

Definition at line 117 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ session_

const tensorflow::Session* DeepCoreSeedGenerator::session_
private

Definition at line 125 of file DeepCoreSeedGenerator.cc.

Referenced by SeedEvaluation().

◆ topoToken_

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

Definition at line 115 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ vertices_

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

Definition at line 109 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ weightfilename_

std::string DeepCoreSeedGenerator::weightfilename_
private

Definition at line 121 of file DeepCoreSeedGenerator.cc.