1 |
bbetchar |
1.1 |
#include "TopQuarkAnalysis/TopRefTuple/interface/Tuple_Electron.h"
|
2 |
bbetchar |
1.3 |
#include "TopQuarkAnalysis/TopRefTuple/interface/fTypes.h"
|
3 |
bbetchar |
1.1 |
|
4 |
bbetchar |
1.2 |
#include <boost/foreach.hpp>
|
5 |
bbetchar |
1.1 |
|
6 |
bbetchar |
1.2 |
#include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
|
7 |
|
|
#include "DataFormats/PatCandidates/interface/Electron.h"
|
8 |
|
|
#include "DataFormats/VertexReco/interface/VertexFwd.h"
|
9 |
|
|
#include "DataFormats/VertexReco/interface/Vertex.h"
|
10 |
|
|
|
11 |
|
|
|
12 |
|
|
|
13 |
|
|
Tuple_Electron::
|
14 |
|
|
Tuple_Electron(const edm::ParameterSet& conf)
|
15 |
|
|
: electronTag( conf.getParameter<edm::InputTag>("electronTag") ),
|
16 |
|
|
vertexTag( conf.getParameter<edm::InputTag>("vertexTag") ),
|
17 |
|
|
electronIDs( conf.getParameter<std::vector<std::string> >("electronIDs") ),
|
18 |
|
|
prefix( conf.getParameter<std::string>("prefix") )
|
19 |
bbetchar |
1.1 |
{
|
20 |
bbetchar |
1.2 |
produces <bool> ( prefix + "HandleValid");
|
21 |
bbetchar |
1.3 |
produces <std::vector<fTypes::dPolarLorentzV> > ( prefix + "P4" );
|
22 |
bbetchar |
1.6 |
produces <std::vector<fTypes::dPolarLorentzV> > ( prefix + "P4PF" );
|
23 |
bbetchar |
1.2 |
produces <std::vector<int> > ( prefix + "Charge");
|
24 |
|
|
|
25 |
bbetchar |
1.3 |
produces <std::vector<float> > ( prefix + "RelIso" );
|
26 |
|
|
produces <std::vector<float> > ( prefix + "ChIso" );
|
27 |
|
|
produces <std::vector<float> > ( prefix + "NhIso" );
|
28 |
|
|
produces <std::vector<float> > ( prefix + "PhIso" );
|
29 |
|
|
produces <std::vector<float> > ( prefix + "PuIso" );
|
30 |
bbetchar |
1.5 |
produces <std::vector<float> > ( prefix + "EAIso" );
|
31 |
bbetchar |
1.2 |
|
32 |
bbetchar |
1.3 |
produces <std::vector<float> > ( prefix + "SuperClusterEta");
|
33 |
bbetchar |
1.2 |
produces <std::vector<bool> > ( prefix + "PassConversionVeto");
|
34 |
bbetchar |
1.3 |
produces <std::vector<float> > ( prefix + "Dxy" );
|
35 |
bbetchar |
1.2 |
produces <std::vector<unsigned> > (prefix + "GsfTrackInnerHits");
|
36 |
|
|
|
37 |
|
|
BOOST_FOREACH(const std::string& idName, electronIDs ) {
|
38 |
bbetchar |
1.3 |
produces <std::vector<float> > ( prefix + idName );
|
39 |
bbetchar |
1.2 |
}
|
40 |
bbetchar |
1.1 |
}
|
41 |
|
|
|
42 |
bbetchar |
1.2 |
void Tuple_Electron::
|
43 |
|
|
produce(edm::Event &event, const edm::EventSetup&) {
|
44 |
bbetchar |
1.3 |
std::auto_ptr<std::vector<fTypes::dPolarLorentzV> > p4 ( new std::vector<fTypes::dPolarLorentzV>() );
|
45 |
bbetchar |
1.6 |
std::auto_ptr<std::vector<fTypes::dPolarLorentzV> > p4pf ( new std::vector<fTypes::dPolarLorentzV>() );
|
46 |
bbetchar |
1.2 |
std::auto_ptr<std::vector<int> > charge ( new std::vector<int>() );
|
47 |
|
|
|
48 |
bbetchar |
1.3 |
std::auto_ptr<std::vector<float> > relIso( new std::vector<float>() );
|
49 |
|
|
std::auto_ptr<std::vector<float> > chIso ( new std::vector<float>() );
|
50 |
|
|
std::auto_ptr<std::vector<float> > nhIso ( new std::vector<float>() );
|
51 |
|
|
std::auto_ptr<std::vector<float> > phIso ( new std::vector<float>() );
|
52 |
|
|
std::auto_ptr<std::vector<float> > puIso ( new std::vector<float>() );
|
53 |
bbetchar |
1.5 |
std::auto_ptr<std::vector<float> > EAIso ( new std::vector<float>() );
|
54 |
bbetchar |
1.2 |
|
55 |
bbetchar |
1.3 |
std::auto_ptr<std::vector<float> > scEta ( new std::vector<float>() );
|
56 |
bbetchar |
1.5 |
std::auto_ptr<std::vector<bool> > passCV ( new std::vector<bool>() );
|
57 |
bbetchar |
1.3 |
std::auto_ptr<std::vector<float> > dxy ( new std::vector<float>() );
|
58 |
bbetchar |
1.2 |
std::auto_ptr<std::vector<unsigned> > nHits ( new std::vector<unsigned>() );
|
59 |
|
|
|
60 |
bbetchar |
1.3 |
std::map<std::string, std::vector<float>* > idMaps;
|
61 |
bbetchar |
1.2 |
BOOST_FOREACH(const std::string& idName, electronIDs ) {
|
62 |
bbetchar |
1.3 |
idMaps[idName] = new std::vector<float>();
|
63 |
bbetchar |
1.2 |
}
|
64 |
|
|
|
65 |
|
|
typedef edm::View<pat::Electron> els_t;
|
66 |
|
|
edm::Handle<els_t> electrons; edm::Handle<reco::VertexCollection> vertices;
|
67 |
|
|
event.getByLabel(electronTag,electrons); event.getByLabel(vertexTag, vertices);
|
68 |
|
|
|
69 |
|
|
const reco::Vertex* vertex = (vertices.isValid() && vertices->size() ) ? &vertices->at(0) : 0;
|
70 |
|
|
|
71 |
|
|
if( electrons.isValid() ) {
|
72 |
|
|
for(els_t::const_iterator el=electrons->begin(); el!=electrons->end(); el++) {
|
73 |
bbetchar |
1.6 |
pat::Electron::LorentzVector p4_ = el->ecalDrivenMomentum();
|
74 |
|
|
p4->push_back(fTypes::dPolarLorentzV(p4_.pt(), p4_.eta(), p4_.phi(), p4_.mass()));
|
75 |
|
|
p4pf->push_back(fTypes::dPolarLorentzV(el->pt(), el->eta(), el->phi(), el->mass()));
|
76 |
bbetchar |
1.2 |
charge->push_back(el->charge());
|
77 |
|
|
chIso->push_back( el->chargedHadronIso() );
|
78 |
|
|
nhIso->push_back( el->neutralHadronIso() );
|
79 |
|
|
phIso->push_back( el->photonIso() );
|
80 |
|
|
puIso->push_back( el->puChargedHadronIso() );
|
81 |
bbetchar |
1.5 |
EAIso->push_back( el->userIsolation("User1Iso") );
|
82 |
|
|
relIso->push_back( (chIso->back() + std::max(float(0), nhIso->back()+phIso->back()-EAIso->back()) ) / el->pt() );
|
83 |
bbetchar |
1.2 |
|
84 |
|
|
scEta->push_back( el->superCluster()->eta() );
|
85 |
|
|
passCV->push_back( el->passConversionVeto() );
|
86 |
bbetchar |
1.4 |
dxy->push_back( vertex? fabs( el->gsfTrack()->dxy(vertex->position()) ) : -9999.9 );
|
87 |
bbetchar |
1.2 |
nHits->push_back( el->gsfTrack()->trackerExpectedHitsInner().numberOfHits() );
|
88 |
|
|
|
89 |
|
|
BOOST_FOREACH(const std::string& idName, electronIDs) {
|
90 |
|
|
idMaps[idName]->push_back( el->electronID(idName) );
|
91 |
|
|
}
|
92 |
|
|
}
|
93 |
|
|
}
|
94 |
bbetchar |
1.1 |
|
95 |
bbetchar |
1.2 |
event.put( std::auto_ptr<bool> ( new bool(electrons.isValid() ) ), prefix + "HandleValid" );
|
96 |
|
|
event.put( p4, prefix+"P4" );
|
97 |
bbetchar |
1.6 |
event.put( p4pf, prefix+"P4PF" );
|
98 |
bbetchar |
1.2 |
event.put( charge, prefix+"Charge");
|
99 |
|
|
event.put( relIso, prefix+"RelIso");
|
100 |
|
|
event.put( chIso, prefix+"ChIso");
|
101 |
|
|
event.put( nhIso, prefix+"NhIso");
|
102 |
|
|
event.put( phIso, prefix+"PhIso");
|
103 |
|
|
event.put( puIso, prefix+"PuIso");
|
104 |
bbetchar |
1.5 |
event.put( EAIso, prefix+"EAIso");
|
105 |
bbetchar |
1.2 |
event.put( scEta, prefix+"SuperClusterEta");
|
106 |
|
|
event.put( passCV,prefix+"PassConversionVeto");
|
107 |
|
|
event.put( dxy, prefix+"Dxy");
|
108 |
|
|
event.put( nHits, prefix+"GsfTrackInnerHits");
|
109 |
|
|
BOOST_FOREACH(const std::string& idName, electronIDs) {
|
110 |
bbetchar |
1.3 |
event.put( std::auto_ptr<std::vector<float> >(idMaps[idName]), prefix + idName);
|
111 |
bbetchar |
1.2 |
}
|
112 |
bbetchar |
1.1 |
}
|