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