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