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 
32 
33 namespace {
34 
35  template <typename T>
36  void getSimTrackId(std::vector<UniqueSimTrackId>& simTrkId,
37  const edm::Handle<edm::DetSetVector<T> >& simLinks,
38  const DetId& detId,
39  uint32_t channel) {
40  auto isearch = simLinks->find(detId);
41  if (isearch != simLinks->end()) {
42  // Loop over DigiSimLink in this det unit
43  edm::DetSet<T> link_detset = (*isearch);
44  for (typename edm::DetSet<T>::const_iterator it = link_detset.data.begin(); it != link_detset.data.end(); ++it) {
45  if (channel == it->channel()) {
46  simTrkId.emplace_back(it->SimTrackId(), it->eventId());
47  }
48  }
49  }
50  }
51 
52 } // namespace
53 
55 public:
56  using OmniClusterCollection = std::vector<OmniClusterRef>;
57 
59  ~ClusterTPAssociationProducer() override = default;
60 
61  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
62 
63 private:
64  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
65 
74 };
75 
77  : sipixelSimLinksToken_(
78  consumes<edm::DetSetVector<PixelDigiSimLink> >(cfg.getParameter<edm::InputTag>("pixelSimLinkSrc"))),
79  sistripSimLinksToken_(
80  consumes<edm::DetSetVector<StripDigiSimLink> >(cfg.getParameter<edm::InputTag>("stripSimLinkSrc"))),
81  siphase2OTSimLinksToken_(
82  consumes<edm::DetSetVector<PixelDigiSimLink> >(cfg.getParameter<edm::InputTag>("phase2OTSimLinkSrc"))),
83  pixelClustersToken_(
84  consumes<edmNew::DetSetVector<SiPixelCluster> >(cfg.getParameter<edm::InputTag>("pixelClusterSrc"))),
85  stripClustersToken_(
86  consumes<edmNew::DetSetVector<SiStripCluster> >(cfg.getParameter<edm::InputTag>("stripClusterSrc"))),
87  phase2OTClustersToken_(consumes<edmNew::DetSetVector<Phase2TrackerCluster1D> >(
88  cfg.getParameter<edm::InputTag>("phase2OTClusterSrc"))),
89  trackingParticleToken_(
90  consumes<TrackingParticleCollection>(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))),
91  throwOnMissingCollections_(cfg.getParameter<bool>("throwOnMissingCollections")) {
92  produces<ClusterTPAssociation>();
93 }
94 
97  desc.add<edm::InputTag>("simTrackSrc", edm::InputTag("g4SimHits"));
98  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis"));
99  desc.add<edm::InputTag>("stripSimLinkSrc", edm::InputTag("simSiStripDigis"));
100  desc.add<edm::InputTag>("phase2OTSimLinkSrc", edm::InputTag("simSiPixelDigis", "Tracker"));
101  desc.add<edm::InputTag>("pixelClusterSrc", edm::InputTag("siPixelClusters"));
102  desc.add<edm::InputTag>("stripClusterSrc", edm::InputTag("siStripClusters"));
103  desc.add<edm::InputTag>("phase2OTClusterSrc", edm::InputTag("siPhase2Clusters"));
104  desc.add<edm::InputTag>("trackingParticleSrc", edm::InputTag("mix", "MergedTrackTruth"));
105  desc.add<bool>("throwOnMissingCollections", true);
106  descriptions.add("tpClusterProducerDefault", desc);
107 }
108 
110  // Pixel Cluster
112  bool foundPixelClusters = iEvent.getByToken(pixelClustersToken_, pixelClusters);
113 
114  // Strip Cluster
116  bool foundStripClusters = iEvent.getByToken(stripClustersToken_, stripClusters);
117 
118  // Phase2 Cluster
120  bool foundPhase2OTClusters = iEvent.getByToken(phase2OTClustersToken_, phase2OTClusters);
121 
122  // Pixel DigiSimLink
124  auto pixelSimLinksFound = iEvent.getByToken(sipixelSimLinksToken_, sipixelSimLinks);
125  if (not throwOnMissingCollections_ and foundPixelClusters and not pixelSimLinksFound) {
126  auto clusterTPList = std::make_unique<ClusterTPAssociation>();
127  iEvent.put(std::move(clusterTPList));
128  return;
129  }
130 
131  // SiStrip DigiSimLink
133  auto stripSimLinksFound = iEvent.getByToken(sistripSimLinksToken_, sistripSimLinks);
134  if (not throwOnMissingCollections_ and foundStripClusters and not stripSimLinksFound) {
135  auto clusterTPList = std::make_unique<ClusterTPAssociation>();
136  iEvent.put(std::move(clusterTPList));
137  return;
138  }
139 
140  // Phase2 OT DigiSimLink
142  auto phase2OTSimLinksFound = iEvent.getByToken(siphase2OTSimLinksToken_, siphase2OTSimLinks);
143  if (not throwOnMissingCollections_ and foundPhase2OTClusters and not phase2OTSimLinksFound) {
144  auto clusterTPList = std::make_unique<ClusterTPAssociation>();
145  iEvent.put(std::move(clusterTPList));
146  return;
147  }
148 
149  // TrackingParticle
151  auto tpFound = iEvent.getByToken(trackingParticleToken_, TPCollectionH);
152  if (not throwOnMissingCollections_ and not tpFound) {
153  auto clusterTPList = std::make_unique<ClusterTPAssociation>();
154  iEvent.put(std::move(clusterTPList));
155  return;
156  }
157 
158  auto clusterTPList = std::make_unique<ClusterTPAssociation>(TPCollectionH);
159 
160  // prepare temporary map between SimTrackId and TrackingParticle index
161  std::unordered_map<UniqueSimTrackId, TrackingParticleRef, UniqueSimTrackIdHash> mapping;
162  auto const& tpColl = *TPCollectionH.product();
163  for (TrackingParticleCollection::size_type itp = 0; itp < tpColl.size(); ++itp) {
164  TrackingParticleRef trackingParticleRef(TPCollectionH, itp);
165  auto const& trackingParticle = tpColl[itp];
166  // SimTracks inside TrackingParticle
167  EncodedEventId eid(trackingParticle.eventId());
168  //size_t index = 0;
169  for (auto const& trk : trackingParticle.g4Tracks()) {
170  UniqueSimTrackId trkid(trk.trackId(), eid);
171  //std::cout << "creating map for id: " << trkid.first << " with tp: " << trackingParticle.key() << std::endl;
172  mapping.insert(std::make_pair(trkid, trackingParticleRef));
173  }
174  }
175 
176  std::unordered_set<UniqueSimTrackId, UniqueSimTrackIdHash> simTkIds;
177  std::vector<UniqueSimTrackId> trkid;
178  if (foundPixelClusters) {
179  // Pixel Clusters
180  clusterTPList->addKeyID(pixelClusters.id());
182  iter != pixelClusters->end();
183  ++iter) {
184  uint32_t detid = iter->id();
185  DetId detId(detid);
186  edmNew::DetSet<SiPixelCluster> link_pixel = (*iter);
187  for (edmNew::DetSet<SiPixelCluster>::const_iterator di = link_pixel.begin(); di != link_pixel.end(); ++di) {
188  const SiPixelCluster& cluster = (*di);
190 
191  simTkIds.clear();
192  for (int irow = cluster.minPixelRow(); irow <= cluster.maxPixelRow(); ++irow) {
193  for (int icol = cluster.minPixelCol(); icol <= cluster.maxPixelCol(); ++icol) {
194  uint32_t channel = PixelChannelIdentifier::pixelToChannel(irow, icol);
195  trkid.clear();
196  getSimTrackId<PixelDigiSimLink>(trkid, sipixelSimLinks, detId, channel);
197  simTkIds.insert(trkid.begin(), trkid.end());
198  }
199  }
200  for (auto iset = simTkIds.begin(); iset != simTkIds.end(); iset++) {
201  auto ipos = mapping.find(*iset);
202  if (ipos != mapping.end()) {
203  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
204  clusterTPList->emplace_back(OmniClusterRef(c_ref), ipos->second);
205  }
206  }
207  }
208  }
209  }
210 
211  if (foundStripClusters) {
212  // Strip Clusters
213  clusterTPList->addKeyID(stripClusters.id());
215  eter = stripClusters->end(false);
216  iter != eter;
217  ++iter) {
218  if (!(*iter).isValid())
219  continue;
220  uint32_t detid = iter->id();
221  DetId detId(detid);
222  edmNew::DetSet<SiStripCluster> link_strip = (*iter);
223  for (edmNew::DetSet<SiStripCluster>::const_iterator di = link_strip.begin(); di != link_strip.end(); di++) {
224  const SiStripCluster& cluster = (*di);
226 
227  simTkIds.clear();
228  int first = cluster.firstStrip();
229  int last = first + cluster.amplitudes().size();
230 
231  for (int istr = first; istr < last; ++istr) {
232  trkid.clear();
233  getSimTrackId<StripDigiSimLink>(trkid, sistripSimLinks, detId, istr);
234  simTkIds.insert(trkid.begin(), trkid.end());
235  }
236  for (auto iset = simTkIds.begin(); iset != simTkIds.end(); iset++) {
237  auto ipos = mapping.find(*iset);
238  if (ipos != mapping.end()) {
239  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
240  clusterTPList->emplace_back(OmniClusterRef(c_ref), ipos->second);
241  }
242  }
243  }
244  }
245  }
246 
247  if (foundPhase2OTClusters) {
248  // Phase2 Clusters
249  clusterTPList->addKeyID(phase2OTClusters.id());
250  if (phase2OTClusters.isValid()) {
252  eter = phase2OTClusters->end(false);
253  iter != eter;
254  ++iter) {
255  if (!(*iter).isValid())
256  continue;
257  uint32_t detid = iter->id();
258  DetId detId(detid);
259  edmNew::DetSet<Phase2TrackerCluster1D> link_phase2 = (*iter);
260  for (edmNew::DetSet<Phase2TrackerCluster1D>::const_iterator di = link_phase2.begin(); di != link_phase2.end();
261  di++) {
262  const Phase2TrackerCluster1D& cluster = (*di);
265 
266  simTkIds.clear();
267 
268  for (unsigned int istr(0); istr < cluster.size(); ++istr) {
269  uint32_t channel = Phase2TrackerDigi::pixelToChannel(cluster.firstRow() + istr, cluster.column());
270  trkid.clear();
271  getSimTrackId<PixelDigiSimLink>(trkid, siphase2OTSimLinks, detId, channel);
272  simTkIds.insert(trkid.begin(), trkid.end());
273  }
274 
275  for (auto iset = simTkIds.begin(); iset != simTkIds.end(); iset++) {
276  auto ipos = mapping.find(*iset);
277  if (ipos != mapping.end()) {
278  clusterTPList->emplace_back(OmniClusterRef(c_ref), ipos->second);
279  }
280  }
281  }
282  }
283  }
284  }
285  clusterTPList->sortAndUnique();
286  iEvent.put(std::move(clusterTPList));
287 }
288 
291 
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)
uint16_t firstStrip() const
std::pair< uint32_t, EncodedEventId > UniqueSimTrackId
std::vector< OmniClusterRef > OmniClusterCollection
edm::EDGetTokenT< edmNew::DetSetVector< Phase2TrackerCluster1D > > phase2OTClustersToken_
int maxPixelRow() const
T const * product() const
Definition: Handle.h:70
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_
int minPixelRow() const
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
int minPixelCol() const
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
ClusterTPAssociationProducer(const edm::ParameterSet &)
decltype(auto) emplace_back(Args &&... args)
Definition: DetSet.h:68
SiStripCluster const & amplitudes() const
auto size() const
static PackedDigiType pixelToChannel(unsigned int row, unsigned int col)
unsigned int firstRow() const
int maxPixelCol() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleToken_
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
unsigned int column() const
edm::EDGetTokenT< edm::DetSetVector< PixelDigiSimLink > > sipixelSimLinksToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~ClusterTPAssociationProducer() override=default
collection_type data
Definition: DetSet.h:80
Pixel cluster – collection of neighboring pixels above threshold.
HLT enums.
iterator end()
Definition: DetSetNew.h:53
static int pixelToChannel(int row, int col)
std::vector< TrackingParticle > TrackingParticleCollection
collection_type::const_iterator const_iterator
Definition: DetSet.h:31
edm::EDGetTokenT< edm::DetSetVector< StripDigiSimLink > > sistripSimLinksToken_
def move(src, dest)
Definition: eostools.py:511
iterator begin()
Definition: DetSetNew.h:51