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  barrelEtaSearchRoad_{static_cast<float>(ps.getParameter<double>("barrelEtaSearchRoad"))},
26  barrelPhiSearchRoad_{static_cast<float>(ps.getParameter<double>("barrelPhiSearchRoad"))},
27  endcapEtaSearchRoad_{static_cast<float>(ps.getParameter<double>("endcapEtaSearchRoad"))},
28  endcapPhiSearchRoad_{static_cast<float>(ps.getParameter<double>("endcapPhiSearchRoad"))},
29  seedTransverseEnergyThreshold_{static_cast<float>(ps.getParameter<double>("seedTransverseEnergyThreshold"))},
30  doBarrel_{ps.getParameter<bool>("doBarrel")},
31  doEndcaps_{ps.getParameter<bool>("doEndcaps")},
32  totalE{0.},
34 {
35 
36  if(doEndcaps_) {
38  consumes<reco::BasicClusterCollection>(ps.getParameter<edm::InputTag>("endcapClusterTag"));
39  }
40  if(doBarrel_) {
42  consumes<reco::BasicClusterCollection>(ps.getParameter<edm::InputTag>("barrelClusterTag"));
43  }
44 
45  const edm::ParameterSet bremRecoveryPset = ps.getParameter<edm::ParameterSet>("bremRecoveryPset");
46  bool dynamicPhiRoad = ps.getParameter<bool>("dynamicPhiRoad");
47 
48  bremAlgo_p = std::make_unique<Multi5x5BremRecoveryClusterAlgo>(bremRecoveryPset, barrelEtaSearchRoad_, barrelPhiSearchRoad_,
51 
52  if(doEndcaps_) {
53  endcapPutToken_ = produces< reco::SuperClusterCollection >(ps.getParameter<std::string>("endcapSuperclusterCollection"));
54  }
55  if(doBarrel_) {
56  barrelPutToken_ = produces< reco::SuperClusterCollection >(ps.getParameter<std::string>("barrelSuperclusterCollection"));
57  }
58 }
59 
60 
61 void
63  double averEnergy = 0.;
64  std::ostringstream str;
65  str << "Multi5x5SuperClusterProducer::endJob()\n"
66  << " total # reconstructed super clusters: " << noSuperClusters << "\n"
67  << " total energy of all clusters: " << totalE << "\n";
68  if(noSuperClusters>0) {
69  averEnergy = totalE / noSuperClusters;
70  str << " average SuperCluster energy = " << averEnergy << "\n";
71  }
72  edm::LogInfo("Multi5x5SuperClusterProducerInfo") << str.str() << "\n";
73 
74 }
75 
76 
78 {
79  if(doEndcaps_)
81 
82  if(doBarrel_)
84 }
85 
86 
91 {
92  // get the cluster collection out and turn it to a BasicClusterRefVector:
93  reco::CaloClusterPtrVector clusterPtrVector_p=getClusterPtrVector(evt, clustersToken);
94 
95  // run the brem recovery and get the SC collection
96  reco::SuperClusterCollection superclusters_ap(bremAlgo_p->makeSuperClusters(clusterPtrVector_p));
97 
98  // count the total energy and the number of superclusters
99  for (auto const& sc: superclusters_ap)
100  {
101  totalE += sc.energy();
102  noSuperClusters++;
103  }
104 
105  // put the SC collection in the event
106  evt.emplace(putToken,std::move(superclusters_ap));
107 
108 }
109 
112  const edm::EDGetTokenT<reco::BasicClusterCollection>& clustersToken) const
113 {
114  reco::CaloClusterPtrVector clusterPtrVector_p;
116  evt.getByToken(clustersToken, bccHandle);
117 
118  const reco::BasicClusterCollection *clusterCollection_p = bccHandle.product();
119  clusterPtrVector_p.reserve(clusterCollection_p->size());
120  for (unsigned int i = 0; i < clusterCollection_p->size(); i++)
121  {
122  clusterPtrVector_p.push_back(reco::CaloClusterPtr(bccHandle, i));
123  }
124  return clusterPtrVector_p;
125 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::BasicClusterCollection > eeClustersToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void produce(edm::Event &, const edm::EventSetup &) override
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
std::unique_ptr< Multi5x5BremRecoveryClusterAlgo > bremAlgo_p
void produceSuperclustersForECALPart(edm::Event &evt, const edm::EDGetTokenT< reco::BasicClusterCollection > &clustersToken, const edm::EDPutTokenT< reco::SuperClusterCollection > &putToken)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
Multi5x5SuperClusterProducer(const edm::ParameterSet &ps)
edm::EDPutTokenT< reco::SuperClusterCollection > endcapPutToken_
float energy() const
Energy. Note this is taken from the first SimTrack only.
Definition: SimCluster.h:104
edm::EDPutTokenT< reco::SuperClusterCollection > barrelPutToken_
edm::EDGetTokenT< reco::BasicClusterCollection > ebClustersToken_
T const * product() const
Definition: Handle.h:74
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:413
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
void reserve(size_type n)
Reserve space for RefVector.
Definition: PtrVectorBase.h:88
reco::CaloClusterPtrVector getClusterPtrVector(edm::Event &evt, const edm::EDGetTokenT< reco::BasicClusterCollection > &clustersToken) const
#define str(s)
def move(src, dest)
Definition: eostools.py:511