CMS 3D CMS Logo

ClusterChargeMasker.cc
Go to the documentation of this file.
5 
8 
9 
14 
16 
18 
21 
22 
23 namespace {
24 
25  class ClusterChargeMasker : public edm::stream::EDProducer<> {
26  public:
27  ClusterChargeMasker(const edm::ParameterSet& iConfig);
28  ~ClusterChargeMasker(){}
29  void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override ;
30  private:
31 
32 
35 
36 
37  bool mergeOld_;
38  float minGoodStripCharge_;
39 
42 
43  edm::EDGetTokenT<PixelMaskContainer> oldPxlMaskToken_;
44  edm::EDGetTokenT<StripMaskContainer> oldStrMaskToken_;
45 
46 
47 
48  };
49 
50 
51  ClusterChargeMasker::ClusterChargeMasker(const edm::ParameterSet& iConfig) :
52  mergeOld_(iConfig.exists("oldClusterRemovalInfo")),
53  minGoodStripCharge_(clusterChargeCut(iConfig))
54  {
55  produces<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > >();
56  produces<edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > >();
57 
58 
59  pixelClusters_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(iConfig.getParameter<edm::InputTag>("pixelClusters"));
60  stripClusters_ = consumes<edmNew::DetSetVector<SiStripCluster> >(iConfig.getParameter<edm::InputTag>("stripClusters"));
61 
62  if (mergeOld_) {
63  oldPxlMaskToken_ = consumes<PixelMaskContainer>(iConfig.getParameter<edm::InputTag>("oldClusterRemovalInfo"));
64  oldStrMaskToken_ = consumes<StripMaskContainer>(iConfig.getParameter<edm::InputTag>("oldClusterRemovalInfo"));
65  }
66 
67  }
68 
69 
70  void
71  ClusterChargeMasker::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
72  {
73 
74 
76  iEvent.getByToken(pixelClusters_, pixelClusters);
78  iEvent.getByToken(stripClusters_, stripClusters);
79 
80 
81  std::vector<bool> collectedStrips;
82  std::vector<bool> collectedPixels;
83 
84  if(mergeOld_) {
87  iEvent.getByToken(oldPxlMaskToken_ ,oldPxlMask);
88  iEvent.getByToken(oldStrMaskToken_ ,oldStrMask);
89  LogDebug("ClusterChargeMasker")<<"to merge in, "<<oldStrMask->size()<<" strp and "<<oldPxlMask->size()<<" pxl";
90  oldStrMask->copyMaskTo(collectedStrips);
91  oldPxlMask->copyMaskTo(collectedPixels);
92  assert(stripClusters->dataSize()>=collectedStrips.size());
93  collectedStrips.resize(stripClusters->dataSize(), false);
94  }else {
95  collectedStrips.resize(stripClusters->dataSize(), false);
96  collectedPixels.resize(pixelClusters->dataSize(), false);
97  }
98 
99  auto const & clusters = stripClusters->data();
100  for (auto const & item : stripClusters->ids()) {
101 
102  if (!item.isValid()) continue; // not umpacked (hlt only)
103 
104  auto detid = item.id;
105 
106  for (int i = item.offset; i<item.offset+int(item.size); ++i) {
107  auto clusCharge = siStripClusterTools::chargePerCM(detid,clusters[i]);
108  if(clusCharge < minGoodStripCharge_) collectedStrips[i] = true;
109  }
110 
111  }
112 
113  LogDebug("ClusterChargeMasker")<<"total strip to skip: "<<std::count(collectedStrips.begin(),collectedStrips.end(),true);
114  // std::cout << "ClusterChargeMasker " <<"total strip to skip: "<<std::count(collectedStrips.begin(),collectedStrips.end(),true)
115  // << " for CCC " << minGoodStripCharge_ <<std::endl;
116  iEvent.put(std::make_unique<StripMaskContainer>(edm::RefProd<edmNew::DetSetVector<SiStripCluster> >(stripClusters),collectedStrips));
117 
118  LogDebug("ClusterChargeMasker")<<"total pxl to skip: "<<std::count(collectedPixels.begin(),collectedPixels.end(),true);
119  iEvent.put(std::make_unique<PixelMaskContainer>(edm::RefProd<edmNew::DetSetVector<SiPixelCluster> >(pixelClusters),collectedPixels));
120 
121 
122 
123  }
124 
125 
126 }
127 
128 
131 DEFINE_FWK_MODULE(ClusterChargeMasker);
132 
#define LogDebug(id)
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
ProductID id() const
Definition: HandleBase.cc:15
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
virtual void produce(Event &, EventSetup const &)=0
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
float chargePerCM(DetId detid, Iter a, Iter b)
float clusterChargeCut(const edm::ParameterSet &conf, const char *name="clusterChargeCut")
int iEvent
Definition: GenABIO.cc:230