CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HiSuperClusterProducer.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
16 
17 // Reconstruction Classes
22 
23 // Class header file
25 
26 
28 {
29  // The verbosity level
30  std::string verbosityString = ps.getParameter<std::string>("VerbosityLevel");
31  if (verbosityString == "DEBUG") verbosity = HiBremRecoveryClusterAlgo::pDEBUG;
32  else if (verbosityString == "WARNING") verbosity = HiBremRecoveryClusterAlgo::pWARNING;
33  else if (verbosityString == "INFO") verbosity = HiBremRecoveryClusterAlgo::pINFO;
35 
36  endcapClusterProducer_ = ps.getParameter<std::string>("endcapClusterProducer");
37  barrelClusterProducer_ = ps.getParameter<std::string>("barrelClusterProducer");
38 
39  endcapClusterCollection_ = ps.getParameter<std::string>("endcapClusterCollection");
40  barrelClusterCollection_ = ps.getParameter<std::string>("barrelClusterCollection");
41 
42  endcapSuperclusterCollection_ = ps.getParameter<std::string>("endcapSuperclusterCollection");
43  barrelSuperclusterCollection_ = ps.getParameter<std::string>("barrelSuperclusterCollection");
44 
45  doBarrel_ = ps.getParameter<bool>("doBarrel");
46  doEndcaps_ = ps.getParameter<bool>("doEndcaps");
47 
48 
49  barrelEtaSearchRoad_ = ps.getParameter<double>("barrelEtaSearchRoad");
50  barrelPhiSearchRoad_ = ps.getParameter<double>("barrelPhiSearchRoad");
51  endcapEtaSearchRoad_ = ps.getParameter<double>("endcapEtaSearchRoad");
52  endcapPhiSearchRoad_ = ps.getParameter<double>("endcapPhiSearchRoad");
53  seedTransverseEnergyThreshold_ = ps.getParameter<double>("seedTransverseEnergyThreshold");
54  barrelBCEnergyThreshold_ = ps.getParameter<double>("barrelBCEnergyThreshold");
55  endcapBCEnergyThreshold_ = ps.getParameter<double>("endcapBCEnergyThreshold");
56 
57  if (verbosityString == "INFO") {
58  std::cout <<"Barrel BC Energy threshold = "<<barrelBCEnergyThreshold_<<std::endl;
59  std::cout <<"Endcap BC Energy threshold = "<<endcapBCEnergyThreshold_<<std::endl;
60  }
61 
62  bremAlgo_p = new HiBremRecoveryClusterAlgo(barrelEtaSearchRoad_, barrelPhiSearchRoad_,
63  endcapEtaSearchRoad_, endcapPhiSearchRoad_,
64  seedTransverseEnergyThreshold_,
65  barrelBCEnergyThreshold_,
66  endcapBCEnergyThreshold_,
67  verbosity);
68 
69  produces< reco::SuperClusterCollection >(endcapSuperclusterCollection_);
70  produces< reco::SuperClusterCollection >(barrelSuperclusterCollection_);
71 
72  totalE = 0;
73  noSuperClusters = 0;
74  nEvt_ = 0;
75 }
76 
77 
79 {
80  delete bremAlgo_p;
81 }
82 
84  double averEnergy = 0.;
85  std::ostringstream str;
86  str << "HiSuperClusterProducer::endJob()\n"
87  << " total # reconstructed super clusters: " << noSuperClusters << "\n"
88  << " total energy of all clusters: " << totalE << "\n";
89  if(noSuperClusters>0) {
90  averEnergy = totalE / noSuperClusters;
91  str << " average SuperCluster energy = " << averEnergy << "\n";
92  }
93  edm::LogInfo("HiSuperClusterProducerInfo") << str.str() << "\n";
94 
95 }
96 
97 
99 {
100  if(doEndcaps_)
102 
103  if(doBarrel_)
105 
106  nEvt_++;
107 }
108 
109 
111  std::string clusterProducer,
112  std::string clusterCollection,
113  std::string superclusterCollection)
114 {
115  // get the cluster collection out and turn it to a BasicClusterRefVector:
116  reco::CaloClusterPtrVector *clusterPtrVector_p = new reco::CaloClusterPtrVector;
117  getClusterPtrVector(evt, clusterProducer, clusterCollection, clusterPtrVector_p);
118 
119  // run the brem recovery and get the SC collection
120  std::auto_ptr<reco::SuperClusterCollection>
121  superclusters_ap(new reco::SuperClusterCollection(bremAlgo_p->makeSuperClusters(*clusterPtrVector_p)));
122 
123  // count the total energy and the number of superclusters
124  reco::SuperClusterCollection::iterator it;
125  for (it = superclusters_ap->begin(); it != superclusters_ap->end(); it++)
126  {
127  totalE += it->energy();
128  noSuperClusters++;
129  }
130 
131  // put the SC collection in the event
132  evt.put(superclusters_ap, superclusterCollection);
133 
134  delete clusterPtrVector_p;
135 }
136 
137 
138 void HiSuperClusterProducer::getClusterPtrVector(edm::Event& evt, std::string clusterProducer_, std::string clusterCollection_, reco::CaloClusterPtrVector *clusterPtrVector_p)
139 {
141 
142  evt.getByLabel(clusterProducer_, clusterCollection_, bccHandle);
143  if (!(bccHandle.isValid()))
144  {
145  edm::LogError("HiSuperClusterProducerError") << "could not get a handle on the BasicCluster Collection!";
146  edm::LogError("HiSuperClusterProducerError") << "Error! can't get the product " << clusterCollection_.c_str();
147  clusterPtrVector_p = 0;
148  }
149 
150  const reco::BasicClusterCollection *clusterCollection_p = bccHandle.product();
151  for (unsigned int i = 0; i < clusterCollection_p->size(); i++)
152  {
153  clusterPtrVector_p->push_back(reco::CaloClusterPtr(bccHandle, i));
154  }
155 }
156 
157 
158 
159 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void getClusterPtrVector(edm::Event &evt, std::string clusterProducer_, std::string clusterCollection_, reco::CaloClusterPtrVector *)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:138
reco::SuperClusterCollection makeSuperClusters(reco::CaloClusterPtrVector &clusters)
void produceSuperclustersForECALPart(edm::Event &evt, std::string clusterProducer, std::string clusterCollection, std::string superclusterColection)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::PtrVector< CaloCluster > CaloClusterPtrVector
HiSuperClusterProducer(const edm::ParameterSet &ps)
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
HiBremRecoveryClusterAlgo::VerbosityLevel verbosity
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
T const * product() const
Definition: Handle.h:81
virtual void produce(edm::Event &, const edm::EventSetup &)
HiBremRecoveryClusterAlgo * bremAlgo_p
tuple cout
Definition: gather_cfg.py:121