CMS 3D CMS Logo

HitExtractorPIX.cc
Go to the documentation of this file.
1 #include "HitExtractorPIX.h"
5 
8 
11 
12 #include <iostream>
13 using namespace ctfseeding;
14 using namespace std;
15 
17  int idLayer,
18  const std::string& hitProducer,
20  : theHitProducer(iC.consumes<SiPixelRecHitCollection>(hitProducer)), theSide(side), theIdLayer(idLayer) {}
21 
24 }
25 
27  const edm::Event& ev,
28  const edm::EventSetup& es) const {
30 
32  es.get<TrackerTopologyRcd>().get(httopo);
33  const TrackerTopology& ttopo = *httopo;
34 
36  ev.getByToken(theHitProducer, pixelHits);
38  range2SeedingHits(*pixelHits, result, ttopo.pxbDetIdLayerComparator(theIdLayer));
39  } else {
40  range2SeedingHits(*pixelHits, result, ttopo.pxfDetIdDiskComparator(static_cast<unsigned int>(theSide), theIdLayer));
41  }
42 
43  if (skipClusters) {
44  LogDebug("HitExtractorPIX") << "getting : " << result.size() << " pixel hits.";
45  //std::cout<<" skipping"<<std::endl;
46  edm::Handle<SkipClustersCollection> pixelClusterMask;
47  ev.getByToken(theSkipClusters, pixelClusterMask);
48  unsigned int skipped = 0;
49  for (unsigned int iH = 0; iH != result.size(); ++iH) {
50  if (result[iH]->isValid()) { // can be NOT valid???
51  auto const& concrete = (SiPixelRecHit const&)(*result[iH]);
52  assert(pixelClusterMask->refProd().id() == concrete.cluster().id());
53  if (pixelClusterMask->mask(concrete.cluster().key())) {
54  //too much debug LogDebug("HitExtractorPIX")<<"skipping a pixel hit on: "<< result[iH]->hit()->geographicalId().rawId()<<" key: "<<find(f->begin(),f->end(),concrete->cluster())->key();
55  skipped++;
56  result[iH].reset();
57  }
58  }
59  }
60  LogDebug("HitExtractorPIX") << "skipped :" << skipped << " pixel clusters";
61  // std::cout << "HitExtractorPIX " <<"skipped :"<<skipped<<" pixel clusters out of " << result.size() << std::endl;
62  if (skipped > 0) {
63  auto last = std::remove_if(result.begin(), result.end(), [](HitPointer const& p) { return p.empty(); });
64  result.resize(last - result.begin());
65  }
66  }
67  LogDebug("HitExtractorPIX") << "giving :" << result.size() << " rechits out";
68  // std::cout << "HitExtractorPIX "<<"giving :"<<result.size()<<" rechits out" << std::endl;
69  return result;
70 }
#define LogDebug(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void range2SeedingHits(DSTV const &dstv, HitExtractor::Hits &v, std::pair< A, B > const &sel)
Definition: HitExtractor.h:54
std::vector< HitPointer > Hits
Definition: HitExtractor.h:28
std::pair< DetId, SameLayerComparator > pxfDetIdDiskComparator(uint32_t side, uint32_t disk) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
TrackerDetSide
Definition: TrackerDetSide.h:4
void useSkipClusters_(const edm::InputTag &m, edm::ConsumesCollector &iC) override
bool ev
HitExtractor::Hits hits(const TkTransientTrackingRecHitBuilder &ttrhBuilder, const edm::Event &, const edm::EventSetup &) const override
edm::EDGetTokenT< SiPixelRecHitCollection > theHitProducer
bool mask(unsigned int iIndex) const
Definition: ContainerMask.h:43
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
const edm::RefProd< T > & refProd() const
Definition: ContainerMask.h:55
T get() const
Definition: EventSetup.h:73
HitExtractorPIX(TrackerDetSide side, int idLayer, const std::string &hitProducer, edm::ConsumesCollector &iC)
edm::EDGetTokenT< SkipClustersCollection > theSkipClusters
ProductIndex id() const
Definition: ProductID.h:35
std::pair< DetId, SameLayerComparator > pxbDetIdLayerComparator(uint32_t layer) const
Our base class.
Definition: SiPixelRecHit.h:23