CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTEcalPFClusterIsolationProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <vector>
3 #include <memory>
4 
6 
7 // Framework
16 
19 
22 
24 
26 
28 
30 
31  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(config.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
32  pfClusterProducer_ = consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducer"));
33 
34  drMax_ = config.getParameter<double>("drMax");
35  drVetoBarrel_ = config.getParameter<double>("drVetoBarrel");
36  drVetoEndcap_ = config.getParameter<double>("drVetoEndcap");
37  etaStripBarrel_ = config.getParameter<double>("etaStripBarrel");
38  etaStripEndcap_ = config.getParameter<double>("etaStripEndcap");
39  energyBarrel_ = config.getParameter<double>("energyBarrel");
40  energyEndcap_ = config.getParameter<double>("energyEndcap");
41 
42  doRhoCorrection_ = config.getParameter<bool>("doRhoCorrection");
43  if (doRhoCorrection_)
44  rhoProducer_ = consumes<double>(config.getParameter<edm::InputTag>("rhoProducer"));
45 
46  rhoMax_ = config.getParameter<double>("rhoMax");
47  rhoScale_ = config.getParameter<double>("rhoScale");
48  effectiveAreaBarrel_ = config.getParameter<double>("effectiveAreaBarrel");
49  effectiveAreaEndcap_ = config.getParameter<double>("effectiveAreaEndcap");
50 
51  produces < reco::RecoEcalCandidateIsolationMap >();
52 
53 }
54 
56 {}
57 
60  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
61  desc.add<edm::InputTag>("pfClusterProducer", edm::InputTag("hltParticleFlowClusterECAL"));
62  desc.add<edm::InputTag>("rhoProducer", edm::InputTag("fixedGridRhoFastjetAllCalo"));
63  desc.add<bool>("doRhoCorrection", false);
64  desc.add<double>("rhoMax", 9.9999999E7);
65  desc.add<double>("rhoScale", 1.0);
66  desc.add<double>("effectiveAreaBarrel", 0.101);
67  desc.add<double>("effectiveAreaEndcap", 0.046);
68  desc.add<double>("drMax", 0.3);
69  desc.add<double>("drVetoBarrel", 0.0);
70  desc.add<double>("drVetoEndcap", 0.0);
71  desc.add<double>("etaStripBarrel", 0.0);
72  desc.add<double>("etaStripEndcap", 0.0);
73  desc.add<double>("energyBarrel", 0.0);
74  desc.add<double>("energyEndcap", 0.0);
75  descriptions.add(("hltEgammaHLTEcalPFClusterIsolationProducer"), desc);
76 }
77 
79 
80  edm::Handle<double> rhoHandle;
81  double rho = 0.0;
82  if (doRhoCorrection_) {
83  iEvent.getByToken(rhoProducer_, rhoHandle);
84  rho = *(rhoHandle.product());
85  }
86 
87  if (rho > rhoMax_)
88  rho = rhoMax_;
89 
90  rho = rho*rhoScale_;
91 
94 
95  iEvent.getByToken(recoEcalCandidateProducer_,recoecalcandHandle);
96  iEvent.getByToken(pfClusterProducer_, clusterHandle);
97 
99 
100  float dRVeto = -1.;
101  float etaStrip = -1;
102 
103  for (unsigned int iReco = 0; iReco < recoecalcandHandle->size(); iReco++) {
104  reco::RecoEcalCandidateRef candRef(recoecalcandHandle, iReco);
105 
106  if (fabs(candRef->eta()) < 1.479) {
107  dRVeto = drVetoBarrel_;
108  etaStrip = etaStripBarrel_;
109  } else {
110  dRVeto = drVetoEndcap_;
111  etaStrip = etaStripEndcap_;
112  }
113 
114  float sum = 0;
115  for (size_t i=0; i<clusterHandle->size(); i++) {
116  reco::PFClusterRef pfclu(clusterHandle, i);
117 
118  if (fabs(candRef->eta()) < 1.479) {
119  if (fabs(pfclu->pt()) < energyBarrel_)
120  continue;
121  } else {
122  if (fabs(pfclu->energy()) < energyEndcap_)
123  continue;
124  }
125 
126  float dEta = fabs(candRef->eta() - pfclu->eta());
127  if(dEta < etaStrip) continue;
128 
129  float dR = deltaR(candRef->eta(), candRef->phi(), pfclu->eta(), pfclu->phi());
130  if(dR > drMax_ || dR < dRVeto) continue;
131 
132  // Exclude clusters that are part of the candidate
133  bool isCandCluster = false;
134  for (reco::CaloCluster_iterator it = candRef->superCluster()->clustersBegin(); it != candRef->superCluster()->clustersEnd(); ++it) {
135  if ((*it)->seed() == pfclu->seed()) {
136  isCandCluster = true;
137  break;
138  }
139  }
140  if(isCandCluster) continue;
141 
142  sum += pfclu->pt();
143  }
144 
145  if (doRhoCorrection_) {
146  if (fabs(candRef->eta()) < 1.479)
147  sum = sum - rho*effectiveAreaBarrel_;
148  else
149  sum = sum - rho*effectiveAreaEndcap_;
150  }
151 
152  recoEcalCandMap.insert(candRef, sum);
153  }
154 
155  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> mapForEvent(new reco::RecoEcalCandidateIsolationMap(recoEcalCandMap));
156  iEvent.put(mapForEvent);
157 
158 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
Definition: DDAxes.h:10
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::EDGetTokenT< reco::PFClusterCollection > pfClusterProducer_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void insert(const key_type &k, const data_type &v)
insert an association
virtual void produce(edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T const * product() const
Definition: Handle.h:81