CMS 3D CMS Logo

PackedCandidateTrackChi2Producer.cc
Go to the documentation of this file.
11 
12 namespace pat {
13 
16 
17  public:
19  : candidateToken_(consumes<pat::PackedCandidateCollection>(iConfig.getParameter<edm::InputTag>("candidates"))),
20  trackToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("trackCollection"))),
21  doLostTracks_(iConfig.getParameter<bool>("doLostTracks")) {
22  if (doLostTracks_) {
23  track2LostTrackToken_ = consumes<edm::Association<pat::PackedCandidateCollection>>(
24  iConfig.getParameter<edm::InputTag>("candidates"));
25  } else {
27  consumes<edm::Association<reco::PFCandidateCollection>>(iConfig.getParameter<edm::InputTag>("candidates"));
28  }
29 
30  produces<FloatMap>();
31  }
32 
33  void produce(edm::Event&, const edm::EventSetup&) override;
34 
36 
37  private:
42  const bool doLostTracks_;
43  static const uint8_t roundingPrecision = 8;
44  };
45 
46 } // namespace pat
47 
49  auto const candidates = iEvent.getHandle(candidateToken_);
50 
51  const edm::Association<reco::PFCandidateCollection>* candidate2PF = nullptr;
52  if (!doLostTracks_) {
53  candidate2PF = &iEvent.get(candidate2PFToken_);
54  }
55 
56  const edm::Association<pat::PackedCandidateCollection>* tracks2LT = nullptr;
58  if (doLostTracks_) {
59  tracks2LT = &iEvent.get(track2LostTrackToken_);
60  iEvent.getByToken(trackToken_, trks);
61  }
62 
63  const auto nCand = candidates->size();
64  std::vector<float> trkChi2Map(nCand, 0);
65 
66  if (doLostTracks_) { //for Lost tracks we don't have references to PFCands, so we must loop over tracks and check keys...
67  for (size_t i = 0; i < trks->size(); i++) {
68  const auto& trk = reco::TrackRef(trks, i);
69  const auto& lostTrack = (*tracks2LT)[trk];
70  if (lostTrack.isNonnull()) {
71  const float nChi2 = trk->normalizedChi2();
72  trkChi2Map.at(lostTrack.key()) = MiniFloatConverter::reduceMantissaToNbitsRounding<roundingPrecision>(nChi2);
73  }
74  }
75  } else { //for the regular PackedPFCands we have direct references...
76  for (size_t i = 0; i < nCand; i++) {
77  const auto& cand = pat::PackedCandidateRef(candidates, i);
78 
79  // ignore neutral candidates or without track
80  if (cand->charge() == 0 || !cand->hasTrackDetails())
81  continue;
82 
83  const auto& candTrack = (*candidate2PF)[cand]->trackRef();
84 
85  if (candTrack.isNonnull()) {
86  const float nChi2 = candTrack->normalizedChi2();
87 
88  trkChi2Map.at(i) = MiniFloatConverter::reduceMantissaToNbitsRounding<roundingPrecision>(nChi2);
89  }
90  }
91  }
92 
93  // fill the value maps
94  std::unique_ptr<FloatMap> valueMap = std::make_unique<FloatMap>();
95  FloatMap::Filler filler(*valueMap);
96  filler.insert(candidates, trkChi2Map.begin(), trkChi2Map.end());
97  filler.fill();
98  iEvent.put(std::move(valueMap), "");
99 }
100 
101 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
104  desc.add<edm::InputTag>("candidates", edm::InputTag("packedPFCandidates"))
105  ->setComment("packed candidate input collection");
106  desc.add<edm::InputTag>("trackCollection", edm::InputTag("generalTracks"))->setComment("track input collection");
107  desc.add<bool>("doLostTracks", false);
108  descriptions.add("packedPFCandidateTrackChi2", desc);
109 }
110 
112 using namespace pat;
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< pat::PackedCandidate > PackedCandidateCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const edm::EDGetTokenT< pat::PackedCandidateCollection > candidateToken_
edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > track2LostTrackToken_
edm::EDGetTokenT< edm::Association< reco::PFCandidateCollection > > candidate2PFToken_
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &)
edm::Ref< pat::PackedCandidateCollection > PackedCandidateRef
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
const edm::EDGetTokenT< reco::TrackCollection > trackToken_
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &, const edm::EventSetup &) override
PackedCandidateTrackChi2Producer(const edm::ParameterSet &iConfig)