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  TrackerDetSide side, int idLayer, const std::string & hitProducer, edm::ConsumesCollector& iC)
18  : theHitProducer(iC.consumes<SiPixelRecHitCollection>(hitProducer)), theSide(side), theIdLayer(idLayer)
19 { }
20 
23 }
24 
26 {
28 
30  es.get<TrackerTopologyRcd>().get(httopo);
31  const TrackerTopology& ttopo = *httopo;
32 
34  ev.getByToken( theHitProducer, pixelHits);
36  range2SeedingHits( *pixelHits, result, ttopo.pxbDetIdLayerComparator(theIdLayer));
37  } else {
38  range2SeedingHits( *pixelHits, result, ttopo.pxfDetIdDiskComparator(static_cast<unsigned int>(theSide),theIdLayer));
39  }
40 
41 
42  if (skipClusters){
43  LogDebug("HitExtractorPIX")<<"getting : "<<result.size()<<" pixel hits.";
44  //std::cout<<" skipping"<<std::endl;
45  edm::Handle<SkipClustersCollection> pixelClusterMask;
46  ev.getByToken(theSkipClusters,pixelClusterMask);
47  unsigned int skipped=0;
48  for (unsigned int iH=0;iH!=result.size();++iH){
49  if (result[iH]->isValid()){ // can be NOT valid???
50  auto const & concrete = (SiPixelRecHit const&)(*result[iH]);
51  assert(pixelClusterMask->refProd().id() == concrete.cluster().id());
52  if(pixelClusterMask->mask(concrete.cluster().key())) {
53  //too much debug LogDebug("HitExtractorPIX")<<"skipping a pixel hit on: "<< result[iH]->hit()->geographicalId().rawId()<<" key: "<<find(f->begin(),f->end(),concrete->cluster())->key();
54  skipped++;
55  result[iH].reset();
56  }
57  }
58  }
59  LogDebug("HitExtractorPIX")<<"skipped :"<<skipped<<" pixel clusters";
60  // std::cout << "HitExtractorPIX " <<"skipped :"<<skipped<<" pixel clusters out of " << result.size() << std::endl;
61  if (skipped>0) {
62  auto last = std::remove_if(result.begin(),result.end(),[]( HitPointer const & p) {return p.empty();});
63  result.resize(last-result.begin());
64  }
65  }
66  LogDebug("HitExtractorPIX")<<"giving :"<<result.size()<<" rechits out";
67  // std::cout << "HitExtractorPIX "<<"giving :"<<result.size()<<" rechits out" << std::endl;
68  return result;
69 }
#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:48
std::vector< HitPointer > Hits
Definition: HitExtractor.h:24
std::pair< DetId, SameLayerComparator > pxfDetIdDiskComparator(uint32_t side, uint32_t disk) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
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:44
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
const edm::RefProd< T > & refProd() const
Definition: ContainerMask.h:57
T get() const
Definition: EventSetup.h:71
HitExtractorPIX(TrackerDetSide side, int idLayer, const std::string &hitProducer, edm::ConsumesCollector &iC)
edm::EDGetTokenT< SkipClustersCollection > theSkipClusters
ProductIndex id() const
Definition: ProductID.h:38
std::pair< DetId, SameLayerComparator > pxbDetIdLayerComparator(uint32_t layer) const
Our base class.
Definition: SiPixelRecHit.h:23