CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTBcHcalIsolationProducersRegional.cc
Go to the documentation of this file.
1 /* \class EgammaHLTBcHcalIsolationProducersRegional
2  *
3  * \author Matteo Sani (UCSD)
4  *
5  */
6 
9 
12 
16 
21 
25 
27 
28  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(config.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
29  caloTowerProducer_ = consumes<CaloTowerCollection>(config.getParameter<edm::InputTag>("caloTowerProducer"));
30  rhoProducer_ = consumes<double>(config.getParameter<edm::InputTag>("rhoProducer"));
31 
32  doRhoCorrection_ = config.getParameter<bool>("doRhoCorrection");
33  rhoMax_ = config.getParameter<double>("rhoMax");
34  rhoScale_ = config.getParameter<double>("rhoScale");
35 
36  etMin_ = config.getParameter<double>("etMin");
37  innerCone_ = config.getParameter<double>("innerCone");
38  outerCone_ = config.getParameter<double>("outerCone");
39  depth_ = config.getParameter<int>("depth");
40  doEtSum_ = config.getParameter<bool>("doEtSum");
41  effectiveAreaBarrel_ = config.getParameter<double>("effectiveAreaBarrel");
42  effectiveAreaEndcap_ = config.getParameter<double>("effectiveAreaEndcap");
43 
44  useSingleTower_ = config.getParameter<bool>("useSingleTower");
45 
47  hcalCfg_.useTowers = true;
50 
52 
53  produces <reco::RecoEcalCandidateIsolationMap>();
54 }
55 
57  delete hcalHelper_;
58 }
59 
61 
63 
64  desc.add<edm::InputTag>(("recoEcalCandidateProducer"), edm::InputTag("hltRecoEcalCandidate"));
65  desc.add<edm::InputTag>(("caloTowerProducer"), edm::InputTag("hltTowerMakerForAll"));
66  desc.add<edm::InputTag>(("rhoProducer"), edm::InputTag("fixedGridRhoFastjetAllCalo"));
67  desc.add<bool>(("doRhoCorrection"), false);
68  desc.add<double>(("rhoMax"), 999999.);
69  desc.add<double>(("rhoScale"), 1.0);
70  desc.add<double>(("etMin"), -1.0);
71  desc.add<double>(("innerCone"), 0);
72  desc.add<double>(("outerCone"), 0.15);
73  desc.add<int>(("depth"), -1);
74  desc.add<bool>(("doEtSum"), false);
75  desc.add<double>(("effectiveAreaBarrel"), 0.021);
76  desc.add<double>(("effectiveAreaEndcap"), 0.040);
77  desc.add<bool>(("useSingleTower"), false);
78  descriptions.add(("hltEgammaHLTBcHcalIsolationProducersRegional"), desc);
79 }
80 
82 
83  // Get the HLT filtered objects
85  iEvent.getByToken(recoEcalCandidateProducer_, recoEcalCandHandle);
86 
87  edm::Handle<CaloTowerCollection> caloTowersHandle;
88  iEvent.getByToken(caloTowerProducer_, caloTowersHandle);
89 
90  edm::Handle<double> rhoHandle;
91  double rho = 0.0;
92 
93  if (doRhoCorrection_) {
94  iEvent.getByToken(rhoProducer_, rhoHandle);
95  rho = *(rhoHandle.product());
96  }
97 
98  if (rho > rhoMax_)
99  rho = rhoMax_;
100 
101  rho = rho*rhoScale_;
102 
103  hcalHelper_->checkSetup(iSetup);
104  hcalHelper_->readEvent(iEvent);
105 
107 
108  for(unsigned int iRecoEcalCand=0; iRecoEcalCand <recoEcalCandHandle->size(); iRecoEcalCand++) {
109 
110  reco::RecoEcalCandidateRef recoEcalCandRef(recoEcalCandHandle, iRecoEcalCand);
111 
112  float isol = 0;
113 
114  std::vector<CaloTowerDetId> towersBehindCluster;
115 
116  if (useSingleTower_)
117  towersBehindCluster = hcalHelper_->hcalTowersBehindClusters(*(recoEcalCandRef->superCluster()));
118 
119  if (doEtSum_) { //calculate hcal isolation excluding the towers behind the cluster which will be used for H for H/E
120  EgammaTowerIsolation isolAlgo(outerCone_, innerCone_, etMin_, depth_, caloTowersHandle.product());
121  if (useSingleTower_)
122  isol = isolAlgo.getTowerEtSum(&(*recoEcalCandRef), &(towersBehindCluster)); // towersBehindCluster are excluded from the isolation sum
123  else
124  isol = isolAlgo.getTowerEtSum(&(*recoEcalCandRef));
125 
126  if (doRhoCorrection_) {
127  if (fabs(recoEcalCandRef->superCluster()->eta()) < 1.442)
128  isol = isol - rho*effectiveAreaBarrel_;
129  else
130  isol = isol - rho*effectiveAreaEndcap_;
131  }
132  } else { //calcuate H for H/E
133  if (useSingleTower_)
134  isol = hcalHelper_->hcalESumDepth1BehindClusters(towersBehindCluster) + hcalHelper_->hcalESumDepth2BehindClusters(towersBehindCluster);
135  else
136  isol = hcalHelper_->hcalESum(*(recoEcalCandRef->superCluster()));
137  }
138 
139  isoMap.insert(recoEcalCandRef, isol);
140  }
141 
142  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> isolMap(new reco::RecoEcalCandidateIsolationMap(isoMap));
143  iEvent.put(isolMap);
144 }
T getParameter(std::string const &) const
double hcalESum(const reco::SuperCluster &, const std::vector< CaloTowerDetId > *excludeTowers=0)
void readEvent(const edm::Event &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: DDAxes.h:10
double hcalESumDepth2BehindClusters(const std::vector< CaloTowerDetId > &towers)
std::vector< CaloTowerDetId > hcalTowersBehindClusters(const reco::SuperCluster &sc)
void checkSetup(const edm::EventSetup &)
double hcalESumDepth1BehindClusters(const std::vector< CaloTowerDetId > &towers)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void insert(const key_type &k, const data_type &v)
insert an association
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T const * product() const
Definition: Handle.h:81
double getTowerEtSum(const reco::Candidate *cand, const std::vector< CaloTowerDetId > *detIdToExclude=0) const
edm::EDGetTokenT< CaloTowerCollection > hcalTowers