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

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_
 
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 161 of file DeepCoreSeedGenerator.cc.

162  : vertices_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
164  consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("pixelClusters"))),
165  cores_(consumes<edm::View<reco::Candidate>>(iConfig.getParameter<edm::InputTag>("cores"))),
168  topoToken_(esConsumes<TrackerTopology, TrackerTopologyRcd>()),
169  ptMin_(iConfig.getParameter<double>("ptMin")),
170  deltaR_(iConfig.getParameter<double>("deltaR")),
171  chargeFracMin_(iConfig.getParameter<double>("chargeFractionMin")),
172  centralMIPCharge_(iConfig.getParameter<double>("centralMIPCharge")),
173  weightfilename_(iConfig.getParameter<edm::FileInPath>("weightFile").fullPath()),
174  inputTensorName_(iConfig.getParameter<std::vector<std::string>>("inputTensorName")),
175  outputTensorName_(iConfig.getParameter<std::vector<std::string>>("outputTensorName")),
176  probThr_(iConfig.getParameter<double>("probThr")),
178 
179 {
180  produces<TrajectorySeedCollection>();
181  produces<reco::TrackCollection>();
182 }
Session * createSession(SessionOptions &sessionOptions)
Definition: TensorFlow.cc:85
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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 edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
tensorflow::Session * session_

◆ ~DeepCoreSeedGenerator()

DeepCoreSeedGenerator::~DeepCoreSeedGenerator ( )
override

Definition at line 184 of file DeepCoreSeedGenerator.cc.

184 {}

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

483  {
484  double minDist = 0.0;
485  GeomDet* output = (GeomDet*)nullptr;
486  for (const auto& detset : clusters) {
487  auto aClusterID = detset.id();
488  if (DetId(aClusterID).subdetId() != 1)
489  continue;
490  const GeomDet* det = geometry_->idToDet(aClusterID);
491  int lay = tTopo->layer(det->geographicalId());
492  if (lay != llay)
493  continue;
494  std::pair<bool, Basic3DVector<float>> interPair =
495  findIntersection(jetDir, (reco::Candidate::Point)jetVertex.position(), det);
496  if (interPair.first == false)
497  continue;
498  Basic3DVector<float> inter = interPair.second;
499  auto localInter = det->specificSurface().toLocal((GlobalPoint)inter);
500  if ((minDist == 0.0 || std::abs(localInter.x()) < minDist) && std::abs(localInter.y()) < 3.35) {
501  minDist = std::abs(localInter.x());
502  output = (GeomDet*)det;
503  }
504  } //detset
505  return output;
506 }
const Point & position() const
position
Definition: Vertex.h:127
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
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 546 of file DeepCoreSeedGenerator.cc.

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

546  {
548  desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
549  desc.add<edm::InputTag>("pixelClusters", edm::InputTag("siPixelClustersPreSplitting"));
550  desc.add<edm::InputTag>("cores", edm::InputTag("jetsForCoreTracking"));
551  desc.add<double>("ptMin", 100);
552  desc.add<double>("deltaR", 0.25); // the current training makes use of 0.1
553  desc.add<double>("chargeFractionMin", 18000.0);
554  desc.add<double>("centralMIPCharge", 2);
555  desc.add<std::string>("pixelCPE", "PixelCPEGeneric");
556  desc.add<edm::FileInPath>(
557  "weightFile",
558  edm::FileInPath("RecoTracker/TkSeedGenerator/data/DeepCore/DeepCoreSeedGenerator_TrainedModel_barrel_2017.pb"));
559  desc.add<std::vector<std::string>>("inputTensorName", {"input_1", "input_2", "input_3"});
560  desc.add<std::vector<std::string>>("outputTensorName", {"output_node0", "output_node1"});
561  desc.add<double>("probThr", 0.85);
562  descriptions.add("deepCoreSeedGenerator", desc);
563 }
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 427 of file DeepCoreSeedGenerator.cc.

References funct::abs(), SiPixelCluster::Pixel::adc, f, HLT_2022v12_cff::flip, dqmMemoryStats::float, mps_fire::i, jetDimX, jetDimY, phase1PixelTopology::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().

