CMS 3D CMS Logo

L1AnalysisRecoElectron.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 
6 {
7 }
8 
9 
11 {
12 }
13 
15  const edm::EventSetup& setup,
16  //const edm::Handle<edm::View<reco::GsfElectron>>& electrons,
18  std::vector<edm::Handle<edm::ValueMap<bool> > > eleVIDDecisionHandles,
19  const unsigned& maxElectron)
20 
21 {
22 
23  recoElectron_.nElectrons=0;
24 
25  for(reco::GsfElectronCollection::const_iterator el=electrons->begin();
26  el!=electrons->end() && recoElectron_.nElectrons < maxElectron;
27  ++el) {
28 
29  recoElectron_.e.push_back(el->energy());
30  recoElectron_.pt.push_back(el->pt());
31  recoElectron_.et.push_back(el->et());
32  recoElectron_.eta.push_back(el->eta());
33  recoElectron_.phi.push_back(el->phi());
34  recoElectron_.eta_SC.push_back((el->superClusterPosition()).eta());
35  recoElectron_.phi_SC.push_back((el->superClusterPosition()).phi());
36  recoElectron_.e_ECAL.push_back(el->ecalEnergy());
37  recoElectron_.e_SC.push_back(el->superCluster()->energy());
38  recoElectron_.charge.push_back(el->charge());
39 
40  edm::Ref<reco::GsfElectronCollection> electronEdmRef(electrons,recoElectron_.nElectrons);
41 
42  recoElectron_.isVetoElectron.push_back( (*(eleVIDDecisionHandles[0]))[electronEdmRef] );
43  recoElectron_.isLooseElectron.push_back( (*(eleVIDDecisionHandles[1]))[electronEdmRef] );
44  recoElectron_.isMediumElectron.push_back( (*(eleVIDDecisionHandles[2]))[electronEdmRef] );
45  recoElectron_.isTightElectron.push_back( (*(eleVIDDecisionHandles[3]))[electronEdmRef] );
46 
47  double iso = (el->pfIsolationVariables().sumChargedHadronPt + max(
48  el->pfIsolationVariables().sumNeutralHadronEt +
49  el->pfIsolationVariables().sumPhotonEt -
50  0.5 * el->pfIsolationVariables().sumPUPt, 0.0)) / el->pt();
51 
52  recoElectron_.iso.push_back(iso);
53  // cout<<ConversionTools::hasMatchedConversion(*el,conversions,theBeamSpot->position())<<endl;
54 
55  recoElectron_.nElectrons++;
56 
57  }
58  // for(reco::GsfElectronCollection::const_iterator it=electrons->begin();
59  // it!=electrons->end() && recoElectron_.nElectrons < maxElectron;
60  // ++it) {
61 
62  // recoElectron_.e.push_back(it->energy());
63  // recoElectron_.pt.push_back(it->pt());
64  // recoElectron_.et.push_back(it->et());
65  // recoElectron_.eta.push_back(it->eta());
66  // recoElectron_.phi.push_back(it->phi());
67 
68  // // cout<<it->superCluster().position().eta()<<endl;
69  // isVetoElectronCustom(*it, vertices, conversions, theBeamSpot, Rho);
70 
71  // // cout<<ConversionTools::hasMatchedConversion(*it,conversions,theBeamSpot->position())<<endl;
72 
73  // recoElectron_.nElectrons++;
74 
75  // }
76 }
77 
78 
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
void SetElectron(const edm::Event &event, const edm::EventSetup &setup, const edm::Handle< reco::GsfElectronCollection > electrons, const std::vector< edm::Handle< edm::ValueMap< bool > > > eleVIDDecisionHandles, const unsigned &maxElectron)
Definition: event.py:1