CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTPFNeutralIsolationProducer.cc
Go to the documentation of this file.
1 
8 #include <iostream>
9 #include <vector>
10 #include <memory>
11 
13 
14 // Framework
20 
23 
25 
27 
28  pfCandidateProducer_ = consumes<reco::PFCandidateCollection>(config.getParameter<edm::InputTag>("pfCandidatesProducer"));
29 
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  etaStripBarrel_ = config.getParameter<double>("etaStripBarrel");
36  etaStripEndcap_ = config.getParameter<double>("etaStripEndcap");
37  energyBarrel_ = config.getParameter<double>("energyBarrel");
38  energyEndcap_ = config.getParameter<double>("energyEndcap");
39  pfToUse_ = config.getParameter<int>("pfCandidateType");
40 
41  doRhoCorrection_ = config.getParameter<bool>("doRhoCorrection");
42  if (doRhoCorrection_)
43  rhoProducer_ = consumes<double>(config.getParameter<edm::InputTag>("rhoProducer"));
44 
45  rhoMax_ = config.getParameter<double>("rhoMax");
46  rhoScale_ = config.getParameter<double>("rhoScale");
47  effectiveAreaBarrel_ = config.getParameter<double>("effectiveAreaBarrel");
48  effectiveAreaEndcap_ = config.getParameter<double>("effectiveAreaEndcap");
49 
50  if(useSCRefs_) {
51  produces < reco::RecoEcalCandidateIsolationMap >();
52  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(config.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
53  } else {
54  produces < reco::ElectronIsolationMap >();
55  electronProducer_ = consumes<reco::ElectronCollection>(config.getParameter<edm::InputTag>("electronProducer"));
56  }
57 }
58 
61  desc.add<edm::InputTag>("electronProducer", edm::InputTag("hltEle27WP80PixelMatchElectronsL1SeededPF"));
62  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
63  desc.add<edm::InputTag>("pfCandidatesProducer", edm::InputTag("hltParticleFlowReg"));
64  desc.add<edm::InputTag>("rhoProducer", edm::InputTag("fixedGridRhoFastjetAllCalo"));
65  desc.add<bool>("doRhoCorrection", false);
66  desc.add<double>("rhoMax", 9.9999999E7);
67  desc.add<double>("rhoScale", 1.0);
68  desc.add<double>("effectiveAreaBarrel", 0.101);
69  desc.add<double>("effectiveAreaEndcap", 0.046);
70  desc.add<bool>("useSCRefs", false);
71  desc.add<double>("drMax", 0.3);
72  desc.add<double>("drVetoBarrel", 0.0);
73  desc.add<double>("drVetoEndcap", 0.0);
74  desc.add<double>("etaStripBarrel", 0.0);
75  desc.add<double>("etaStripEndcap", 0.0);
76  desc.add<double>("energyBarrel", 0.0);
77  desc.add<double>("energyEndcap", 0.0);
78  desc.add<int>("pfCandidateType", 5);
79  descriptions.add(("hltEgammaHLTPFNeutralIsolationProducer"), desc);
80 }
81 
83 
84  edm::Handle<double> rhoHandle;
85  double rho = 0.0;
86  if (doRhoCorrection_) {
87  iEvent.getByToken(rhoProducer_, rhoHandle);
88  rho = *(rhoHandle.product());
89  }
90 
91  if (rho > rhoMax_)
92  rho = rhoMax_;
93 
94  rho = rho*rhoScale_;
95 
99 
100  iEvent.getByToken(pfCandidateProducer_, pfHandle);
101  const reco::PFCandidateCollection* forIsolation = pfHandle.product();
102 
104  reco::RecoEcalCandidateIsolationMap recoEcalCandMap;
105 
106  if(useSCRefs_) {
107 
108  iEvent.getByToken(recoEcalCandidateProducer_,recoecalcandHandle);
109 
110  float dRVeto = -1.;
111  float etaStrip = -1;
112 
113  for (unsigned int iReco = 0; iReco < recoecalcandHandle->size(); iReco++) {
114  reco::RecoEcalCandidateRef candRef(recoecalcandHandle, iReco);
115 
116  if (fabs(candRef->eta()) < 1.479) {
117  dRVeto = drVetoBarrel_;
118  etaStrip = etaStripBarrel_;
119  } else {
120  dRVeto = drVetoEndcap_;
121  etaStrip = etaStripEndcap_;
122  }
123 
124  float sum = 0;
125 
126  // Loop over the PFCandidates
127  for(unsigned i=0; i<forIsolation->size(); i++) {
128  const reco::PFCandidate& pfc = (*forIsolation)[i];
129 
130  //require that the PFCandidate is a neutral hadron
131  if (pfc.particleId() == pfToUse_) {
132 
133  if (fabs(candRef->eta()) < 1.479) {
134  if (fabs(pfc.pt()) < energyBarrel_)
135  continue;
136  } else {
137  if (fabs(pfc.energy()) < energyEndcap_)
138  continue;
139  }
140 
141  // Shift the RecoEcalCandidate direction vector according to the PF vertex
142  math::XYZPoint pfvtx = pfc.vertex();
143  math::XYZVector candDirectionWrtVtx(candRef->superCluster()->x() - pfvtx.x(),
144  candRef->superCluster()->y() - pfvtx.y(),
145  candRef->superCluster()->z() - pfvtx.z());
146 
147  float dEta = fabs(candDirectionWrtVtx.Eta() - pfc.momentum().Eta());
148  if(dEta < etaStrip) continue;
149 
150  float dR = deltaR(candDirectionWrtVtx.Eta(), candDirectionWrtVtx.Phi(), pfc.momentum().Eta(), pfc.momentum().Phi());
151  if(dR > drMax_ || dR < dRVeto) continue;
152 
153  sum += pfc.pt();
154  }
155  }
156 
157  if (doRhoCorrection_) {
158  if (fabs(candRef->eta()) < 1.479)
159  sum = sum - rho*effectiveAreaBarrel_;
160  else
161  sum = sum - rho*effectiveAreaEndcap_;
162  }
163 
164  recoEcalCandMap.insert(candRef, sum);
165  }
166 
167  } else {
168 
169  iEvent.getByToken(electronProducer_,electronHandle);
170 
171  float dRVeto = -1.;
172  float etaStrip = -1;
173 
174  for(unsigned int iEl=0; iEl<electronHandle->size(); iEl++) {
175  reco::ElectronRef eleRef(electronHandle, iEl);
176 
177  if (fabs(eleRef->eta()) < 1.479) {
178  dRVeto = drVetoBarrel_;
179  etaStrip = etaStripBarrel_;
180  } else {
181  dRVeto = drVetoEndcap_;
182  etaStrip = etaStripEndcap_;
183  }
184 
185  float sum = 0;
186 
187  // Loop over the PFCandidates
188  for(unsigned i=0; i<forIsolation->size(); i++) {
189  const reco::PFCandidate& pfc = (*forIsolation)[i];
190 
191  //require that the PFCandidate is a neutral hadron
192  if (pfc.particleId() == pfToUse_) {
193 
194  if (fabs(eleRef->eta()) < 1.479) {
195  if (fabs(pfc.pt()) < energyBarrel_)
196  continue;
197  } else {
198  if (fabs(pfc.energy()) < energyEndcap_)
199  continue;
200  }
201 
202  float dEta = fabs(eleRef->eta() - pfc.momentum().Eta());
203  if(dEta < etaStrip) continue;
204 
205  float dR = deltaR(eleRef->eta(), eleRef->phi(), pfc.momentum().Eta(), pfc.momentum().Phi());
206  if(dR > drMax_ || dR < dRVeto) continue;
207 
208  sum += pfc.pt();
209  }
210  }
211 
212  if (doRhoCorrection_) {
213  if (fabs(eleRef->superCluster()->eta()) < 1.479)
214  sum = sum - rho*effectiveAreaBarrel_;
215  else
216  sum = sum - rho*effectiveAreaEndcap_;
217  }
218 
219  eleMap.insert(eleRef, sum);
220  }
221 
222  }
223 
224  if(useSCRefs_){
225  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> mapForEvent(new reco::RecoEcalCandidateIsolationMap(recoEcalCandMap));
226  iEvent.put(mapForEvent);
227  }else{
228  std::auto_ptr<reco::ElectronIsolationMap> mapForEvent(new reco::ElectronIsolationMap(eleMap));
229  iEvent.put(mapForEvent);
230  }
231 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
Definition: DDAxes.h:10
virtual Vector momentum() const
spatial momentum vector
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
virtual double pt() const
transverse momentum
virtual double energy() const
energy
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::PFCandidateCollection > pfCandidateProducer_
edm::EDGetTokenT< reco::ElectronCollection > electronProducer_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
virtual const Point & vertex() const
vertex position (overwritten by PF...)
Definition: PFCandidate.cc:647
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: Handle.h:81
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
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
virtual ParticleType particleId() const
Definition: PFCandidate.h:373