Go to the documentation of this file.00001 #include <memory>
00002 #include <vector>
00003 #include <iostream>
00004 #include <fstream>
00005 #include <utility>
00006
00007 #include "FWCore/Utilities/interface/InputTag.h"
00008
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 #include "FWCore/Utilities/interface/EDMException.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013
00014 #include "DataFormats/Common/interface/Handle.h"
00015 #include "DataFormats/Candidate/interface/Candidate.h"
00016 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00017 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00018 #include "DataFormats/SiPixelDetId/interface/PixelChannelIdentifier.h"
00019
00020 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00021 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00022 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00023 #include "SimDataFormats/TrackerDigiSimLink/interface/StripDigiSimLink.h"
00024 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"
00025
00026 #include "SimTracker/TrackerHitAssociation/interface/ClusterTPAssociationProducer.h"
00027
00028 ClusterTPAssociationProducer::ClusterTPAssociationProducer(const edm::ParameterSet & cfg)
00029 : _verbose(cfg.getParameter<bool>("verbose")),
00030 _pixelSimLinkSrc(cfg.getParameter<edm::InputTag>("pixelSimLinkSrc")),
00031 _stripSimLinkSrc(cfg.getParameter<edm::InputTag>("stripSimLinkSrc")),
00032 _pixelClusterSrc(cfg.getParameter<edm::InputTag>("pixelClusterSrc")),
00033 _stripClusterSrc(cfg.getParameter<edm::InputTag>("stripClusterSrc")),
00034 _trackingParticleSrc(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))
00035 {
00036 produces<ClusterTPAssociationList>();
00037 }
00038
00039 ClusterTPAssociationProducer::~ClusterTPAssociationProducer() {
00040 }
00041
00042 void ClusterTPAssociationProducer::produce(edm::Event& iEvent, const edm::EventSetup& es) {
00043 std::auto_ptr<ClusterTPAssociationList> clusterTPList(new ClusterTPAssociationList);
00044
00045
00046 edm::Handle<edm::DetSetVector<PixelDigiSimLink> > sipixelSimLinks;
00047 iEvent.getByLabel(_pixelSimLinkSrc, sipixelSimLinks);
00048
00049
00050 edm::Handle<edm::DetSetVector<StripDigiSimLink> > sistripSimLinks;
00051 iEvent.getByLabel(_stripSimLinkSrc, sistripSimLinks);
00052
00053
00054 edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00055 iEvent.getByLabel(_pixelClusterSrc, pixelClusters);
00056
00057
00058 edm::Handle<edmNew::DetSetVector<SiStripCluster> > stripClusters;
00059 iEvent.getByLabel(_stripClusterSrc, stripClusters);
00060
00061
00062 edm::Handle<TrackingParticleCollection> TPCollectionH;
00063 iEvent.getByLabel(_trackingParticleSrc, TPCollectionH);
00064
00065
00066 std::map<std::pair<size_t, EncodedEventId>, TrackingParticleRef> mapping;
00067 for (TrackingParticleCollection::size_type itp = 0;
00068 itp < TPCollectionH.product()->size(); ++itp) {
00069 TrackingParticleRef trackingParticle(TPCollectionH, itp);
00070
00071
00072 EncodedEventId eid(trackingParticle->eventId());
00073
00074 for (std::vector<SimTrack>::const_iterator itrk = trackingParticle->g4Track_begin();
00075 itrk != trackingParticle->g4Track_end(); ++itrk) {
00076 std::pair<uint32_t, EncodedEventId> trkid(itrk->trackId(), eid);
00077
00078 mapping.insert(std::make_pair(trkid, trackingParticle));
00079 }
00080 }
00081
00082
00083 for (edmNew::DetSetVector<SiPixelCluster>::const_iterator iter = pixelClusters->begin();
00084 iter != pixelClusters->end(); ++iter)
00085 {
00086 uint32_t detid = iter->id();
00087 DetId detId(detid);
00088 edmNew::DetSet<SiPixelCluster> link_pixel = (*iter);
00089 for (edmNew::DetSet<SiPixelCluster>::const_iterator di = link_pixel.begin();
00090 di != link_pixel.end(); ++di)
00091 {
00092 const SiPixelCluster& cluster = (*di);
00093 edm::Ref<edmNew::DetSetVector<SiPixelCluster>, SiPixelCluster> c_ref =
00094 edmNew::makeRefTo(pixelClusters, di);
00095
00096 std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
00097 for (int irow = cluster.minPixelRow(); irow <= cluster.maxPixelRow(); ++irow) {
00098 for (int icol = cluster.minPixelCol(); icol <= cluster.maxPixelCol(); ++icol) {
00099 uint32_t channel = PixelChannelIdentifier::pixelToChannel(irow, icol);
00100 std::vector<std::pair<uint32_t, EncodedEventId> > trkid(getSimTrackId<PixelDigiSimLink>(sipixelSimLinks, detId, channel));
00101 if (trkid.size()==0) continue;
00102 simTkIds.insert(trkid.begin(),trkid.end());
00103 }
00104 }
00105 for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
00106 iset != simTkIds.end(); iset++) {
00107 auto ipos = mapping.find(*iset);
00108 if (ipos != mapping.end()) {
00109
00110 clusterTPList->push_back(std::make_pair(OmniClusterRef(c_ref), ipos->second));
00111 }
00112 }
00113 }
00114 }
00115
00116
00117 for (edmNew::DetSetVector<SiStripCluster>::const_iterator iter = stripClusters->begin();
00118 iter != stripClusters->end(); ++iter) {
00119 uint32_t detid = iter->id();
00120 DetId detId(detid);
00121 edmNew::DetSet<SiStripCluster> link_strip = (*iter);
00122 for (edmNew::DetSet<SiStripCluster>::const_iterator di = link_strip.begin();
00123 di != link_strip.end(); di++) {
00124 const SiStripCluster& cluster = (*di);
00125 edm::Ref<edmNew::DetSetVector<SiStripCluster>, SiStripCluster> c_ref =
00126 edmNew::makeRefTo(stripClusters, di);
00127
00128 std::set<std::pair<uint32_t, EncodedEventId> > simTkIds;
00129 int first = cluster.firstStrip();
00130 int last = first + cluster.amplitudes().size();
00131
00132 for (int istr = first; istr < last; ++istr) {
00133 std::vector<std::pair<uint32_t, EncodedEventId> > trkid(getSimTrackId<StripDigiSimLink>(sistripSimLinks, detId, istr));
00134 if (trkid.size()==0) continue;
00135 simTkIds.insert(trkid.begin(),trkid.end());
00136 }
00137 for (std::set<std::pair<uint32_t, EncodedEventId> >::const_iterator iset = simTkIds.begin();
00138 iset != simTkIds.end(); iset++) {
00139 auto ipos = mapping.find(*iset);
00140 if (ipos != mapping.end()) {
00141
00142 clusterTPList->push_back(std::make_pair(OmniClusterRef(c_ref), ipos->second));
00143 }
00144 }
00145 }
00146 }
00147 iEvent.put(clusterTPList);
00148 }
00149 template <typename T>
00150 std::vector<std::pair<uint32_t, EncodedEventId> >
00151
00152 ClusterTPAssociationProducer::getSimTrackId(const edm::Handle<edm::DetSetVector<T> >& simLinks,
00153 const DetId& detId, uint32_t channel) const
00154 {
00155
00156 std::vector<std::pair<uint32_t, EncodedEventId> > simTrkId;
00157 auto isearch = simLinks->find(detId);
00158 if (isearch != simLinks->end()) {
00159
00160 edm::DetSet<T> link_detset = (*isearch);
00161 for (typename edm::DetSet<T>::const_iterator it = link_detset.data.begin();
00162 it != link_detset.data.end(); ++it) {
00163 if (channel == it->channel()) {
00164 simTrkId.push_back(std::make_pair(it->SimTrackId(), it->eventId()));
00165 }
00166 }
00167 }
00168 return simTrkId;
00169 }
00170 #include "FWCore/PluginManager/interface/ModuleDef.h"
00171 #include "FWCore/Framework/interface/MakerMacros.h"
00172
00173 DEFINE_FWK_MODULE(ClusterTPAssociationProducer);