CMS 3D CMS Logo

Multi5x5SuperClusterProducer.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
16 
19 
20 // Class header file
22 
23 
25 {
26 
28  consumes<reco::BasicClusterCollection>(ps.getParameter<edm::InputTag>("endcapClusterTag"));
30  consumes<reco::BasicClusterCollection>(ps.getParameter<edm::InputTag>("barrelClusterTag"));
31 
32 
33  endcapSuperclusterCollection_ = ps.getParameter<std::string>("endcapSuperclusterCollection");
34  barrelSuperclusterCollection_ = ps.getParameter<std::string>("barrelSuperclusterCollection");
35 
36  doBarrel_ = ps.getParameter<bool>("doBarrel");
37  doEndcaps_ = ps.getParameter<bool>("doEndcaps");
38 
39  barrelEtaSearchRoad_ = ps.getParameter<double>("barrelEtaSearchRoad");
40  barrelPhiSearchRoad_ = ps.getParameter<double>("barrelPhiSearchRoad");
41  endcapEtaSearchRoad_ = ps.getParameter<double>("endcapEtaSearchRoad");
42  endcapPhiSearchRoad_ = ps.getParameter<double>("endcapPhiSearchRoad");
43  seedTransverseEnergyThreshold_ = ps.getParameter<double>("seedTransverseEnergyThreshold");
44 
45  const edm::ParameterSet bremRecoveryPset = ps.getParameter<edm::ParameterSet>("bremRecoveryPset");
46  bool dynamicPhiRoad = ps.getParameter<bool>("dynamicPhiRoad");
47 
48  bremAlgo_p = new Multi5x5BremRecoveryClusterAlgo(bremRecoveryPset, barrelEtaSearchRoad_, barrelPhiSearchRoad_,
49  endcapEtaSearchRoad_, endcapPhiSearchRoad_,
50  dynamicPhiRoad, seedTransverseEnergyThreshold_);
51 
52  produces< reco::SuperClusterCollection >(endcapSuperclusterCollection_);
53  produces< reco::SuperClusterCollection >(barrelSuperclusterCollection_);
54 
55  totalE = 0;
56  noSuperClusters = 0;
57  nEvt_ = 0;
58 }
59 
60 
62 {
63  delete bremAlgo_p;
64 }
65 
66 void
68  double averEnergy = 0.;
69  std::ostringstream str;
70  str << "Multi5x5SuperClusterProducer::endJob()\n"
71  << " total # reconstructed super clusters: " << noSuperClusters << "\n"
72  << " total energy of all clusters: " << totalE << "\n";
73  if(noSuperClusters>0) {
74  averEnergy = totalE / noSuperClusters;
75  str << " average SuperCluster energy = " << averEnergy << "\n";
76  }
77  edm::LogInfo("Multi5x5SuperClusterProducerInfo") << str.str() << "\n";
78 
79 }
80 
81 
83 {
84  if(doEndcaps_)
86 
87  if(doBarrel_)
89 
90  nEvt_++;
91 }
92 
93 
98 {
99  // get the cluster collection out and turn it to a BasicClusterRefVector:
100  reco::CaloClusterPtrVector *clusterPtrVector_p = new reco::CaloClusterPtrVector;
101  getClusterPtrVector(evt, clustersToken, clusterPtrVector_p);
102 
103  // run the brem recovery and get the SC collection
104  auto superclusters_ap = std::make_unique<reco::SuperClusterCollection>(bremAlgo_p->makeSuperClusters(*clusterPtrVector_p));
105 
106  // count the total energy and the number of superclusters
107  reco::SuperClusterCollection::iterator it;
108  for (it = superclusters_ap->begin(); it != superclusters_ap->end(); it++)
109  {
110  totalE += it->energy();
111  noSuperClusters++;
112  }
113 
114  // put the SC collection in the event
115  evt.put(std::move(superclusters_ap), superclusterCollection);
116 
117  delete clusterPtrVector_p;
118 }
119 
120 
123  reco::CaloClusterPtrVector *clusterPtrVector_p)
124 {
126  evt.getByToken(clustersToken, bccHandle);
127 
128  const reco::BasicClusterCollection *clusterCollection_p = bccHandle.product();
129  for (unsigned int i = 0; i < clusterCollection_p->size(); i++)
130  {
131  clusterPtrVector_p->push_back(reco::CaloClusterPtr(bccHandle, i));
132  }
133 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::BasicClusterCollection > eeClustersToken_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
void produce(edm::Event &, const edm::EventSetup &) override
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
Multi5x5BremRecoveryClusterAlgo * bremAlgo_p
Multi5x5SuperClusterProducer(const edm::ParameterSet &ps)
edm::PtrVector< CaloCluster > CaloClusterPtrVector
edm::EDGetTokenT< reco::BasicClusterCollection > ebClustersToken_
void produceSuperclustersForECALPart(edm::Event &evt, const edm::EDGetTokenT< reco::BasicClusterCollection > &clustersToken, std::string superclusterColection)
T const * product() const
Definition: Handle.h:81
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
void getClusterPtrVector(edm::Event &evt, const edm::EDGetTokenT< reco::BasicClusterCollection > &clustersToken, reco::CaloClusterPtrVector *)
def move(src, dest)
Definition: eostools.py:510
reco::SuperClusterCollection makeSuperClusters(const reco::CaloClusterPtrVector &clusters)