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 
14 #include <string>
15 
16 
17 using namespace edm;
18 using namespace std;
19 using namespace reco;
20 
22  conf_(iConfig)
23 {
24  LogInfo("ElectronSeedMerger")<<"Electron SeedMerger started ";
25 
26 
27  ecalSeedToken_ = consumes<ElectronSeedCollection>(iConfig.getParameter<InputTag>("EcalBasedSeeds"));
28  tkSeedToken_ = consumes<ElectronSeedCollection>(iConfig.getParameter<InputTag>("TkBasedSeeds"));
29 
30  produces<ElectronSeedCollection>();
31 
32 }
33 
34 
36 {
37 
38  // do anything here that needs to be done at desctruction time
39  // (e.g. close files, deallocate resources etc.)
40 
41 }
42 
43 
44 //
45 // member functions
46 //
47 
48 // ------------ method called to produce the data ------------
49 void
51 {
52  //CREATE OUTPUT COLLECTION
53  auto output = std::make_unique<ElectronSeedCollection>();
54 
55  //HANDLE THE INPUT SEED COLLECTIONS
57  iEvent.getByToken(ecalSeedToken_,EcalBasedSeeds);
58  ElectronSeedCollection ESeed = *(EcalBasedSeeds.product());
59 
61  iEvent.getByToken(tkSeedToken_,TkBasedSeeds);
62  ElectronSeedCollection TSeed = *(TkBasedSeeds.product());
63 
64 
65  //VECTOR FOR MATCHED SEEDS
66  vector<bool> TSeedMatched;
67  for (unsigned int it=0;it<TSeed.size();it++){
68  TSeedMatched.push_back(false);
69  }
70 
71 
72  //LOOP OVER THE ECAL SEED COLLECTION
73  ElectronSeedCollection::const_iterator e_beg= ESeed.begin();
74  ElectronSeedCollection::const_iterator e_end= ESeed.end();
75  for (;e_beg!=e_end;++e_beg){
76 
77  ElectronSeed NewSeed=*(e_beg);
78  bool AlreadyMatched =false;
79 
80  //LOOP OVER THE TK SEED COLLECTION
81  for (unsigned int it=0;it<TSeed.size();it++){
82  if (AlreadyMatched) continue;
83 
84  //HITS FOR ECAL SEED
85  TrajectorySeed::const_iterator eh = e_beg->recHits().first;
86  TrajectorySeed::const_iterator eh_end = e_beg->recHits().second;
87 
88  //HITS FOR TK SEED
89  unsigned int hitShared=0;
90  unsigned int hitSeed=0;
91  for (;eh!=eh_end;++eh){
92 
93  if (!eh->isValid()) continue;
94  hitSeed++;
95  bool Shared=false;
96  TrajectorySeed::const_iterator th = TSeed[it].recHits().first;
97  TrajectorySeed::const_iterator th_end = TSeed[it].recHits().second;
98  for (;th!=th_end;++th){
99  if (!th->isValid()) continue;
100  //CHECK THE HIT COMPATIBILITY: put back sharesInput
101  // as soon Egamma solves the bug on the seed collection
102  if (eh->sharesInput(&(*th),TrackingRecHit::all)) Shared = true;
103  // if(eh->geographicalId() == th->geographicalId() &&
104 // (eh->localPosition() - th->localPosition()).mag() < 0.001) Shared=true;
105  }
106  if (Shared) hitShared++;
107  }
108  if (hitShared==hitSeed){
109  AlreadyMatched=true;
110  TSeedMatched[it]=true;
111  NewSeed.setCtfTrack(TSeed[it].ctfTrack());
112  }
113  if ( hitShared == (hitSeed-1)){
114  NewSeed.setCtfTrack(TSeed[it].ctfTrack());
115  }
116  }
117 
118  output->push_back(NewSeed);
119  }
120 
121  //FILL THE COLLECTION WITH UNMATCHED TK-BASED SEED
122  for (unsigned int it=0;it<TSeed.size();it++){
123  if (!TSeedMatched[it]) output->push_back(TSeed[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
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
ElectronSeedMerger(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
void setCtfTrack(const CtfTrackRef &)
Set additional info.
Definition: ElectronSeed.cc:35
virtual void produce(edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:230
recHitContainer::const_iterator const_iterator
edm::EDGetTokenT< reco::ElectronSeedCollection > ecalSeedToken_
SEED COLLECTIONS.
std::vector< ElectronSeed > ElectronSeedCollection
collection of ElectronSeed objects
T const * product() const
Definition: Handle.h:81
range recHits() const
fixed size matrix
HLT enums.
edm::EDGetTokenT< reco::ElectronSeedCollection > tkSeedToken_
def move(src, dest)
Definition: eostools.py:510