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 
22 
24  : _verbose(cfg.getParameter<bool>("verbose")),
25  _pixelSimLinkSrc(cfg.getParameter<edm::InputTag>("pixelSimLinkSrc")),
26  _stripSimLinkSrc(cfg.getParameter<edm::InputTag>("stripSimLinkSrc")),
27  _pixelClusterSrc(cfg.getParameter<edm::InputTag>("pixelClusterSrc")),
28  _stripClusterSrc(cfg.getParameter<edm::InputTag>("stripClusterSrc")),
29  _trackingParticleSrc(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))
30 {
31 
32  sipixelSimLinksToken_ = consumes<edm::DetSetVector<PixelDigiSimLink> >(cfg.getParameter<edm::InputTag>("pixelSimLinkSrc"));
33  sistripSimLinksToken_ = consumes<edm::DetSetVector<StripDigiSimLink> >(cfg.getParameter<edm::InputTag>("stripSimLinkSrc"));
34  pixelClustersToken_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(cfg.getParameter<edm::InputTag>("pixelClusterSrc"));
35  stripClustersToken_ = consumes<edmNew::DetSetVector<SiStripCluster> >(cfg.getParameter<edm::InputTag>("stripClusterSrc"));
36  trackingParticleToken_ = consumes<TrackingParticleCollection>(cfg.getParameter<edm::InputTag>("trackingParticleSrc"));
37 
38  produces<ClusterTPAssociationList>();
39 }
40 
42 }
43 
45  std::auto_ptr<ClusterTPAssociationList> clusterTPList(new ClusterTPAssociationList);
46 
47  // Pixel DigiSimLink
49  // iEvent.getByLabel(_pixelSimLinkSrc, sipixelSimLinks);
50  iEvent.getByToken(sipixelSimLinksToken_,sipixelSimLinks);
51 
52  // SiStrip DigiSimLink
54  iEvent.getByToken(sistripSimLinksToken_,sistripSimLinks);
55 
56  // Pixel Cluster
58  bool foundPixelClusters = iEvent.getByToken(pixelClustersToken_,pixelClusters);
59 
60  // Strip Cluster
62  bool foundStripClusters = iEvent.getByToken(stripClustersToken_,stripClusters);
63 
64 
65  // TrackingParticle
67  iEvent.getByToken(trackingParticleToken_,TPCollectionH);
68 
69  // prepare temporary map between SimTrackId and TrackingParticle index
70  std::map<std::pair<size_t, EncodedEventId>, TrackingParticleRef> mapping;
72  itp < TPCollectionH.product()->size(); ++itp) {
73  TrackingParticleRef trackingParticle(TPCollectionH, itp);
74 
75  // SimTracks inside TrackingParticle
76  EncodedEventId eid(trackingParticle->eventId());
77  //size_t index = 0;
78  for (std::vector<SimTrack>::const_iterator itrk = trackingParticle->g4Track_begin();
79  itrk != trackingParticle->g4Track_end(); ++itrk) {
80  std::pair<uint32_t, EncodedEventId> trkid(itrk->trackId(), eid);
81  //std::cout << "creating map for id: " << trkid.first << " with tp: " << trackingParticle.key() << std::endl;
82  mapping.insert(std::make_pair(trkid, trackingParticle));
83  }
84  }
85 
86  if ( foundPixelClusters ) {
87  // Pixel Clusters
88  for (edmNew::DetSetVector<SiPixelCluster>::const_iterator iter = pixelClusters->begin();
89  iter != pixelClusters->end(); ++iter) {
90  uint32_t detid = iter->id();
91  DetId detId(detid);
92  edmNew::DetSet<SiPixelCluster> link_pixel = (*iter);
94  di != link_pixel.end(); ++di) {
95  const SiPixelCluster& cluster = (*di);
97  edmNew::makeRefTo(pixelClusters, di);
98 
99  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
100  for (int irow = cluster.minPixelRow(); irow <= cluster.maxPixelRow(); ++irow) {
101  for (int icol = cluster.minPixelCol(); icol <= cluster.maxPixelCol(); ++icol) {
102  uint32_t channel = PixelChannelIdentifier::pixelToChannel(irow, icol);
103  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(getSimTrackId<PixelDigiSimLink>(sipixelSimLinks, detId, channel));
104  if (trkid.size()==0) continue;
105  simTkIds.insert(trkid.begin(),trkid.end());
106  }
107  }
108  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
109  iset != simTkIds.end(); iset++) {
110  auto ipos = mapping.find(*iset);
111  if (ipos != mapping.end()) {
112  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
113  clusterTPList->push_back(std::make_pair(OmniClusterRef(c_ref), ipos->second));
114  }
115  }
116  }
117  }
118  }
119 
120  if ( foundStripClusters ) {
121  // Strip Clusters
122  for (edmNew::DetSetVector<SiStripCluster>::const_iterator iter = stripClusters->begin(false), eter = stripClusters->end(false);
123  iter != eter; ++iter) {
124  if (!(*iter).isValid()) continue;
125  uint32_t detid = iter->id();
126  DetId detId(detid);
127  edmNew::DetSet<SiStripCluster> link_strip = (*iter);
129  di != link_strip.end(); di++) {
130  const SiStripCluster& cluster = (*di);
132  edmNew::makeRefTo(stripClusters, di);
133 
134  std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
135  int first = cluster.firstStrip();
136  int last = first + cluster.amplitudes().size();
137 
138  for (int istr = first; istr < last; ++istr) {
139  std::vector<std::pair<uint32_t, EncodedEventId> > trkid(getSimTrackId<StripDigiSimLink>(sistripSimLinks, detId, istr));
140  if (trkid.size()==0) continue;
141  simTkIds.insert(trkid.begin(),trkid.end());
142  }
143  for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
144  iset != simTkIds.end(); iset++) {
145  auto ipos = mapping.find(*iset);
146  if (ipos != mapping.end()) {
147  //std::cout << "cluster in detid: " << detid << " from tp: " << ipos->second.key() << " " << iset->first << std::endl;
148  clusterTPList->push_back(std::make_pair(OmniClusterRef(c_ref), ipos->second));
149  }
150  }
151  }
152  }
153  }
154 
155  iEvent.put(clusterTPList);
156 }
157 
158 template <typename T>
159 std::vector<std::pair<uint32_t, EncodedEventId> >
160 //std::pair<uint32_t, EncodedEventId>
162  const DetId& detId, uint32_t channel) const
163 {
164  //std::pair<uint32_t, EncodedEventId> simTrkId;
165  std::vector<std::pair<uint32_t, EncodedEventId> > simTrkId;
166  auto isearch = simLinks->find(detId);
167  if (isearch != simLinks->end()) {
168  // Loop over DigiSimLink in this det unit
169  edm::DetSet<T> link_detset = (*isearch);
170  for (typename edm::DetSet<T>::const_iterator it = link_detset.data.begin();
171  it != link_detset.data.end(); ++it) {
172  if (channel == it->channel()) {
173  simTrkId.push_back(std::make_pair(it->SimTrackId(), it->eventId()));
174  }
175  }
176  }
177  return simTrkId;
178 }
181 
T getParameter(std::string const &) const
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
tuple cfg
Definition: looper.py:259
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#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
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > pixelClustersToken_
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > stripClustersToken_
int maxPixelRow() const
int iEvent
Definition: GenABIO.cc:230
int minPixelRow() const
ClusterTPAssociationProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
edm::EDGetTokenT< TrackingParticleCollection > trackingParticleToken_
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< edm::DetSetVector< PixelDigiSimLink > > sipixelSimLinksToken_
std::vector< std::pair< OmniClusterRef, TrackingParticleRef > > ClusterTPAssociationList
int maxPixelCol() const
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
edm::EDGetTokenT< edm::DetSetVector< StripDigiSimLink > > sistripSimLinksToken_
const std::vector< uint8_t > & amplitudes() const
virtual void produce(edm::Event &, const edm::EventSetup &)
iterator begin()
Definition: DetSetNew.h:67