CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackMCQuality.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TrackMCQuality
4 // Class: TrackMCQuality
5 //
13 //
14 // Original Author: Jean-Roch Vlimant
15 // Created: Fri Mar 27 15:19:03 CET 2009
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
32 
34 
37 
38 //
39 // class decleration
40 //
41 
43  public:
44  explicit TrackMCQuality(const edm::ParameterSet&);
46 
47  private:
48  virtual void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
49 
50  // ----------member data ---------------------------
51 
55 };
56 
57 //
58 // constants, enums and typedefs
59 //
60 
61 
62 //
63 // static data member definitions
64 //
65 
66 //
67 // constructors and destructor
68 //
70  label_tr(consumes<reco::TrackToTrackingParticleAssociator>(pset.getParameter< edm::InputTag >("label_tr"))),
71  label_tp(consumes<TrackingParticleCollection>(pset.getParameter< edm::InputTag >("label_tp"))),
72  label_associator(consumes<edm::View<reco::Track> >(pset.getParameter< edm::InputTag >("associator")))
73 {
74  produces<reco::TrackCollection>();
75 }
76 
77 
79 {
80 }
81 
82 
83 //
84 // member functions
85 //
86 
87 // ------------ method called to produce the data ------------
88 void
90 {
91 
92  using namespace edm;
94  iEvent.getByToken(label_associator,associator);
95 
98 
101 
102  reco::RecoToSimCollection recSimColl=associator->associateRecoToSim(trackCollection,
103  TPCollection);
104 
105  //then loop the track collection
106  std::auto_ptr<reco::TrackCollection> outTracks(new reco::TrackCollection(trackCollection->size()));
107 
108  for (unsigned int iT=0;iT!=trackCollection->size();++iT){
110  bool matched=false;
111  //find it in the map
112  if (recSimColl.find(track)!=recSimColl.end()){
113  // you can get the data if you want
114  std::vector<std::pair<TrackingParticleRef, double> > tp= recSimColl[track];
115  matched=true;
116  }
117  else{
118  matched=false;
119  }
120 
121  //copy the track into the new container
122  (*outTracks)[iT] = reco::Track(*track);
123  if (matched){
124  (*outTracks)[iT].setQuality(reco::TrackBase::qualitySize); //is not assigned to any quality. use it as a fake/matched flag
125  }
126  }
127 
128  iEvent.put(outTracks);
129 }
130 
131 //define this as a plug-in
std::vector< TrackingParticle > TrackingParticleCollection
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
edm::EDGetTokenT< edm::View< reco::Track > > label_associator
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
TrackMCQuality(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
edm::EDGetTokenT< TrackingParticleCollection > label_tp
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > label_tr