CMS 3D CMS Logo

EgammaIsoESDetIdCollectionProducer.cc
Go to the documentation of this file.
2 
5 
8 
9 
11 {
12 
14  consumes<reco::PFCluster::EEtoPSAssociation>(iConfig.getParameter< edm::InputTag > ("eeClusToESMapLabel"));
15 
17  consumes<reco::PFClusterCollection>(iConfig.getParameter<edm::InputTag>("ecalPFClustersLabel"));
18  elesToken_ =
19  consumes<reco::GsfElectronCollection>(iConfig.getParameter< edm::InputTag >("elesLabel"));
20 
21  phosToken_ =
22  consumes<reco::PhotonCollection>(iConfig.getParameter< edm::InputTag >("phosLabel"));
23 
25  consumes<reco::SuperClusterCollection>(iConfig.getParameter< edm::InputTag >("superClustersLabel"));
26 
27  minSCEt_ = iConfig.getParameter<double>("minSCEt");
28  minEleEt_ = iConfig.getParameter<double>("minEleEt");
29  minPhoEt_ = iConfig.getParameter<double>("minPhoEt");
30 
31 
32  interestingDetIdCollection_ = iConfig.getParameter<std::string>("interestingDetIdCollection");
33 
34  maxDR_ = iConfig.getParameter<double>("maxDR");
35 
36  //register your products
37  produces< DetIdCollection > (interestingDetIdCollection_) ;
38 
39 
40 }
41 
42 
44 {
45 
46 }
47 
48 // ------------ method called to produce the data ------------
49 void
51  const edm::EventSetup& iSetup)
52 {
53 
54  // take BasicClusters
56  iEvent.getByToken(superClustersToken_, superClusters);
57 
59  iEvent.getByToken(elesToken_,eles);
60 
62  iEvent.getByToken(phosToken_,phos);
63 
65  iEvent.getByToken(eeClusToESMapToken_,eeClusToESMap);
66 
68  iEvent.getByToken(ecalPFClustersToken_,ecalPFClusters);
69 
70  //Create empty output collections
71  std::vector<DetId> indexToStore;
72  indexToStore.reserve(100);
73 
74  if(eles.isValid() && eeClusToESMap.isValid() && ecalPFClusters.isValid()){
75  for(auto& ele : *eles){
76 
77  float scEt = ele.superCluster()->energy()*std::sin(ele.superCluster()->position().theta());
78  if(scEt > minEleEt_ || ele.et()> minEleEt_) addDetIds(*ele.superCluster(),*ecalPFClusters,*eeClusToESMap,indexToStore);
79  }
80  }
81  if(phos.isValid() && eeClusToESMap.isValid() && ecalPFClusters.isValid()){
82  for(auto& pho : *phos){
83  float scEt = pho.superCluster()->energy()*std::sin(pho.superCluster()->position().theta());
84  if(scEt > minPhoEt_ || pho.et()> minPhoEt_) addDetIds(*pho.superCluster(),*ecalPFClusters,*eeClusToESMap,indexToStore);
85  }
86  }
87  if(superClusters.isValid() && eeClusToESMap.isValid() && ecalPFClusters.isValid()){
88  for(auto& sc : *superClusters){
89  float scEt = sc.energy()*std::sin(sc.position().theta());
90  if(scEt > minSCEt_) addDetIds(sc,*ecalPFClusters,*eeClusToESMap,indexToStore);
91  }
92  }
93 
94  //unify the vector
95  std::sort(indexToStore.begin(),indexToStore.end());
96  std::unique(indexToStore.begin(),indexToStore.end());
97 
98  auto detIdCollection = std::make_unique<DetIdCollection>(indexToStore);
99 
100  iEvent.put(std::move(detIdCollection), interestingDetIdCollection_ );
101 
102 }
103 
104 
105 
106 
107 //find all clusters within dR2<maxDR2_ of supercluster and then save det id's of hits of all es clusters matched to said ecal clusters
108 void
110 {
111 
112  const float scEta = superClus.eta();
113  // if(std::abs(scEta)+maxDR_<1.5) return; //not possible to have a endcap cluster, let alone one with preshower (eta>1.65) so exit without checking further
114  const float scPhi = superClus.phi();
115 
116  const float maxDR2=maxDR_*maxDR_;
117 
118  for (size_t clusNr=0;clusNr<clusters.size();clusNr++){
119  const reco::PFCluster& clus = clusters[clusNr];
120  if(clus.layer()==PFLayer::ECAL_ENDCAP &&
121  reco::deltaR2(scEta,scPhi,clus.eta(),clus.phi())<maxDR2){
122 
123  auto keyVal = std::make_pair(clusNr,edm::Ptr<reco::PFCluster>());
124  const auto esClusters = std::equal_range(eeClusToESMap.begin(),eeClusToESMap.end(),keyVal,
125  [](const reco::PFCluster::EEtoPSAssociation::value_type& rhs, //roll on c++14, auto & lambda 4 evar!
127  bool{return rhs.first<lhs.first;}
128  );
129  // std::cout <<"cluster "<<clus.eta()<<" had "<<std::distance(esClusters.first,esClusters.second)<<" es clusters"<<std::endl;
130  for(auto esIt = esClusters.first;esIt!=esClusters.second;++esIt){
131  // std::cout <<"es clus "<<esIt->second->hitsAndFractions().size()<<std::endl;
132  for(const auto& hitAndFrac : esIt->second->hitsAndFractions()){
133  detIdsToStore.push_back(hitAndFrac.first);
134  }
135  }
136 
137  }//end of endcap & dR check
138  }//end of cluster loop
139 
140 
141 
142 }
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:128
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
virtual void produce(edm::Event &, const edm::EventSetup &) override
producer
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:47
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
edm::EDGetTokenT< reco::GsfElectronCollection > elesToken_
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:166
edm::EDGetTokenT< reco::PFCluster::EEtoPSAssociation > eeClusToESMapToken_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::PFClusterCollection > ecalPFClustersToken_
def unique(seq, keepstr=True)
Definition: tier0.py:24
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Definition: PFCluster.h:50
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:111
bool isValid() const
Definition: HandleBase.h:74
float theta() const
Momentum polar angle. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:166
virtual void beginRun(edm::Run const &, const edm::EventSetup &) override final
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
edm::EDGetTokenT< reco::PhotonCollection > phosToken_
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
edm::EDGetTokenT< reco::SuperClusterCollection > superClustersToken_
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:169
void addDetIds(const reco::SuperCluster &superClus, reco::PFClusterCollection clusters, const reco::PFCluster::EEtoPSAssociation &eeClusToESMap, std::vector< DetId > &detIdsToStore)
def move(src, dest)
Definition: eostools.py:510
Definition: Run.h:42
EgammaIsoESDetIdCollectionProducer(const edm::ParameterSet &)
ctor