CMS 3D CMS Logo

EgammaPCAHelper.cc
Go to the documentation of this file.
2 
7 
8 // To retrieve HGCalImagingAlgo::maxlayer
9 // Also available as HGCal3DClustering::lastLayerBH and HGCalDepthPreClusterer::lastLayerBH
11 
12 #include <algorithm>
13 #include <iostream>
14 
15 using namespace hgcal;
16 
18  // Thickness correction to dEdx weights
19  // (100um, 200um, 300um silicon)
20  // See RecoLocalCalo.HGCalRecProducers.HGCalRecHit_cfi
21  invThicknessCorrection_({1. / 1.132, 1. / 1.092, 1. / 1.084}),
22  pca_(new TPrincipal(3, "D")) {
23  hitMapOrigin_ = 0;
24  hitMap_ = new std::map<DetId, const HGCRecHit *>();
25  debug_ = false;
26 }
27 
29  if (hitMapOrigin_ == 2) delete hitMap_;
30 }
31 
32 void EGammaPCAHelper::setHitMap( std::map<DetId,const HGCRecHit *> * hitMap) {
33  hitMapOrigin_ = 1;
34  hitMap_ = hitMap ;
35  pcaIteration_ = 0;
36 }
37 
39  recHitTools_ = recHitTools;
40 }
41 
43  const HGCRecHitCollection & rechitsFH,
44  const HGCRecHitCollection & rechitsBH) {
45  hitMap_->clear();
46  for (const auto& hit : rechitsEE) {
47  hitMap_->emplace_hint(hitMap_->end(), hit.detid(), &hit);
48  }
49 
50  for (const auto& hit : rechitsFH) {
51  hitMap_->emplace_hint(hitMap_->end(), hit.detid(), &hit);
52  }
53 
54  for (const auto& hit : rechitsBH) {
55  hitMap_->emplace_hint(hitMap_->end(), hit.detid(), &hit);
56  }
57 
58  pcaIteration_ = 0;
59  hitMapOrigin_ = 2;
60 }
61 
63  theCluster_ = &cluster;
64  std::vector<std::pair<DetId, float>> result;
65  for (reco::HGCalMultiCluster::component_iterator it = cluster.begin(); it != cluster.end();
66  it++) {
67  const std::vector<std::pair<DetId, float>> &hf = (*it)->hitsAndFractions();
68  result.insert(result.end(),hf.begin(),hf.end());
69  }
70  storeRecHits(result);
71 }
72 
74  theCluster_ = &cluster;
75  storeRecHits(cluster.hitsAndFractions());
76 }
77 
78 void EGammaPCAHelper::storeRecHits(const std::vector<std::pair<DetId, float>> &hf) {
79  std::vector<double> pcavars;
80  pcavars.resize(3,0.);
81  theSpots_.clear();
82  pcaIteration_ = 0;
83 
84  sigu_ = 0.;
85  sigv_ = 0.;
86  sigp_ = 0.;
87  sige_ = 0.;
88 
89  unsigned hfsize = hf.size();
90  if (debug_)
91  std::cout << "The seed cluster constains " << hfsize << " hits " << std::endl;
92 
93  if (hfsize == 0) return;
94 
95 
96  for (unsigned int j = 0; j < hfsize; j++) {
97  unsigned int layer = recHitTools_->getLayerWithOffset(hf[j].first);
98 
99  const DetId rh_detid = hf[j].first;
100  std::map<DetId,const HGCRecHit *>::const_iterator itcheck= hitMap_->find(rh_detid);
101  if (itcheck == hitMap_->end()) {
102  edm::LogWarning("EgammaPCAHelper") << " Big problem, unable to find a hit " << rh_detid.rawId() << " " << rh_detid.det() << " " << HGCalDetId(rh_detid) << std::endl;
103  continue;
104  }
105  if (debug_) {
106  std::cout << "DetId " << rh_detid.rawId() << " " << layer << " " << itcheck->second->energy() <<std::endl;
107  std::cout << " Hit " << itcheck->second << " " << itcheck->second->energy() << std::endl;
108  }
109  float fraction = hf[j].second;
110 
111  int thickIndex = recHitTools_->getSiThickIndex(rh_detid);
112  double mip = dEdXWeights_[layer] * 0.001; // convert in GeV
113  if(thickIndex>-1 and thickIndex<3) mip *= invThicknessCorrection_[thickIndex];
114 
115  pcavars[0] = recHitTools_->getPosition(rh_detid).x();
116  pcavars[1] = recHitTools_->getPosition(rh_detid).y();
117  pcavars[2] = recHitTools_->getPosition(rh_detid).z();
118  if (pcavars[2] == 0.)
119  edm::LogWarning("EgammaPCAHelper") << " Problem, hit with z =0 ";
120  else {
121  Spot mySpot(rh_detid,itcheck->second->energy(),pcavars,layer,fraction,mip);
122  theSpots_.push_back(mySpot);
123  }
124  }
125  if (debug_) {
126  std::cout << " Stored " << theSpots_.size() << " hits " << std::endl;
127  }
128 }
129 
130 void EGammaPCAHelper::computePCA(float radius , bool withHalo) {
131  // very important - to reset
132  pca_.reset(new TPrincipal(3, "D"));
133  bool initialCalculation = radius < 0;
134  if (debug_)
135  std::cout << " Initial calculation " << initialCalculation << std::endl;
136  if (initialCalculation && withHalo) {
137  edm::LogWarning("EGammaPCAHelper") << "Warning - in the first iteration, the halo hits are excluded " << std::endl;
138  withHalo=false;
139  }
140 
141  float radius2 = radius*radius;
142  if (! initialCalculation) {
143  math::XYZVector mainAxis(axis_);
144  mainAxis.unit();
145  math::XYZVector phiAxis(barycenter_.x(), barycenter_.y(), 0);
146  math::XYZVector udir(mainAxis.Cross(phiAxis));
147  udir = udir.unit();
149  Point(0., 0., 1.), Point(1., 0., 0.));
150  }
151 
152  std::set<int> layers;
153  for (const auto& spot : theSpots_) {
154  if (spot.layer() > recHitTools_->lastLayerEE()) continue;
155  if (!withHalo && (! spot.isCore() ))
156  continue;
157  if (initialCalculation) {
158  // initial calculation, take only core hits
159  if ( ! spot.isCore() ) continue;
160  layers.insert(spot.layer());
161  for (int i = 0; i < spot.multiplicity(); ++i)
162  pca_->AddRow(spot.row());
163  }
164  else {
165  // use a cylinder, include all hits
166  math::XYZPoint local = trans_(Point( spot.row()[0],spot.row()[1],spot.row()[2]));
167  if (local.Perp2() > radius2) continue;
168  layers.insert(spot.layer());
169  for (int i = 0; i < spot.multiplicity(); ++i)
170  pca_->AddRow(spot.row());
171  }
172  }
173  if (debug_)
174  std::cout << " Nlayers " << layers.size() << std::endl;
175  if (layers.size() < 3) {
176  pcaIteration_ = -1;
177  return;
178  }
179  pca_->MakePrincipals();
180  ++pcaIteration_;
181  const TVectorD& means = *(pca_->GetMeanValues());
182  const TMatrixD& eigens = *(pca_->GetEigenVectors());
183 
184  barycenter_ = math::XYZPoint(means[0], means[1], means[2]);
185  axis_ = math::XYZVector(eigens(0, 0), eigens(1, 0), eigens(2, 0));
186  if (axis_.z() * barycenter_.z() < 0.0) {
187  axis_ = -1. * axis_;
188  }
189 }
190 
191  void EGammaPCAHelper::computeShowerWidth(float radius, bool withHalo){
192  sigu_ = 0.;
193  sigv_ = 0.;
194  sigp_ = 0.;
195  sige_ = 0.;
196  double cyl_ene = 0.;
197 
198  float radius2 = radius * radius;
199  for (const auto& spot : theSpots_) {
200  Point globalPoint(spot.row()[0],spot.row()[1],spot.row()[2]);
201  math::XYZPoint local = trans_(globalPoint);
202  if (local.Perp2() > radius2) continue;
203 
204  // Select halo hits or not
205  if (withHalo && spot.fraction() < 0) continue;
206  if (!withHalo && !(spot.isCore())) continue;
207 
208  sige_ += (globalPoint.eta() - theCluster_->eta()) * (globalPoint.eta() - theCluster_->eta()) * spot.energy();
209  sigp_ += deltaPhi(globalPoint.phi(), theCluster_->phi()) * deltaPhi(globalPoint.phi(), theCluster_->phi()) *
210  spot.energy();
211 
212  sigu_ += local.x() * local.x() * spot.energy();
213  sigv_ += local.y() * local.y() * spot.energy();
214  cyl_ene += spot.energy();
215  }
216 
217  if (cyl_ene > 0.) {
218  const double inv_cyl_ene = 1. / cyl_ene;
219  sigu_ = sigu_ * inv_cyl_ene;
220  sigv_ = sigv_ * inv_cyl_ene;
221  sigp_ = sigp_ * inv_cyl_ene;
222  sige_ = sige_ * inv_cyl_ene;
223  }
224  sigu_ = std::sqrt(sigu_);
225  sigv_ = std::sqrt(sigv_);
226  sigp_ = std::sqrt(sigp_);
227  sige_ = std::sqrt(sige_);
228 }
229 
231  if (pcaIteration_ == 0) {
232  if(debug_)
233  std::cout << " The PCA has not been run yet " << std::endl;
234  return false;
235  } else if (pcaIteration_ == 1) {
236  if (debug_)
237  std::cout << " The PCA has been run only once - careful " << std::endl;
238  return false;
239  } else if (pcaIteration_ == -1){
240  if (debug_)
241  std::cout << " Not enough layers to perform PCA " << std::endl;
242  return false;
243  }
244  return true;
245 }
246 
248  theSpots_.clear();
249  pcaIteration_ = 0;
250  sigu_ = 0.;
251  sigv_ = 0.;
252  sigp_ = 0.;
253  sige_ = 0.;
254 }
255 
257  if (debug_) checkIteration();
258  std::set<int> layers;
259  float radius2 = radius*radius;
260  std::vector<float> energyPerLayer(HGCalImagingAlgo::maxlayer+1, 0.f);
261  math::XYZVector mainAxis(axis_);
262  mainAxis.unit();
263  math::XYZVector phiAxis(barycenter_.x(), barycenter_.y(), 0);
264  math::XYZVector udir(mainAxis.Cross(phiAxis));
265  udir = udir.unit();
267  Point(0., 0., 1.), Point(1., 0., 0.));
268  float energyEE = 0.;
269  float energyFH = 0.;
270  float energyBH = 0.;
271 
272  for (const auto& spot : theSpots_) {
273  if (!withHalo && ! spot.isCore())
274  continue;
275  math::XYZPoint local = trans_(Point( spot.row()[0],spot.row()[1],spot.row()[2]));
276  if (local.Perp2() > radius2) continue;
277  energyPerLayer[spot.layer()] += spot.energy();
278  layers.insert(spot.layer());
279  if (spot.detId().det() == DetId::HGCalEE or spot.subdet() == HGCEE) { energyEE += spot.energy();}
280  else if (spot.detId().det() == DetId::HGCalHSi or spot.subdet() == HGCHEF) { energyFH += spot.energy();}
281  else if (spot.detId().det() == DetId::HGCalHSc or spot.subdet() == HGCHEB) { energyBH += spot.energy();}
282 
283  }
284  return LongDeps(radius,energyPerLayer,energyEE,energyFH,energyBH,layers);
285 }
286 
288  unsigned nSpots = theSpots_.size();
289  float radius2=radius*radius;
290  for ( unsigned i =0; i< nSpots ; ++i) {
291  Spot spot(theSpots_[i]);
292  math::XYZPoint local = trans_(Point( spot.row()[0],spot.row()[1],spot.row()[2]));
293  if (local.Perp2() < radius2 ) {
294  std::cout << i << " " << theSpots_[i].detId().rawId() << " " << theSpots_[i].layer() << " " << theSpots_[i].energy() << " " <<theSpots_[i].isCore() ;
295  std::cout << " " << std::sqrt(local.Perp2()) << std::endl;
296  }
297  }
298 }
299 
301  unsigned int firstLayer = 0;
302  for(unsigned il=1;il<=HGCalImagingAlgo::maxlayer;++il) {
303  if (ld.energyPerLayer()[il] > 0.) {
304  firstLayer = il;
305  break;
306  }
307  }
308  // Make dummy DetId to get abs(z) for layer
309  return recHitTools_->getPositionLayer(firstLayer).z();
310 }
311 
312 float EGammaPCAHelper::clusterDepthCompatibility(const LongDeps & ld, float & measuredDepth, float& expectedDepth, float&expectedSigma) {
313  expectedDepth = -999.;
314  expectedSigma = -999.;
315  measuredDepth = -999.;
316  if (!checkIteration()) return -999.;
317 
318  float z = findZFirstLayer(ld);
319  math::XYZVector dir=axis_.unit();
320  measuredDepth = std::abs((z-std::abs(barycenter_.z()))/dir.z());
321  return showerDepth_.getClusterDepthCompatibility(measuredDepth,ld.energyEE(), expectedDepth,expectedSigma);
322 }
void computeShowerWidth(float radius, bool withHalo=true)
GlobalPoint getPositionLayer(int layer) const
Definition: RecHitTools.cc:145
math::XYZPoint barycenter_
float energyEE() const
Definition: LongDeps.h:20
static const unsigned int maxlayer
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
void computePCA(float radius, bool withHalo=true)
std::map< DetId, const HGCRecHit * > * hitMap_
LongDeps energyPerLayer(float radius, bool withHalo=true)
float findZFirstLayer(const LongDeps &) const
std::unique_ptr< TPrincipal > pca_
component_iterator begin() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
ROOT::Math::Transform3DPJ Transform3D
T y() const
Definition: PV3DBase.h:63
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:197
float getClusterDepthCompatibility(float measuredDepth, float emEnergy, float &expectedDepth, float &expectedSigma) const
Definition: ShowerDepth.cc:6
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:168
void storeRecHits(const reco::CaloCluster &theCluster)
ROOT::Math::Transform3DPJ::Point Point
T sqrt(T t)
Definition: SSEVec.h:18
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
T z() const
Definition: PV3DBase.h:64
void setHitMap(std::map< DetId, const HGCRecHit * > *hitMap)
to set from outside - once per event
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
std::vector< Spot > theSpots_
unsigned int lastLayerEE() const
Definition: RecHitTools.h:58
component_iterator end() const
Definition: DetId.h:18
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:304
float clusterDepthCompatibility(const LongDeps &, float &measuredDepth, float &expectedDepth, float &expectedSigma)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:189
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const std::vector< float > & energyPerLayer() const
Definition: LongDeps.h:23
void printHits(float radius) const
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:133
const double * row() const
Definition: Spot.h:17
std::vector< double > invThicknessCorrection_
void fillHitMap(const HGCRecHitCollection &HGCEERecHits, const HGCRecHitCollection &HGCFHRecHits, const HGCRecHitCollection &HGCBHRecHits)
to compute from inside - once per event
std::vector< double > dEdXWeights_
const reco::CaloCluster * theCluster_
bool checkIteration() const
dbl *** dir
Definition: mlp_gen.cc:35
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:171
T x() const
Definition: PV3DBase.h:62
const hgcal::RecHitTools * recHitTools_
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39