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/MessageLogger/interface/MessageLogger.h"
00012
00013 #include "DataFormats/Common/interface/Handle.h"
00014
00015 #include "SimGeneral/TrackingAnalysis/interface/SimHitTPAssociationProducer.h"
00016
00017 SimHitTPAssociationProducer::SimHitTPAssociationProducer(const edm::ParameterSet & cfg)
00018 : _simHitSrc(cfg.getParameter<std::vector<edm::InputTag> >("simHitSrc")),
00019 _trackingParticleSrc(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))
00020 {
00021 produces<SimHitTPAssociationList>();
00022 }
00023
00024 SimHitTPAssociationProducer::~SimHitTPAssociationProducer() {
00025 }
00026
00027 void SimHitTPAssociationProducer::produce(edm::Event& iEvent, const edm::EventSetup& es) {
00028 std::auto_ptr<SimHitTPAssociationList> simHitTPList(new SimHitTPAssociationList);
00029
00030
00031 edm::Handle<TrackingParticleCollection> TPCollectionH;
00032 iEvent.getByLabel(_trackingParticleSrc, TPCollectionH);
00033
00034
00035 std::map<std::pair<size_t, EncodedEventId>, TrackingParticleRef> mapping;
00036 for (TrackingParticleCollection::size_type itp = 0; itp < TPCollectionH.product()->size(); ++itp) {
00037 TrackingParticleRef trackingParticle(TPCollectionH, itp);
00038
00039 EncodedEventId eid(trackingParticle->eventId());
00040 for (auto itrk = trackingParticle->g4Track_begin(); itrk != trackingParticle->g4Track_end(); ++itrk) {
00041 std::pair<uint32_t, EncodedEventId> trkid(itrk->trackId(), eid);
00042 mapping.insert(std::make_pair(trkid, trackingParticle));
00043 }
00044 }
00045
00046
00047 for (auto psit=_simHitSrc.begin();psit<_simHitSrc.end();++psit) {
00048 edm::Handle<edm::PSimHitContainer> PSimHitCollectionH;
00049 iEvent.getByLabel(*psit, PSimHitCollectionH);
00050 for (unsigned int psimHit = 0;psimHit != PSimHitCollectionH->size();++psimHit) {
00051 TrackPSimHitRef pSimHitRef(PSimHitCollectionH,psimHit);
00052 std::pair<uint32_t, EncodedEventId> simTkIds(pSimHitRef->trackId(),pSimHitRef->eventId());
00053 auto ipos = mapping.find(simTkIds);
00054 if (ipos != mapping.end()) {
00055 simHitTPList->push_back(std::make_pair(ipos->second,pSimHitRef));
00056 }
00057 }
00058 }
00059
00060 std::sort(simHitTPList->begin(),simHitTPList->end(),simHitTPAssociationListGreater);
00061 iEvent.put(simHitTPList);
00062
00063 }
00064
00065 #include "FWCore/PluginManager/interface/ModuleDef.h"
00066 #include "FWCore/Framework/interface/MakerMacros.h"
00067
00068 DEFINE_FWK_MODULE(SimHitTPAssociationProducer);