CMS 3D CMS Logo

HybridClusterProducer.cc
Go to the documentation of this file.
1 // C/C++ headers
2 #include <iostream>
3 #include <vector>
4 #include <memory>
5 
6 // Framework
13 
14 // Reconstruction Classes
21 
22 // Geometry
30 
33 
34 // Class header file
37 
38 
39 
41 {
42 
43 
44  basicclusterCollection_ = ps.getParameter<std::string>("basicclusterCollection");
45  superclusterCollection_ = ps.getParameter<std::string>("superclusterCollection");
46  hitsToken_ =
47  consumes<EcalRecHitCollection>(ps.getParameter<edm::InputTag>("recHitsCollection"));
48 
49  //Setup for core tools objects.
51  ps.getParameter<edm::ParameterSet>("posCalcParameters");
52 
53  posCalculator_ = PositionCalc(posCalcParameters);
54 
55  const std::vector<std::string> flagnames =
56  ps.getParameter<std::vector<std::string> >("RecHitFlagToBeExcluded");
57 
58  const std::vector<int> flagsexcl=
59  StringToEnumValue<EcalRecHit::Flags>(flagnames);
60 
61  const std::vector<std::string> severitynames =
62  ps.getParameter<std::vector<std::string> >("RecHitSeverityToBeExcluded");
63 
64  const std::vector<int> severitiesexcl=
65  StringToEnumValue<EcalSeverityLevel::SeverityLevel>(severitynames);
66 
67  hybrid_p = new HybridClusterAlgo(ps.getParameter<double>("HybridBarrelSeedThr"),
68  ps.getParameter<int>("step"),
69  ps.getParameter<double>("ethresh"),
70  ps.getParameter<double>("eseed"),
71  ps.getParameter<double>("xi"),
72  ps.getParameter<bool>("useEtForXi"),
73  ps.getParameter<double>("ewing"),
74  flagsexcl,
75  posCalculator_,
76  ps.getParameter<bool>("dynamicEThresh"),
77  ps.getParameter<double>("eThreshA"),
78  ps.getParameter<double>("eThreshB"),
79  severitiesexcl,
80  ps.getParameter<bool>("excludeFlagged")
81  );
82  //bremRecoveryPset,
83 
84  // get brem recovery parameters
85  bool dynamicPhiRoad = ps.getParameter<bool>("dynamicPhiRoad");
86  if (dynamicPhiRoad) {
88  hybrid_p->setDynamicPhiRoad(bremRecoveryPset);
89  }
90 
91  produces< reco::BasicClusterCollection >(basicclusterCollection_);
92  produces< reco::SuperClusterCollection >(superclusterCollection_);
93  nEvt_ = 0;
94 }
95 
96 
98 {
99  delete hybrid_p;
100 }
101 
102 
104 {
105  // get the hit collection from the event:
107 
108  evt.getByToken(hitsToken_, rhcHandle);
109  if (!(rhcHandle.isValid())){
110  edm::LogError("MissingProduct") << "could not get a handle on the EcalRecHitCollection!";
111  return;
112 
113  }
114  const EcalRecHitCollection *hit_collection = rhcHandle.product();
115 
116  // get the collection geometry:
117  edm::ESHandle<CaloGeometry> geoHandle;
118  es.get<CaloGeometryRecord>().get(geoHandle);
119  const CaloGeometry& geometry = *geoHandle;
120  const CaloSubdetectorGeometry *geometry_p;
121  std::unique_ptr<const CaloSubdetectorTopology> topology;
122 
124  es.get<EcalSeverityLevelAlgoRcd>().get(sevLv);
125 
126  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
127  topology = std::make_unique<EcalBarrelTopology>(geoHandle);
128 
129  // make the Basic clusters!
130  reco::BasicClusterCollection basicClusters;
131  hybrid_p->makeClusters(hit_collection, geometry_p, basicClusters, sevLv.product(),false,
132  std::vector<RectangularEtaPhiRegion>());
133 
134  LogTrace("EcalClusters") << "Finished clustering - BasicClusterCollection returned to producer..." ;
135 
136  // create a unique_ptr to a BasicClusterCollection, copy the clusters into it and put in the Event:
137  auto basicclusters_p = std::make_unique<reco::BasicClusterCollection>();
138  basicclusters_p->assign(basicClusters.begin(), basicClusters.end());
140 
141  //Basic clusters now in the event.
142  LogTrace("EcalClusters") << "Basic Clusters now put into event." ;
143 
144 
145  //Weird though it is, get the BasicClusters back out of the event. We need the
146  //edm::Ref to these guys to make our superclusters for Hybrid.
147 
148  if (!(bccHandle.isValid())) {
149  edm::LogError("Missing Product") << "could not get a handle on the BasicClusterCollection!" ;
150  return;
151  }
152 
153  reco::BasicClusterCollection clusterCollection = *bccHandle;
154 
155  LogTrace("EcalClusters")<< "Got the BasicClusterCollection" << std::endl;
156 
157  reco::CaloClusterPtrVector clusterPtrVector;
158  for (unsigned int i = 0; i < clusterCollection.size(); i++){
159  clusterPtrVector.push_back(reco::CaloClusterPtr(bccHandle, i));
160  }
161 
162  reco::SuperClusterCollection superClusters = hybrid_p->makeSuperClusters(clusterPtrVector);
163  LogTrace("EcalClusters") << "Found: " << superClusters.size() << " superclusters." ;
164 
165  auto superclusters_p = std::make_unique<reco::SuperClusterCollection>();
166  superclusters_p->assign(superClusters.begin(), superClusters.end());
167 
168  evt.put(std::move(superclusters_p), superclusterCollection_);
169  LogTrace("EcalClusters") << "Hybrid Clusters (Basic/Super) added to the Event! :-)" ;
170 
171 
172  nEvt_++;
173 }
174 
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:49
HybridClusterProducer(const edm::ParameterSet &ps)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
CaloTopology const * topology(0)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
edm::EDGetTokenT< EcalRecHitCollection > hitsToken_
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
reco::SuperClusterCollection makeSuperClusters(const reco::CaloClusterPtrVector &)
bool isValid() const
Definition: HandleBase.h:74
#define LogTrace(id)
T const * product() const
Definition: Handle.h:81
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
T get() const
Definition: EventSetup.h:68
HybridClusterAlgo * hybrid_p
void makeClusters(const EcalRecHitCollection *, const CaloSubdetectorGeometry *geometry, reco::BasicClusterCollection &basicClusters, const EcalSeverityLevelAlgo *sevLv, bool regional=false, const std::vector< RectangularEtaPhiRegion > &regions=std::vector< RectangularEtaPhiRegion >())
T const * product() const
Definition: ESHandle.h:84
def move(src, dest)
Definition: eostools.py:511