CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHcalPFClusterIsolationProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <vector>
3 #include <memory>
4 
6 
7 // Framework
16 
18 
19 template<typename T1>
21 
22  std::string recoCandidateProducerName = "recoCandidateProducer";
23  if ((typeid(HLTHcalPFClusterIsolationProducer<T1>) == typeid(HLTHcalPFClusterIsolationProducer<reco::RecoEcalCandidate>))) recoCandidateProducerName = "recoEcalCandidateProducer";
24 
25  recoCandidateProducer_ = consumes<T1Collection>(config.getParameter<edm::InputTag>(recoCandidateProducerName));
26  pfClusterProducerHCAL_ = consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducerHCAL"));
27  useHF_ = config.getParameter<bool>("useHF");
28  if (useHF_) {
29  pfClusterProducerHFEM_ = consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducerHFEM"));
30  pfClusterProducerHFHAD_ = consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducerHFHAD"));
31  }
32 
33  drMax_ = config.getParameter<double>("drMax");
34  drVetoBarrel_ = config.getParameter<double>("drVetoBarrel");
35  drVetoEndcap_ = config.getParameter<double>("drVetoEndcap");
36  etaStripBarrel_ = config.getParameter<double>("etaStripBarrel");
37  etaStripEndcap_ = config.getParameter<double>("etaStripEndcap");
38  energyBarrel_ = config.getParameter<double>("energyBarrel");
39  energyEndcap_ = config.getParameter<double>("energyEndcap");
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  produces <T1IsolationMap >();
51 }
52 
53 template<typename T1>
55 {}
56 
57 template<typename T1>
59 
60  std::string recoCandidateProducerName = "recoCandidateProducer";
61  if ((typeid(HLTHcalPFClusterIsolationProducer<T1>) == typeid(HLTHcalPFClusterIsolationProducer<reco::RecoEcalCandidate>))) recoCandidateProducerName = "recoEcalCandidateProducer";
62 
64  desc.add<edm::InputTag>(recoCandidateProducerName, edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
65  desc.add<edm::InputTag>("pfClusterProducerHCAL", edm::InputTag("hltParticleFlowClusterHCAL"));
66  desc.ifValue(edm::ParameterDescription<bool>("useHF", false, true),
67  true >> (edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFEM", edm::InputTag("hltParticleFlowClusterHFEM"), true) and
68  edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFHAD", edm::InputTag("hltParticleFlowClusterHFHAD"), true)) or
69  false >> (edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFEM", edm::InputTag(""), true) and
70  edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFHAD", edm::InputTag(""), true)));
71  desc.add<edm::InputTag>("rhoProducer", edm::InputTag("fixedGridRhoFastjetAllCalo"));
72  desc.add<bool>("doRhoCorrection", false);
73  desc.add<double>("rhoMax", 9.9999999E7);
74  desc.add<double>("rhoScale", 1.0);
75  desc.add<double>("effectiveAreaBarrel", 0.101);
76  desc.add<double>("effectiveAreaEndcap", 0.046);
77  desc.add<double>("drMax", 0.3);
78  desc.add<double>("drVetoBarrel", 0.0);
79  desc.add<double>("drVetoEndcap", 0.0);
80  desc.add<double>("etaStripBarrel", 0.0);
81  desc.add<double>("etaStripEndcap", 0.0);
82  desc.add<double>("energyBarrel", 0.0);
83  desc.add<double>("energyEndcap", 0.0);
84  descriptions.add(std::string("hlt")+std::string(typeid(HLTHcalPFClusterIsolationProducer<T1>).name()), desc);
85 }
86 
87 template<typename T1>
89 
90  edm::Handle<double> rhoHandle;
91  double rho = 0.0;
92  if (doRhoCorrection_) {
93  iEvent.getByToken(rhoProducer_, rhoHandle);
94  rho = *(rhoHandle.product());
95  }
96 
97  if (rho > rhoMax_)
98  rho = rhoMax_;
99 
100  rho = rho*rhoScale_;
101 
102  edm::Handle<T1Collection> recoCandHandle;
103  edm::Handle<reco::PFClusterCollection> clusterHcalHandle;
104  edm::Handle<reco::PFClusterCollection> clusterHfemHandle;
105  edm::Handle<reco::PFClusterCollection> clusterHfhadHandle;
106 
107  iEvent.getByToken(recoCandidateProducer_,recoCandHandle);
108  iEvent.getByToken(pfClusterProducerHCAL_, clusterHcalHandle);
109  const reco::PFClusterCollection* forIsolationHcal = clusterHcalHandle.product();
110 
111  if (useHF_) {
112  iEvent.getByToken(pfClusterProducerHFEM_, clusterHfemHandle);
113  iEvent.getByToken(pfClusterProducerHFHAD_, clusterHfhadHandle);
114  }
115 
116  T1IsolationMap recoCandMap;
117 
118  float dRVeto = -1.;
119  float etaStrip = -1;
120 
121  for (unsigned int iReco = 0; iReco < recoCandHandle->size(); iReco++) {
122  T1Ref candRef(recoCandHandle, iReco);
123 
124  if (fabs(candRef->eta()) < 1.479) {
125  dRVeto = drVetoBarrel_;
126  etaStrip = etaStripBarrel_;
127  } else {
128  dRVeto = drVetoEndcap_;
129  etaStrip = etaStripEndcap_;
130  }
131 
132  float sum = 0;
133 
134  // Loop over the 3 types of PFClusters
135 
136  for(unsigned i=0; i<forIsolationHcal->size(); i++) {
137  const reco::PFCluster& pfclu = (*forIsolationHcal)[i];
138 
139  if (fabs(candRef->eta()) < 1.479) {
140  if (fabs(pfclu.pt()) < energyBarrel_)
141  continue;
142  } else {
143  if (fabs(pfclu.energy()) < energyEndcap_)
144  continue;
145  }
146 
147  float dEta = fabs(candRef->eta() - pfclu.eta());
148  if(dEta < etaStrip) continue;
149 
150  float dR = deltaR(candRef->eta(), candRef->phi(), pfclu.eta(), pfclu.phi());
151  if(dR > drMax_ || dR < dRVeto) continue;
152 
153  sum += pfclu.pt();
154  }
155 
156  if (useHF_) {
157  const reco::PFClusterCollection* forIsolationHfem = clusterHfemHandle.product();
158  const reco::PFClusterCollection* forIsolationHfhad = clusterHfhadHandle.product();
159 
160  for(unsigned i=0; i<forIsolationHfem->size(); i++) {
161  const reco::PFCluster& pfclu = (*forIsolationHfem)[i];
162 
163  if (fabs(candRef->eta()) < 1.479) {
164  if (fabs(pfclu.pt()) < energyBarrel_)
165  continue;
166  } else {
167  if (fabs(pfclu.energy()) < energyEndcap_)
168  continue;
169  }
170 
171  float dEta = fabs(candRef->eta() - pfclu.eta());
172  if(dEta < etaStrip) continue;
173 
174  float dR = deltaR(candRef->eta(), candRef->phi(), pfclu.eta(), pfclu.phi());
175  if(dR > drMax_ || dR < dRVeto) continue;
176 
177  sum += pfclu.pt();
178  }
179 
180  for(unsigned i=0; i<forIsolationHfhad->size(); i++) {
181  const reco::PFCluster& pfclu = (*forIsolationHfhad)[i];
182 
183  if (fabs(candRef->eta()) < 1.479) {
184  if (fabs(pfclu.pt()) < energyBarrel_)
185  continue;
186  } else {
187  if (fabs(pfclu.energy()) < energyEndcap_)
188  continue;
189  }
190 
191  float dEta = fabs(candRef->eta() - pfclu.eta());
192  if(dEta < etaStrip) continue;
193 
194  float dR = deltaR(candRef->eta(), candRef->phi(), pfclu.eta(), pfclu.phi());
195  if(dR > drMax_ || dR < dRVeto) continue;
196 
197  sum += pfclu.pt();
198  }
199  }
200 
201  if (doRhoCorrection_) {
202  if (fabs(candRef->eta()) < 1.479)
203  sum = sum - rho*effectiveAreaBarrel_;
204  else
205  sum = sum - rho*effectiveAreaEndcap_;
206  }
207 
208  recoCandMap.insert(candRef, sum);
209  }
210 
211  std::auto_ptr<T1IsolationMap> mapForEvent(new T1IsolationMap(recoCandMap));
212  iEvent.put(mapForEvent);
213 }
214 
217 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventIDconst &, edm::Timestampconst & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Definition: DDAxes.h:10
double pt() const
transverse momentum, massless approximation
Definition: PFCluster.h:135
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:163
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
virtual void produce(edm::Event &, const edm::EventSetup &)
HLTHcalPFClusterIsolationProducer< reco::RecoEcalCandidate > EgammaHLTHcalPFClusterIsolationProducer
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double energy() const
cluster energy
Definition: PFCluster.h:79
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
ParameterDescriptionNode * ifValue(ParameterDescription< T > const &switchParameter, std::auto_ptr< ParameterDescriptionCases< T > > cases)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T const * product() const
Definition: Handle.h:81
void insert(const key_type &k, const data_type &v)
insert an association
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:166
HLTHcalPFClusterIsolationProducer(const edm::ParameterSet &)
HLTHcalPFClusterIsolationProducer< reco::RecoChargedCandidate > MuonHLTHcalPFClusterIsolationProducer