CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ClusterTPAssociationProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <vector>
3 #include <iostream>
4 #include <fstream>
5 #include <utility>
6 
8 
13 
19 
25 
27 
29  : _verbose(cfg.getParameter<bool>("verbose")),
30  _pixelSimLinkSrc(cfg.getParameter<edm::InputTag>("pixelSimLinkSrc")),
31  _stripSimLinkSrc(cfg.getParameter<edm::InputTag>("stripSimLinkSrc")),
32  _pixelClusterSrc(cfg.getParameter<edm::InputTag>("pixelClusterSrc")),
33  _stripClusterSrc(cfg.getParameter<edm::InputTag>("stripClusterSrc")),
34  _trackingParticleSrc(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))
35 {
36  produces<ClusterTPAssociationList>();
37 }
38 
40 }
41 
43  std::auto_ptr<ClusterTPAssociationList> clusterTPList(new ClusterTPAssociationList);
44 
45  // Pixel DigiSimLink
47  iEvent.getByLabel(_pixelSimLinkSrc, sipixelSimLinks);
48 
49  // SiStrip DigiSimLink
51  iEvent.getByLabel(_stripSimLinkSrc, sistripSimLinks);
52 
53  // Pixel Cluster
55  iEvent.getByLabel(_pixelClusterSrc, pixelClusters);
56 
57  // Strip Cluster
59  iEvent.getByLabel(_stripClusterSrc, stripClusters);
60 
61  // TrackingParticle
63  iEvent.getByLabel(_trackingParticleSrc, TPCollectionH);
64 
65  // prepare temporary map between SimTrackId and TrackingParticle index
66  std::map<std::pair<size_t, EncodedEventId>, TrackingParticleRef> mapping;
68  itp < TPCollectionH.product()->size(); ++itp) {
69  TrackingParticleRef trackingParticle(TPCollectionH, itp);
70 
71  // SimTracks inside TrackingParticle
72  EncodedEventId eid(trackingParticle->eventId());
73  //size_t index = 0;
74  for (std::vector<SimTrack>::const_iterator itrk = trackingParticle->g4Track_begin();
75  itrk != trackingParticle->g4Track_end(); ++itrk) {
76  std::pair<uint32_t, EncodedEventId> trkid(itrk->trackId(), eid);
77  //std::cout << "creating map for id: " << trkid.first << " with tp: " << trackingParticle.key() << std::endl;
78  mapping.insert(std::make_pair(trkid, trackingParticle));
79  }
80  }
81 
82  // Pixel Clusters
83  for (edmNew::DetSetVector<SiPixelCluster>::const_iterator iter = pixelClusters->begin();
84  iter != pixelClusters->end(); ++iter)
85  {
86  uint32_t detid = iter->id();
87  DetId detId(detid);
88  edmNew::DetSet<SiPixelCluster> link_pixel = (*iter);
90  di != link_pixel.end(); ++di)
91  {
92  const SiPixelCluster& cluster = (*di);
94  edmNew::makeRefTo(pixelClusters, di);
95 
96  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
97  for (int irow = cluster.minPixelRow(); irow <= cluster.maxPixelRow(); ++irow) {
98  for (int icol = cluster.minPixelCol(); icol <= cluster.maxPixelCol(); ++icol) {
99  uint32_t channel = PixelChannelIdentifier::pixelToChannel(irow, icol);
100  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(getSimTrackId<PixelDigiSimLink>(sipixelSimLinks, detId, channel));
101  if (trkid.size()==0) continue;
102  simTkIds.insert(trkid.begin(),trkid.end());
103  }
104  }
105  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
106  iset != simTkIds.end(); iset++) {
107  auto ipos = mapping.find(*iset);
108  if (ipos != mapping.end()) {
109  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
110  clusterTPList->push_back(std::make_pair(OmniClusterRef(c_ref), ipos->second));
111  }
112  }
113  }
114  }
115 
116  // Strip Clusters
117  for (edmNew::DetSetVector<SiStripCluster>::const_iterator iter = stripClusters->begin();
118  iter != stripClusters->end(); ++iter) {
119  uint32_t detid = iter->id();
120  DetId detId(detid);
121  edmNew::DetSet<SiStripCluster> link_strip = (*iter);
123  di != link_strip.end(); di++) {
124  const SiStripCluster& cluster = (*di);
126  edmNew::makeRefTo(stripClusters, di);
127 
128  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
129  int first = cluster.firstStrip();
130  int last = first + cluster.amplitudes().size();
131 
132  for (int istr = first; istr < last; ++istr) {
133  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(getSimTrackId<StripDigiSimLink>(sistripSimLinks, detId, istr));
134  if (trkid.size()==0) continue;
135  simTkIds.insert(trkid.begin(),trkid.end());
136  }
137  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
138  iset != simTkIds.end(); iset++) {
139  auto ipos = mapping.find(*iset);
140  if (ipos != mapping.end()) {
141  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
142  clusterTPList->push_back(std::make_pair(OmniClusterRef(c_ref), ipos->second));
143  }
144  }
145  }
146  }
147  iEvent.put(clusterTPList);
148 }
149 template <typename T>
150 std::vector<std::pair<uint32_t, EncodedEventId> >
151 //std::pair<uint32_t, EncodedEventId>
153  const DetId& detId, uint32_t channel) const
154 {
155  //std::pair<uint32_t, EncodedEventId> simTrkId;
156  std::vector<std::pair<uint32_t, EncodedEventId> > simTrkId;
157  auto isearch = simLinks->find(detId);
158  if (isearch != simLinks->end()) {
159  // Loop over DigiSimLink in this det unit
160  edm::DetSet<T> link_detset = (*isearch);
161  for (typename edm::DetSet<T>::const_iterator it = link_detset.data.begin();
162  it != link_detset.data.end(); ++it) {
163  if (channel == it->channel()) {
164  simTrkId.push_back(std::make_pair(it->SimTrackId(), it->eventId()));
165  }
166  }
167  }
168  return simTrkId;
169 }
172 
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
void push_back(const T &t)
Definition: DetSet.h:68
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
uint16_t firstStrip() const
data_type const * const_iterator
Definition: DetSetNew.h:30
uint16_t size_type
int maxPixelRow() const
int iEvent
Definition: GenABIO.cc:243
int minPixelRow() const
ClusterTPAssociationProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
bool first
Definition: L1TdeRCT.cc:79
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
Definition: DetId.h:18
std::vector< std::pair< OmniClusterRef, TrackingParticleRef > > ClusterTPAssociationList
int maxPixelCol() const
T const * product() const
Definition: Handle.h:81
collection_type data
Definition: DetSet.h:78
Pixel cluster – collection of neighboring pixels above threshold.
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
const std::vector< uint8_t > & amplitudes() const
virtual void produce(edm::Event &, const edm::EventSetup &)
iterator begin()
Definition: DetSetNew.h:67