CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTEcalPFClusterIsolationProducer.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 
27 
28 template<typename T1>
30  pfClusterProducer_ (consumes<reco::PFClusterCollection>(config.getParameter<edm::InputTag>("pfClusterProducer"))),
31  rhoProducer_ (consumes<double>(config.getParameter<edm::InputTag>("rhoProducer"))),
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  doRhoCorrection_ (config.getParameter<bool>("doRhoCorrection")),
40  rhoMax_ (config.getParameter<double>("rhoMax")),
41  rhoScale_ (config.getParameter<double>("rhoScale")),
42  effectiveAreas_ (config.getParameter<std::vector<double> >("effectiveAreas")),
43  absEtaLowEdges_ (config.getParameter<std::vector<double> >("absEtaLowEdges")) {
44 
45  if (doRhoCorrection_) {
46  if (absEtaLowEdges_.size() != effectiveAreas_.size())
47  throw cms::Exception("IncompatibleVects") << "absEtaLowEdges and effectiveAreas should be of the same size. \n";
48 
49  if (absEtaLowEdges_.at(0) != 0.0)
50  throw cms::Exception("IncompleteCoverage") << "absEtaLowEdges should start from 0. \n";
51 
52  for (unsigned int aIt = 0; aIt < absEtaLowEdges_.size() - 1; aIt++) {
53  if ( !(absEtaLowEdges_.at( aIt ) < absEtaLowEdges_.at( aIt + 1 )) )
54  throw cms::Exception("ImproperBinning") << "absEtaLowEdges entries should be in increasing order. \n";
55  }
56  }
57 
58  std::string recoCandidateProducerName = "recoCandidateProducer";
59  if ((typeid(HLTEcalPFClusterIsolationProducer<T1>) == typeid(HLTEcalPFClusterIsolationProducer<reco::RecoEcalCandidate>))) recoCandidateProducerName = "recoEcalCandidateProducer";
60 
61  recoCandidateProducer_ = consumes<T1Collection>(config.getParameter<edm::InputTag>(recoCandidateProducerName));
62  produces <T1IsolationMap>();
63 
64 }
65 
66 template<typename T1>
68 {}
69 
70 template<typename T1>
72 
73  std::string recoCandidateProducerName = "recoCandidateProducer";
74  if ((typeid(HLTEcalPFClusterIsolationProducer<T1>) == typeid(HLTEcalPFClusterIsolationProducer<reco::RecoEcalCandidate>))) recoCandidateProducerName = "recoEcalCandidateProducer";
75 
77  desc.add<edm::InputTag>(recoCandidateProducerName, edm::InputTag("hltL1SeededRecoEcalCandidatePF"));
78  desc.add<edm::InputTag>("pfClusterProducer", edm::InputTag("hltParticleFlowClusterECAL"));
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<std::vector<double> >("effectiveAreas", {0.29, 0.21}); // 2016 post-ichep sinEle default
91  desc.add<std::vector<double> >("absEtaLowEdges", {0.0, 1.479}); // Barrel, Endcap
93 }
94 
95 template<typename T1>
97 
98  edm::Handle<double> rhoHandle;
99  double rho = 0.0;
100  if (doRhoCorrection_) {
101  iEvent.getByToken(rhoProducer_, rhoHandle);
102  rho = *(rhoHandle.product());
103  }
104 
105  if (rho > rhoMax_)
106  rho = rhoMax_;
107 
108  rho = rho*rhoScale_;
109 
110  edm::Handle<T1Collection> recoCandHandle;
112 
113  iEvent.getByToken(recoCandidateProducer_,recoCandHandle);
114  iEvent.getByToken(pfClusterProducer_, clusterHandle);
115 
116  EcalPFClusterIsolation<T1> isoAlgo(drMax_, drVetoBarrel_, drVetoEndcap_, etaStripBarrel_, etaStripEndcap_, energyBarrel_, energyEndcap_);
117  T1IsolationMap recoCandMap(recoCandHandle);
118 
119  for (unsigned int iReco = 0; iReco < recoCandHandle->size(); iReco++) {
120  T1Ref candRef(recoCandHandle, iReco);
121 
122  float sum = isoAlgo.getSum(candRef, clusterHandle);
123 
124  if (doRhoCorrection_) {
125  int iEA = -1;
126  auto cEta = std::abs(candRef->eta());
127  for (int bIt = absEtaLowEdges_.size() - 1; bIt > -1; bIt--) {
128  if ( cEta > absEtaLowEdges_.at(bIt) ) {
129  iEA = bIt;
130  break;
131  }
132  }
133 
134  sum = sum - rho*effectiveAreas_.at(iEA);
135  }
136 
137  recoCandMap.insert(candRef, sum);
138  }
139 
140  std::auto_ptr<T1IsolationMap> mapForEvent(new T1IsolationMap(recoCandMap));
141  iEvent.put(mapForEvent);
142 }
143 
146 
T getParameter(std::string const &) const
std::string defaultModuleLabel()
edm::EDGetTokenT< T1Collection > recoCandidateProducer_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int iEvent
Definition: GenABIO.cc:230
void produce(edm::Event &, const edm::EventSetup &) override
HLTEcalPFClusterIsolationProducer< reco::RecoChargedCandidate > MuonHLTEcalPFClusterIsolationProducer
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLTEcalPFClusterIsolationProducer(const edm::ParameterSet &)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
T const * product() const
Definition: Handle.h:81
void insert(const key_type &k, const data_type &v)
insert an association
double getSum(T1Ref, edm::Handle< std::vector< reco::PFCluster > >)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
HLTEcalPFClusterIsolationProducer< reco::RecoEcalCandidate > EgammaHLTEcalPFClusterIsolationProducer