CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CleanAndMergeProducer.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
21 
22 // Geometry
30 
31 // Class header file
35 
36 
37 /*
38 CleanAndMergeProducer:
39 ^^^^^^^^^^^^^^^^^^^^^^
40 
41 Takes as input the cleaned and the uncleaned collection of SC
42 and produces an uncleaned collection of SC without the duplicates
43 and a collection of references to the cleaned collection of SC
44 
45 25 June 2010
46 Nikolaos Rompotis and Chris Seez - Imperial College London
47 many thanks to Shahram Rahatlou and Federico Ferri
48 
49 
50 */
51 
52 
54 {
55 
56  // get the parameters
57  // the cleaned collection:
58  cleanScToken_ =
59  consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("cleanScInputTag"));
61  consumes<reco::SuperClusterCollection>(ps.getParameter<edm::InputTag>("uncleanScInputTag"));
62 
63  // the names of the products to be produced:
64  bcCollection_ = ps.getParameter<std::string>("bcCollection");
65  scCollection_ = ps.getParameter<std::string>("scCollection");
66  refScCollection_ = ps.getParameter<std::string>("refScCollection");
67 
68  std::map<std::string,double> providedParameters;
69  providedParameters.insert(std::make_pair("LogWeighted",ps.getParameter<bool>("posCalc_logweight")));
70  providedParameters.insert(std::make_pair("T0_barl",ps.getParameter<double>("posCalc_t0")));
71  providedParameters.insert(std::make_pair("W0",ps.getParameter<double>("posCalc_w0")));
72  providedParameters.insert(std::make_pair("X0",ps.getParameter<double>("posCalc_x0")));
73 
74 
75  // the products:
76  produces< reco::BasicClusterCollection >(bcCollection_);
77  produces< reco::SuperClusterCollection >(scCollection_);
78  produces< reco::SuperClusterRefVector >(refScCollection_);
79 
80 }
81 
83 
85  const edm::EventSetup& es)
86 {
87  // get the input collections
88  // ______________________________________________________________________________________
89  //
90  // cluster collections:
91 
94 
95  evt.getByToken(cleanScToken_, pCleanSC);
96  evt.getByToken(uncleanScToken_, pUncleanSC);
97 
98  //
99  // collections are all taken now _____________________________________________________
100  //
101  //
102  // the collections to be produced:
103  reco::BasicClusterCollection basicClusters;
104  reco::SuperClusterCollection superClusters;
106 
107  //
108  // run over the uncleaned SC and check how many of them are matched to the cleaned ones
109  // if you find a matched one, create a reference to the cleaned collection and store it
110  // if you find an unmatched one, then keep all its basic clusters in the basic Clusters
111  // vector
112  int uncleanSize = pUncleanSC->size();
113  int cleanSize = pCleanSC->size();
114 
115  LogTrace("EcalCleaning") << "Size of Clean Collection: " << cleanSize
116  << ", uncleanSize: " << uncleanSize << std::endl;
117  //
118  // keep whether the SC in unique in the uncleaned collection
119  std::vector<int> isUncleanOnly; // 1 if unique in the uncleaned
120  std::vector<int> basicClusterOwner; // contains the index of the SC that owns that BS
121  std::vector<int> isSeed; // if this basic cluster is a seed it is 1
122  for (int isc =0; isc< uncleanSize; ++isc) {
123  reco::SuperClusterRef unscRef( pUncleanSC, isc );
124  const std::vector< std::pair<DetId, float> > & uhits = unscRef->hitsAndFractions();
125  int uhitsSize = uhits.size();
126  bool foundTheSame = false;
127  for (int jsc=0; jsc < cleanSize; ++jsc) { // loop over the cleaned SC
128  reco::SuperClusterRef cscRef( pCleanSC, jsc );
129  const std::vector< std::pair<DetId, float> > & chits = cscRef->hitsAndFractions();
130  int chitsSize = chits.size();
131  foundTheSame = true;
132  if (unscRef->seed()->seed() == cscRef->seed()->seed() && chitsSize == uhitsSize) {
133  // if the clusters are exactly the same then because the clustering
134  // algorithm works in a deterministic way, the order of the rechits
135  // will be the same
136  for (int i=0; i< chitsSize; ++i) {
137  if (uhits[i].first != chits[i].first ) { foundTheSame=false; break;}
138  }
139  if (foundTheSame) { // ok you have found it!
140  // make the reference
141  //scRefs->push_back( edm::Ref<reco::SuperClusterCollection>(pCleanSC, jsc) );
142  scRefs->push_back( cscRef );
143  isUncleanOnly.push_back(0);
144  break;
145  }
146  }
147  }
148  if (not foundTheSame) {
149  // mark it as unique in the uncleaned
150  isUncleanOnly.push_back(1);
151  // keep all its basic clusters
152  for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
153  // the basic clusters
154  basicClusters.push_back(**bciter);
155  basicClusterOwner.push_back(isc);
156  }
157  }
158  }
159  int bcSize = basicClusters.size();
160 
161  LogDebug("EcalCleaning") << "Found cleaned SC: " << cleanSize << " uncleaned SC: "
162  << uncleanSize << " from which " << scRefs->size()
163  << " will become refs to the cleaned collection" ;
164 
165 
166  // now you have the collection of basic clusters of the SC to be remain in the
167  // in the clean collection, export them to the event
168  // you will need to reread them later in order to make correctly the refs to the SC
169  std::auto_ptr< reco::BasicClusterCollection> basicClusters_p(new reco::BasicClusterCollection);
170  basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
172  evt.put(basicClusters_p, bcCollection_);
173  if (!(bccHandle.isValid())) {
174  edm::LogWarning("MissingInput")<<"could not get a handle on the BasicClusterCollection!" << std::endl;
175  return;
176  }
177 
178  LogDebug("EcalCleaning")<< "Got the BasicClusters from the event again";
179  // now you have to create again your superclusters
180  // you run over the uncleaned SC, but now you know which of them are
181  // the ones that are needed and which are their basic clusters
182  for (int isc=0; isc< uncleanSize; ++isc) {
183  if (isUncleanOnly[isc]==1) { // look for sc that are unique in the unclean collection
184  // make the basic cluster collection
185  reco::CaloClusterPtrVector clusterPtrVector;
186  reco::CaloClusterPtr seed; // the seed is the basic cluster with the highest energy
187  double energy = -1;
188  for (int jbc=0; jbc< bcSize; ++jbc) {
189  if (basicClusterOwner[jbc]==isc) {
190  reco::CaloClusterPtr currentClu = reco::CaloClusterPtr(bccHandle, jbc);
191  clusterPtrVector.push_back(currentClu);
192  if (energy< currentClu->energy()) {
193  energy = currentClu->energy(); seed = currentClu;
194  }
195  }
196  }
197  reco::SuperClusterRef unscRef( pUncleanSC, isc );
198  reco::SuperCluster newSC(unscRef->energy(), unscRef->position(), seed, clusterPtrVector );
199  superClusters.push_back(newSC);
200  }
201 
202  }
203  // export the collection of references to the clean collection
204  std::auto_ptr< reco::SuperClusterRefVector > scRefs_p( scRefs );
205  evt.put(scRefs_p, refScCollection_);
206 
207  // the collection of basic clusters is already in the event
208  // the collection of uncleaned SC
209  std::auto_ptr< reco::SuperClusterCollection > superClusters_p(new reco::SuperClusterCollection);
210  superClusters_p->assign(superClusters.begin(), superClusters.end());
211  evt.put(superClusters_p, scCollection_);
212  LogDebug("EcalCleaning")<< "Hybrid Clusters (Basic/Super) added to the Event! :-)";
213 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:141
edm::EDGetTokenT< reco::SuperClusterCollection > uncleanScToken_
edm::Ptr< CaloCluster > CaloClusterPtr
edm::EDGetTokenT< reco::SuperClusterCollection > cleanScToken_
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
virtual void produce(edm::Event &, const edm::EventSetup &)
#define LogTrace(id)
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:62
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
edm::RefVector< SuperClusterCollection > SuperClusterRefVector
vector of references to objects in the same colletion of SuperCluster objects
CleanAndMergeProducer(const edm::ParameterSet &ps)