CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTPFChargedIsolationProducer.cc
Go to the documentation of this file.
1 
9 
15 
17 
18 //#include "DataFormats/TrackReco/interface/TrackFwd.h"
19 //#include "DataFormats/TrackReco/interface/Track.h"
20 //#include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
21 
23 
25 
26  pfCandidateProducer_ = consumes<reco::PFCandidateCollection>(config.getParameter<edm::InputTag>("pfCandidatesProducer"));
27  beamSpotProducer_ = consumes<reco::BeamSpot>(config.getParameter<edm::InputTag>("beamSpotProducer"));
28 
29  useGsfTrack_ = config.getParameter<bool>("useGsfTrack");
30  useSCRefs_ = config.getParameter<bool>("useSCRefs");
31 
32  drMax_ = config.getParameter<double>("drMax");
33  drVetoBarrel_ = config.getParameter<double>("drVetoBarrel");
34  drVetoEndcap_ = config.getParameter<double>("drVetoEndcap");
35  ptMin_ = config.getParameter<double>("ptMin");
36  dzMax_ = config.getParameter<double>("dzMax");
37  dxyMax_ = config.getParameter<double>("dxyMax");
38  pfToUse_ = config.getParameter<int>("pfCandidateType");
39 
40  if(useSCRefs_) {
41  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(config.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
42  produces < reco::RecoEcalCandidateIsolationMap >();
43  } else {
44  electronProducer_ = consumes<reco::ElectronCollection>(config.getParameter<edm::InputTag>("electronProducer"));
45  produces < reco::ElectronIsolationMap >();
46  }
47 }
48 
51  desc.add<edm::InputTag>("electronProducer", edm::InputTag("hltEle27WP80PixelMatchElectronsL1SeededPF"));
52  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
53  desc.add<edm::InputTag>("pfCandidatesProducer", edm::InputTag("hltParticleFlowReg"));
54  desc.add<edm::InputTag>("beamSpotProducer", edm::InputTag("hltOnlineBeamSpot"));
55  desc.add<bool>("useGsfTrack", false);
56  desc.add<bool>("useSCRefs", false);
57  desc.add<double>("drMax", 0.3);
58  desc.add<double>("drVetoBarrel", 0.02);
59  desc.add<double>("drVetoEndcap", 0.02);
60  desc.add<double>("ptMin", 0.0);
61  desc.add<double>("dzMax", 0.2);
62  desc.add<double>("dxyMax", 0.1);
63  desc.add<int>("pfCandidateType", 1);
64  descriptions.add(("hltEgammaHLTPFChargedIsolationProducer"), desc);
65 }
66 
68 
72  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
73 
74  iEvent.getByToken(pfCandidateProducer_, pfHandle);
75  const reco::PFCandidateCollection* forIsolation = pfHandle.product();
76 
79 
80  if(useSCRefs_) {
81 
82  iEvent.getByToken(recoEcalCandidateProducer_, recoEcalCandHandle);
83  iEvent.getByToken(beamSpotProducer_, recoBeamSpotHandle);
84  const reco::BeamSpot::Point& beamSpotPosition = recoBeamSpotHandle->position();
85 
86  float dRveto = -1;
87 
88  for(unsigned int iReco=0; iReco<recoEcalCandHandle->size(); iReco++) {
89  reco::RecoEcalCandidateRef candRef(recoEcalCandHandle, iReco);
90 
91  if (fabs(candRef->eta())<1.479)
92  dRveto = drVetoBarrel_;
93  else
94  dRveto = drVetoEndcap_;
95 
96  // Shift the RecoEcalCandidate direction vector according to the vertex
97  math::XYZVector candDirectionWrtVtx(candRef->superCluster()->x() - beamSpotPosition.x(),
98  candRef->superCluster()->y() - beamSpotPosition.y(),
99  candRef->superCluster()->z() - beamSpotPosition.z());
100 
101  float sum = 0;
102 
103  // Loop over the PFCandidates
104  for(unsigned i=0; i<forIsolation->size(); i++) {
105  const reco::PFCandidate& pfc = (*forIsolation)[i];
106 
107  //require that the PFCandidate is a charged hadron
108  if (pfc.particleId() == pfToUse_) {
109 
110  if(pfc.pt() < ptMin_) continue;
111 
112  float dz = fabs(pfc.trackRef()->dz(beamSpotPosition));
113  if(dz > dzMax_) continue;
114 
115  float dxy = fabs(pfc.trackRef()->dxy(beamSpotPosition));
116  if(fabs(dxy) > dxyMax_) continue;
117 
118  float dR = deltaR(candDirectionWrtVtx.Eta(), candDirectionWrtVtx.Phi(), pfc.momentum().Eta(), pfc.momentum().Phi());
119  if(dR > drMax_ || dR < dRveto) continue;
120 
121  sum += pfc.pt();
122  }
123  }
124 
125  recoEcalCandMap.insert(candRef, sum);
126  }
127 
128  } else {
129 
130  iEvent.getByToken(electronProducer_,electronHandle);
131 
132  float dRveto = -1;
133 
134  for(unsigned int iEl=0; iEl<electronHandle->size(); iEl++) {
135  reco::ElectronRef eleRef(electronHandle, iEl);
136  //const reco::Track* eleTrk = useGsfTrack_ ? &*eleRef->gsfTrack() : &*eleRef->track();
137 
138  if (fabs(eleRef->eta())<1.479)
139  dRveto = drVetoBarrel_;
140  else
141  dRveto = drVetoEndcap_;
142 
143  float sum = 0;
144 
145  // Loop over the PFCandidates
146  for(unsigned i=0; i<forIsolation->size(); i++) {
147  const reco::PFCandidate& pfc = (*forIsolation)[i];
148 
149  //require that the PFCandidate is a charged hadron
150  if (pfc.particleId() == pfToUse_) {
151 
152  if(pfc.pt() < ptMin_) continue;
153 
154  float dz = fabs(pfc.trackRef()->dz(eleRef->vertex()));
155  if(dz > dzMax_) continue;
156 
157  float dxy = fabs(pfc.trackRef()->dxy(eleRef->vertex()));
158  if(fabs(dxy) > dxyMax_) continue;
159 
160  float dR = deltaR(eleRef->eta(), eleRef->phi(), pfc.momentum().Eta(), pfc.momentum().Phi());
161  if(dR > drMax_ || dR < dRveto) continue;
162 
163  sum += pfc.pt();
164  }
165  }
166 
167  eleMap.insert(eleRef, sum);
168  }
169 
170  }
171 
172  if(useSCRefs_){
173  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> mapForEvent(new reco::RecoEcalCandidateIsolationMap(recoEcalCandMap));
174  iEvent.put(mapForEvent);
175  }else{
176  std::auto_ptr<reco::ElectronIsolationMap> mapForEvent(new reco::ElectronIsolationMap(eleMap));
177  iEvent.put(mapForEvent);
178  }
179 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual float pt() const
transverse momentum
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual Vector momentum() const
spatial momentum vector
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
edm::EDGetTokenT< reco::BeamSpot > beamSpotProducer_
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:429
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateProducer_
virtual void produce(edm::Event &, const edm::EventSetup &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::EDGetTokenT< reco::ElectronCollection > electronProducer_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
void insert(const key_type &k, const data_type &v)
insert an association
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T const * product() const
Definition: Handle.h:81
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
virtual ParticleType particleId() const
Definition: PFCandidate.h:355