CMS 3D CMS Logo

AlCaECALRecHitReducer.cc
Go to the documentation of this file.
4 
9 
10 //#define ALLrecHits
11 //#define DEBUG
12 
13 //#define QUICK -> if commented loop over the recHits of the SC and add them to the list of recHits to be saved
14 // comment it if you want a faster module but be sure the window is large enough
15 
20 
23  ebRecHitsToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("ebRecHitsLabel"));
24  eeRecHitsToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("eeRecHitsLabel"));
25  // esRecHitsLabel_ = iConfig.getParameter< edm::InputTag > ("esRecHitsLabel");
26 
27  std::vector<edm::InputTag> srcLabels = iConfig.getParameter<std::vector<edm::InputTag> >("srcLabels");
28  for (auto inputTag = srcLabels.begin(); inputTag != srcLabels.end(); ++inputTag) {
30  }
31 
32  //eleViewToken_ = consumes<edm::View <reco::RecoCandidate> > (iConfig.getParameter< edm::InputTag > ("electronLabel"));
33  photonToken_ = consumes<reco::PhotonCollection>(iConfig.getParameter<edm::InputTag>("photonLabel"));
35  consumes<reco::SuperClusterCollection>(iConfig.getParameter<edm::InputTag>("EESuperClusterCollection"));
36 
37  caloTopologyToken_ = esConsumes<CaloTopology, CaloTopologyRecord>();
38 
39  minEta_highEtaSC_ = iConfig.getParameter<double>("minEta_highEtaSC");
40 
41  alcaBarrelHitsCollection_ = iConfig.getParameter<std::string>("alcaBarrelHitCollection");
42  alcaEndcapHitsCollection_ = iConfig.getParameter<std::string>("alcaEndcapHitCollection");
43  alcaCaloClusterCollection_ = iConfig.getParameter<std::string>("alcaCaloClusterCollection");
44 
45  // alcaPreshowerHitsCollection_ = iConfig.getParameter<std::string>("alcaPreshowerHitCollection");
46 
47  etaSize_ = iConfig.getParameter<int>("etaSize");
48  phiSize_ = iConfig.getParameter<int>("phiSize");
49  // FIXME: minimum size of etaSize_ and phiSize_
50  if (phiSize_ % 2 == 0 || etaSize_ % 2 == 0)
51  edm::LogError("AlCaECALRecHitReducerError") << "Size of eta/phi should be odd numbers";
52 
53  // esNstrips_ = iConfig.getParameter<int> ("esNstrips");
54  // esNcolumns_ = iConfig.getParameter<int> ("esNcolumns");
55 
56  //register your products
57  produces<EBRecHitCollection>(alcaBarrelHitsCollection_);
58  produces<EERecHitCollection>(alcaEndcapHitsCollection_);
59  produces<reco::CaloClusterCollection>(alcaCaloClusterCollection_);
60  // produces< ESRecHitCollection > (alcaPreshowerHitsCollection_) ;
61 }
62 
64 
65 // ------------ method called to produce the data ------------
67  using namespace edm;
68  //using namespace std;
69 
70  const CaloTopology* caloTopology = &iSetup.getData(caloTopologyToken_);
71 
72  // Get Photons
74  iEvent.getByToken(photonToken_, phoHandle);
75 
76  // get RecHits
77  Handle<EBRecHitCollection> barrelRecHitsHandle;
78  iEvent.getByToken(ebRecHitsToken_, barrelRecHitsHandle);
79  const EBRecHitCollection* barrelHitsCollection = barrelRecHitsHandle.product();
80 
81  // get RecHits
82  Handle<EERecHitCollection> endcapRecHitsHandle;
83  iEvent.getByToken(eeRecHitsToken_, endcapRecHitsHandle);
84  const EERecHitCollection* endcapHitsCollection = endcapRecHitsHandle.product();
85 
86  // // get ES RecHits
87  // Handle<ESRecHitCollection> preshowerRecHitsHandle;
88  // iEvent.getByToken(esRecHitsToken_,preshowerRecHitsHandle);
89 
90  // const ESRecHitCollection * preshowerHitsCollection = 0 ;
91  // if (preshowerIsFull)
92  // preshowerHitsCollection = preshowerRecHitsHandle.product () ;
93 
94  // // make a vector to store the used ES rechits:
95  // set<ESDetId> used_strips;
96  // used_strips.clear();
97 
98  // for Z->ele+SC
100  iEvent.getByToken(EESuperClusterToken_, EESCHandle);
101 
102  //Create empty output collections
103  auto miniEBRecHitCollection = std::make_unique<EBRecHitCollection>();
104  auto miniEERecHitCollection = std::make_unique<EERecHitCollection>();
105  // auto miniESRecHitCollection = std::make_unique<ESRecHitCollection>();
106 
107  std::set<DetId> reducedRecHit_EBmap;
108  std::set<DetId> reducedRecHit_EEmap;
109 
110  // std::set< edm::Ref<reco::CaloCluster> > reducedCaloClusters_map;
111 
112  auto reducedCaloClusterCollection = std::make_unique<reco::CaloClusterCollection>();
113 
114  //Photons:
115 #ifdef shervin
116  for (reco::PhotonCollection::const_iterator phoIt = phoHandle->begin(); phoIt != phoHandle->end(); phoIt++) {
117  const reco::SuperCluster& sc = *(phoIt->superCluster());
118 
119  if (phoIt->isEB()) {
120  AddMiniRecHitCollection(sc, reducedRecHit_EBmap, caloTopology);
121  } else { // endcap
122  AddMiniRecHitCollection(sc, reducedRecHit_EEmap, caloTopology);
123  } // end of endcap
124 
128  for (; it != itend; ++it) {
129  reco::CaloCluster caloClus(**it);
130  reducedCaloClusterCollection->push_back(caloClus);
131  }
132  }
133 #endif
134 
136  for (auto iToken = eleViewTokens_.begin(); iToken != eleViewTokens_.end(); iToken++) {
137  iEvent.getByToken(*iToken, eleViewHandle);
138 
139  //Electrons:
140  for (auto eleIt = eleViewHandle->begin(); eleIt != eleViewHandle->end(); eleIt++) {
141  const reco::SuperCluster& sc = *(eleIt->superCluster());
142 
143  if (fabs(sc.eta()) < 1.479) {
144  AddMiniRecHitCollection(sc, reducedRecHit_EBmap, caloTopology);
145  } else { // endcap
146  AddMiniRecHitCollection(sc, reducedRecHit_EEmap, caloTopology);
147  } // end of endcap
148 
151  for (; it != itend; ++it) {
152  reco::CaloCluster caloClus(**it);
153  reducedCaloClusterCollection->push_back(caloClus);
154  }
155  }
156  }
157 
158  //saving recHits for highEta SCs for highEta studies
159  for (reco::SuperClusterCollection::const_iterator SC_iter = EESCHandle->begin(); SC_iter != EESCHandle->end();
160  SC_iter++) {
161  if (fabs(SC_iter->eta()) < minEta_highEtaSC_)
162  continue;
163  AddMiniRecHitCollection(*SC_iter, reducedRecHit_EEmap, caloTopology);
164 
165  const reco::SuperCluster& sc = *(SC_iter);
168  for (; it != itend; ++it) {
169  reco::CaloCluster caloClus(**it);
170  reducedCaloClusterCollection->push_back(caloClus);
171  }
172  }
173 
174  //------------------------------ fill the alcareco reduced recHit collection
175  for (std::set<DetId>::const_iterator itr = reducedRecHit_EBmap.begin(); itr != reducedRecHit_EBmap.end(); itr++) {
176  if (barrelHitsCollection->find(*itr) != barrelHitsCollection->end())
177  miniEBRecHitCollection->push_back(*(barrelHitsCollection->find(*itr)));
178  }
179 
180  for (std::set<DetId>::const_iterator itr = reducedRecHit_EEmap.begin(); itr != reducedRecHit_EEmap.end(); itr++) {
181  if (endcapHitsCollection->find(*itr) != endcapHitsCollection->end())
182  miniEERecHitCollection->push_back(*(endcapHitsCollection->find(*itr)));
183  }
184 
185  //--------------------------------------- Put selected information in the event
186  iEvent.put(std::move(miniEBRecHitCollection), alcaBarrelHitsCollection_);
187  iEvent.put(std::move(miniEERecHitCollection), alcaEndcapHitsCollection_);
188  // iEvent.put(std::move(miniESRecHitCollection),alcaPreshowerHitsCollection_ );
189  iEvent.put(std::move(reducedCaloClusterCollection), alcaCaloClusterCollection_);
190 }
191 
193  std::set<DetId>& reducedRecHitMap,
194  const CaloTopology* caloTopology) const {
195  DetId seed = (sc.seed()->seed());
196  int phiSize = phiSize_, etaSize = etaSize_;
197  if (seed.subdetId() != EcalBarrel) { // if not EB, take a square window
199  phiSize = etaSize;
200  }
201 
202  std::vector<DetId> recHit_window = caloTopology->getWindow(seed, phiSize, etaSize);
203  for (unsigned int i = 0; i < recHit_window.size(); i++) {
204  reducedRecHitMap.insert(recHit_window[i]);
205  }
206 
207  const std::vector<std::pair<DetId, float> >& scHits = sc.hitsAndFractions();
208  for (std::vector<std::pair<DetId, float> >::const_iterator scHit_itr = scHits.begin(); scHit_itr != scHits.end();
209  scHit_itr++) {
210  // the map fills just one time (avoiding double insert of recHits)
211  reducedRecHitMap.insert(scHit_itr->first);
212  }
213 
214  return;
215 }
216 
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
producer
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:88
T const * product() const
Definition: Handle.h:70
std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
Get the neighbors of the given cell in a window of given size.
Definition: CaloTopology.cc:64
Log< level::Error, false > LogError
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:91
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::ESGetToken< CaloTopology, CaloTopologyRecord > caloTopologyToken_
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::SuperClusterCollection > EESuperClusterToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const_iterator end() const
Definition: DetId.h:17
AlCaECALRecHitReducer(const edm::ParameterSet &)
ctor
std::string alcaCaloClusterCollection_
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79
edm::EDGetTokenT< reco::PhotonCollection > photonToken_
iterator find(key_type k)
HLT enums.
std::vector< edm::EDGetTokenT< edm::View< reco::RecoCandidate > > > eleViewTokens_
edm::EDGetTokenT< EcalRecHitCollection > eeRecHitsToken_
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:181
void AddMiniRecHitCollection(const reco::SuperCluster &sc, std::set< DetId > &reducedRecHitMap, const CaloTopology *caloTopology) const
edm::EDGetTokenT< EcalRecHitCollection > ebRecHitsToken_
def move(src, dest)
Definition: eostools.py:511