CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SuperClusterProducer.cc
Go to the documentation of this file.
1 // C/C++ headers
2 #include <iostream>
3 #include <vector>
4 #include <memory>
5 #include <sstream>
6 
7 // Framework
13 
14 // Reconstruction Classes
19 
20 // Class header file
22 
23 
25 {
26  // The verbosity level
27  std::string verbosityString = ps.getParameter<std::string>("VerbosityLevel");
28  if (verbosityString == "DEBUG") verbosity = BremRecoveryClusterAlgo::pDEBUG;
29  else if (verbosityString == "WARNING") verbosity = BremRecoveryClusterAlgo::pWARNING;
30  else if (verbosityString == "INFO") verbosity = BremRecoveryClusterAlgo::pINFO;
32 
33  endcapClusterProducer_ = ps.getParameter<std::string>("endcapClusterProducer");
34  barrelClusterProducer_ = ps.getParameter<std::string>("barrelClusterProducer");
35 
36  endcapClusterCollection_ = ps.getParameter<std::string>("endcapClusterCollection");
37  barrelClusterCollection_ = ps.getParameter<std::string>("barrelClusterCollection");
38 
39  endcapSuperclusterCollection_ = ps.getParameter<std::string>("endcapSuperclusterCollection");
40  barrelSuperclusterCollection_ = ps.getParameter<std::string>("barrelSuperclusterCollection");
41 
42  doBarrel_ = ps.getParameter<bool>("doBarrel");
43  doEndcaps_ = ps.getParameter<bool>("doEndcaps");
44 
45 
46  barrelEtaSearchRoad_ = ps.getParameter<double>("barrelEtaSearchRoad");
47  barrelPhiSearchRoad_ = ps.getParameter<double>("barrelPhiSearchRoad");
48  endcapEtaSearchRoad_ = ps.getParameter<double>("endcapEtaSearchRoad");
49  endcapPhiSearchRoad_ = ps.getParameter<double>("endcapPhiSearchRoad");
50  seedTransverseEnergyThreshold_ = ps.getParameter<double>("seedTransverseEnergyThreshold");
51 
52  bremAlgo_p = new BremRecoveryClusterAlgo(barrelEtaSearchRoad_, barrelPhiSearchRoad_,
53  endcapEtaSearchRoad_, endcapPhiSearchRoad_,
54  seedTransverseEnergyThreshold_, verbosity);
55 
56  produces< reco::SuperClusterCollection >(endcapSuperclusterCollection_);
57  produces< reco::SuperClusterCollection >(barrelSuperclusterCollection_);
58 
59  totalE = 0;
60  noSuperClusters = 0;
61  nEvt_ = 0;
62 }
63 
64 
66 {
67  delete bremAlgo_p;
68 }
69 
70 void
72  double averEnergy = 0.;
73  std::ostringstream str;
74  str << "SuperClusterProducer::endJob()\n"
75  << " total # reconstructed super clusters: " << noSuperClusters << "\n"
76  << " total energy of all clusters: " << totalE << "\n";
77  if(noSuperClusters>0) {
78  averEnergy = totalE / noSuperClusters;
79  str << " average SuperCluster energy = " << averEnergy << "\n";
80  }
81  edm::LogInfo("SuperClusterProducerInfo") << str.str() << "\n";
82 
83 }
84 
85 
87 {
88  if(doEndcaps_)
90 
91  if(doBarrel_)
93 
94  nEvt_++;
95 }
96 
97 
99  std::string clusterProducer,
100  std::string clusterCollection,
101  std::string superclusterCollection)
102 {
103  // get the cluster collection out and turn it to a BasicClusterRefVector:
104  reco::CaloClusterPtrVector *clusterPtrVector_p = new reco::CaloClusterPtrVector;
105  getClusterPtrVector(evt, clusterProducer, clusterCollection, clusterPtrVector_p);
106 
107  // run the brem recovery and get the SC collection
108  std::auto_ptr<reco::SuperClusterCollection>
109  superclusters_ap(new reco::SuperClusterCollection(bremAlgo_p->makeSuperClusters(*clusterPtrVector_p)));
110 
111  // count the total energy and the number of superclusters
112  reco::SuperClusterCollection::iterator it;
113  for (it = superclusters_ap->begin(); it != superclusters_ap->end(); it++)
114  {
115  totalE += it->energy();
116  noSuperClusters++;
117  }
118 
119  // put the SC collection in the event
120  evt.put(superclusters_ap, superclusterCollection);
121 
122  delete clusterPtrVector_p;
123 }
124 
125 
126 void SuperClusterProducer::getClusterPtrVector(edm::Event& evt, std::string clusterProducer_, std::string clusterCollection_, reco::CaloClusterPtrVector *clusterPtrVector_p)
127 {
129 
130  evt.getByLabel(clusterProducer_, clusterCollection_, bccHandle);
131  if (!(bccHandle.isValid()))
132  {
133  edm::LogError("SuperClusterProducerError") << "could not get a handle on the BasicCluster Collection!";
134  edm::LogError("SuperClusterProducerError") << "Error! can't get the product " << clusterCollection_.c_str();
135  clusterPtrVector_p = 0;
136  }
137 
138  const reco::BasicClusterCollection *clusterCollection_p = bccHandle.product();
139  for (unsigned int i = 0; i < clusterCollection_p->size(); i++)
140  {
141  clusterPtrVector_p->push_back(reco::CaloClusterPtr(bccHandle, i));
142  }
143 }
144 
145 
146 
147 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::string barrelClusterProducer_
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:137
void getClusterPtrVector(edm::Event &evt, std::string clusterProducer_, std::string clusterCollection_, reco::CaloClusterPtrVector *)
reco::SuperClusterCollection makeSuperClusters(reco::CaloClusterPtrVector &clusters)
BremRecoveryClusterAlgo::VerbosityLevel verbosity
std::string barrelSuperclusterCollection_
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
BremRecoveryClusterAlgo * bremAlgo_p
void produceSuperclustersForECALPart(edm::Event &evt, std::string clusterProducer, std::string clusterCollection, std::string superclusterColection)
edm::PtrVector< CaloCluster > CaloClusterPtrVector
bool isValid() const
Definition: HandleBase.h:76
SuperClusterProducer(const edm::ParameterSet &ps)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual void produce(edm::Event &, const edm::EventSetup &)
std::string endcapClusterProducer_
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
T const * product() const
Definition: Handle.h:74
std::string barrelClusterCollection_
std::string endcapClusterCollection_
std::string endcapSuperclusterCollection_