CMS 3D CMS Logo

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

#include <Calibration/TkAlCaRecoProducers/plugins/NearbyPixelClustersProducer.cc>

Inheritance diagram for NearbyPixelClustersProducer:
edm::stream::EDProducer<>

Public Member Functions

 NearbyPixelClustersProducer (const edm::ParameterSet &)
 
 ~NearbyPixelClustersProducer () override=default
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 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)
 

Private Member Functions

bool detidIsOnPixel (const DetId &detid)
 
const std::vector< edmNew::DetSet< SiPixelCluster >::const_iterator > findAllNearbyClusters (const SiPixelClusterCollectionNew::const_iterator &clusterSet, const uint32_t rawId, const std::vector< LocalPoint > &vLocalPos)
 
const std::vector< edmNew::DetSet< SiPixelCluster >::const_iterator > findAllNearbyClusters (const SiPixelClusterCollectionNew &clusterSet, const uint32_t rawId, const std::vector< LocalPoint > &vLocalPos)
 
const trajCrossings_t findAllTrajectoriesCrossings (const edm::Handle< TrajTrackAssociationCollection > &trajTrackCollectionHandle)
 
TrajectoryStateOnSurface getTrajectoryStateOnSurface (const TrajectoryMeasurement &measurement)
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const edm::ESGetToken< SiPixelQuality, SiPixelQualityFromDbRcdbadModuleToken_
 
edm::EDPutTokenT< SiPixelClusterCollectionNewclusterPutToken_
 
edm::EDGetTokenT< SiPixelClusterCollectionNewclustersToken_
 
const bool dumpWholeDetId_
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordgeomEsToken_
 
const PixelClusterParameterEstimatorpixelCPE_
 
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecordpixelCPEEsToken_
 
const bool throwBadComponents_
 
const TrackerGeometrytrackerGeometry_
 
edm::EDGetTokenT< TrajTrackAssociationCollectiontrajTrackCollectionToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
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: Class to produce the collection of SiPixelClusters closest, hit by hit, to the trajectory measurements of a given track

Implementation: Implementation of this class is heavily endebted to https://github.com/jkarancs/PhaseIPixelNtuplizer/blob/master/plugins/PhaseIPixelNtuplizer.h

Definition at line 63 of file NearbyPixelClustersProducer.cc.

Constructor & Destructor Documentation

◆ NearbyPixelClustersProducer()

NearbyPixelClustersProducer::NearbyPixelClustersProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 112 of file NearbyPixelClustersProducer.cc.

References dumpWholeDetId_.

