CMS 3D CMS Logo

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