CMS 3D CMS Logo

EgammaHLTMulti5x5ClusterProducer.cc
Go to the documentation of this file.
1 // C/C++ headers
2 #include <iostream>
3 #include <vector>
4 #include <memory>
5 
10 
11 // Reconstruction Classes
16 
17 // Geometry
24 
25 // Level 1 Trigger
30 
31 // EgammaCoreTools
34 
35 // Class header file
37 
39 
40  doBarrel_ = ps.getParameter<bool>("doBarrel");
41  doEndcaps_ = ps.getParameter<bool>("doEndcaps");
42  doIsolated_ = ps.getParameter<bool>("doIsolated");
43 
44  // Parameters to identify the hit collections
45  barrelHitCollection_ = ps.getParameter<edm::InputTag>("barrelHitProducer");
46  endcapHitCollection_ = ps.getParameter<edm::InputTag>("endcapHitProducer");
47  barrelHitToken_ = consumes<EcalRecHitCollection>( barrelHitCollection_);
48  endcapHitToken_ = consumes<EcalRecHitCollection>( endcapHitCollection_);
49 
50  // The names of the produced cluster collections
51  barrelClusterCollection_ = ps.getParameter<std::string>("barrelClusterCollection");
52  endcapClusterCollection_ = ps.getParameter<std::string>("endcapClusterCollection");
53 
54  // Multi5x5 algorithm parameters
55  double barrelSeedThreshold = ps.getParameter<double>("Multi5x5BarrelSeedThr");
56  double endcapSeedThreshold = ps.getParameter<double>("Multi5x5EndcapSeedThr");
57 
58  // L1 matching parameters
59  l1TagIsolated_ = consumes<l1extra::L1EmParticleCollection>(ps.getParameter< edm::InputTag > ("l1TagIsolated"));
60  l1TagNonIsolated_ = consumes<l1extra::L1EmParticleCollection>(ps.getParameter< edm::InputTag > ("l1TagNonIsolated"));
61  l1LowerThr_ = ps.getParameter<double> ("l1LowerThr");
62  l1UpperThr_ = ps.getParameter<double> ("l1UpperThr");
63  l1LowerThrIgnoreIsolation_ = ps.getParameter<double> ("l1LowerThrIgnoreIsolation");
64 
65  regionEtaMargin_ = ps.getParameter<double>("regionEtaMargin");
66  regionPhiMargin_ = ps.getParameter<double>("regionPhiMargin");
67 
68  // Parameters for the position calculation:
69  posCalculator_ = PositionCalc( ps.getParameter<edm::ParameterSet>("posCalcParameters") );
70 
71  const std::vector<std::string> flagnames =
72  ps.getParameter<std::vector<std::string> >("RecHitFlagToBeExcluded");
73 
74  // exclude recHit flags from seeding
75  std::vector<int> v_chstatus = StringToEnumValue<EcalRecHit::Flags>(flagnames);
76 
77  // Produces a collection of barrel and a collection of endcap clusters
78  produces< reco::BasicClusterCollection >(endcapClusterCollection_);
79  produces< reco::BasicClusterCollection >(barrelClusterCollection_);
80 
81  Multi5x5_p = new Multi5x5ClusterAlgo(barrelSeedThreshold, endcapSeedThreshold, v_chstatus, posCalculator_);
82  nEvt_ = 0;
83 }
84 
86  delete Multi5x5_p;
87 }
88 
90 
92  desc.add<bool>(("doBarrel"), false);
93  desc.add<bool>(("doEndcaps"), true);
94  desc.add<bool>(("doIsolated"), true);
95  desc.add<std::string>("VerbosityLevel" ,"ERROR");
96 
97  edm::ParameterSetDescription posCalcPSET;
98  posCalcPSET.add<double>("T0_barl", 7.4);
99  posCalcPSET.add<double>("T0_endc", 3.1);
100  posCalcPSET.add<double>("T0_endcPresh", 1.2);
101  posCalcPSET.add<double>("W0", 4.2);
102  posCalcPSET.add<double>("X0", 0.89);
103  posCalcPSET.add<bool>("LogWeighted", true);
104  desc.add<edm::ParameterSetDescription>("posCalcParameters", posCalcPSET);
105 
106  desc.add<edm::InputTag>(("barrelHitProducer"), edm::InputTag("hltEcalRegionalEgammaRecHit", "EcalRecHitsEB"));
107  desc.add<edm::InputTag>(("endcapHitProducer"), edm::InputTag("hltEcalRegionalEgammaRecHit", "EcalRecHitsEE"));
108  desc.add<std::string>(("barrelClusterCollection"), "notused");
109  desc.add<std::string>(("endcapClusterCollection"), "multi5x5EndcapBasicClusters");
110  desc.add<double>(("Multi5x5BarrelSeedThr"), 0.5);
111  desc.add<double>(("Multi5x5EndcapSeedThr"), 0.5);
112  desc.add<edm::InputTag>(("l1TagIsolated"), edm::InputTag("hltL1extraParticles","Isolated"));
113  desc.add<edm::InputTag>(("l1TagNonIsolated"), edm::InputTag("hltL1extraParticles","NonIsolated"));
114  desc.add<double>(("l1LowerThr"), 5.0);
115  desc.add<double>(("l1UpperThr"), 9999.);
116  desc.add<double>(("l1LowerThrIgnoreIsolation"), 999.0);
117  desc.add<double>(("regionEtaMargin"), 0.3);
118  desc.add<double>(("regionPhiMargin"), 0.4);
119 
120  desc.add<std::vector<std::string> >(("RecHitFlagToBeExcluded"), std::vector<std::string>());
121  descriptions.add(("hltEgammaHLTMulti5x5ClusterProducer"), desc);
122 }
123 
125 
126  //Get the L1 EM Particle Collection
128  if(doIsolated_)
129  evt.getByToken(l1TagIsolated_, emIsolColl);
130 
131  //Get the L1 EM Particle Collection
133  evt.getByToken(l1TagNonIsolated_, emNonIsolColl);
134 
135  // Get the CaloGeometry
136  edm::ESHandle<L1CaloGeometry> l1CaloGeom ;
137  es.get<L1CaloGeometryRecord>().get(l1CaloGeom) ;
138 
139  std::vector<EcalEtaPhiRegion> barrelRegions;
140  std::vector<EcalEtaPhiRegion> endcapRegions;
141 
142  if(doIsolated_) {
143  for( l1extra::L1EmParticleCollection::const_iterator emItr = emIsolColl->begin(); emItr != emIsolColl->end() ;++emItr ){
144 
145  if (emItr->et() > l1LowerThr_ && emItr->et() < l1UpperThr_) {
146 
147  // Access the GCT hardware object corresponding to the L1Extra EM object.
148  int etaIndex = emItr->gctEmCand()->etaIndex() ;
149 
150 
151  int phiIndex = emItr->gctEmCand()->phiIndex() ;
152  // Use the L1CaloGeometry to find the eta, phi bin boundaries.
153  double etaLow = l1CaloGeom->etaBinLowEdge( etaIndex ) ;
154  double etaHigh = l1CaloGeom->etaBinHighEdge( etaIndex ) ;
155  double phiLow = l1CaloGeom->emJetPhiBinLowEdge( phiIndex ) ;
156  double phiHigh = l1CaloGeom->emJetPhiBinHighEdge( phiIndex ) ;
157 
158  //Attention isForward does not work
159  int isforw=0;
160  int isbarl=0;
161  if((float)(etaHigh)>1.479 || (float)(etaLow)<-1.479) isforw=1;
162  if(((float)(etaLow)>-1.479 && (float)(etaLow)<1.479) ||
163  ((float)(etaHigh)>-1.479 && (float)(etaHigh)<1.479)) isbarl=1;
164 
165  //std::cout<<"Multi5x5 etaindex "<<etaIndex<<" low hig : "<<etaLow<<" "<<etaHigh<<" phi low hig" <<phiLow<<" " << phiHigh<<" isforw "<<emItr->gctEmCand()->regionId().isForward()<<" isforwnew" <<isforw<< std::endl;
166 
167  etaLow -= regionEtaMargin_;
168  etaHigh += regionEtaMargin_;
169  phiLow -= regionPhiMargin_;
170  phiHigh += regionPhiMargin_;
171 
172  //if (emItr->gctEmCand()->regionId().isForward()) {
173  if (isforw) {
174  if (etaHigh>-1.479 && etaHigh<1.479) etaHigh=-1.479;
175  if ( etaLow>-1.479 && etaLow<1.479) etaLow=1.479;
176  EcalEtaPhiRegion region(etaLow,etaHigh,phiLow,phiHigh);
177  endcapRegions.push_back(region);
178  }
179  if (isbarl) {
180  if (etaHigh>1.479) etaHigh=1.479;
181  if (etaLow<-1.479) etaLow=-1.479;
182  EcalEtaPhiRegion region(etaLow,etaHigh,phiLow,phiHigh);
183  barrelRegions.push_back(region);
184  }
185  EcalEtaPhiRegion region(etaLow,etaHigh,phiLow,phiHigh);
186 
187  }
188  }
189  }
190 
191 
193  for( l1extra::L1EmParticleCollection::const_iterator emItr = emNonIsolColl->begin(); emItr != emNonIsolColl->end() ;++emItr ){
194 
195  if(doIsolated_&&emItr->et()<l1LowerThrIgnoreIsolation_) continue;
196 
197  if (emItr->et() > l1LowerThr_ && emItr->et() < l1UpperThr_) {
198 
199  // Access the GCT hardware object corresponding to the L1Extra EM object.
200  int etaIndex = emItr->gctEmCand()->etaIndex() ;
201 
202 
203  int phiIndex = emItr->gctEmCand()->phiIndex() ;
204  // Use the L1CaloGeometry to find the eta, phi bin boundaries.
205  double etaLow = l1CaloGeom->etaBinLowEdge( etaIndex ) ;
206  double etaHigh = l1CaloGeom->etaBinHighEdge( etaIndex ) ;
207  double phiLow = l1CaloGeom->emJetPhiBinLowEdge( phiIndex ) ;
208  double phiHigh = l1CaloGeom->emJetPhiBinHighEdge( phiIndex ) ;
209 
210 
211  int isforw=0;
212  int isbarl=0;
213  if((float)(etaHigh)>1.479 || (float)(etaLow)<-1.479) isforw=1;
214  if(((float)(etaLow)>-1.479 && (float)(etaLow)<1.479) ||
215  ((float)(etaHigh)>-1.479 && (float)(etaHigh)<1.479)) isbarl=1;
216 
217  //std::cout<<"Multi5x5 etaindex "<<etaIndex<<" low hig : "<<etaLow<<" "<<etaHigh<<" phi low hig" <<phiLow<<" " << phiHigh<<" isforw "<<emItr->gctEmCand()->regionId().isForward()<<" isforwnew" <<isforw<< std::endl;
218 
219  etaLow -= regionEtaMargin_;
220  etaHigh += regionEtaMargin_;
221  phiLow -= regionPhiMargin_;
222  phiHigh += regionPhiMargin_;
223 
224  //if (emItr->gctEmCand()->regionId().isForward()) {
225  if (isforw) {
226  if (etaHigh>-1.479 && etaHigh<1.479) etaHigh=-1.479;
227  if ( etaLow>-1.479 && etaLow<1.479) etaLow=1.479;
228  EcalEtaPhiRegion region(etaLow,etaHigh,phiLow,phiHigh);
229  endcapRegions.push_back(region);
230  }
231  if (isbarl) {
232  if (etaHigh>1.479) etaHigh=1.479;
233  if (etaLow<-1.479) etaLow=-1.479;
234  EcalEtaPhiRegion region(etaLow,etaHigh,phiLow,phiHigh);
235  barrelRegions.push_back(region);
236  }
237 
238  }
239  }
240  }
241 
242 
243  if (doEndcaps_) {
245  }
246  if (doBarrel_) {
248  }
249  nEvt_++;
250 }
251 
252 
255 
257  evt.getByToken(hitToken, rhcHandle);
258 
259  if (!(rhcHandle.isValid()))
260  {
261  std::cout << "could not get a handle on the EcalRecHitCollection!" << std::endl;
262  edm::LogError("EgammaHLTMulti5x5ClusterProducerError") << "Error! can't get the product ";
263  return 0;
264  }
265  return rhcHandle.product();
266 }
267 
268 
271  const std::string& clusterCollection,
272  const std::vector<EcalEtaPhiRegion>& regions,
274 
275  // get the hit collection from the event:
276  const EcalRecHitCollection *hitCollection_p = getCollection(evt, hitToken);
277 
278  // get the geometry and topology from the event setup:
279  edm::ESHandle<CaloGeometry> geoHandle;
280  es.get<CaloGeometryRecord>().get(geoHandle);
281 
282  const CaloSubdetectorGeometry *geometry_p;
283  CaloSubdetectorTopology *topology_p;
284 
285  if (detector == reco::CaloID::DET_ECAL_BARREL)
286  {
287  geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
288  topology_p = new EcalBarrelTopology(geoHandle);
289  }
290  else
291  {
292  geometry_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalEndcap);
293  topology_p = new EcalEndcapTopology(geoHandle);
294  }
295 
296 
297  const CaloSubdetectorGeometry *geometryES_p;
298  geometryES_p = geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalPreshower);
299 
300  // Run the clusterization algorithm:
302  clusters = Multi5x5_p->makeClusters(hitCollection_p, geometry_p, topology_p, geometryES_p, detector, true, regions);
303 
304  // create an unique_ptr to a BasicClusterCollection, copy the barrel clusters into it and put in the Event:
305  auto clusters_p = std::make_unique<reco::BasicClusterCollection>();
306  clusters_p->assign(clusters.begin(), clusters.end());
308  if (detector == reco::CaloID::DET_ECAL_BARREL)
309  bccHandle = evt.put(std::move(clusters_p), barrelClusterCollection_);
310  else
311  bccHandle = evt.put(std::move(clusters_p), endcapClusterCollection_);
312 
313  delete topology_p;
314 }
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:45
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
void clusterizeECALPart(edm::Event &evt, const edm::EventSetup &es, edm::EDGetTokenT< EcalRecHitCollection > &hitToken, const std::string &clusterCollection, const std::vector< EcalEtaPhiRegion > &regions, const reco::CaloID::Detectors detector)
std::vector< reco::BasicCluster > makeClusters(const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry, const CaloSubdetectorTopology *topology_p, const CaloSubdetectorGeometry *geometryES_p, reco::CaloID::Detectors detector, bool regional=false, const std::vector< EcalEtaPhiRegion > &regions=std::vector< EcalEtaPhiRegion >())
double etaBinHighEdge(unsigned int etaIndex, bool central=true) const
double etaBinLowEdge(unsigned int etaIndex, bool central=true) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
edm::EDGetTokenT< l1extra::L1EmParticleCollection > l1TagNonIsolated_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< EcalRecHitCollection > barrelHitToken_
edm::EDGetTokenT< EcalRecHitCollection > endcapHitToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:74
const EcalRecHitCollection * getCollection(edm::Event &evt, edm::EDGetTokenT< EcalRecHitCollection > &hitToken)
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:56
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< BasicCluster > BasicClusterCollection
collection of BasicCluster objects
void produce(edm::Event &, const edm::EventSetup &) override
double emJetPhiBinLowEdge(unsigned int phiIndex) const
def move(src, dest)
Definition: eostools.py:510
double emJetPhiBinHighEdge(unsigned int phiIndex) const
EgammaHLTMulti5x5ClusterProducer(const edm::ParameterSet &ps)
edm::EDGetTokenT< l1extra::L1EmParticleCollection > l1TagIsolated_