CMS 3D CMS Logo

ElectronSeedMerger.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PFTracking
4 // Class: ElectronSeedMerger
5 //
6 // Original Author: Michele Pioppi
7 
17 
19 public:
20  explicit ElectronSeedMerger(const edm::ParameterSet&);
21 
22 private:
23  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
24 
27 };
28 
29 using namespace edm;
30 using namespace std;
31 using namespace reco;
32 
34  : ecalSeedToken_{consumes<ElectronSeedCollection>(iConfig.getParameter<InputTag>("EcalBasedSeeds"))} {
35  edm::InputTag tkSeedLabel_ = iConfig.getParameter<InputTag>("TkBasedSeeds");
36  if (!tkSeedLabel_.label().empty())
37  tkSeedToken_ = consumes<ElectronSeedCollection>(tkSeedLabel_);
38 
39  produces<ElectronSeedCollection>();
40 }
41 
42 // ------------ method called to produce the data ------------
44  //HANDLE THE INPUT SEED COLLECTIONS
45  auto const& eSeeds = iEvent.get(ecalSeedToken_);
46 
47  ElectronSeedCollection tSeedsEmpty;
48  auto const& tSeeds = tkSeedToken_.isUninitialized() ? tSeedsEmpty : iEvent.get(tkSeedToken_);
49 
50  //CREATE OUTPUT COLLECTION
51  auto output = std::make_unique<ElectronSeedCollection>();
52  output->reserve(eSeeds.size() + tSeeds.size());
53 
54  //VECTOR FOR MATCHED SEEDS
55  vector<bool> tSeedsMatched(tSeeds.size(), false);
56 
57  //LOOP OVER THE ECAL SEED COLLECTION
58  for (auto newSeed : eSeeds) { // make a copy
59 
60  //LOOP OVER THE TK SEED COLLECTION
61  int it = -1;
62  for (auto const& tSeed : tSeeds) {
63  it++;
64 
65  //HITS FOR TK SEED
66  unsigned int hitShared = 0;
67  unsigned int hitSeed = 0;
68  for (auto const& eh : newSeed.recHits()) {
69  if (!eh.isValid())
70  continue;
71  hitSeed++;
72 
73  for (auto const& th : tSeed.recHits()) {
74  if (!th.isValid())
75  continue;
76  //CHECK THE HIT COMPATIBILITY
77  if (eh.sharesInput(&th, TrackingRecHit::all)) {
78  hitShared++;
79  break;
80  }
81  }
82  }
83  if (hitShared == hitSeed) {
84  tSeedsMatched[it] = true;
85  newSeed.setCtfTrack(tSeed.ctfTrack());
86  break;
87  }
88  if (hitShared == (hitSeed - 1)) {
89  newSeed.setCtfTrack(tSeed.ctfTrack());
90  } else if ((hitShared > 0 || tSeed.nHits() == 0) && !newSeed.isTrackerDriven()) {
91  //try to find hits in the full track
92  unsigned int hitSharedOnTrack = 0;
93  for (auto const& eh : newSeed.recHits()) {
94  if (!eh.isValid())
95  continue;
96  for (auto const* th : tSeed.ctfTrack()->recHits()) {
97  if (!th->isValid())
98  continue;
99  // hits on tracks are not matched : use ::some
100  if (eh.sharesInput(th, TrackingRecHit::some)) {
101  hitSharedOnTrack++;
102  break;
103  }
104  }
105  }
106  if (hitSharedOnTrack == hitSeed) {
107  tSeedsMatched[it] = true;
108  newSeed.setCtfTrack(tSeed.ctfTrack());
109  break;
110  }
111  if (hitSharedOnTrack == (hitSeed - 1)) {
112  newSeed.setCtfTrack(tSeed.ctfTrack());
113  }
114  }
115  }
116 
117  output->push_back(newSeed);
118  }
119 
120  //FILL THE COLLECTION WITH UNMATCHED TK-BASED SEED
121  for (unsigned int it = 0; it < tSeeds.size(); it++) {
122  if (!tSeedsMatched[it])
123  output->push_back(tSeeds[it]);
124  }
125 
126  //PUT THE MERGED COLLECTION IN THE EVENT
127  iEvent.put(std::move(output));
128 }
129 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
ElectronSeedMerger(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:99
std::string const & label() const
Definition: InputTag.h:36
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
int iEvent
Definition: GenABIO.cc:224
std::vector< ElectronSeed > ElectronSeedCollection
collection of ElectronSeed objects
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::ElectronSeedCollection > tkSeedToken_
def move(src, dest)
Definition: eostools.py:511
const edm::EDGetTokenT< reco::ElectronSeedCollection > ecalSeedToken_