CMS 3D CMS Logo

PreshowerPhiClusterAlgo.cc
Go to the documentation of this file.
1 #include <vector>
2 #include <map>
3 #include <iostream>
4 
12 
13 #include "TH1.h"
14 
16  HitsID *used_strips,
17  RecHitsMap *the_rechitsMap_p,
18  const CaloSubdetectorGeometry* geometry_p,
19  double deltaEta, double minDeltaPhi, double maxDeltaPhi)
20 {
21 
22  rechits_map = the_rechitsMap_p;
23  used_s = used_strips;
24 
25  int plane = strip.plane();
26 
27  // create null-cluster
28  std::vector< std::pair<DetId,float> > dummy;
29  Point posi(0,0,0);
30  reco::PreshowerCluster nullcluster = reco::PreshowerCluster(0., posi, dummy, plane);
31  if (strip == ESDetId(0)) return nullcluster; // works in case of no intersected strip found (e.g. in the Barrel)
32 
33  auto refCell = geometry_p->getGeometry(strip);
34  const GlobalPoint& refpos = refCell->getPosition();
35  double refEta = refpos.eta();
36  double refPhi = refpos.phi();
37 
38  // Collection of cluster strips
39  EcalRecHitCollection clusterRecHits;
40 
41  double x_pos = 0;
42  double y_pos = 0;
43  double z_pos = 0;
44 
45  RecHitsMap::iterator strip_it;
46  for (strip_it = rechits_map->begin(); strip_it != rechits_map->end(); ++strip_it) {
47 
48  if (!goodStrip(strip_it)) continue;
49 
50  ESDetId mystrip = (strip_it->first == DetId(0)) ? ESDetId(0) : ESDetId(strip_it->first);
51  if (mystrip.plane() != plane) continue;
52 
53  auto thisCell = geometry_p->getGeometry(strip_it->first);
54  const GlobalPoint& position = thisCell->getPosition();
55 
56  if (fabs(position.eta() - refEta) < deltaEta) {
57 
58  //std::cout<<"all strips : "<<mystrip.plane()<<" "<<position.phi()<<" "<<reco::deltaPhi(position.phi(), refPhi)<<std::endl;
59 
60  if (reco::deltaPhi(position.phi(), refPhi) > 0 && reco::deltaPhi(position.phi(), refPhi) > maxDeltaPhi) continue;
61  if (reco::deltaPhi(position.phi(), refPhi) < 0 && reco::deltaPhi(position.phi(), refPhi) < minDeltaPhi) continue;
62 
63  //std::cout<<mystrip.zside()<<" "<<mystrip.plane()<<" "<<mystrip.six()<<" "<<mystrip.siy()<<" "<<mystrip.strip()<<" "<<position.eta()<<" "<<position.phi()<<" "<<strip_it->second.energy()<<" "<<strip_it->second.recoFlag()<<" "<<refEta<<" "<<refPhi<<" "<<reco::deltaPhi(position.phi(), refPhi)<<std::endl;
64 
65  clusterRecHits.push_back(strip_it->second);
66  used_s->insert(strip_it->first);
67 
68  x_pos += strip_it->second.energy() * position.x();
69  y_pos += strip_it->second.energy() * position.y();
70  z_pos += strip_it->second.energy() * position.z();
71  }
72 
73  }
74 
76  double Eclust = 0;
77 
78  std::vector<std::pair<DetId,float > > usedHits;
79  for (it=clusterRecHits.begin(); it != clusterRecHits.end(); it++) {
80  Eclust += it->energy();
81  usedHits.push_back(std::pair<DetId,float > (it->id(), 1.));
82  }
83 
84  if (Eclust > 0.) {
85  x_pos /= Eclust;
86  y_pos /= Eclust;
87  z_pos /= Eclust;
88  }
89  Point pos(x_pos,y_pos,z_pos);
90 
91  reco::PreshowerCluster cluster = reco::PreshowerCluster(Eclust, pos, usedHits, plane);
92 
93  return cluster;
94 }
95 
96 // returns true if the candidate strip fulfills the requirements to be added to the cluster:
97 bool PreshowerPhiClusterAlgo::goodStrip(RecHitsMap::iterator candidate_it) {
98 
99  if ( used_s->find(candidate_it->first) != used_s->end())
100  LogTrace("PreShowerPhiClusterAlgo") << " This strip is in use";
101  if (candidate_it == rechits_map->end() )
102  LogTrace("PreShowerPhiClusterAlgo") << " No such a strip in rechits_map";
103  if (candidate_it->second.energy() <= esStripEnergyCut_)
104  LogTrace("PreShowerPhiClusterAlgo") << "Strip energy" << candidate_it->second.energy() <<"is below threshold";
105 
106  // crystal should not be included...
107  if ( (used_s->find(candidate_it->first) != used_s->end()) || // ...if it already belongs to a cluster
108  (candidate_it == rechits_map->end() ) || // ...if it corresponds to a hit
109  (candidate_it->second.energy() <= esStripEnergyCut_ ) ) // ...if it has a negative or zero energy
110  {
111  return false;
112  }
113  return true;
114 }
115 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
reco::PreshowerCluster makeOneCluster(ESDetId strip, HitsID *used_strips, RecHitsMap *rechits_map, const CaloSubdetectorGeometry *geometry_p, double deltaEta, double minDeltaPhi, double maxDeltaPhi)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
void push_back(T const &t)
static const double deltaEta
Definition: CaloConstants.h:8
bool goodStrip(RecHitsMap::iterator candidate_it)
T z() const
Definition: PV3DBase.h:64
#define LogTrace(id)
std::vector< EcalRecHit >::iterator iterator
const_iterator end() const
Definition: DetId.h:18
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
T eta() const
Definition: PV3DBase.h:76
std::map< DetId, EcalRecHit > RecHitsMap
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
static int position[264][3]
Definition: ReadPGInfo.cc:509
int plane() const
Definition: ESDetId.h:41
T x() const
Definition: PV3DBase.h:62
const_iterator begin() const