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