113  : throwBadComponents_(iConfig.getParameter<bool>("throwBadComponents")),
114  dumpWholeDetId_(iConfig.getParameter<bool>("dumpWholeDetIds")),
116  pixelCPEEsToken_(esConsumes(edm::ESInputTag("", "PixelCPEGeneric"))),
118  clustersToken_(consumes<SiPixelClusterCollectionNew>(iConfig.getParameter<edm::InputTag>("clusterCollection"))),
120  consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("trajectoryInput"))),
121  clusterPutToken_(produces<SiPixelClusterCollectionNew>()) {
122  if (dumpWholeDetId_) {
123  edm::LogPrint("NearbyPixelClustersProducer") << "WARNING: selected to dump the whole DetId's worth of clusters.\n "
124  "This will have consequences on the event size!"
125  << std::endl;
126  }
127 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDPutTokenT< SiPixelClusterCollectionNew > clusterPutToken_
edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackCollectionToken_
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > pixelCPEEsToken_
Log< level::Warning, true > LogPrint
const edm::ESGetToken< SiPixelQuality, SiPixelQualityFromDbRcd > badModuleToken_
edm::EDGetTokenT< SiPixelClusterCollectionNew > clustersToken_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_

◆ ~NearbyPixelClustersProducer()

NearbyPixelClustersProducer::~NearbyPixelClustersProducer ( )
overridedefault

Member Function Documentation

◆ detidIsOnPixel()

bool NearbyPixelClustersProducer::detidIsOnPixel ( const DetId detid)
private

Definition at line 383 of file NearbyPixelClustersProducer.cc.

References DetId::det(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, DetId::subdetId(), and DetId::Tracker.

Referenced by findAllTrajectoriesCrossings().

385 {
386  if (detid.det() != DetId::Tracker)
387  return false;
389  return true;
391  return true;
392  return false;
393 }
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48

◆ fillDescriptions()

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

Definition at line 421 of file NearbyPixelClustersProducer.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and HLT_2022v12_cff::InputTag.

421  {
423  desc.setComment(
424  "Produces the collection of SiPixelClusters closest, hit by hit, to the trajectory measurements of a given "
425  "track");
426  desc.add<bool>("throwBadComponents", false)
427  ->setComment(
428  "do not consider modules flagged as bad components. Careful, it changes the efficiency denominator!");
429  desc.add<bool>("dumpWholeDetIds", false)
430  ->setComment("put in the event all the pixel cluster on the impacted module, by default only the closest");
431  ;
432  desc.add<edm::InputTag>("clusterCollection", edm::InputTag("siPixelClusters"));
433  desc.add<edm::InputTag>("trajectoryInput", edm::InputTag("myRefitter"));
434  descriptions.addWithDefaultLabel(desc);
435 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ findAllNearbyClusters() [1/2]

const std::vector< edmNew::DetSet< SiPixelCluster >::const_iterator > NearbyPixelClustersProducer::findAllNearbyClusters ( const SiPixelClusterCollectionNew::const_iterator clusterSet,
const uint32_t  rawId,
const std::vector< LocalPoint > &  vLocalPos 
)
private

Definition at line 247 of file NearbyPixelClustersProducer.cc.

References submitPVResolutionJobs::count, dumpWholeDetId_, edmNew::DetSet< T >::end(), PixelClusterParameterEstimator::getParameters(), TrackerGeometry::idToDetUnit(), LogDebug, slimmedTrackExtras_cff::outputClusters, submitPVValidationJobs::params, pixelCPE_, mathSSE::sqrt(), trackerGeometry_, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by produce().

252 {
253  std::vector<edmNew::DetSet<SiPixelCluster>::const_iterator> outputClusters;
254 
255  static constexpr unsigned int k_maxClustersInDet = 1024;
256 
257  // something funny is going on here ...
258  if ((*clusterSet).size() > k_maxClustersInDet) {
259  edm::LogWarning("NearbyPixelClustersProducer")
260  << __func__ << "() number of clusters in det " << rawId /*(*clusterSet).id()*/ << " is " << (*clusterSet).size()
261  << ", which is larger than maximum (1024).\n Something funny with the data is going on!" << std::endl;
262  return outputClusters;
263  }
264 
265  const PixelGeomDetUnit* pixdet = (const PixelGeomDetUnit*)trackerGeometry_->idToDetUnit(rawId);
266  edmNew::DetSet<SiPixelCluster>::const_iterator itCluster = clusterSet->begin();
267 
268  // just copy straight the whole set of clusters in the detid
269  if (dumpWholeDetId_) {
270  for (; itCluster != clusterSet->end(); ++itCluster) {
271  outputClusters.push_back(itCluster);
272  }
273  return outputClusters;
274  }
275 
276  int count = 0;
277  for (const auto& localPos : vLocalPos) {
278  count++;
279  //trajectory crossing local coordinates
280  const auto& traj_lx = localPos.x();
281  const auto& traj_ly = localPos.y();
282 
283  float minD = 10000.;
285 
286  //std::cout<< rawId << " count: " << count << " n. clusters: " << (*clusterSet).size() << std::endl;
287  LogDebug("NearbyPixelClustersProducer")
288  << __func__ << rawId << " count: " << count << " n. clusters: " << (*clusterSet).size() << std::endl;
289 
290  for (; itCluster != clusterSet->end(); ++itCluster) {
291  LocalPoint lp(itCluster->x(), itCluster->y(), 0.);
292  const auto& params = pixelCPE_->getParameters(*itCluster, *pixdet);
293  lp = std::get<0>(params);
294 
295  float D = sqrt((lp.x() - traj_lx) * (lp.x() - traj_lx) + (lp.y() - traj_ly) * (lp.y() - traj_ly));
296  if (D < minD) {
297  closest = itCluster;
298  minD = D;
299  }
300  } // loop on cluster sets
301 
302  if (closest) {
303  outputClusters.push_back(closest);
304  }
305  } // loop on the trajectory crossings
306 
307  return outputClusters;
308 }
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
const PixelClusterParameterEstimator * pixelCPE_
data_type const * const_iterator
Definition: DetSetNew.h:31
virtual ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const =0
T sqrt(T t)
Definition: SSEVec.h:19
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
iterator end()
Definition: DetSetNew.h:56
Log< level::Warning, false > LogWarning
#define LogDebug(id)

◆ findAllNearbyClusters() [2/2]

const std::vector< edmNew::DetSet< SiPixelCluster >::const_iterator > NearbyPixelClustersProducer::findAllNearbyClusters ( const SiPixelClusterCollectionNew clusterSet,
const uint32_t  rawId,
const std::vector< LocalPoint > &  vLocalPos 
)
private

Definition at line 312 of file NearbyPixelClustersProducer.cc.

References ALCARECOSiPixelCalSingleMuonDQM_cff::clusterCollection, submitPVResolutionJobs::count, dumpWholeDetId_, edmNew::DetSet< T >::end(), PixelClusterParameterEstimator::getParameters(), TrackerGeometry::idToDetUnit(), LogDebug, slimmedTrackExtras_cff::outputClusters, submitPVValidationJobs::params, PixelSubdetector::PixelBarrel, pixelCPE_, PixelSubdetector::PixelEndcap, mathSSE::sqrt(), GeomDetEnumerators::subDetId, trackerGeometry_, PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

317 {
318  std::vector<edmNew::DetSet<SiPixelCluster>::const_iterator> outputClusters;
319 
320  int count = 0;
321  for (const auto& localPos : vLocalPos) {
322  count++;
323 
324  //trajectory crossing local coordinates
325  const auto& traj_lx = localPos.x();
326  const auto& traj_ly = localPos.y();
327 
328  float minD = 10000.;
329 
331  for (; itClusterSet != clusterCollection.end(); itClusterSet++) {
332  DetId detId(itClusterSet->id());
333  if (detId.rawId() != rawId)
334  continue;
335 
336  unsigned int subDetId = detId.subdetId();
338  edm::LogError("NearByPixelClusterProducer")
339  << "ERROR: not a pixel cluster!!!" << std::endl; // should not happen
340  continue;
341  }
342 
343  edmNew::DetSet<SiPixelCluster>::const_iterator itCluster = itClusterSet->begin();
344 
345  // just copy straight the whole set of clusters in the detid
346  if (dumpWholeDetId_ && count == 1) {
347  for (; itCluster != itClusterSet->end(); ++itCluster) {
348  outputClusters.push_back(itCluster);
349  }
350  return outputClusters;
351  }
352 
353  //std::cout<< rawId << " count: " << count << " n. clusters: " << (*clusterSet).size() << std::endl;
354  LogDebug("NearbyPixelClustersProducer")
355  << __func__ << rawId << " count: " << count << " n. clusters: " << (*itClusterSet).size() << std::endl;
356 
357  const PixelGeomDetUnit* pixdet = (const PixelGeomDetUnit*)trackerGeometry_->idToDetUnit(rawId);
358 
360 
361  for (; itCluster != itClusterSet->end(); ++itCluster) {
362  LocalPoint lp(itCluster->x(), itCluster->y(), 0.);
363  const auto& params = pixelCPE_->getParameters(*itCluster, *pixdet);
364  lp = std::get<0>(params);
365 
366  float D = sqrt((lp.x() - traj_lx) * (lp.x() - traj_lx) + (lp.y() - traj_ly) * (lp.y() - traj_ly));
367  if (D < minD) {
368  closest = itCluster;
369  minD = D;
370  }
371  } // loop on cluster sets
372 
373  if (closest) {
374  outputClusters.push_back(closest);
375  }
376  } // loop on all clusters
377  } // loop on the trajectory crossings
378 
379  return outputClusters;
380 }
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
const PixelClusterParameterEstimator * pixelCPE_
constexpr unsigned int subDetId[21]
data_type const * const_iterator
Definition: DetSetNew.h:31
Log< level::Error, false > LogError
virtual ReturnType getParameters(const SiPixelCluster &cl, const GeomDetUnit &det) const =0
T sqrt(T t)
Definition: SSEVec.h:19
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
iterator end()
Definition: DetSetNew.h:56
#define LogDebug(id)

◆ findAllTrajectoriesCrossings()

const trajCrossings_t NearbyPixelClustersProducer::findAllTrajectoriesCrossings ( const edm::Handle< TrajTrackAssociationCollection > &  trajTrackCollectionHandle)
private

Definition at line 199 of file NearbyPixelClustersProducer.cc.

References detidIsOnPixel(), spr::find(), getTrajectoryStateOnSurface(), TrajectoryStateOnSurface::isValid(), edm::Ref< C, T, F >::key(), TrajectoryStateOnSurface::localPosition(), DetId::rawId(), and rpcPointValidation_cfi::recHit.

Referenced by produce().

202 {
203  trajCrossings_t crossings;
204 
205  std::vector<uint32_t> treatedIds;
206 
207  for (const auto& pair : *trajTrackCollectionHandle) {
208  const edm::Ref<std::vector<Trajectory>> traj = pair.key;
209 
210  for (const TrajectoryMeasurement& measurement : traj->measurements()) {
211  //Check if the measurement infos can be read
212  if (!measurement.updatedState().isValid())
213  continue;
214 
215  const TransientTrackingRecHit::ConstRecHitPointer& recHit = measurement.recHit();
216 
217  // Only looking for pixel hits
218  DetId recHitDetid = recHit->geographicalId();
219  const auto& rawId = recHitDetid.rawId();
220 
221  if (!this->detidIsOnPixel(recHitDetid))
222  continue;
223 
224  // Skipping hits with undeterminable positions
225  TrajectoryStateOnSurface trajStateOnSurface = this->getTrajectoryStateOnSurface(measurement);
226 
227  if (!(trajStateOnSurface.isValid()))
228  continue;
229 
230  // Position measurements
231  // Looking for valid and missing hits
232  LocalPoint localPosition = trajStateOnSurface.localPosition();
233 
234  if (std::find(treatedIds.begin(), treatedIds.end(), rawId) != treatedIds.end()) {
235  crossings.at(rawId).push_back(localPosition);
236  } else {
237  crossings.insert(std::pair<uint32_t, std::vector<LocalPoint>>(rawId, {localPosition}));
238  treatedIds.push_back(rawId);
239  }
240  } // loop on measurements in trajectory
241  } // loop on trajectories
242 
243  return crossings;
244 }
std::map< uint32_t, std::vector< LocalPoint > > trajCrossings_t
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
key_type key() const
Accessor for product key.
Definition: Ref.h:250
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
TrajectoryStateOnSurface getTrajectoryStateOnSurface(const TrajectoryMeasurement &measurement)
Definition: DetId.h:17
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57

◆ getTrajectoryStateOnSurface()

TrajectoryStateOnSurface NearbyPixelClustersProducer::getTrajectoryStateOnSurface ( const TrajectoryMeasurement measurement)
private

Definition at line 396 of file NearbyPixelClustersProducer.cc.

References TrajectoryMeasurement::backwardPredictedState(), and TrajectoryMeasurement::forwardPredictedState().

Referenced by findAllTrajectoriesCrossings().

399 {
400  const static TrajectoryStateCombiner trajStateCombiner;
401 
402  const auto& forwardPredictedState = measurement.forwardPredictedState();
403  const auto& backwardPredictedState = measurement.backwardPredictedState();
404 
405  if (forwardPredictedState.isValid() && backwardPredictedState.isValid())
406  return trajStateCombiner(forwardPredictedState, backwardPredictedState);
407 
408  else if (backwardPredictedState.isValid())
409  return backwardPredictedState;
410 
411  else if (forwardPredictedState.isValid())
412  return forwardPredictedState;
413 
414  edm::LogError("NearbyPixelClustersProducer") << "Error saving traj. measurement data."
415  << " Trajectory state on surface cannot be determined." << std::endl;
416 
417  return TrajectoryStateOnSurface();
418 }
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
Log< level::Error, false > LogError
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)

◆ produce()

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

Definition at line 134 of file NearbyPixelClustersProducer.cc.

References edmNew::DetSetVector< T >::FastFiller::abort(), badModuleToken_, ALCARECOSiPixelCalSingleMuonDQM_cff::clusterCollection, clusterPutToken_, clustersToken_, edmNew::DetSetVector< T >::FastFiller::empty(), findAllNearbyClusters(), findAllTrajectoriesCrossings(), geomEsToken_, edm::EventSetup::getData(), iEvent, LogDebug, eostools::move(), slimmedTrackExtras_cff::outputClusters, pixelCPE_, pixelCPEEsToken_, edmNew::DetSetVector< T >::FastFiller::push_back(), throwBadComponents_, trackerGeometry_, and trajTrackCollectionToken_.

134  {
135  using namespace edm;
136 
137  auto outputClusters = std::make_unique<SiPixelClusterCollectionNew>();
138 
139  // get the Tracker geometry from event setup
141 
142  // get the Pixel CPE from event setup
144 
145  const auto& SiPixelBadModule_ = &iSetup.getData(badModuleToken_);
146 
147  // get cluster collection
148  const auto& clusterCollectionHandle = iEvent.getHandle(clustersToken_);
149  const SiPixelClusterCollectionNew& clusterCollection = *clusterCollectionHandle;
150 
151  // get Traj-Track Collection
152  const auto& trajTrackCollectionHandle = iEvent.getHandle(trajTrackCollectionToken_);
153  if (!trajTrackCollectionHandle.isValid())
154  return;
155 
156  // find all trajectory crossings in the event
157  const auto& allCrossings = this->findAllTrajectoriesCrossings(trajTrackCollectionHandle);
158 
159  LogDebug("NearbyPixelClustersProducer") << allCrossings.size() << std::endl;
160 
161  // now find all nearby clusters
162  for (const auto& [id, vLocalPos] : allCrossings) {
163  // prepare the filler
165 
166  // retrieve the clusters of the right detId
167  const auto& clustersOnDet = clusterCollection.find(DetId(id));
168 
169  if (clustersOnDet == clusterCollection.end())
170  continue;
171 
172  // if the cluster DetSet is not valid, move on
173  if (!(*clustersOnDet).isValid())
174  continue;
175 
176  // if the module is bad continue
177  if (throwBadComponents_ && SiPixelBadModule_->IsModuleBad(id))
178  continue;
179 
180  // find all the clusters to put into the event
181  const auto& clustersToPut = this->findAllNearbyClusters(clustersOnDet, id, vLocalPos);
182 
183  // find all the clusters to put into the event (different interface)
184  //const auto& clustersToPut = this->findAllNearbyClusters(clusterCollection, id, vLocalPos);
185 
186  for (const auto& cluster : clustersToPut) {
187  spc.push_back(*cluster);
188  }
189 
190  if (spc.empty())
191  spc.abort();
192 
193  } // loop on trajectory crossings
194 
196 }
const std::vector< edmNew::DetSet< SiPixelCluster >::const_iterator > findAllNearbyClusters(const SiPixelClusterCollectionNew::const_iterator &clusterSet, const uint32_t rawId, const std::vector< LocalPoint > &vLocalPos)
edm::EDPutTokenT< SiPixelClusterCollectionNew > clusterPutToken_
const PixelClusterParameterEstimator * pixelCPE_
edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackCollectionToken_
const trajCrossings_t findAllTrajectoriesCrossings(const edm::Handle< TrajTrackAssociationCollection > &trajTrackCollectionHandle)
int iEvent
Definition: GenABIO.cc:224
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const edm::ESGetToken< PixelClusterParameterEstimator, TkPixelCPERecord > pixelCPEEsToken_
const edm::ESGetToken< SiPixelQuality, SiPixelQualityFromDbRcd > badModuleToken_
Definition: DetId.h:17
HLT enums.
edm::EDGetTokenT< SiPixelClusterCollectionNew > clustersToken_
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomEsToken_

Member Data Documentation

◆ badModuleToken_

const edm::ESGetToken<SiPixelQuality, SiPixelQualityFromDbRcd> NearbyPixelClustersProducer::badModuleToken_
private

Definition at line 95 of file NearbyPixelClustersProducer.cc.

Referenced by produce().

◆ clusterPutToken_

edm::EDPutTokenT<SiPixelClusterCollectionNew> NearbyPixelClustersProducer::clusterPutToken_
private

Definition at line 102 of file NearbyPixelClustersProducer.cc.

Referenced by produce().

◆ clustersToken_

edm::EDGetTokenT<SiPixelClusterCollectionNew> NearbyPixelClustersProducer::clustersToken_
private

Definition at line 98 of file NearbyPixelClustersProducer.cc.

Referenced by produce().

◆ dumpWholeDetId_

const bool NearbyPixelClustersProducer::dumpWholeDetId_
private

◆ geomEsToken_

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> NearbyPixelClustersProducer::geomEsToken_
private

Definition at line 93 of file NearbyPixelClustersProducer.cc.

Referenced by produce().

◆ pixelCPE_

const PixelClusterParameterEstimator* NearbyPixelClustersProducer::pixelCPE_
private

Definition at line 106 of file NearbyPixelClustersProducer.cc.

Referenced by findAllNearbyClusters(), and produce().

◆ pixelCPEEsToken_

const edm::ESGetToken<PixelClusterParameterEstimator, TkPixelCPERecord> NearbyPixelClustersProducer::pixelCPEEsToken_
private

Definition at line 94 of file NearbyPixelClustersProducer.cc.

Referenced by produce().

◆ throwBadComponents_

const bool NearbyPixelClustersProducer::throwBadComponents_
private

Definition at line 89 of file NearbyPixelClustersProducer.cc.

Referenced by produce().

◆ trackerGeometry_

const TrackerGeometry* NearbyPixelClustersProducer::trackerGeometry_
private

Definition at line 105 of file NearbyPixelClustersProducer.cc.

Referenced by findAllNearbyClusters(), and produce().

◆ trajTrackCollectionToken_

edm::EDGetTokenT<TrajTrackAssociationCollection> NearbyPixelClustersProducer::trajTrackCollectionToken_
private

Definition at line 99 of file NearbyPixelClustersProducer.cc.

Referenced by produce().