CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SimHitTPAssociationProducer.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 
12 
14 
16 
18  : _simHitSrc(cfg.getParameter<std::vector<edm::InputTag> >("simHitSrc")),
19  _trackingParticleSrc(cfg.getParameter<edm::InputTag>("trackingParticleSrc"))
20 {
21  produces<SimHitTPAssociationList>();
22 }
23 
25 }
26 
28  std::auto_ptr<SimHitTPAssociationList> simHitTPList(new SimHitTPAssociationList);
29 
30  // TrackingParticle
32  iEvent.getByLabel(_trackingParticleSrc, TPCollectionH);
33 
34  // prepare temporary map between SimTrackId and TrackingParticle index
35  std::map<std::pair<size_t, EncodedEventId>, TrackingParticleRef> mapping;
36  for (TrackingParticleCollection::size_type itp = 0; itp < TPCollectionH.product()->size(); ++itp) {
37  TrackingParticleRef trackingParticle(TPCollectionH, itp);
38  // SimTracks inside TrackingParticle
39  EncodedEventId eid(trackingParticle->eventId());
40  for (auto itrk = trackingParticle->g4Track_begin(); itrk != trackingParticle->g4Track_end(); ++itrk) {
41  std::pair<uint32_t, EncodedEventId> trkid(itrk->trackId(), eid);
42  mapping.insert(std::make_pair(trkid, trackingParticle));
43  }
44  }
45 
46  // PSimHits
47  for (auto psit=_simHitSrc.begin();psit<_simHitSrc.end();++psit) {
48  edm::Handle<edm::PSimHitContainer> PSimHitCollectionH;
49  iEvent.getByLabel(*psit, PSimHitCollectionH);
50  for (unsigned int psimHit = 0;psimHit != PSimHitCollectionH->size();++psimHit) {
51  TrackPSimHitRef pSimHitRef(PSimHitCollectionH,psimHit);
52  std::pair<uint32_t, EncodedEventId> simTkIds(pSimHitRef->trackId(),pSimHitRef->eventId());
53  auto ipos = mapping.find(simTkIds);
54  if (ipos != mapping.end()) {
55  simHitTPList->push_back(std::make_pair(ipos->second,pSimHitRef));
56  }
57  }
58  }
59 
60  std::sort(simHitTPList->begin(),simHitTPList->end(),simHitTPAssociationListGreater);
61  iEvent.put(simHitTPList);
62 
63 }
64 
67 
SimHitTPAssociationProducer(const edm::ParameterSet &)
static bool simHitTPAssociationListGreater(SimHitTPPair i, SimHitTPPair j)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< edm::InputTag > _simHitSrc
uint16_t size_type
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
T const * product() const
Definition: Handle.h:74
std::vector< SimHitTPPair > SimHitTPAssociationList
virtual void produce(edm::Event &, const edm::EventSetup &)