CMS 3D CMS Logo

EgammaHLTEcalIsolationProducersRegional.cc
Go to the documentation of this file.
1 
8 #include <iostream>
9 #include <vector>
10 #include <memory>
11 
16 
18  : conf_(config)
19  , recoEcalCandidateProducer_ (consumes<reco::RecoEcalCandidateCollection>(conf_.getParameter<edm::InputTag>("recoEcalCandidateProducer")))
20  , bcBarrelProducer_ (consumes<reco::BasicClusterCollection>(conf_.getParameter<edm::InputTag>("bcBarrelProducer")))
21  , bcEndcapProducer_ (consumes<reco::BasicClusterCollection>(conf_.getParameter<edm::InputTag>("bcEndcapProducer")))
22  , scIslandBarrelProducer_ (consumes<reco::SuperClusterCollection>(conf_.getParameter<edm::InputTag>("scIslandBarrelProducer")))
23  , scIslandEndcapProducer_ (consumes<reco::SuperClusterCollection>(conf_.getParameter<edm::InputTag>("scIslandEndcapProducer")))
24  , egEcalIsoEtMin_ (conf_.getParameter<double>("egEcalIsoEtMin"))
25  , egEcalIsoConeSize_ (conf_.getParameter<double>("egEcalIsoConeSize"))
26  , algoType_ (conf_.getParameter<int>("SCAlgoType"))
27  , test_ (new EgammaHLTEcalIsolation(egEcalIsoEtMin_,egEcalIsoConeSize_,algoType_))
28 {
29 
30  //register your products
31  produces < reco::RecoEcalCandidateIsolationMap >();
32 }
33 
35  delete test_;
36 }
37 
40  desc.add<edm::InputTag>("bcBarrelProducer", edm::InputTag(""));
41  desc.add<edm::InputTag>("bcEndcapProducer", edm::InputTag(""));
42  desc.add<edm::InputTag>("scIslandBarrelProducer", edm::InputTag(""));
43  desc.add<edm::InputTag>("scIslandEndcapProducer", edm::InputTag(""));
44  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag(""));
45  desc.add<double>("egEcalIsoEtMin", 0.);
46  desc.add<double>("egEcalIsoConeSize", 0.3);
47  desc.add<int>("SCAlgoType", 1);
48  descriptions.add("hltEgammaHLTEcalIsolationProducersRegional", desc);
49 }
50 
51 
53 
54  // Get the basic cluster collection in the Barrel
56  iEvent.getByToken(bcBarrelProducer_, bcBarrelHandle);
57  const reco::BasicClusterCollection* clusterBarrelCollection = bcBarrelHandle.product();
58  // Get the basic cluster collection in the endcap
60  iEvent.getByToken(bcEndcapProducer_, bcEndcapHandle);
61  const reco::BasicClusterCollection* clusterEndcapCollection = (bcEndcapHandle.product());
62  // Get the Barrel Super Cluster collection
64  iEvent.getByToken(scIslandBarrelProducer_,scBarrelHandle);
65  const reco::SuperClusterCollection* scBarrelCollection = (scBarrelHandle.product());
66  // Get the Endcap Super Cluster collection
68  iEvent.getByToken(scIslandEndcapProducer_,scEndcapHandle);
69  const reco::SuperClusterCollection* scEndcapCollection = (scEndcapHandle.product());
70  // Get the RecoEcalCandidate Collection
72  iEvent.getByToken(recoEcalCandidateProducer_,recoecalcandHandle);
73 
74  std::vector<const reco::BasicCluster*> clusterCollection;
75  for (reco::BasicClusterCollection::const_iterator ibc = clusterBarrelCollection->begin();
76  ibc < clusterBarrelCollection->end(); ibc++ ){clusterCollection.push_back(&(*ibc));}
77  for (reco::BasicClusterCollection::const_iterator iec = clusterEndcapCollection->begin();
78  iec < clusterEndcapCollection->end(); iec++ ){clusterCollection.push_back(&(*iec));}
79  std::vector<const reco::SuperCluster*> scCollection;
80  for (reco::SuperClusterCollection::const_iterator ibsc = scBarrelCollection->begin();
81  ibsc < scBarrelCollection->end(); ibsc++ ){scCollection.push_back(&(*ibsc));}
82  for (reco::SuperClusterCollection::const_iterator iesc = scEndcapCollection->begin();
83  iesc < scEndcapCollection->end(); iesc++ ){scCollection.push_back(&(*iesc));}
84 
86 
87 
88 
89  for (reco::RecoEcalCandidateCollection::const_iterator iRecoEcalCand= recoecalcandHandle->begin(); iRecoEcalCand!=recoecalcandHandle->end(); iRecoEcalCand++) {
90 
91 
92  reco::RecoEcalCandidateRef recoecalcandref(reco::RecoEcalCandidateRef(recoecalcandHandle,iRecoEcalCand -recoecalcandHandle ->begin()));
93 
94 
95  const reco::RecoCandidate *tempiRecoEcalCand = &(*recoecalcandref);
96  float isol = test_->isolPtSum(tempiRecoEcalCand,scCollection, clusterCollection);
97 
98  isoMap.insert(recoecalcandref, isol);
99 
100  }
101 
102  iEvent.put(std::make_unique<reco::RecoEcalCandidateIsolationMap>(isoMap));
103 
104 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
const edm::EDGetTokenT< reco::BasicClusterCollection > bcBarrelProducer_
Definition: config.py:1
const edm::EDGetTokenT< reco::SuperClusterCollection > scIslandEndcapProducer_
int iEvent
Definition: GenABIO.cc:224
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
void produce(edm::StreamID sid, edm::Event &, const edm::EventSetup &) const override
const edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
const edm::EDGetTokenT< reco::SuperClusterCollection > scIslandBarrelProducer_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDGetTokenT< reco::BasicClusterCollection > bcEndcapProducer_
float isolPtSum(const reco::RecoCandidate *recocandidate, const std::vector< const reco::SuperCluster * > &sclusters, const std::vector< const reco::BasicCluster * > &bclusters) const
T const * product() const
Definition: Handle.h:74
void insert(const key_type &k, const data_type &v)
insert an association
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
std::vector< RecoEcalCandidate > RecoEcalCandidateCollection
collectin of RecoEcalCandidate objects
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.