CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaHLTClusterShapeProducer.cc
Go to the documentation of this file.
1 
10 
12 
17 
19 
20  // use configuration file to setup input/output collection names
21  recoEcalCandidateProducer_ = consumes<reco::RecoEcalCandidateCollection>(conf_.getParameter<edm::InputTag>("recoEcalCandidateProducer"));
22 
23  ecalRechitEBToken_ = consumes<EcalRecHitCollection>(conf_.getParameter< edm::InputTag > ("ecalRechitEB"));
24  ecalRechitEEToken_ = consumes<EcalRecHitCollection>(conf_.getParameter< edm::InputTag > ("ecalRechitEE"));
25 
26  EtaOrIeta_ = conf_.getParameter< bool > ("isIeta");
27 
28  //register your products
29  produces < reco::RecoEcalCandidateIsolationMap >();
30  produces < reco::RecoEcalCandidateIsolationMap >("sigmaIEtaIEta5x5");
31 }
32 
34 {}
35 
37 
39  desc.add<edm::InputTag>(("recoEcalCandidateProducer"), edm::InputTag("hltL1SeededRecoEcalCandidate"));
40  desc.add< edm::InputTag >(("ecalRechitEB"), edm::InputTag("hltEcalRegionalEgammaRecHit","EcalRecHitsEB"));
41  desc.add< edm::InputTag >(("ecalRechitEE"), edm::InputTag("hltEcalRegionalEgammaRecHit","EcalRecHitsEE"));
42  desc.add< bool >(("isIeta"), true);
43  descriptions.add(("hltEgammaHLTClusterShapeProducer"), desc);
44 }
45 
47 
48  // Get the HLT filtered objects
50  iEvent.getByToken(recoEcalCandidateProducer_,recoecalcandHandle);
51 
52  EcalClusterLazyTools lazyTools( iEvent, iSetup, ecalRechitEBToken_, ecalRechitEEToken_ );
54 
57 
58 
59  for(unsigned int iRecoEcalCand = 0; iRecoEcalCand<recoecalcandHandle->size(); iRecoEcalCand++) {
60 
61  reco::RecoEcalCandidateRef recoecalcandref(recoecalcandHandle, iRecoEcalCand);
62 
63  std::vector<float> vCov ;
64  double sigmaee;
65  if (EtaOrIeta_) {
66  vCov = lazyTools.localCovariances( *(recoecalcandref->superCluster()->seed()) );
67  sigmaee = sqrt(vCov[0]);
68  } else {
69  vCov = lazyTools.covariances( *(recoecalcandref->superCluster()->seed()) );
70  sigmaee = sqrt(vCov[0]);
71  double EtaSC = recoecalcandref->eta();
72  if (EtaSC > 1.479) sigmaee = sigmaee - 0.02*(EtaSC - 2.3);
73  }
74 
75  double sigmaee5x5 = sqrt(lazyTools5x5.localCovariances(*(recoecalcandref->superCluster()->seed()) )[0]);
76  clshMap.insert(recoecalcandref, sigmaee);
77  clsh5x5Map.insert(recoecalcandref,sigmaee5x5);
78 
79 
80  }
81 
82 
83 
84  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> clushMap(new reco::RecoEcalCandidateIsolationMap(clshMap));
85  std::auto_ptr<reco::RecoEcalCandidateIsolationMap> clush5x5Map(new reco::RecoEcalCandidateIsolationMap(clsh5x5Map));
86  iEvent.put(clushMap);
87  iEvent.put(clush5x5Map,"sigmaIEtaIEta5x5");
88 }
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::RecoEcalCandidateCollection > recoEcalCandidateProducer_
EgammaHLTClusterShapeProducer(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< EcalRecHitCollection > ecalRechitEEToken_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
T sqrt(T t)
Definition: SSEVec.h:48
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void insert(const key_type &k, const data_type &v)
insert an association
std::vector< float > localCovariances(const reco::BasicCluster &cluster, float w0=4.7)
edm::EDGetTokenT< EcalRecHitCollection > ecalRechitEBToken_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
virtual void produce(edm::Event &, const edm::EventSetup &)