CMS 3D CMS Logo

UnifiedSCCollectionProducer.cc
Go to the documentation of this file.
1 /*
2 UnifiedSCCollectionProducer:
3 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^
4 
5 Takes as input the cleaned and the uncleaned collection of SC
6 and produces two collections of SC: one with the clean SC, but flagged
7 such that with the algoID value one can identify the SC that are also
8 in the unclean collection and a collection with the unclean only SC.
9 This collection has the algoID enumeration of the SC altered
10 such that:
11 flags = 0 (cleanedOnly) cluster is only in the cleaned collection
12 flags = 100 (common) cluster is common in both collections
13 flags = 200 (uncleanedOnly) cluster is only in the uncleaned collection
14 
15 In that way the user can get hold of objects from the
16 - cleaned collection only if they choose flags < 200
17 - uncleaned collection only if they choose flags >= 100
18 
19 18 Aug 2010
20 Nikolaos Rompotis and Chris Seez - Imperial College London
21 many thanks to David Wardrope, Shahram Rahatlou and Federico Ferri
22 */
23 
43 
44 #include <iostream>
45 #include <memory>
46 #include <vector>
47 
49 public:
51 
52  void produce(edm::Event&, const edm::EventSetup&) override;
53 
54 private:
55  // the clean collection
58  // the uncleaned collection
61 
62  // the names of the products to be produced:
67 };
68 
71 
75  // get the parameters
76  // the cleaned collection:
77  cleanBcCollection_ = consumes<BasicClusterCollection>(ps.getParameter<edm::InputTag>("cleanBcCollection"));
78  cleanScCollection_ = consumes<SuperClusterCollection>(ps.getParameter<edm::InputTag>("cleanScCollection"));
79 
80  // the uncleaned collection
81  uncleanBcCollection_ = consumes<BasicClusterCollection>(ps.getParameter<edm::InputTag>("uncleanBcCollection"));
82  uncleanScCollection_ = consumes<SuperClusterCollection>(ps.getParameter<edm::InputTag>("uncleanScCollection"));
83 
84  // the names of the products to be produced:
85  //
86  // the clean collection: this is as it was before, but labeled
87  bcCollection_ = ps.getParameter<std::string>("bcCollection");
88  scCollection_ = ps.getParameter<std::string>("scCollection");
89  // the unclean only collection: SC unique to the unclean collection
90  bcCollectionUncleanOnly_ = ps.getParameter<std::string>("bcCollectionUncleanOnly");
91  scCollectionUncleanOnly_ = ps.getParameter<std::string>("scCollectionUncleanOnly");
92  // the products:
93  produces<reco::BasicClusterCollection>(bcCollection_);
94  produces<reco::SuperClusterCollection>(scCollection_);
95  produces<reco::BasicClusterCollection>(bcCollectionUncleanOnly_);
96  produces<reco::SuperClusterCollection>(scCollectionUncleanOnly_);
97 }
98 
100  edm::LogInfo("UnifiedSC") << ">>>>> Entering UnifiedSCCollectionProducer <<<<<";
101  // get the input collections
102  // __________________________________________________________________________
103  //
104  // cluster collections:
107  //
110 
111  evt.getByToken(cleanScCollection_, pCleanSC);
112  evt.getByToken(cleanBcCollection_, pCleanBC);
113  evt.getByToken(uncleanBcCollection_, pUncleanBC);
114  evt.getByToken(uncleanScCollection_, pUncleanSC);
115 
116  // the collections to be produced ___________________________________________
117  reco::BasicClusterCollection basicClusters;
119  //
120  reco::BasicClusterCollection basicClustersUncleanOnly;
121  reco::SuperClusterCollection superClustersUncleanOnly;
122  //
123  // run over the uncleaned SC and check how many of them are matched to
124  // the cleaned ones
125  // if you find a matched one, then keep the info that it is matched
126  // along with which clean SC was matched + its basic clusters
127  // if you find an unmatched one, keep the info and store its basic clusters
128  //
129  //
130  int uncleanSize = pUncleanSC->size();
131  int cleanSize = pCleanSC->size();
132 
133  LogTrace("UnifiedSC") << "Size of Clean Collection: " << cleanSize << ", uncleanSize: " << uncleanSize;
134 
135  // keep the indices
136  std::vector<int> inUncleanOnlyInd; // counting the unclean
137  std::vector<int> inCleanInd; // counting the unclean
138  std::vector<int> inCleanOnlyInd; // counting the clean
139  std::vector<DetId> scUncleanSeedDetId; // counting the unclean
140  std::vector<DetId> scCleanSeedDetId; // counting the clean
141  // ontains the index of the SC that owns that BS
142  // first basic cluster index, second: 0 for unclean and 1 for clean
143  std::vector<std::pair<int, int> > basicClusterOwner;
144  std::vector<std::pair<int, int> > basicClusterOwnerUncleanOnly;
145  // if this basic cluster is a seed it is 1
146  std::vector<int> uncleanBasicClusterIsSeed;
147 
148  // loop over unclean SC _____________________________________________________
149  for (int isc = 0; isc < uncleanSize; ++isc) {
150  reco::SuperClusterRef unscRef(pUncleanSC, isc);
151  const std::vector<std::pair<DetId, float> >& uhits = unscRef->hitsAndFractions();
152  int uhitsSize = uhits.size();
153  bool foundTheSame = false;
154  for (int jsc = 0; jsc < cleanSize; ++jsc) { // loop over the cleaned SC
155  reco::SuperClusterRef cscRef(pCleanSC, jsc);
156  const std::vector<std::pair<DetId, float> >& chits = cscRef->hitsAndFractions();
157  int chitsSize = chits.size();
158  foundTheSame = false;
159  if (unscRef->seed()->seed() == cscRef->seed()->seed() && chitsSize == uhitsSize) {
160  // if the clusters are exactly the same then because the clustering
161  // algorithm works in a deterministic way, the order of the rechits
162  // will be the same
163  foundTheSame = true;
164  for (int i = 0; i < chitsSize; ++i) {
165  if (uhits[i].first != chits[i].first) {
166  foundTheSame = false;
167  break;
168  }
169  }
170  }
171  if (foundTheSame) { // ok you have found it:
172  // this supercluster belongs to both collections
173  inUncleanOnlyInd.push_back(0);
174  inCleanInd.push_back(jsc); // keeps the index of the clean SC
175  scUncleanSeedDetId.push_back(unscRef->seed()->seed());
176  //
177  // keep its basic clusters:
178  for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
179  // the basic clusters
180  basicClusters.push_back(**bciter);
181  // index of the unclean SC
182  basicClusterOwner.push_back(std::make_pair(isc, 0));
183  }
184  break; // break the loop over unclean sc
185  }
186  }
187  if (not foundTheSame) { // this SC is only in the unclean collection
188  // mark it as unique in the uncleaned
189  inUncleanOnlyInd.push_back(1);
190  scUncleanSeedDetId.push_back(unscRef->seed()->seed());
191  // keep all its basic clusters
192  for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
193  // the basic clusters
194  basicClustersUncleanOnly.push_back(**bciter);
195  basicClusterOwnerUncleanOnly.push_back(std::make_pair(isc, 0));
196  }
197  }
198  } // loop over the unclean SC _______________________________________________
199  //
200  int inCleanSize = inCleanInd.size();
201  //
202  // loop over the clean SC, check that are not in common with the unclean
203  // ones and then store their SC as before ___________________________________
204  for (int jsc = 0; jsc < cleanSize; ++jsc) {
205  // check whether this index is already in the common collection
206  bool takenAlready = false;
207  for (int j = 0; j < inCleanSize; ++j) {
208  if (jsc == inCleanInd[j]) {
209  takenAlready = true;
210  break;
211  }
212  }
213  if (takenAlready) {
214  inCleanOnlyInd.push_back(0);
215  scCleanSeedDetId.push_back(DetId(0));
216  continue;
217  }
218  inCleanOnlyInd.push_back(1);
219  reco::SuperClusterRef cscRef(pCleanSC, jsc);
220  scCleanSeedDetId.push_back(cscRef->seed()->seed());
221  for (reco::CaloCluster_iterator bciter = cscRef->clustersBegin(); bciter != cscRef->clustersEnd(); ++bciter) {
222  // the basic clusters
223  basicClusters.push_back(**bciter);
224  basicClusterOwner.push_back(std::make_pair(jsc, 1));
225  }
226  } // end loop over clean SC _________________________________________________
227  //
228  //
229 
230  // Final check: in the endcap BC may exist that are not associated to SC,
231  // we need to recover them as well (e.g. multi5x5 algo)
232  // This is should be optimized (SA, 20110621)
233 
234  // loop on original clean BC collection and see if the BC is missing from the new one
235  for (reco::BasicClusterCollection::const_iterator bc = pCleanBC->begin(); bc != pCleanBC->end(); ++bc) {
236  bool foundTheSame = false;
237  for (reco::BasicClusterCollection::const_iterator cleanonly_bc = basicClusters.begin();
238  cleanonly_bc != basicClusters.end();
239  ++cleanonly_bc) {
240  const std::vector<std::pair<DetId, float> >& chits = bc->hitsAndFractions();
241  int chitsSize = chits.size();
242 
243  const std::vector<std::pair<DetId, float> >& uhits = cleanonly_bc->hitsAndFractions();
244  int uhitsSize = uhits.size();
245 
246  if (cleanonly_bc->seed() == bc->seed() && chitsSize == uhitsSize) {
247  foundTheSame = true;
248  for (int i = 0; i < chitsSize; ++i) {
249  if (uhits[i].first != chits[i].first) {
250  foundTheSame = false;
251  break;
252  }
253  }
254  }
255 
256  } // loop on new clean BC collection
257 
258  // clean basic cluster is not associated to SC and does not belong to the
259  // new collection, add it
260  if (!foundTheSame) {
261  basicClusters.push_back(*bc);
262  LogTrace("UnifiedSC") << "found BC to add that was not associated to any SC";
263  }
264 
265  } // loop on original clean BC collection
266 
267  // at this point we have the basic cluster collection ready
268  // Up to index basicClusterOwner.size() we have the BC owned by a SC
269  // The remaining are BCs not owned by a SC
270 
271  int bcSize = (int)basicClusterOwner.size();
272  int bcSizeUncleanOnly = (int)basicClustersUncleanOnly.size();
273 
274  LogTrace("UnifiedSC") << "Found cleaned SC: " << cleanSize << " uncleaned SC: " << uncleanSize;
275  //
276  // export the clusters to the event from the clean clusters
277  auto basicClusters_p = std::make_unique<reco::BasicClusterCollection>();
278  basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
280  if (!(bccHandle.isValid())) {
281  edm::LogWarning("MissingInput") << "could not handle the new BasicClusters!";
282  return;
283  }
284 
285  LogTrace("UnifiedSC") << "Got the BasicClusters from the event again";
286  //
287  // export the clusters to the event: from the unclean only clusters
288  auto basicClustersUncleanOnly_p = std::make_unique<reco::BasicClusterCollection>();
289  basicClustersUncleanOnly_p->assign(basicClustersUncleanOnly.begin(), basicClustersUncleanOnly.end());
291  evt.put(std::move(basicClustersUncleanOnly_p), bcCollectionUncleanOnly_);
292  if (!(bccHandleUncleanOnly.isValid())) {
293  edm::LogWarning("MissingInput") << "could not handle the new BasicClusters (Unclean Only)!";
294  return;
295  }
296  LogTrace("UnifiedSC") << "Got the BasicClusters from the event again (Unclean Only)";
297  //
298 
299  // now we can build the SC collection
300  //
301  // start again from the unclean collection
302  // all the unclean SC will become members of the new collection
303  // with different algoIDs ___________________________________________________
304  for (int isc = 0; isc < uncleanSize; ++isc) {
305  reco::CaloClusterPtrVector clusterPtrVector;
306  // the seed is the basic cluster with the highest energy
308  if (inUncleanOnlyInd[isc] == 1) { // unclean SC Unique in Unclean
309  for (int jbc = 0; jbc < bcSizeUncleanOnly; ++jbc) {
310  std::pair<int, int> theBcOwner = basicClusterOwnerUncleanOnly[jbc];
311  if (theBcOwner.first == isc && theBcOwner.second == 0) {
312  reco::CaloClusterPtr currentClu = reco::CaloClusterPtr(bccHandleUncleanOnly, jbc);
313  clusterPtrVector.push_back(currentClu);
314  if (scUncleanSeedDetId[isc] == currentClu->seed()) {
315  seed = currentClu;
316  }
317  }
318  }
319 
320  } else { // unclean SC common in clean and unclean
321  for (int jbc = 0; jbc < bcSize; ++jbc) {
322  std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
323  if (theBcOwner.first == isc && theBcOwner.second == 0) {
324  reco::CaloClusterPtr currentClu = reco::CaloClusterPtr(bccHandle, jbc);
325  clusterPtrVector.push_back(currentClu);
326  if (scUncleanSeedDetId[isc] == currentClu->seed()) {
327  seed = currentClu;
328  }
329  }
330  }
331  }
332  //std::cout << "before getting the uncl" << std::endl;
333  reco::SuperClusterRef unscRef(pUncleanSC, isc);
334  reco::SuperCluster newSC(unscRef->energy(), unscRef->position(), seed, clusterPtrVector);
335  // now set the algoID for this SC again
336  if (inUncleanOnlyInd[isc] == 1) {
337  // set up the quality to unclean only .............
339  superClustersUncleanOnly.push_back(newSC);
340  } else {
341  // set up the quality to common .............
342  newSC.setFlags(reco::CaloCluster::common);
343  superClusters.push_back(newSC);
344  }
345  // now you can store your SC
346 
347  } // end loop over unclean SC _______________________________________________
348  // flags numbering scheme
349  // flags = 0 = cleanedOnly cluster is only in the cleaned collection
350  // flags = 100 = common cluster is common in both collections
351  // flags = 200 = uncleanedOnly cluster is only in the uncleaned collection
352 
353  // now loop over the clean SC and do the same but now you have to avoid the
354  // the duplicated ones ______________________________________________________
355  for (int jsc = 0; jsc < cleanSize; ++jsc) {
356  //std::cout << "working in cl #" << jsc << std::endl;
357  // check that the SC is not in the unclean collection
358  if (inCleanOnlyInd[jsc] == 0)
359  continue;
360  reco::CaloClusterPtrVector clusterPtrVector;
361  // the seed is the basic cluster with the highest energy
363  for (int jbc = 0; jbc < bcSize; ++jbc) {
364  std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
365  if (theBcOwner.first == jsc && theBcOwner.second == 1) {
366  reco::CaloClusterPtr currentClu = reco::CaloClusterPtr(bccHandle, jbc);
367  clusterPtrVector.push_back(currentClu);
368  if (scCleanSeedDetId[jsc] == currentClu->seed()) {
369  seed = currentClu;
370  }
371  }
372  }
373  reco::SuperClusterRef cscRef(pCleanSC, jsc);
374  reco::SuperCluster newSC(cscRef->energy(), cscRef->position(), seed, clusterPtrVector);
376 
377  // add it to the collection:
378  superClusters.push_back(newSC);
379 
380  } // end loop over clean SC _________________________________________________
381 
382  LogTrace("UnifiedSC") << "New SC collection was created";
383 
384  auto superClusters_p = std::make_unique<reco::SuperClusterCollection>();
385  superClusters_p->assign(superClusters.begin(), superClusters.end());
386 
387  evt.put(std::move(superClusters_p), scCollection_);
388 
389  LogTrace("UnifiedSC") << "Clusters (Basic/Super) added to the Event! :-)";
390 
391  auto superClustersUncleanOnly_p = std::make_unique<reco::SuperClusterCollection>();
392  superClustersUncleanOnly_p->assign(superClustersUncleanOnly.begin(), superClustersUncleanOnly.end());
393 
394  evt.put(std::move(superClustersUncleanOnly_p), scCollectionUncleanOnly_);
395 
396  // ----- debugging ----------
397  // print the new collection SC quantities
398 
399  // print out the clean collection SC
400  LogTrace("UnifiedSC") << "Clean Collection SC ";
401  for (int i = 0; i < cleanSize; ++i) {
402  reco::SuperClusterRef cscRef(pCleanSC, i);
403  LogTrace("UnifiedSC") << " >>> clean #" << i << "; Energy: " << cscRef->energy() << " eta: " << cscRef->eta()
404  << " sc seed detid: " << cscRef->seed()->seed().rawId() << std::endl;
405  }
406  // the unclean SC
407  LogTrace("UnifiedSC") << "Unclean Collection SC ";
408  for (int i = 0; i < uncleanSize; ++i) {
409  reco::SuperClusterRef uscRef(pUncleanSC, i);
410  LogTrace("UnifiedSC") << " >>> unclean #" << i << "; Energy: " << uscRef->energy() << " eta: " << uscRef->eta()
411  << " sc seed detid: " << uscRef->seed()->seed().rawId();
412  }
413  // the new collection
414  LogTrace("UnifiedSC") << "The new SC clean collection with size " << superClusters.size() << std::endl;
415 
416  int new_unclean = 0, new_clean = 0;
417  for (int i = 0; i < (int)superClusters.size(); ++i) {
418  const reco::SuperCluster& nsc = superClusters[i];
419  LogTrace("UnifiedSC") << "SC was got" << std::endl
420  << " ---> energy: " << nsc.energy() << std::endl
421  << " ---> eta: " << nsc.eta() << std::endl
422  << " ---> inClean: " << nsc.isInClean() << std::endl
423  << " ---> id: " << nsc.seed()->seed().rawId() << std::endl
424  << " >>> newSC #" << i << "; Energy: " << nsc.energy() << " eta: " << nsc.eta()
425  << " isClean=" << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
426  << " sc seed detid: " << nsc.seed()->seed().rawId();
427 
428  if (nsc.isInUnclean())
429  ++new_unclean;
430  if (nsc.isInClean())
431  ++new_clean;
432  }
433  LogTrace("UnifiedSC") << "The new SC unclean only collection with size " << superClustersUncleanOnly.size();
434  for (int i = 0; i < (int)superClustersUncleanOnly.size(); ++i) {
435  const reco::SuperCluster nsc = superClustersUncleanOnly[i];
436  LogTrace("UnifiedSC") << " >>> newSC #" << i << "; Energy: " << nsc.energy() << " eta: " << nsc.eta()
437  << " isClean=" << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
438  << " sc seed detid: " << nsc.seed()->seed().rawId();
439  if (nsc.isInUnclean())
440  ++new_unclean;
441  if (nsc.isInClean())
442  ++new_clean;
443  }
444  if ((new_unclean != uncleanSize) || (new_clean != cleanSize)) {
445  LogTrace("UnifiedSC") << ">>>>!!!!!! MISMATCH: new unclean/ old unclean= " << new_unclean << " / " << uncleanSize
446  << ", new clean/ old clean" << new_clean << " / " << cleanSize;
447  }
448 }
void produce(edm::Event &, const edm::EventSetup &) override
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::EDGetTokenT< reco::SuperClusterCollection > uncleanScCollection_
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
edm::Ptr< CaloCluster > CaloClusterPtr
#define LogTrace(id)
void setFlags(uint32_t flags)
Definition: CaloCluster.h:193
bool isInUnclean() const
Definition: CaloCluster.h:198
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
bool isInClean() const
Definition: CaloCluster.h:197
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::BasicClusterCollection > uncleanBcCollection_
edm::EDGetTokenT< reco::SuperClusterCollection > cleanScCollection_
edm::EDGetTokenT< reco::BasicClusterCollection > cleanBcCollection_
double energy() const
cluster energy
Definition: CaloCluster.h:148
Log< level::Info, false > LogInfo
Definition: DetId.h:17
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:79
UnifiedSCCollectionProducer(const edm::ParameterSet &ps)
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:180
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511