431  {
432  int flip = pixelFlipper(det); // 1=not flip, -1=flip
433 
434  for (int i = 0; i < cluster.size(); i++) {
435  SiPixelCluster::Pixel pix = cluster.pixel(i);
436  std::pair<int, int> pixInter = local2Pixel(inter.x(), inter.y(), det);
437  int nx = pix.x - pixInter.first;
438  int ny = pix.y - pixInter.second;
439  nx = flip * nx;
440 
441  if (abs(nx) < jetDimX / 2 && abs(ny) < jetDimY / 2) {
442  nx = nx + jetDimX / 2;
443  ny = ny + jetDimY / 2;
444  //14000 = normalization of ACD counts used in the NN
445  input_tensors[2].second.tensor<float, 4>()(0, nx, ny, layer - 1) += (pix.adc) / (14000.f);
446  }
447  }
448 }
Pixel pixel(int i) const
int pixelFlipper(const GeomDet *)
constexpr std::array< uint8_t, layerIndexSize > layer
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 389 of file DeepCoreSeedGenerator.cc.

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

Referenced by DetectorSelector(), and produce().

391  {
393  Basic3DVector<float>(dir.x(), dir.y(), dir.z()));
394 
395  std::pair<bool, Basic3DVector<float>> pos = vertexPlane.position(det->specificSurface());
396 
397  return pos;
398 }
const Plane & specificSurface() const
Same as surface(), kept for backward compatibility.
Definition: GeomDet.h:40

◆ globalEndJob()

void DeepCoreSeedGenerator::globalEndJob ( DeepCoreCache cache)
static

Definition at line 543 of file DeepCoreSeedGenerator.cc.

References utilities::cache().

543 { delete cache->graph_def; }
def cache(function)
Definition: utilities.py:3

◆ initializeGlobalCache()

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

Definition at line 535 of file DeepCoreSeedGenerator.cc.

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

