CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalHaloAlgo.cc
Go to the documentation of this file.
3 
4 /*
5  [class]: EcalHaloAlgo
6  [authors]: R. Remington, The University of Florida
7  [description]: See EcalHaloAlgo.h
8  [date]: October 15, 2009
9 */
10 namespace {
11  constexpr float c_cm_per_ns = 29.9792458;
12 };
13 using namespace std;
14 using namespace reco;
15 using namespace edm;
16 
17 bool CompareTime(const EcalRecHit* x, const EcalRecHit* y){ return x->time() < y->time(); }
18 
20 {
21  RoundnessCut = 0 ;
22  AngleCut = 0 ;
23  EBRecHitEnergyThreshold = 0.;
24  EERecHitEnergyThreshold = 0.;
25  ESRecHitEnergyThreshold = 0.;
26  SumEnergyThreshold = 0.;
27  NHitsThreshold =0;
28 
29  geo = 0;
30 }
31 
33 {
34 
35  EcalHaloData TheEcalHaloData;
36 
37  // Store energy sum of rechits as a function of iPhi (iphi goes from 1 to 72)
38  float SumE[361];
39  // Store number of rechits as a function of iPhi
40  int NumHits[361];
41  // Store minimum time of rechit as a function of iPhi
42  float MinTimeHits[361];
43  // Store maximum time of rechit as a function of iPhi
44  float MaxTimeHits[361];
45 
46  // initialize
47  for(int i = 0 ; i < 361 ; i++ )
48  {
49  SumE[i] = 0.;
50  NumHits[i] = 0 ;
51  MinTimeHits[i] = 9999.;
52  MaxTimeHits[i] = -9999.;
53  }
54 
55  // Loop over EB RecHits
56  for(EBRecHitCollection::const_iterator hit = TheEBRecHits->begin() ; hit != TheEBRecHits->end() ; hit++ )
57  {
58  // Arbitrary threshold to kill noise (needs to be optimized with data)
59  if (hit->energy() < EBRecHitEnergyThreshold ) continue;
60 
61  // Get Det Id of the rechit
62  DetId id = DetId(hit->id());
63  const CaloSubdetectorGeometry* TheSubGeometry = 0;
64  const CaloCellGeometry* cell = 0 ;
65 
66  // Get EB geometry
67  TheSubGeometry = TheCaloGeometry.getSubdetectorGeometry(DetId::Ecal, 1);
68  EBDetId EcalID(id.rawId());
69  if( TheSubGeometry )
70  cell = TheSubGeometry->getGeometry(id);
71 
72  if(cell)
73  {
74  // GlobalPoint globalpos = cell->getPosition();
75  // float r = TMath::Sqrt ( globalpos.y()*globalpos.y() + globalpos.x()*globalpos.x());
76  int iPhi = EcalID.iphi();
77 
78  if( iPhi < 361 ) // just to be safe
79  {
80  //iPhi = (iPhi-1)/5 +1; // convert ecal iphi to phiwedge iphi (e.g. there are 5 crystal per phi wedge, as in calotowers )
81  SumE[iPhi] += hit->energy();
82  NumHits[iPhi] ++;
83 
84  float time = hit->time();
85  MinTimeHits[iPhi] = time < MinTimeHits[iPhi] ? time : MinTimeHits[iPhi];
86  MaxTimeHits[iPhi] = time > MaxTimeHits[iPhi] ? time : MaxTimeHits[iPhi];
87  }
88  }
89  }
90 
91  //for( int iPhi = 1 ; iPhi < 73; iPhi++ )
92  for( int iPhi = 1 ; iPhi < 361; iPhi++ )
93  {
94  if( SumE[iPhi] >= SumEnergyThreshold && NumHits[iPhi] > NHitsThreshold )
95  {
96  // Build PhiWedge and store to EcalHaloData if energy or #hits pass thresholds
97  PhiWedge wedge(SumE[iPhi], iPhi, NumHits[iPhi], MinTimeHits[iPhi], MaxTimeHits[iPhi]);
98 
99  // Loop over rechits again to calculate direction based on timing info
100 
101  // Loop over EB RecHits
102  std::vector<const EcalRecHit*> Hits;
103  for(EBRecHitCollection::const_iterator hit = TheEBRecHits->begin() ; hit != TheEBRecHits->end() ; hit++ )
104  {
105  if (hit->energy() < EBRecHitEnergyThreshold ) continue;
106 
107  // Get Det Id of the rechit
108  DetId id = DetId(hit->id());
109  EBDetId EcalID(id.rawId());
110  int Hit_iPhi = EcalID.iphi();
111  //Hit_iPhi = (Hit_iPhi-1)/5 +1; // convert ecal iphi to phiwedge iphi
112  if( Hit_iPhi != iPhi ) continue;
113  Hits.push_back( &(*hit) );
114 
115  }
116  std::sort( Hits.begin() , Hits.end(), CompareTime);
117  float MinusToPlus = 0.;
118  float PlusToMinus = 0.;
119  for( unsigned int i = 0 ; i < Hits.size() ; i++ )
120  {
121  DetId id_i = DetId(Hits[i]->id());
122  EBDetId EcalID_i(id_i.rawId());
123  int ieta_i = EcalID_i.ieta();
124  for( unsigned int j = (i+1) ; j < Hits.size() ; j++ )
125  {
126  DetId id_j = DetId(Hits[j]->id() );
127  EBDetId EcalID_j(id_j.rawId());
128  int ieta_j = EcalID_j.ieta();
129  if( ieta_i > ieta_j ) PlusToMinus += TMath::Abs(ieta_j - ieta_i );
130  else MinusToPlus += TMath::Abs(ieta_j -ieta_i) ;
131  }
132  }
133 
134  float PlusZOriginConfidence = (PlusToMinus+MinusToPlus) ? PlusToMinus / (PlusToMinus+MinusToPlus) : -1.;
135  wedge.SetPlusZOriginConfidence(PlusZOriginConfidence);
136  TheEcalHaloData.GetPhiWedges().push_back(wedge);
137  }
138  }
139 
140  std::vector<float> vShowerShapes_Roundness;
141  std::vector<float> vShowerShapes_Angle ;
142  if(TheSuperClusters.isValid()){
143  for(reco::SuperClusterCollection::const_iterator cluster = TheSuperClusters->begin() ; cluster != TheSuperClusters->end() ; cluster++ )
144  {
145  if( abs(cluster->eta()) <= 1.48 )
146  {
147  vector<float> shapes = EcalClusterTools::roundnessBarrelSuperClusters( *cluster, (*TheEBRecHits.product()));
148  float roundness = shapes[0];
149  float angle = shapes[1];
150 
151  // Check if supercluster belongs to photon and passes the cuts on roundness and angle, if so store the reference to it
152  if( (roundness >=0 && roundness < GetRoundnessCut()) && angle >= 0 && angle < GetAngleCut() )
153  {
154  edm::Ref<SuperClusterCollection> TheClusterRef( TheSuperClusters, cluster - TheSuperClusters->begin());
155  bool BelongsToPhoton = false;
156  if( ThePhotons.isValid() )
157  {
158  for(reco::PhotonCollection::const_iterator iPhoton = ThePhotons->begin() ; iPhoton != ThePhotons->end() ; iPhoton++ )
159  {
160  if(iPhoton->isEB())
161  if ( TheClusterRef == iPhoton->superCluster() )
162  {
163  BelongsToPhoton = true;
164  break;
165  }
166  }
167  }
168  //Only store refs to suspicious EB SuperClusters which belong to Photons
169  //Showershape variables are more discriminating for these cases
170  if( BelongsToPhoton )
171  {
172  TheEcalHaloData.GetSuperClusters().push_back( TheClusterRef ) ;
173  }
174  }
175  vShowerShapes_Roundness.push_back(shapes[0]);
176  vShowerShapes_Angle.push_back(shapes[1]);
177  }
178  else
179  {
180  vShowerShapes_Roundness.push_back(-1.);
181  vShowerShapes_Angle.push_back(-1.);
182  }
183  }
184 
185 
186  edm::ValueMap<float>::Filler TheRoundnessFiller( TheEcalHaloData.GetShowerShapesRoundness() );
187  TheRoundnessFiller.insert( TheSuperClusters, vShowerShapes_Roundness.begin(), vShowerShapes_Roundness.end() );
188  TheRoundnessFiller.fill();
189 
190  edm::ValueMap<float>::Filler TheAngleFiller( TheEcalHaloData.GetShowerShapesAngle() );
191  TheAngleFiller.insert( TheSuperClusters, vShowerShapes_Angle.begin() , vShowerShapes_Angle.end() );
192  TheAngleFiller.fill();
193  }
194 
195 
196 
197  geo = 0;
199  TheSetup.get<CaloGeometryRecord>().get(pGeo);
200  geo = pGeo.product();
201 
202  //Halo cluster building:
203  //Various clusters are built, depending on the subdetector.
204  //In barrel, one looks for deposits narrow in phi.
205  //In endcaps, one looks for localized deposits (dr condition in EE where r =sqrt(dphi*dphi+deta*deta)
206  //H/E condition is also applied in EB.
207  //The halo cluster building step targets a large efficiency (ideally >99%) for beam halo deposits.
208  //These clusters are used as input for the halo pattern finding methods in EcalHaloAlgo and for the CSC-calo matching methods in GlobalHaloAlgo.
209 
210  //Et threshold hardcoded for now. Might one to get it from config
211  std::vector<HaloClusterCandidateECAL> haloclustercands_EB;
212  haloclustercands_EB= GetHaloClusterCandidateEB(TheEBRecHits , TheHBHERecHits, 5);
213 
214  std::vector<HaloClusterCandidateECAL> haloclustercands_EE;
215  haloclustercands_EE= GetHaloClusterCandidateEE(TheEERecHits , TheHBHERecHits, 10);
216 
217  TheEcalHaloData.setHaloClusterCandidatesEB(haloclustercands_EB);
218  TheEcalHaloData.setHaloClusterCandidatesEE(haloclustercands_EE);
219 
220  return TheEcalHaloData;
221 }
222 
223 
224 
225 
226 std::vector<HaloClusterCandidateECAL> EcalHaloAlgo::GetHaloClusterCandidateEB(edm::Handle<EcalRecHitCollection>& ecalrechitcoll, edm::Handle<HBHERecHitCollection>& hbherechitcoll,float et_thresh_seedrh){
227 
228  std::vector<HaloClusterCandidateECAL> TheHaloClusterCandsEB;
229  reco::Vertex::Point vtx(0,0,0);
230 
231  for(size_t ihit = 0; ihit<ecalrechitcoll->size(); ++ ihit){
232  HaloClusterCandidateECAL clustercand;
233 
234  const EcalRecHit & rechit = (*ecalrechitcoll)[ ihit ];
235  math::XYZPoint rhpos = getPosition(rechit.id(),vtx);
236  //Et condition
237 
238  double rhet = rechit.energy() * sqrt(rhpos.perp2()/rhpos.mag2());
239  if(rhet<et_thresh_seedrh) continue;
240  double eta = rhpos.eta();
241  double phi = rhpos.phi();
242 
243  bool isiso = true;
244  double etcluster(0);
245  int nbcrystalsameeta(0);
246  double timediscriminator(0);
247  double etstrip_iphiseedplus1(0), etstrip_iphiseedminus1(0);
248 
249  //Building the cluster
251  for(size_t jhit = 0; jhit<ecalrechitcoll->size(); ++ jhit){
252  const EcalRecHit & rechitj = (*ecalrechitcoll)[ jhit ];
253  EcalRecHitRef rhRef(ecalrechitcoll,jhit);
254  math::XYZPoint rhposj = getPosition(rechitj.id(),vtx);
255 
256  double etaj = rhposj.eta();
257  double phij = rhposj.phi();
258 
259  double deta = eta - etaj;
260  double dphi = deltaPhi(phi,phij);
261  if(std::abs(deta)>0.2) continue;//This means +/-11 crystals in eta
262  if(std::abs(dphi)>0.08) continue;//This means +/-4 crystals in phi
263 
264  double rhetj = rechitj.energy()* sqrt(rhposj.perp2()/rhposj.mag2());
265  //Rechits with et between 1 and 2 GeV are saved in the rh list but not used in the calculation of the halocluster variables
266  if(rhetj<1) continue;
267  bhrhcandidates.push_back(rhRef);
268  if(rhetj<2) continue;
269 
270  if(std::abs(dphi)>0.03){isiso=false;break;}//The strip should be isolated
271  if(std::abs(dphi)<0.01) nbcrystalsameeta++;
272  if(dphi>0.01) etstrip_iphiseedplus1+=rhetj;
273  if(dphi<-0.01) etstrip_iphiseedminus1+=rhetj;
274  etcluster+=rhetj;
275  //Timing discriminator
276  //We assign a weight to the rechit defined as:
277  //Log10(Et)*f(T,R,Z)
278  //where f(T,R,Z) is the separation curve between halo-like and IP-like times.
279  //The time difference between a deposit from a outgoing IT halo and a deposit coming from a particle emitted at the IP is given by:
280  //dt= ( - sqrt(R^2+z^2) + |z| )/c
281  //Here we take R to be 130 cm.
282  //For EB, the function was parametrized as a function of ieta instead of Z.
283  double rhtj = rechitj.time();
284  EBDetId detj = rechitj.id();
285  int rhietaj= detj.ieta();
286  timediscriminator+= std::log10( rhetj )* ( rhtj +0.5*(sqrt(16900+9*rhietaj*rhietaj)-3*std::abs(rhietaj))/c_cm_per_ns );
287 
288  }
289  //Isolation condition
290  if(!isiso) continue;
291 
292  //Calculate H/E
293  double hoe(0);
294  for(size_t jhit = 0; jhit<hbherechitcoll->size(); ++ jhit){
295  const HBHERecHit & rechitj = (*hbherechitcoll)[ jhit ];
296  math::XYZPoint rhposj = getPosition(rechitj.id(),vtx);
297  double rhetj = rechitj.energy()* sqrt(rhposj.perp2()/rhposj.mag2());
298  if(rhetj<2) continue;
299  double etaj = rhposj.eta();
300  double phij = rhposj.phi();
301  double deta = eta - etaj;
302  double dphi = deltaPhi(phi,phij);
303  if(std::abs(deta)>0.2) continue;
304  if(std::abs(dphi)>0.2) continue;
305  hoe+=rhetj/etcluster;
306  }
307  //H/E condition
308  if(hoe>0.1) continue;
309 
310 
311  clustercand.setClusterEt(etcluster);
312  clustercand.setSeedEt(rhet);
313  clustercand.setSeedEta(eta);
314  clustercand.setSeedPhi(phi);
315  clustercand.setSeedZ(rhpos.Z());
316  clustercand.setSeedR(sqrt(rhpos.perp2()));
317  clustercand.setSeedTime(rechit.time());
318  clustercand.setHoverE(hoe);
319  clustercand.setNbofCrystalsInEta(nbcrystalsameeta);
320  clustercand.setEtStripIPhiSeedPlus1(etstrip_iphiseedplus1);
321  clustercand.setEtStripIPhiSeedMinus1(etstrip_iphiseedminus1);
322  clustercand.setTimeDiscriminator(timediscriminator);
323  clustercand.setBeamHaloRecHitsCandidates(bhrhcandidates);
324 
325 
326  bool isbeamhalofrompattern = EBClusterShapeandTimeStudy(clustercand,false);
327  clustercand.setIsHaloFromPattern(isbeamhalofrompattern);
328 
329  bool isbeamhalofrompattern_hlt = EBClusterShapeandTimeStudy(clustercand,true);
330  clustercand.setIsHaloFromPattern_HLT(isbeamhalofrompattern_hlt);
331 
332 
333  TheHaloClusterCandsEB.push_back(clustercand);
334  }
335 
336  return TheHaloClusterCandsEB;
337 }
338 
339 
340 
341 
342 std::vector<HaloClusterCandidateECAL> EcalHaloAlgo::GetHaloClusterCandidateEE(edm::Handle<EcalRecHitCollection>& ecalrechitcoll, edm::Handle<HBHERecHitCollection>& hbherechitcoll,float et_thresh_seedrh){
343 
344  std::vector<HaloClusterCandidateECAL> TheHaloClusterCandsEE;
345 
346  reco::Vertex::Point vtx(0,0,0);
347 
348  for(size_t ihit = 0; ihit<ecalrechitcoll->size(); ++ ihit){
349  HaloClusterCandidateECAL clustercand;
350 
351  const EcalRecHit & rechit = (*ecalrechitcoll)[ ihit ];
352  math::XYZPoint rhpos = getPosition(rechit.id(),vtx);
353  //Et condition
354  double rhet = rechit.energy()* sqrt(rhpos.perp2()/rhpos.mag2());
355  if(rhet<et_thresh_seedrh) continue;
356  double eta = rhpos.eta();
357  double phi = rhpos.phi();
358  double rhr = sqrt(rhpos.perp2());
359 
360  bool isiso = true;
361  double etcluster(0);
362  double timediscriminator(0);
363  int clustersize(0);
364  int nbcrystalssmallt(0);
365  int nbcrystalshight(0);
366  //Building the cluster
368  for(size_t jhit = 0; jhit<ecalrechitcoll->size(); ++ jhit){
369  const EcalRecHit & rechitj = (*ecalrechitcoll)[ jhit ];
370  EcalRecHitRef rhRef(ecalrechitcoll,jhit);
371  math::XYZPoint rhposj = getPosition(rechitj.id(),vtx);
372 
373  //Ask the hits to be in the same endcap
374  if(rhposj.z()*rhpos.z()<0)continue;
375 
376  double etaj = rhposj.eta();
377  double phij = rhposj.phi();
378  double dr = sqrt((eta-etaj)*(eta-etaj)+deltaPhi(phi,phij)*deltaPhi(phi,phij));
379 
380  //Outer cone
381  if(dr>0.3) continue;
382 
383  double rhetj = rechitj.energy()* sqrt(rhposj.perp2()/rhposj.mag2());
384  //Rechits with et between 1 and 2 GeV are saved in the rh list but not used in the calculation of the halocluster variables
385  if(rhetj<1) continue;
386  bhrhcandidates.push_back(rhRef);
387  if(rhetj<2) continue;
388 
389  //Isolation between outer and inner cone
390  if(dr>0.05){isiso=false;break;}//The deposit should be isolated
391 
392  etcluster+=rhetj;
393 
394  //Timing infos:
395  //Here we target both IT and OT beam halo
396  double rhtj=rechitj.time();
397 
398  //Discriminating variables for OT beam halo:
399  if(rhtj>1) nbcrystalshight++;
400  if(rhtj<0) nbcrystalssmallt++;
401  //Timing test (likelihood ratio), only for seeds with large R (100 cm) and for crystals with et>5,
402  //This targets IT beam halo (t around - 1ns)
403  if(rhtj>5){
404  double corrt_j = rhtj + sqrt(rhposj.x()*rhposj.x()+rhposj.y()*rhposj.y() + 320.*320.)/c_cm_per_ns - 320./c_cm_per_ns;
405  //BH is modeled by a Gaussian peaking at 0.
406  //Collisions is modeled by a Gaussian peaking at 0.3
407  //The width is similar and taken to be 0.4
408  timediscriminator+= 0.5*(pow( (corrt_j-0.3)/0.4,2)-pow( (corrt_j-0.)/0.4,2));
409  clustersize++;
410  }
411 
412  }
413  //Isolation condition
414  if(!isiso) continue;
415 
416  //Calculate H2/E
417  //Only second hcal layer is considered as it can happen that a shower initiated in EE reaches HCAL first layer
418  double h2oe(0);
419  for(size_t jhit = 0; jhit<hbherechitcoll->size(); ++ jhit){
420  const HBHERecHit & rechitj = (*hbherechitcoll)[ jhit ];
421  math::XYZPoint rhposj = getPosition(rechitj.id(),vtx);
422 
423  //Ask the hits to be in the same endcap
424  if(rhposj.z()*rhpos.z()<0)continue;
425  //Selects only second HCAL layer
426  if(std::abs(rhposj.z())<425) continue;
427 
428  double rhetj = rechitj.energy()* sqrt(rhposj.perp2()/rhposj.mag2());
429  if(rhetj<2) continue;
430 
431  double phij = rhposj.phi();
432  if(std::abs(deltaPhi(phi,phij))>0.4 ) continue;
433 
434  double rhrj = sqrt(rhposj.perp2());
435  if(std::abs(rhr-rhrj)>50) continue;
436 
437  h2oe+=rhetj/etcluster;
438  }
439  //H/E condition
440  if(h2oe>0.1) continue;
441 
442 
443  clustercand.setClusterEt(etcluster);
444  clustercand.setSeedEt(rhet);
445  clustercand.setSeedEta(eta);
446  clustercand.setSeedPhi(phi);
447  clustercand.setSeedZ(rhpos.Z());
448  clustercand.setSeedR(sqrt(rhpos.perp2()));
449  clustercand.setSeedTime(rechit.time());
450  clustercand.setH2overE(h2oe);
451  clustercand.setNbEarlyCrystals(nbcrystalssmallt);
452  clustercand.setNbLateCrystals(nbcrystalshight);
453  clustercand.setClusterSize(clustersize);
454  clustercand.setTimeDiscriminator(timediscriminator);
455  clustercand.setBeamHaloRecHitsCandidates(bhrhcandidates);
456 
457  bool isbeamhalofrompattern = EEClusterShapeandTimeStudy_ITBH(clustercand,false) || EEClusterShapeandTimeStudy_OTBH(clustercand,false);
458  clustercand.setIsHaloFromPattern(isbeamhalofrompattern);
459 
460  bool isbeamhalofrompattern_hlt = EEClusterShapeandTimeStudy_ITBH(clustercand,true) || EEClusterShapeandTimeStudy_OTBH(clustercand,true);
461  clustercand.setIsHaloFromPattern_HLT(isbeamhalofrompattern_hlt);
462 
463  TheHaloClusterCandsEE.push_back(clustercand);
464  }
465 
466  return TheHaloClusterCandsEE;
467 }
468 
469 
470 
471 
472 
474  //Conditions on the central strip size in eta.
475  //For low size, extra conditions on seed et, isolation and cluster timing
476  //The time condition only targets IT beam halo.
477  //EB rechits from OT beam halos are typically too late (around 5 ns or more) and seem therefore already cleaned by the reconstruction.
478 
479  if(hcand.getSeedEt()<5)return false;
480  if(hcand.getNbofCrystalsInEta()<4) return false;
481  if(hcand.getNbofCrystalsInEta()==4&&hcand.getSeedEt()<10) return false;
482  if(hcand.getNbofCrystalsInEta()==4 && hcand.getEtStripIPhiSeedPlus1()>0.1 &&hcand.getEtStripIPhiSeedMinus1()>0.1 ) return false;
483  if(hcand.getNbofCrystalsInEta()<=5 && hcand.getTimeDiscriminator()>=0.)return false;
484 
485  //For HLT, only use conditions without timing and tighten seed et condition
486  if(ishlt &&hcand.getNbofCrystalsInEta()<=5)return false;
487  if(ishlt && hcand.getSeedEt()<10)return false;
488 
489  hcand.setIsHaloFromPattern(true);
490 
491  return true;
492 }
493 
494 
495 
497  //Separate conditions targeting IT and OT beam halos
498  //For OT beam halos, just require enough crystals with large T
499  if(hcand.getSeedEt()<20)return false;
500  if(hcand.getSeedTime()<0.5)return false;
501  if(hcand.getNbLateCrystals()-hcand.getNbEarlyCrystals() <2)return false;
502 
503  //The use of time information does not allow this method to work at HLT
504  if(ishlt)return false;
505 
506  hcand.setIsHaloFromPattern(true);
507 
508  return true;
509 }
510 
512  //Separate conditions targeting IT and OT beam halos
513  //For IT beam halos, fakes from collisions are higher => require the cluster size to be small.
514  //Only halos with R>100 cm are considered here.
515  //For lower values, the time difference with particles from collisions is too small
516  //IT outgoing beam halos that interact in EE at low R is probably the most difficult category to deal with:
517  //Their signature is very close to the one of photon from collisions (similar cluster shape and timing)
518  if(hcand.getSeedEt()<20)return false;
519  if(hcand.getSeedR()<100)return false;
520  if(hcand.getTimeDiscriminator()<1) return false;
521  if(hcand.getClusterSize()<2) return false;
522  if(hcand.getClusterSize()>4) return false;
523 
524  //The use of time information does not allow this method to work at HLT
525  if(ishlt)return false;
526 
527  hcand.setIsHaloFromPattern(true);
528 
529  return true;
530 }
531 
532 
534 
535  const GlobalPoint& pos=geo->getPosition(id);
536  math::XYZPoint posV(pos.x() - vtx.x(),pos.y() - vtx.y(),pos.z() - vtx.z());
537  return posV;
538 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
int i
Definition: DBlmapReader.cc:9
edm::ValueMap< float > & GetShowerShapesRoundness()
Definition: EcalHaloData.h:43
const std::vector< PhiWedge > & GetPhiWedges() const
Definition: EcalHaloData.h:34
edm::RefVector< reco::SuperClusterCollection > & GetSuperClusters()
Definition: EcalHaloData.h:38
edm::ValueMap< float > & GetShowerShapesAngle()
Definition: EcalHaloData.h:46
HcalDetId id() const
get the id
Definition: HBHERecHit.h:23
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
std::vector< EcalRecHit >::const_iterator const_iterator
T y() const
Definition: PV3DBase.h:63
float time() const
Definition: EcalRecHit.h:70
#define constexpr
bool EEClusterShapeandTimeStudy_OTBH(reco::HaloClusterCandidateECAL hcand, bool ishlt)
int iphi() const
get the crystal iphi
Definition: EBDetId.h:53
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
float energy() const
Definition: CaloRecHit.h:17
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
T Abs(T a)
Definition: MathUtil.h:49
void setHaloClusterCandidatesEE(const std::vector< HaloClusterCandidateECAL > &x)
Definition: EcalHaloData.h:56
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
void SetPlusZOriginConfidence(float x)
Definition: PhiWedge.h:67
float energy() const
Definition: EcalRecHit.h:68
math::XYZPoint Point
point in the space
Definition: Vertex.h:39
void setHaloClusterCandidatesEB(const std::vector< HaloClusterCandidateECAL > &x)
Definition: EcalHaloData.h:55
std::vector< reco::HaloClusterCandidateECAL > GetHaloClusterCandidateEB(edm::Handle< EcalRecHitCollection > &ecalrechitcoll, edm::Handle< HBHERecHitCollection > &hbherechitcoll, float et_thresh_seedrh)
int ieta() const
get the crystal ieta
Definition: EBDetId.h:51
bool isValid() const
Definition: HandleBase.h:75
unsigned int id
std::vector< reco::HaloClusterCandidateECAL > GetHaloClusterCandidateEE(edm::Handle< EcalRecHitCollection > &ecalrechitcoll, edm::Handle< HBHERecHitCollection > &hbherechitcoll, float et_thresh_seedrh)
Definition: DetId.h:18
DetId id() const
get the id
Definition: EcalRecHit.h:76
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool EEClusterShapeandTimeStudy_ITBH(reco::HaloClusterCandidateECAL hcand, bool ishlt)
reco::EcalHaloData Calculate(const CaloGeometry &TheCaloGeometry, edm::Handle< reco::PhotonCollection > &ThePhotons, edm::Handle< reco::SuperClusterCollection > &TheSuperClusters, edm::Handle< EBRecHitCollection > &TheEBRecHits, edm::Handle< EERecHitCollection > &TheEERecHits, edm::Handle< ESRecHitCollection > &TheESRecHits, edm::Handle< HBHERecHitCollection > &TheHBHERecHits, const edm::EventSetup &TheSetup)
Definition: EcalHaloAlgo.cc:32
Geom::Phi< T > phi() const
math::XYZPoint getPosition(const DetId &id, reco::Vertex::Point vtx)
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:62
bool CompareTime(const EcalRecHit *x, const EcalRecHit *y)
Definition: EcalHaloAlgo.cc:17
void setBeamHaloRecHitsCandidates(edm::RefVector< EcalRecHitCollection > x)
bool EBClusterShapeandTimeStudy(reco::HaloClusterCandidateECAL hcand, bool ishlt)
T x() const
Definition: PV3DBase.h:62
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11