CMS 3D CMS Logo

ClusterTPAssociationProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <vector>
3 #include <utility>
4 
13 
23 
31 
33 public:
34  typedef std::vector<OmniClusterRef> OmniClusterCollection;
35 
38 
39  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
40 
41 private:
42  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
43 
44  template <typename T>
45  std::vector<std::pair<uint32_t, EncodedEventId> > getSimTrackId(const edm::Handle<edm::DetSetVector<T> >& simLinks,
46  const DetId& detId,
47  uint32_t channel) const;
48 
56 };
57 
60  consumes<edm::DetSetVector<PixelDigiSimLink> >(cfg.getParameter<edm::InputTag>("pixelSimLinkSrc"))),
62  consumes<edm::DetSetVector<StripDigiSimLink> >(cfg.getParameter<edm::InputTag>("stripSimLinkSrc"))),
64  consumes<edm::DetSetVector<PixelDigiSimLink> >(cfg.getParameter<edm::InputTag>("phase2OTSimLinkSrc"))),
66  consumes<edmNew::DetSetVector<SiPixelCluster> >(cfg.getParameter<edm::InputTag>("pixelClusterSrc"))),
68  consumes<edmNew::DetSetVector<SiStripCluster> >(cfg.getParameter<edm::InputTag>("stripClusterSrc"))),
70  cfg.getParameter<edm::InputTag>("phase2OTClusterSrc"))),
72  consumes<TrackingParticleCollection>(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))) {
73  produces<ClusterTPAssociation>();
74 }
75 
77 
80  desc.add<edm::InputTag>("simTrackSrc", edm::InputTag("g4SimHits"));
81  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis"));
82  desc.add<edm::InputTag>("stripSimLinkSrc", edm::InputTag("simSiStripDigis"));
83  desc.add<edm::InputTag>("phase2OTSimLinkSrc", edm::InputTag("simSiPixelDigis", "Tracker"));
84  desc.add<edm::InputTag>("pixelClusterSrc", edm::InputTag("siPixelClusters"));
85  desc.add<edm::InputTag>("stripClusterSrc", edm::InputTag("siStripClusters"));
86  desc.add<edm::InputTag>("phase2OTClusterSrc", edm::InputTag("siPhase2Clusters"));
87  desc.add<edm::InputTag>("trackingParticleSrc", edm::InputTag("mix", "MergedTrackTruth"));
88  descriptions.add("tpClusterProducerDefault", desc);
89 }
90 
92  // Pixel DigiSimLink
94  // iEvent.getByLabel(_pixelSimLinkSrc, sipixelSimLinks);
95  iEvent.getByToken(sipixelSimLinksToken_, sipixelSimLinks);
96 
97  // SiStrip DigiSimLink
99  iEvent.getByToken(sistripSimLinksToken_, sistripSimLinks);
100 
101  // Phase2 OT DigiSimLink
103  iEvent.getByToken(siphase2OTSimLinksToken_, siphase2OTSimLinks);
104 
105  // Pixel Cluster
107  bool foundPixelClusters = iEvent.getByToken(pixelClustersToken_, pixelClusters);
108 
109  // Strip Cluster
111  bool foundStripClusters = iEvent.getByToken(stripClustersToken_, stripClusters);
112 
113  // Phase2 Cluster
115  bool foundPhase2OTClusters = iEvent.getByToken(phase2OTClustersToken_, phase2OTClusters);
116 
117  // TrackingParticle
119  iEvent.getByToken(trackingParticleToken_, TPCollectionH);
120 
121  auto clusterTPList = std::make_unique<ClusterTPAssociation>(TPCollectionH);
122 
123  // prepare temporary map between SimTrackId and TrackingParticle index
124  std::map<std::pair<size_t, EncodedEventId>, TrackingParticleRef> mapping;
125  for (TrackingParticleCollection::size_type itp = 0; itp < TPCollectionH.product()->size(); ++itp) {
126  TrackingParticleRef trackingParticle(TPCollectionH, itp);
127 
128  // SimTracks inside TrackingParticle
129  EncodedEventId eid(trackingParticle->eventId());
130  //size_t index = 0;
131  for (std::vector<SimTrack>::const_iterator itrk = trackingParticle->g4Track_begin();
132  itrk != trackingParticle->g4Track_end();
133  ++itrk) {
134  std::pair<uint32_t, EncodedEventId> trkid(itrk->trackId(), eid);
135  //std::cout << "creating map for id: " << trkid.first << " with tp: " << trackingParticle.key() << std::endl;
136  mapping.insert(std::make_pair(trkid, trackingParticle));
137  }
138  }
139 
140  if (foundPixelClusters) {
141  // Pixel Clusters
142  clusterTPList->addKeyID(pixelClusters.id());
143  for (edmNew::DetSetVector<SiPixelCluster>::const_iterator iter = pixelClusters->begin();
144  iter != pixelClusters->end();
145  ++iter) {
146  uint32_t detid = iter->id();
147  DetId detId(detid);
148  edmNew::DetSet<SiPixelCluster> link_pixel = (*iter);
149  for (edmNew::DetSet<SiPixelCluster>::const_iterator di = link_pixel.begin(); di != link_pixel.end(); ++di) {
150  const SiPixelCluster& cluster = (*di);
152 
153  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
154  for (int irow = cluster.minPixelRow(); irow <= cluster.maxPixelRow(); ++irow) {
155  for (int icol = cluster.minPixelCol(); icol <= cluster.maxPixelCol(); ++icol) {
156  uint32_t channel = PixelChannelIdentifier::pixelToChannel(irow, icol);
157  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(
158  getSimTrackId<PixelDigiSimLink>(sipixelSimLinks, detId, channel));
159  if (trkid.empty())
160  continue;
161  simTkIds.insert(trkid.begin(), trkid.end());
162  }
163  }
164  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
165  iset != simTkIds.end();
166  iset++) {
167  auto ipos = mapping.find(*iset);
168  if (ipos != mapping.end()) {
169  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
170  clusterTPList->emplace_back(OmniClusterRef(c_ref), ipos->second);
171  }
172  }
173  }
174  }
175  }
176 
177  if (foundStripClusters) {
178  // Strip Clusters
179  clusterTPList->addKeyID(stripClusters.id());
180  for (edmNew::DetSetVector<SiStripCluster>::const_iterator iter = stripClusters->begin(false),
181  eter = stripClusters->end(false);
182  iter != eter;
183  ++iter) {
184  if (!(*iter).isValid())
185  continue;
186  uint32_t detid = iter->id();
187  DetId detId(detid);
188  edmNew::DetSet<SiStripCluster> link_strip = (*iter);
189  for (edmNew::DetSet<SiStripCluster>::const_iterator di = link_strip.begin(); di != link_strip.end(); di++) {
190  const SiStripCluster& cluster = (*di);
192 
193  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
194  int first = cluster.firstStrip();
195  int last = first + cluster.amplitudes().size();
196 
197  for (int istr = first; istr < last; ++istr) {
198  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(
199  getSimTrackId<StripDigiSimLink>(sistripSimLinks, detId, istr));
200  if (trkid.empty())
201  continue;
202  simTkIds.insert(trkid.begin(), trkid.end());
203  }
204  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
205  iset != simTkIds.end();
206  iset++) {
207  auto ipos = mapping.find(*iset);
208  if (ipos != mapping.end()) {
209  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
210  clusterTPList->emplace_back(OmniClusterRef(c_ref), ipos->second);
211  }
212  }
213  }
214  }
215  }
216 
217  if (foundPhase2OTClusters) {
218  // Phase2 Clusters
219  clusterTPList->addKeyID(phase2OTClusters.id());
220  if (phase2OTClusters.isValid()) {
221  for (edmNew::DetSetVector<Phase2TrackerCluster1D>::const_iterator iter = phase2OTClusters->begin(false),
222  eter = phase2OTClusters->end(false);
223  iter != eter;
224  ++iter) {
225  if (!(*iter).isValid())
226  continue;
227  uint32_t detid = iter->id();
228  DetId detId(detid);
229  edmNew::DetSet<Phase2TrackerCluster1D> link_phase2 = (*iter);
230  for (edmNew::DetSet<Phase2TrackerCluster1D>::const_iterator di = link_phase2.begin(); di != link_phase2.end();
231  di++) {
232  const Phase2TrackerCluster1D& cluster = (*di);
234  edmNew::makeRefTo(phase2OTClusters, di);
235 
236  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
237 
238  for (unsigned int istr(0); istr < cluster.size(); ++istr) {
239  uint32_t channel = Phase2TrackerDigi::pixelToChannel(cluster.firstRow() + istr, cluster.column());
240  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(
241  getSimTrackId<PixelDigiSimLink>(siphase2OTSimLinks, detId, channel));
242  if (trkid.empty())
243  continue;
244  simTkIds.insert(trkid.begin(), trkid.end());
245  }
246 
247  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
248  iset != simTkIds.end();
249  iset++) {
250  auto ipos = mapping.find(*iset);
251  if (ipos != mapping.end()) {
252  clusterTPList->emplace_back(OmniClusterRef(c_ref), ipos->second);
253  }
254  }
255  }
256  }
257  }
258  }
259  clusterTPList->sortAndUnique();
260  iEvent.put(std::move(clusterTPList));
261 }
262 
263 template <typename T>
264 std::vector<std::pair<uint32_t, EncodedEventId> >
265 //std::pair<uint32_t, EncodedEventId>
267  const DetId& detId,
268  uint32_t channel) const {
269  //std::pair<uint32_t, EncodedEventId> simTrkId;
270  std::vector<std::pair<uint32_t, EncodedEventId> > simTrkId;
271  auto isearch = simLinks->find(detId);
272  if (isearch != simLinks->end()) {
273  // Loop over DigiSimLink in this det unit
274  edm::DetSet<T> link_detset = (*isearch);
275  for (typename edm::DetSet<T>::const_iterator it = link_detset.data.begin(); it != link_detset.data.end(); ++it) {
276  if (channel == it->channel()) {
277  simTrkId.push_back(std::make_pair(it->SimTrackId(), it->eventId()));
278  }
279  }
280  }
281  return simTrkId;
282 }
285 
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
edm::Ref< typename HandleT::element_type, typename HandleT::element_type::value_type::value_type > makeRefTo(const HandleT &iHandle, typename HandleT::element_type::value_type::const_iterator itIter)
int minPixelCol() const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
void push_back(const T &t)
Definition: DetSet.h:67
std::vector< TrackingParticle > TrackingParticleCollection
ProductID id() const
Definition: HandleBase.cc:13
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< edmNew::DetSetVector< Phase2TrackerCluster1D > > phase2OTClustersToken_
uint16_t firstStrip() const
data_type const * const_iterator
Definition: DetSetNew.h:31
uint16_t size_type
edm::EDGetTokenT< edm::DetSetVector< PixelDigiSimLink > > siphase2OTSimLinksToken_
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClustersToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
int maxPixelRow() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
unsigned int firstRow() const
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int minPixelRow() const
ClusterTPAssociationProducer(const edm::ParameterSet &)
std::vector< OmniClusterRef > OmniClusterCollection
unsigned int column() const
static PackedDigiType pixelToChannel(unsigned int row, unsigned int col)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleToken_
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
T const * product() const
Definition: Handle.h:69
edm::EDGetTokenT< edm::DetSetVector< PixelDigiSimLink > > sipixelSimLinksToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int maxPixelCol() const
collection_type data
Definition: DetSet.h:81
Pixel cluster – collection of neighboring pixels above threshold.
HLT enums.
iterator end()
Definition: DetSetNew.h:56
static int pixelToChannel(int row, int col)
std::vector< std::pair< uint32_t, EncodedEventId > > getSimTrackId(const edm::Handle< edm::DetSetVector< T > > &simLinks, const DetId &detId, uint32_t channel) const
collection_type::const_iterator const_iterator
Definition: DetSet.h:32
edm::EDGetTokenT< edm::DetSetVector< StripDigiSimLink > > sistripSimLinksToken_
const std::vector< uint8_t > & amplitudes() const
def move(src, dest)
Definition: eostools.py:511
iterator begin()
Definition: DetSetNew.h:54