CMS 3D CMS Logo

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