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 handle clean super clusters";
103  return;
104  }
105 
106  evt.getByLabel(uncleanScCollection_, pUncleanSC);
107  if (!(pUncleanSC.isValid()))
108  {
109 
110  edm::LogError("MissingInput")<< "could not handle unclean super clusters!" ;
111  return;
112  }
113 
114  // the collections to be produced ___________________________________________
115  reco::BasicClusterCollection basicClusters;
116  reco::SuperClusterCollection superClusters;
117  //
118  reco::BasicClusterCollection basicClustersUncleanOnly;
119  reco::SuperClusterCollection superClustersUncleanOnly;
120  //
121  // run over the uncleaned SC and check how many of them are matched to
122  // the cleaned ones
123  // if you find a matched one, then keep the info that it is matched
124  // along with which clean SC was matched + its basic clusters
125  // if you find an unmatched one, keep the info and store its basic clusters
126  //
127  //
128  int uncleanSize = pUncleanSC->size();
129  int cleanSize = pCleanSC->size();
130 
131  LogDebug("UnifiedSC") << "Size of Clean Collection: " << cleanSize
132  << ", uncleanSize: " << uncleanSize;
133 
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  for (int i=0; i< chitsSize; ++i) {
164  if (uhits[i].first != chits[i].first ) { break;}
165  }
166  foundTheSame = true;
167  }
168  if (foundTheSame) { // ok you have found it:
169  // this supercluster belongs to both collections
170  inUncleanOnlyInd.push_back(0);
171  inCleanInd.push_back(jsc); // keeps the index of the clean SC
172  scUncleanSeedDetId.push_back(unscRef->seed()->seed());
173  //
174  // keep its basic clusters:
175  for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
176  // the basic clusters
177  basicClusters.push_back(**bciter);
178  // index of the unclean SC
179  basicClusterOwner.push_back( std::make_pair(isc,0) );
180  }
181  break; // break the loop over unclean sc
182  }
183  }
184  if (not foundTheSame) { // this SC is only in the unclean collection
185  // mark it as unique in the uncleaned
186  inUncleanOnlyInd.push_back(1);
187  scUncleanSeedDetId.push_back(unscRef->seed()->seed());
188  // keep all its basic clusters
189  for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
190  // the basic clusters
191  basicClustersUncleanOnly.push_back(**bciter);
192  basicClusterOwnerUncleanOnly.push_back( std::make_pair(isc,0) );
193  }
194  }
195  } // loop over the unclean SC _______________________________________________
196  //
197  int inCleanSize = inCleanInd.size();
198  //
199  // loop over the clean SC, check that are not in common with the unclean
200  // ones and then store their SC as before ___________________________________
201  for (int jsc =0; jsc< cleanSize; ++jsc) {
202  // check whether this index is already in the common collection
203  bool takenAlready = false;
204  for (int j=0; j< inCleanSize; ++j) {
205  if (jsc == inCleanInd[j]) { takenAlready = true ;break;}
206  }
207  if (takenAlready) {
208  inCleanOnlyInd.push_back(0);
209  scCleanSeedDetId.push_back(DetId(0));
210  continue;
211  }
212  inCleanOnlyInd.push_back(1);
213  reco::SuperClusterRef cscRef( pCleanSC, jsc );
214  scCleanSeedDetId.push_back(cscRef->seed()->seed());
215  for (reco::CaloCluster_iterator bciter = cscRef->clustersBegin(); bciter != cscRef->clustersEnd(); ++bciter) {
216  // the basic clusters
217  basicClusters.push_back(**bciter);
218  basicClusterOwner.push_back( std::make_pair(jsc,1) );
219  }
220  } // end loop over clean SC _________________________________________________
221  //
222  //
223  // at this point we have the basic cluster collection ready
224  //
225  int bcSize = (int) basicClusters.size();
226  int bcSizeUncleanOnly = (int) basicClustersUncleanOnly.size();
227 
228  LogDebug("UnifiedSC") << "Found cleaned SC: " << cleanSize
229  << " uncleaned SC: " << uncleanSize ;
230  //
231  // export the clusters to the event from the clean clusters
232  std::auto_ptr< reco::BasicClusterCollection>
233  basicClusters_p(new reco::BasicClusterCollection);
234  basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
236  evt.put(basicClusters_p, bcCollection_);
237  if (!(bccHandle.isValid())) {
238 
239  edm::LogWarning("MissingInput")<< "could not handle the new BasicClusters!";
240  return;
241  }
242  reco::BasicClusterCollection basicClustersProd = *bccHandle;
243 
244  LogDebug("UnifiedSC")<< "Got the BasicClusters from the event again" ;
245  //
246  // export the clusters to the event: from the unclean only clusters
247  std::auto_ptr< reco::BasicClusterCollection>
248  basicClustersUncleanOnly_p(new reco::BasicClusterCollection);
249  basicClustersUncleanOnly_p->assign(basicClustersUncleanOnly.begin(),
250  basicClustersUncleanOnly.end());
252  evt.put(basicClustersUncleanOnly_p, bcCollectionUncleanOnly_);
253  if (!(bccHandleUncleanOnly.isValid())) {
254 
255  edm::LogWarning("MissingInput")<< "could not handle the new BasicClusters (Unclean Only)!" ;
256  return;
257  }
258  reco::BasicClusterCollection basicClustersUncleanOnlyProd = *bccHandleUncleanOnly;
259  LogDebug("UnifiedSC")<< "Got the BasicClusters from the event again (Unclean Only)" ;
260  //
261 
262  // now we can build the SC collection
263  //
264  // start again from the unclean collection
265  // all the unclean SC will become members of the new collection
266  // with different algoIDs ___________________________________________________
267  for (int isc=0; isc< uncleanSize; ++isc) {
268  reco::CaloClusterPtrVector clusterPtrVector;
269  // the seed is the basic cluster with the highest energy
270  reco::CaloClusterPtr seed;
271  if (inUncleanOnlyInd[isc] == 1) { // unclean SC Unique in Unclean
272  for (int jbc=0; jbc< bcSizeUncleanOnly; ++jbc) {
273  std::pair<int, int> theBcOwner = basicClusterOwnerUncleanOnly[jbc];
274  if (theBcOwner.first == isc && theBcOwner.second == 0) {
275  reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandleUncleanOnly,jbc);
276  clusterPtrVector.push_back(currentClu);
277  if (scUncleanSeedDetId[isc] == currentClu->seed()) {
278  seed = currentClu;
279  }
280  }
281  }
282 
283  }
284  else { // unclean SC common in clean and unclean
285  for (int jbc=0; jbc< bcSize; ++jbc) {
286  std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
287  if (theBcOwner.first == isc && theBcOwner.second == 0) {
288  reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
289  clusterPtrVector.push_back(currentClu);
290  if (scUncleanSeedDetId[isc] == currentClu->seed()) {
291  seed = currentClu;
292  }
293  }
294  }
295  }
296  //std::cout << "before getting the uncl" << std::endl;
297  reco::SuperClusterRef unscRef( pUncleanSC, isc );
298  reco::SuperCluster newSC(unscRef->energy(), unscRef->position(),
299  seed, clusterPtrVector );
300  // now set the algoID for this SC again
301  if (inUncleanOnlyInd[isc] == 1) {
302  // set up the quality to unclean only .............
304  superClustersUncleanOnly.push_back(newSC);
305  }
306  else {
307  // set up the quality to common .............
308  newSC.setFlags(reco::CaloCluster::common);
309  superClusters.push_back(newSC);
310  }
311  // now you can store your SC
312 
313  } // end loop over unclean SC _______________________________________________
314  // flags numbering scheme
315  // flags = 0 = cleanedOnly cluster is only in the cleaned collection
316  // flags = 100 = common cluster is common in both collections
317  // flags = 200 = uncleanedOnly cluster is only in the uncleaned collection
318 
319  // now loop over the clean SC and do the same but now you have to avoid the
320  // the duplicated ones ______________________________________________________
321  for (int jsc=0; jsc< cleanSize; ++jsc) {
322  //std::cout << "working in cl #" << jsc << std::endl;
323  // check that the SC is not in the unclean collection
324  if (inCleanOnlyInd[jsc] == 0) continue;
325  reco::CaloClusterPtrVector clusterPtrVector;
326  // the seed is the basic cluster with the highest energy
327  reco::CaloClusterPtr seed;
328  for (int jbc=0; jbc< bcSize; ++jbc) {
329  std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
330  if (theBcOwner.first == jsc && theBcOwner.second == 1) {
331  reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
332  clusterPtrVector.push_back(currentClu);
333  if (scCleanSeedDetId[jsc] == currentClu->seed()) {
334  seed = currentClu;
335  }
336  }
337  }
338  reco::SuperClusterRef cscRef( pCleanSC, jsc );
339  reco::SuperCluster newSC(cscRef->energy(), cscRef->position(),
340  seed, clusterPtrVector );
342 
343  // add it to the collection:
344  superClusters.push_back(newSC);
345 
346  } // end loop over clean SC _________________________________________________
347 
348 
349  LogDebug("UnifiedSC")<< "New SC collection was created";
350 
351  std::auto_ptr< reco::SuperClusterCollection>
352  superClusters_p(new reco::SuperClusterCollection);
353  superClusters_p->assign(superClusters.begin(), superClusters.end());
354 
355  evt.put(superClusters_p, scCollection_);
356 
357  LogDebug("UnifiedSC") << "Clusters (Basic/Super) added to the Event! :-)";
358 
359  std::auto_ptr< reco::SuperClusterCollection>
360  superClustersUncleanOnly_p(new reco::SuperClusterCollection);
361  superClustersUncleanOnly_p->assign(superClustersUncleanOnly.begin(),
362  superClustersUncleanOnly.end());
363 
364  evt.put(superClustersUncleanOnly_p, scCollectionUncleanOnly_);
365 
366  // ----- debugging ----------
367  // print the new collection SC quantities
368 
369  // print out the clean collection SC
370  LogDebug("UnifiedSC") << "Clean Collection SC ";
371  for (int i=0; i < cleanSize; ++i) {
372  reco::SuperClusterRef cscRef( pCleanSC, i );
373  LogDebug("UnifiedSC") << " >>> clean #" << i << "; Energy: " << cscRef->energy()
374  << " eta: " << cscRef->eta()
375  << " sc seed detid: " << cscRef->seed()->seed().rawId()
376  << std::endl;
377  }
378  // the unclean SC
379  LogDebug("UnifiedSC") << "Unclean Collection SC ";
380  for (int i=0; i < uncleanSize; ++i) {
381  reco::SuperClusterRef uscRef( pUncleanSC, i );
382  LogDebug("UnifiedSC") << " >>> unclean #" << i << "; Energy: " << uscRef->energy()
383  << " eta: " << uscRef->eta()
384  << " sc seed detid: " << uscRef->seed()->seed().rawId();
385  }
386  // the new collection
387  LogDebug("UnifiedSC")<< "The new SC clean collection with size "<< superClusters.size() << std::endl;
388 
389  int new_unclean = 0, new_clean=0;
390  for (int i=0; i < (int) superClusters.size(); ++i) {
391  const reco::SuperCluster & nsc = superClusters[i];
392  LogDebug("UnifiedSC") << "SC was got" << std::endl
393  << " ---> energy: " << nsc.energy() << std::endl
394  << " ---> eta: " << nsc.eta() << std::endl
395  << " ---> inClean: " << nsc.isInClean() << std::endl
396  << " ---> id: " << nsc.seed()->seed().rawId() << std::endl
397  << " >>> newSC #" << i << "; Energy: " << nsc.energy()
398  << " eta: " << nsc.eta() << " isClean="
399  << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
400  << " sc seed detid: " << nsc.seed()->seed().rawId();
401 
402  if (nsc.isInUnclean()) ++new_unclean;
403  if (nsc.isInClean()) ++new_clean;
404  }
405  LogDebug("UnifiedSC")<< "The new SC unclean only collection with size "<< superClustersUncleanOnly.size();
406  for (int i=0; i < (int) superClustersUncleanOnly.size(); ++i) {
407  const reco::SuperCluster nsc = superClustersUncleanOnly[i];
408  LogDebug ("UnifiedSC") << " >>> newSC #" << i << "; Energy: " << nsc.energy()
409  << " eta: " << nsc.eta() << " isClean="
410  << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
411  << " sc seed detid: " << nsc.seed()->seed().rawId();
412  if (nsc.isInUnclean()) ++new_unclean;
413  if (nsc.isInClean()) ++new_clean;
414  }
415  if ( (new_unclean != uncleanSize) || (new_clean != cleanSize) ) {
416  LogDebug("UnifiedSC") << ">>>>!!!!!! MISMATCH: new unclean/ old unclean= "
417  << new_unclean << " / " << uncleanSize
418  << ", new clean/ old clean" << new_clean << " / " << cleanSize;
419  }
420 
421 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool isInUnclean() const
Definition: CaloCluster.h:167
virtual void produce(edm::Event &, const edm::EventSetup &)
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:135
edm::Ptr< CaloCluster > CaloClusterPtr
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:149
void setFlags(uint32_t flags)
Definition: CaloCluster.h:162
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
int j
Definition: DBlmapReader.cc:9
double energy() const
cluster energy
Definition: CaloCluster.h:109
bool first
Definition: L1TdeRCT.cc:79
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:359
bool isInClean() const
Definition: CaloCluster.h:166
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:61