CMS 3D CMS Logo

TrackMCQuality.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: TrackMCQuality
4 // Class: TrackMCQuality
5 //
14 //
15 // Original Author: Jean-Roch Vlimant
16 // Created: Fri Mar 27 15:19:03 CET 2009
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
31 
33 
36 
37 //
38 // class decleration
39 //
40 
41 class TrackMCQuality final : public edm::global::EDProducer<> {
42 public:
43  explicit TrackMCQuality(const edm::ParameterSet &);
44  ~TrackMCQuality() override;
45 
46 private:
47  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
48 
49  // ----------member data ---------------------------
50 
54 
55  using Product = std::vector<float>;
56 };
57 
58 //
59 // constants, enums and typedefs
60 //
61 
62 //
63 // static data member definitions
64 //
65 
66 //
67 // constructors and destructor
68 //
70  : label_assoc(consumes<reco::TrackToTrackingParticleAssociator>(pset.getParameter<edm::InputTag>("associator"))),
71  label_tp(consumes<TrackingParticleCollection>(pset.getParameter<edm::InputTag>("trackingParticles"))),
72  label_tr(consumes<edm::View<reco::Track>>(pset.getParameter<edm::InputTag>("tracks"))) {
73  produces<Product>();
74 }
75 
77 
78 //
79 // member functions
80 //
81 
82 // ------------ method called to produce the data ------------
84  using namespace edm;
86  iEvent.getByToken(label_assoc, associator);
87 
89  iEvent.getByToken(label_tp, TPCollection);
90 
92  iEvent.getByToken(label_tr, trackCollection);
93 
94  reco::RecoToSimCollection recSimColl = associator->associateRecoToSim(trackCollection, TPCollection);
95 
96  // then loop the track collection
97  std::unique_ptr<Product> product(new Product(trackCollection->size(), 0));
98 
99  for (unsigned int iT = 0; iT != trackCollection->size(); ++iT) {
100  auto &prod = (*product)[iT];
101 
103 
104  // find it in the map
105  if (recSimColl.find(track) == recSimColl.end())
106  continue;
107 
108  auto const &tp = recSimColl[track];
109 
110  if (tp.empty())
111  continue; // can it be?
112  // nSimHits = tp[0].first->numberOfTrackerHits();
113  prod = tp[0].second;
114  // if (tp[0].first->charge() != track->charge()) isChargeMatched = false;
115  if ((tp[0].first->eventId().event() != 0) || (tp[0].first->eventId().bunchCrossing() != 0))
116  prod = -prod;
117  }
118 
119  iEvent.put(std::move(product));
120 }
121 
122 // define this as a plug-in
edm::EDGetTokenT< edm::View< reco::Track > > label_tr
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const_iterator find(const key_type &k) const
find element with specified reference key
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const_iterator end() const
last iterator over the map (read only)
TrackMCQuality(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
std::vector< float > Product
edm::EDGetTokenT< reco::TrackToTrackingParticleAssociator > label_assoc
~TrackMCQuality() override
fixed size matrix
HLT enums.
std::vector< TrackingParticle > TrackingParticleCollection
edm::EDGetTokenT< TrackingParticleCollection > label_tp
def move(src, dest)
Definition: eostools.py:511