535  {
536  // this method is supposed to create, initialize and return a DeepCoreCache instance
537  std::unique_ptr<DeepCoreCache> cache = std::make_unique<DeepCoreCache>();
538  std::string graphPath = iConfig.getParameter<edm::FileInPath>("weightFile").fullPath();
540  return cache;
541 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
GraphDef * loadGraphDef(const std::string &pbFile)
Definition: TensorFlow.cc:68
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 400 of file DeepCoreSeedGenerator.cc.

References MillePedeFileConverter_cfg::out.

Referenced by fillPixelMatrix(), and produce().

400  {
401  LocalPoint locXY(locX, locY);
402  float pixX = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().pixel(locXY).first;
403  float pixY = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().pixel(locXY).second;
404  std::pair<int, int> out(pixX, pixY);
405  return out;
406 }

◆ pixel2Local()

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

Definition at line 408 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

408  {
409  float locX = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().localX(pixX);
410  float locY = (dynamic_cast<const PixelGeomDetUnit*>(det))->specificTopology().localY(pixY);
411  LocalPoint locXY(locX, locY);
412  return locXY;
413 }

◆ pixelFlipper()

int DeepCoreSeedGenerator::pixelFlipper ( const GeomDet det)
private

Definition at line 415 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().

415  {
416  int out = 1;
417  LocalVector locZdir(0, 0, 1);
418  GlobalVector globZdir = det->specificSurface().toGlobal(locZdir);
419  const GlobalPoint& globDetCenter = det->position();
420  float direction =
421  globZdir.x() * globDetCenter.x() + globZdir.y() * globDetCenter.y() + globZdir.z() * globDetCenter.z();
422  if (direction < 0)
423  out = -1;
424  return out;
425 }
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 186 of file DeepCoreSeedGenerator.cc.

References funct::abs(), alongMomentum, HLT_FULL_cff::cores, cores_, DetectorSelector(), MillePedeFileConverter_cfg::e, PV3DBase< T, PVType, FrameType >::eta(), JetChargeProducer_cfi::exp, fillPixelMatrix(), findIntersection(), first, HLT_2022v12_cff::flip, dqmMemoryStats::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, convertSQLitetoXML_cfg::output, HLT_2022v12_cff::output_names, outputTensorName_, PV3DBase< T, PVType, FrameType >::phi(), pitchX_, pitchY_, pixel2Local(), PixelSubdetector::PixelBarrel, pixelClusters_, HLT_2022v12_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_, 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().

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

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

Referenced by produce().

454  {
455  std::vector<tensorflow::Tensor> outputs;
456  tensorflow::run(session_, input_tensors, output_names, &outputs);
457  auto matrix_output_par = outputs.at(0).tensor<float, 5>();
458  auto matrix_output_prob = outputs.at(1).tensor<float, 5>();
459 
460  std::pair<double[jetDimX][jetDimY][Nover][Npar], double[jetDimX][jetDimY][Nover]> output_combined;
461 
462  for (int x = 0; x < jetDimX; x++) {
463  for (int y = 0; y < jetDimY; y++) {
464  for (int trk = 0; trk < Nover; trk++) {
465  output_combined.second[x][y][trk] =
466  matrix_output_prob(0, x, y, trk, 0); //outputs.at(1).matrix<double>()(0,x,y,trk);
467 
468  for (int p = 0; p < Npar; p++) {
469  output_combined.first[x][y][trk][p] =
470  matrix_output_par(0, x, y, trk, p); //outputs.at(0).matrix<double>()(0,x,y,trk,p);
471  }
472  }
473  }
474  }
475  return output_combined;
476 }
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:213
static constexpr int jetDimY
static constexpr int Nover
static constexpr int jetDimX
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 508 of file DeepCoreSeedGenerator.cc.

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

Referenced by produce().

514  {
515  std::vector<GlobalVector> clustDirs;
516  for (const auto& detset_int : clusters) {
517  const GeomDet* det_int = geometry_->idToDet(detset_int.id());
518  int lay = tTopo->layer(det_int->geographicalId());
519  if (lay != layer)
520  continue; //NB: saved bigClusters on all the layers!!
521  auto detUnit = geometry_->idToDetUnit(detset_int.id());
522  for (const auto& aCluster : detset_int) {
523  GlobalPoint clustPos = det_int->surface().toGlobal(pixelCPE->localParametersV(aCluster, (*detUnit))[0].first);
524  GlobalPoint vertexPos(jetVertex.position().x(), jetVertex.position().y(), jetVertex.position().z());
525  GlobalVector clusterDir = clustPos - vertexPos;
526  GlobalVector jetDir(jet.px(), jet.py(), jet.pz());
527  if (Geom::deltaR(jetDir, clusterDir) < deltaR_) {
528  clustDirs.emplace_back(clusterDir);
529  }
530  }
531  }
532  return clustDirs;
533 }
const Point & position() const
position
Definition: Vertex.h:127
const GeomDet * idToDet(DetId) const override
unsigned int layer(const DetId &id) const
constexpr std::array< uint8_t, layerIndexSize > layer
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 119 of file DeepCoreSeedGenerator.cc.

◆ chargeFracMin_

double DeepCoreSeedGenerator::chargeFracMin_
private

Definition at line 118 of file DeepCoreSeedGenerator.cc.

◆ cores_

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

Definition at line 111 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ deltaR_

double DeepCoreSeedGenerator::deltaR_
private

Definition at line 117 of file DeepCoreSeedGenerator.cc.

Referenced by splittedClusterDirections().

◆ geometry_

const GlobalTrackingGeometry* DeepCoreSeedGenerator::geometry_ = nullptr
private

Definition at line 106 of file DeepCoreSeedGenerator.cc.

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

◆ geometryToken_

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

Definition at line 112 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ inputPixelClusters_

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

Definition at line 110 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ inputTensorName_

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

Definition at line 121 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ jetDimX

constexpr int DeepCoreSeedGenerator::jetDimX = 30
static

Definition at line 95 of file DeepCoreSeedGenerator.cc.

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

◆ jetDimY

constexpr int DeepCoreSeedGenerator::jetDimY = 30
static

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

◆ Nlayer

constexpr int DeepCoreSeedGenerator::Nlayer = 4
static

Definition at line 97 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ Nover

constexpr int DeepCoreSeedGenerator::Nover = 3
static

Definition at line 98 of file DeepCoreSeedGenerator.cc.

Referenced by produce(), and SeedEvaluation().

◆ Npar

constexpr int DeepCoreSeedGenerator::Npar = 5
static

Definition at line 99 of file DeepCoreSeedGenerator.cc.

Referenced by SeedEvaluation().

◆ outputTensorName_

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

Definition at line 122 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 109 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ pixelCPEToken_

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

Definition at line 113 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ probThr_

double DeepCoreSeedGenerator::probThr_
private

Definition at line 123 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ propagatorName_

std::string DeepCoreSeedGenerator::propagatorName_
private

Definition at line 105 of file DeepCoreSeedGenerator.cc.

◆ ptMin_

double DeepCoreSeedGenerator::ptMin_
private

Definition at line 116 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ session_

tensorflow::Session* DeepCoreSeedGenerator::session_
private

Definition at line 124 of file DeepCoreSeedGenerator.cc.

Referenced by SeedEvaluation().

◆ topoToken_

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

Definition at line 114 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ vertices_

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

Definition at line 108 of file DeepCoreSeedGenerator.cc.

Referenced by produce().

◆ weightfilename_

std::string DeepCoreSeedGenerator::weightfilename_
private

Definition at line 120 of file DeepCoreSeedGenerator.cc.