CMS 3D CMS Logo

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