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 
19 
20 template<typename T1>
22  pfClusterProducerHCAL_ ( consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducerHCAL"))),
23  rhoProducer_ ( consumes<double>(config.getParameter<edm::InputTag>("rhoProducer"))),
24  pfClusterProducerHFEM_ ( consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducerHFEM"))),
25  pfClusterProducerHFHAD_ ( consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducerHFHAD"))),
26  useHF_ ( config.getParameter<bool>("useHF")),
27  drMax_ ( config.getParameter<double>("drMax")),
28  drVetoBarrel_ ( config.getParameter<double>("drVetoBarrel")),
29  drVetoEndcap_ ( config.getParameter<double>("drVetoEndcap")),
30  etaStripBarrel_ ( config.getParameter<double>("etaStripBarrel")),
31  etaStripEndcap_ ( config.getParameter<double>("etaStripEndcap")),
32  energyBarrel_ ( config.getParameter<double>("energyBarrel")),
33  energyEndcap_ ( config.getParameter<double>("energyEndcap")),
34  useEt_ ( config.getParameter<bool>("useEt")),
35  doRhoCorrection_ ( config.getParameter<bool>("doRhoCorrection")),
36  rhoMax_ ( config.getParameter<double>("rhoMax")),
37  rhoScale_ ( config.getParameter<double>("rhoScale")),
38  effectiveAreas_ ( config.getParameter<std::vector<double> >("effectiveAreas")) ,
39  absEtaLowEdges_ ( config.getParameter<std::vector<double> >("absEtaLowEdges")) {
40 
41  if (doRhoCorrection_) {
42  if (absEtaLowEdges_.size() != effectiveAreas_.size())
43  throw cms::Exception("IncompatibleVects") << "absEtaLowEdges and effectiveAreas should be of the same size. \n";
44 
45  if (absEtaLowEdges_.at(0) != 0.0)
46  throw cms::Exception("IncompleteCoverage") << "absEtaLowEdges should start from 0. \n";
47 
48  for (unsigned int aIt = 0; aIt < absEtaLowEdges_.size() - 1; aIt++) {
49  if ( !(absEtaLowEdges_.at( aIt ) < absEtaLowEdges_.at( aIt + 1 )) )
50  throw cms::Exception("ImproperBinning") << "absEtaLowEdges entries should be in increasing order. \n";
51  }
52  }
53 
54  std::string recoCandidateProducerName = "recoCandidateProducer";
55  if ((typeid(HLTHcalPFClusterIsolationProducer<T1>) == typeid(HLTHcalPFClusterIsolationProducer<reco::RecoEcalCandidate>))) recoCandidateProducerName = "recoEcalCandidateProducer";
56  recoCandidateProducer_ = consumes<T1Collection>(config.getParameter<edm::InputTag>(recoCandidateProducerName));
57 
58  produces <T1IsolationMap >();
59 }
60 
61 template<typename T1>
63 {}
64 
65 template<typename T1>
67 
68  std::string recoCandidateProducerName = "recoCandidateProducer";
69  if ((typeid(HLTHcalPFClusterIsolationProducer<T1>) == typeid(HLTHcalPFClusterIsolationProducer<reco::RecoEcalCandidate>))) recoCandidateProducerName = "recoEcalCandidateProducer";
70 
72  desc.add<edm::InputTag>(recoCandidateProducerName, edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
73  desc.add<edm::InputTag>("pfClusterProducerHCAL", edm::InputTag("hltParticleFlowClusterHCAL"));
74  desc.ifValue(edm::ParameterDescription<bool>("useHF", false, true),
75  true >> (edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFEM", edm::InputTag("hltParticleFlowClusterHFEM"), true) and
76  edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFHAD", edm::InputTag("hltParticleFlowClusterHFHAD"), true)) or
77  false >> (edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFEM", edm::InputTag(""), true) and
78  edm::ParameterDescription<edm::InputTag>("pfClusterProducerHFHAD", edm::InputTag(""), true)));
79  desc.add<edm::InputTag>("rhoProducer", edm::InputTag("fixedGridRhoFastjetAllCalo"));
80  desc.add<bool>("doRhoCorrection", false);
81  desc.add<double>("rhoMax", 9.9999999E7);
82  desc.add<double>("rhoScale", 1.0);
83  desc.add<double>("drMax", 0.3);
84  desc.add<double>("drVetoBarrel", 0.0);
85  desc.add<double>("drVetoEndcap", 0.0);
86  desc.add<double>("etaStripBarrel", 0.0);
87  desc.add<double>("etaStripEndcap", 0.0);
88  desc.add<double>("energyBarrel", 0.0);
89  desc.add<double>("energyEndcap", 0.0);
90  desc.add<bool>("useEt", true);
91  desc.add<std::vector<double> >("effectiveAreas", {0.2, 0.25}); // 2016 post-ichep sinEle default
92  desc.add<std::vector<double> >("absEtaLowEdges", {0.0, 1.479}); // Barrel, Endcap
94 }
95 
96 template<typename T1>
98 
99  edm::Handle<double> rhoHandle;
100  double rho = 0.0;
101  if (doRhoCorrection_) {
102  iEvent.getByToken(rhoProducer_, rhoHandle);
103  rho = *(rhoHandle.product());
104  }
105 
106  if (rho > rhoMax_)
107  rho = rhoMax_;
108 
109  rho = rho*rhoScale_;
110 
111  edm::Handle<T1Collection> recoCandHandle;
112 
113  std::vector<edm::Handle<reco::PFClusterCollection>> clusterHandles;
114  edm::Handle<reco::PFClusterCollection> clusterHcalHandle;
115  edm::Handle<reco::PFClusterCollection> clusterHfemHandle;
116  edm::Handle<reco::PFClusterCollection> clusterHfhadHandle;
117 
118  iEvent.getByToken(recoCandidateProducer_,recoCandHandle);
119  iEvent.getByToken(pfClusterProducerHCAL_, clusterHcalHandle);
120  //const reco::PFClusterCollection* forIsolationHcal = clusterHcalHandle.product();
121  clusterHandles.push_back(clusterHcalHandle);
122 
123  if (useHF_) {
124  iEvent.getByToken(pfClusterProducerHFEM_, clusterHfemHandle);
125  clusterHandles.push_back(clusterHfemHandle);
126  iEvent.getByToken(pfClusterProducerHFHAD_, clusterHfhadHandle);
127  clusterHandles.push_back(clusterHfhadHandle);
128  }
129 
130  T1IsolationMap recoCandMap(recoCandHandle);
131  HcalPFClusterIsolation<T1> isoAlgo(drMax_, drVetoBarrel_, drVetoEndcap_, etaStripBarrel_, etaStripEndcap_, energyBarrel_, energyEndcap_, useEt_);
132 
133  for (unsigned int iReco = 0; iReco < recoCandHandle->size(); iReco++) {
134  T1Ref candRef(recoCandHandle, iReco);
135 
136  float sum = isoAlgo.getSum(candRef, clusterHandles);
137 
138  if (doRhoCorrection_) {
139  int iEA = -1;
140  auto cEta = std::abs(candRef->eta());
141  for (int bIt = absEtaLowEdges_.size() - 1; bIt > -1; bIt--) {
142  if ( cEta > absEtaLowEdges_.at(bIt) ) {
143  iEA = bIt;
144  break;
145  }
146  }
147 
148  sum = sum - rho*effectiveAreas_.at(iEA);
149  }
150 
151  recoCandMap.insert(candRef, sum);
152  }
153 
154  std::auto_ptr<T1IsolationMap> mapForEvent(new T1IsolationMap(recoCandMap));
155  iEvent.put(mapForEvent);
156 }
157 
160 
T getParameter(std::string const &) const
std::string defaultModuleLabel()
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
void produce(edm::StreamID sid, edm::Event &, const edm::EventSetup &) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HLTHcalPFClusterIsolationProducer< reco::RecoEcalCandidate > EgammaHLTHcalPFClusterIsolationProducer
ParameterDescriptionBase * add(U const &iLabel, T const &value)
double getSum(const T1Ref candRef, const std::vector< edm::Handle< reco::PFClusterCollection >> &clusterHandles)
edm::EDGetTokenT< T1Collection > recoCandidateProducer_
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
HLTHcalPFClusterIsolationProducer(const edm::ParameterSet &)
HLTHcalPFClusterIsolationProducer< reco::RecoChargedCandidate > MuonHLTHcalPFClusterIsolationProducer