CMS 3D CMS Logo

EgammaHLTHybridClusterProducer.cc
Go to the documentation of this file.
1 // C/C++ headers
2 #include <iostream>
3 #include <vector>
4 #include <memory>
5 
6 // Framework
10 
13 
14 // Reconstruction Classes
25 
26 // Geometry
34 
35 // Level 1 Trigger
38 
39 // EgammaCoreTools
41 
42 // Class header file
44 
45 
47  : basicclusterCollection_ (ps.getParameter<std::string>("basicclusterCollection"))
48  , superclusterCollection_ (ps.getParameter<std::string>("superclusterCollection"))
49  , hittoken_ (consumes<EcalRecHitCollection>(hitcollection_))
50  , hitcollection_ (ps.getParameter<edm::InputTag>("ecalhitcollection"))
51 
52  // L1 matching parameters
53  , l1TagIsolated_ (consumes<l1extra::L1EmParticleCollection>(ps.getParameter< edm::InputTag > ("l1TagIsolated")))
54  , l1TagNonIsolated_ (consumes<l1extra::L1EmParticleCollection>(ps.getParameter< edm::InputTag > ("l1TagNonIsolated")))
55 
56  , doIsolated_ (ps.getParameter<bool>("doIsolated"))
57 
58  , l1LowerThr_ (ps.getParameter<double> ("l1LowerThr"))
59  , l1UpperThr_ (ps.getParameter<double> ("l1UpperThr"))
60  , l1LowerThrIgnoreIsolation_ (ps.getParameter<double> ("l1LowerThrIgnoreIsolation"))
61 
62  , regionEtaMargin_ (ps.getParameter<double>("regionEtaMargin"))
63  , regionPhiMargin_ (ps.getParameter<double>("regionPhiMargin"))
64 
65  // Parameters for the position calculation:
66  , posCalculator_ (PositionCalc( ps.getParameter<edm::ParameterSet>("posCalcParameters") ))
67 
68  , hybrid_p (new HybridClusterAlgo(
69  ps.getParameter<double>("HybridBarrelSeedThr"),
70  ps.getParameter<int>("step"),
71  ps.getParameter<double>("ethresh"),
72  ps.getParameter<double>("eseed"),
73  ps.getParameter<double>("xi"),
74  ps.getParameter<bool>("useEtForXi"),
75  ps.getParameter<double>("ewing"),
77  ps.getParameter<std::vector<std::string> >("RecHitFlagToBeExcluded")),
78  posCalculator_,
79  ps.getParameter<bool>("dynamicEThresh"),
80  ps.getParameter<double>("eThreshA"),
81  ps.getParameter<double>("eThreshB"),
83  ps.getParameter<std::vector<std::string> >("RecHitSeverityToBeExcluded")),
84  ps.getParameter<bool>("excludeFlagged")))
85 {
86  if (ps.getParameter<bool>("dynamicPhiRoad")) {
88  hybrid_p->setDynamicPhiRoad(bremRecoveryPset);
89  }
90 
91  produces< reco::BasicClusterCollection >(basicclusterCollection_);
92  produces< reco::SuperClusterCollection >(superclusterCollection_);
93 }
94 
95 
97 {
98  delete hybrid_p;
99 }
100 
102 
104  desc.add<std::string>("debugLevel" , "INFO");
105  desc.add<std::string>("basicclusterCollection", "");
106  desc.add<std::string>("superclusterCollection", "");
107  desc.add<edm::InputTag>("ecalhitcollection", edm::InputTag("ecalRecHit","EcalRecHitsEB"));
108  desc.add<edm::InputTag>("l1TagIsolated", edm::InputTag("l1extraParticles","Isolated"));
109  desc.add<edm::InputTag>("l1TagNonIsolated", edm::InputTag("l1extraParticles","NonIsolated"));
110  desc.add<bool>("doIsolated", true);
111  desc.add<double>("l1LowerThr", 0);
112  desc.add<double>("l1UpperThr", 9999.0);
113  desc.add<double>("l1LowerThrIgnoreIsolation", 999.0);
114  desc.add<double>("regionEtaMargin", 0.14);
115  desc.add<double>("regionPhiMargin", 0.4);
116 
117  edm::ParameterSetDescription posCalcPSET;
118  posCalcPSET.add<double>("T0_barl", 7.4);
119  posCalcPSET.add<double>("T0_endc", 3.1);
120  posCalcPSET.add<double>("T0_endcPresh", 1.2);
121  posCalcPSET.add<double>("W0", 4.2);
122  posCalcPSET.add<double>("X0", 0.89);
123  posCalcPSET.add<bool>("LogWeighted", true);
124  desc.add<edm::ParameterSetDescription>("posCalcParameters", posCalcPSET);
125 
126  desc.add<std::vector<std::string>>("RecHitFlagToBeExcluded", std::vector<std::string>());
127  desc.add<std::vector<std::string> >("RecHitSeverityToBeExcluded", std::vector<std::string>());
128  desc.add<double>("severityRecHitThreshold", 4.0);
129  desc.add<double>("HybridBarrelSeedThr", 1.0);
130  desc.add<int>("step", 10);
131  desc.add<double>("ethresh", 0.1);
132  desc.add<double>("eseed", 0.35);
133  desc.add<double>("xi", 0);
134  desc.add<bool>("useEtForXi", true);
135  desc.add<double>("ewing", 1.0);
136  desc.add<bool>("dynamicEThresh", false);
137  desc.add<double>("eThreshA", 0.003);
138  desc.add<double>("eThreshB", 0.1);
139  desc.add<bool>("excludeFlagged", false);
140  desc.add<bool>("dynamicPhiRoad", false);
141  //desc.add<edm::ParameterSet>("bremRecoveryPset", edm::ParameterSet());
142 
143  descriptions.add("hltEgammaHLTHybridClusterProducer", desc);
144 }
145 
146 
148 
149  // get the hit collection from the event:
151  evt.getByToken(hittoken_, rhcHandle);
152 
153  if (!(rhcHandle.isValid()))
154  {
155  edm::LogError("ProductNotFound")<< "could not get a handle on the EcalRecHitCollection!" << std::endl;
156  return;
157  }
158  const EcalRecHitCollection *hit_collection = rhcHandle.product();
159 
160  // get the collection geometry:
161  edm::ESHandle<CaloGeometry> geoHandle;
162  es.get<CaloGeometryRecord>().get(geoHandle);
163  const CaloGeometry& geometry = *geoHandle;
164  const CaloSubdetectorGeometry *geometry_p;
165  std::unique_ptr<const CaloSubdetectorTopology> topology;
166 
167  //edm::ESHandle<EcalChannelStatus> chStatus;
168  //es.get<EcalChannelStatusRcd>().get(chStatus);
169  //const EcalChannelStatus* theEcalChStatus = chStatus.product();
170 
172  es.get<EcalSeverityLevelAlgoRcd>().get(sevlv);
173  const EcalSeverityLevelAlgo* sevLevel = sevlv.product();
174 
175  if(hitcollection_.instance() == "EcalRecHitsEB") {
176  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
177  topology = std::make_unique<EcalBarrelTopology>(*geoHandle);
178  } else if(hitcollection_.instance() == "EcalRecHitsEE") {
179  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalEndcap);
180  topology = std::make_unique<EcalEndcapTopology>(*geoHandle);
181  } else if(hitcollection_.instance() == "EcalRecHitsPS") {
182  geometry_p = geometry.getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
183  topology = std::make_unique<EcalPreshowerTopology>();
184  } else throw(std::runtime_error("\n\nHybrid Cluster Producer encountered invalied ecalhitcollection type.\n\n"));
185 
186  //Get the L1 EM Particle Collection
187  //Get the L1 EM Particle Collection
189  if(doIsolated_)
190  evt.getByToken(l1TagIsolated_, emIsolColl);
191 
192  //Get the L1 EM Particle Collection
194  evt.getByToken(l1TagNonIsolated_, emNonIsolColl);
195 
196  // Get the CaloGeometry
197  edm::ESHandle<L1CaloGeometry> l1CaloGeom ;
198  es.get<L1CaloGeometryRecord>().get(l1CaloGeom) ;
199 
200  std::vector<RectangularEtaPhiRegion> regions;
201 
202  if(doIsolated_) {
203  for( l1extra::L1EmParticleCollection::const_iterator emItr = emIsolColl->begin(); emItr != emIsolColl->end() ;++emItr ){
204 
205  if (emItr->et() > l1LowerThr_ && emItr->et() < l1UpperThr_) {
206 
207  // Access the GCT hardware object corresponding to the L1Extra EM object.
208  int etaIndex = emItr->gctEmCand()->etaIndex() ;
209  int phiIndex = emItr->gctEmCand()->phiIndex() ;
210  // Use the L1CaloGeometry to find the eta, phi bin boundaries.
211  double etaLow = l1CaloGeom->etaBinLowEdge( etaIndex ) ;
212  double etaHigh = l1CaloGeom->etaBinHighEdge( etaIndex ) ;
213  double phiLow = l1CaloGeom->emJetPhiBinLowEdge( phiIndex ) ;
214  double phiHigh = l1CaloGeom->emJetPhiBinHighEdge( phiIndex ) ;
215 
216  int isbarl=0;
217  //Part of the region is in the barel if either the upper or lower
218  //edge of the region is within the barrel
219  if(((float)(etaLow)>-1.479 && (float)(etaLow)<1.479) ||
220  ((float)(etaHigh)>-1.479 && (float)(etaHigh)<1.479)) isbarl=1;
221 
222 
223  etaLow -= regionEtaMargin_;
224  etaHigh += regionEtaMargin_;
225  phiLow -= regionPhiMargin_;
226  phiHigh += regionPhiMargin_;
227 
228  if (etaHigh>1.479) etaHigh=1.479;
229  if (etaLow<-1.479) etaLow=-1.479;
230 
231  if(isbarl) regions.push_back(RectangularEtaPhiRegion(etaLow,etaHigh,phiLow,phiHigh));
232 
233  }
234  }
235  }
236 
238  for( l1extra::L1EmParticleCollection::const_iterator emItr = emNonIsolColl->begin(); emItr != emNonIsolColl->end() ;++emItr ){
239 
240  if(doIsolated_&&emItr->et()<l1LowerThrIgnoreIsolation_) continue;
241 
242  if (emItr->et() > l1LowerThr_ && emItr->et() < l1UpperThr_) {
243 
244  // Access the GCT hardware object corresponding to the L1Extra EM object.
245  int etaIndex = emItr->gctEmCand()->etaIndex() ;
246  int phiIndex = emItr->gctEmCand()->phiIndex() ;
247  // Use the L1CaloGeometry to find the eta, phi bin boundaries.
248  double etaLow = l1CaloGeom->etaBinLowEdge( etaIndex ) ;
249  double etaHigh = l1CaloGeom->etaBinHighEdge( etaIndex ) ;
250  double phiLow = l1CaloGeom->emJetPhiBinLowEdge( phiIndex ) ;
251  double phiHigh = l1CaloGeom->emJetPhiBinHighEdge( phiIndex ) ;
252 
253  int isbarl=0;
254  //Part of the region is in the barel if either the upper or lower
255  //edge of the region is within the barrel
256  if(((float)(etaLow)>-1.479 && (float)(etaLow)<1.479) ||
257  ((float)(etaHigh)>-1.479 && (float)(etaHigh)<1.479)) isbarl=1;
258 
259 
260  etaLow -= regionEtaMargin_;
261  etaHigh += regionEtaMargin_;
262  phiLow -= regionPhiMargin_;
263  phiHigh += regionPhiMargin_;
264 
265  if (etaHigh>1.479) etaHigh=1.479;
266  if (etaLow<-1.479) etaLow=-1.479;
267 
268  if(isbarl) regions.push_back(RectangularEtaPhiRegion(etaLow,etaHigh,phiLow,phiHigh));
269 
270  }
271  }
272  }
273 
274  // make the Basic clusters!
275  reco::BasicClusterCollection basicClusters;
276  hybrid_p->makeClusters(hit_collection, geometry_p, basicClusters, sevLevel, true, regions);
277 
278  // create an unique_ptr to a BasicClusterCollection, copy the clusters into it and put in the Event:
279  auto basicclusters_p = std::make_unique<reco::BasicClusterCollection>();
280  basicclusters_p->assign(basicClusters.begin(), basicClusters.end());
281  edm::OrphanHandle<reco::BasicClusterCollection> bccHandle = evt.put(std::move(basicclusters_p),
283  if (!(bccHandle.isValid())) {
284  return;
285  }
286  reco::BasicClusterCollection clusterCollection = *bccHandle;
287 
288  reco::CaloClusterPtrVector clusterRefVector;
289  for (unsigned int i = 0; i < clusterCollection.size(); i++){
290  clusterRefVector.push_back(reco::CaloClusterPtr(bccHandle, i));
291  }
292 
293  reco::SuperClusterCollection superClusters = hybrid_p->makeSuperClusters(clusterRefVector);
294 
295  auto superclusters_p = std::make_unique<reco::SuperClusterCollection>();
296  superclusters_p->assign(superClusters.begin(), superClusters.end());
297  evt.put(std::move(superclusters_p), superclusterCollection_);
298 }
299 
300 
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
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
double etaBinHighEdge(unsigned int etaIndex, bool central=true) const
CaloTopology const * topology(0)
double etaBinLowEdge(unsigned int etaIndex, bool central=true) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
const edm::EDGetTokenT< l1extra::L1EmParticleCollection > l1TagIsolated_
void produce(edm::Event &, const edm::EventSetup &) override
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
reco::SuperClusterCollection makeSuperClusters(const reco::CaloClusterPtrVector &)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
void setDynamicPhiRoad(const edm::ParameterSet &bremRecoveryPset)
T const * product() const
Definition: Handle.h:74
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
const edm::EDGetTokenT< l1extra::L1EmParticleCollection > l1TagNonIsolated_
int StringToEnumValue(std::string const &enumConstName)
HLT enums.
T get() const
Definition: EventSetup.h:71
double emJetPhiBinLowEdge(unsigned int phiIndex) const
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:86
const edm::EDGetTokenT< EcalRecHitCollection > hittoken_
std::string const & instance() const
Definition: InputTag.h:37
EgammaHLTHybridClusterProducer(const edm::ParameterSet &ps)
def move(src, dest)
Definition: eostools.py:511
double emJetPhiBinHighEdge(unsigned int phiIndex) const