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 
19  // use configuration file to setup input/output collection names
20 
21  bcBarrelProducer_ = consumes<reco::BasicClusterCollection>(conf_.getParameter<edm::InputTag>("bcBarrelProducer"));
22  bcEndcapProducer_ = consumes<reco::BasicClusterCollection>(conf_.getParameter<edm::InputTag>("bcEndcapProducer"));
23  scIslandBarrelProducer_ = consumes<reco::SuperClusterCollection>(conf_.getParameter<edm::InputTag>("scIslandBarrelProducer"));
24  scIslandEndcapProducer_ = consumes<reco::SuperClusterCollection>(conf_.getParameter<edm::InputTag>("scIslandEndcapProducer"));
25  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(conf_.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
26 
27  egEcalIsoEtMin_ = conf_.getParameter<double>("egEcalIsoEtMin");
28  egEcalIsoConeSize_ = conf_.getParameter<double>("egEcalIsoConeSize");
29  algoType_ = conf_.getParameter<int>("SCAlgoType");
30  test_ = new EgammaHLTEcalIsolation(egEcalIsoEtMin_,egEcalIsoConeSize_,algoType_);
31 
32  //register your products
33  produces < reco::RecoEcalCandidateIsolationMap >();
34 }
35 
37  delete test_;
38 }
39 
42  desc.add<edm::InputTag>("bcBarrelProducer", edm::InputTag(""));
43  desc.add<edm::InputTag>("bcEndcapProducer", edm::InputTag(""));
44  desc.add<edm::InputTag>("scIslandBarrelProducer", edm::InputTag(""));
45  desc.add<edm::InputTag>("scIslandEndcapProducer", edm::InputTag(""));
46  desc.add<edm::InputTag>("recoEcalCandidateProducer", edm::InputTag(""));
47  desc.add<double>("egEcalIsoEtMin", 0.);
48  desc.add<double>("egEcalIsoConeSize", 0.3);
49  desc.add<int>("SCAlgoType", 1);
50  descriptions.add("hltEgammaHLTEcalIsolationProducersRegional", desc);
51 }
52 
53 
55 
56  // Get the basic cluster collection in the Barrel
58  iEvent.getByToken(bcBarrelProducer_, bcBarrelHandle);
59  const reco::BasicClusterCollection* clusterBarrelCollection = bcBarrelHandle.product();
60  // Get the basic cluster collection in the endcap
62  iEvent.getByToken(bcEndcapProducer_, bcEndcapHandle);
63  const reco::BasicClusterCollection* clusterEndcapCollection = (bcEndcapHandle.product());
64  // Get the Barrel Super Cluster collection
66  iEvent.getByToken(scIslandBarrelProducer_,scBarrelHandle);
67  const reco::SuperClusterCollection* scBarrelCollection = (scBarrelHandle.product());
68  // Get the Endcap Super Cluster collection
70  iEvent.getByToken(scIslandEndcapProducer_,scEndcapHandle);
71  const reco::SuperClusterCollection* scEndcapCollection = (scEndcapHandle.product());
72  // Get the RecoEcalCandidate Collection
74  iEvent.getByToken(recoEcalCandidateProducer_,recoecalcandHandle);
75 
76  std::vector<const reco::BasicCluster*> clusterCollection;
77  for (reco::BasicClusterCollection::const_iterator ibc = clusterBarrelCollection->begin();
78  ibc < clusterBarrelCollection->end(); ibc++ ){clusterCollection.push_back(&(*ibc));}
79  for (reco::BasicClusterCollection::const_iterator iec = clusterEndcapCollection->begin();
80  iec < clusterEndcapCollection->end(); iec++ ){clusterCollection.push_back(&(*iec));}
81  std::vector<const reco::SuperCluster*> scCollection;
82  for (reco::SuperClusterCollection::const_iterator ibsc = scBarrelCollection->begin();
83  ibsc < scBarrelCollection->end(); ibsc++ ){scCollection.push_back(&(*ibsc));}
84  for (reco::SuperClusterCollection::const_iterator iesc = scEndcapCollection->begin();
85  iesc < scEndcapCollection->end(); iesc++ ){scCollection.push_back(&(*iesc));}
86 
88 
89 
90 
91  for (reco::RecoEcalCandidateCollection::const_iterator iRecoEcalCand= recoecalcandHandle->begin(); iRecoEcalCand!=recoecalcandHandle->end(); iRecoEcalCand++) {
92 
93 
94  reco::RecoEcalCandidateRef recoecalcandref(reco::RecoEcalCandidateRef(recoecalcandHandle,iRecoEcalCand -recoecalcandHandle ->begin()));
95 
96 
97  const reco::RecoCandidate *tempiRecoEcalCand = &(*recoecalcandref);
98  float isol = test_->isolPtSum(tempiRecoEcalCand,scCollection, clusterCollection);
99 
100  isoMap.insert(recoecalcandref, isol);
101 
102  }
103 
104  iEvent.put(std::make_unique<reco::RecoEcalCandidateIsolationMap>(isoMap));
105 
106 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
Definition: config.py:1
void produce(edm::Event &, const edm::EventSetup &) override
int iEvent
Definition: GenABIO.cc:230
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
float isolPtSum(const reco::RecoCandidate *recocandidate, const std::vector< const reco::SuperCluster * > &sclusters, const std::vector< const reco::BasicCluster * > &bclusters)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::BasicClusterCollection > bcEndcapProducer_
T const * product() const
Definition: Handle.h:81
void insert(const key_type &k, const data_type &v)
insert an association
edm::EDGetTokenT< reco::SuperClusterCollection > scIslandBarrelProducer_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
#define begin
Definition: vmac.h:30
edm::EDGetTokenT< reco::SuperClusterCollection > scIslandEndcapProducer_
edm::EDGetTokenT< reco::BasicClusterCollection > bcBarrelProducer_