CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
hgcal::EGammaPCAHelper Class Reference

#include <EgammaPCAHelper.h>

Public Types

typedef ROOT::Math::Transform3D::Point Point
 
typedef ROOT::Math::Transform3D Transform3D
 

Public Member Functions

const math::XYZVectoraxis () const
 
const math::XYZPointbarycenter () const
 
void clear ()
 
float clusterDepthCompatibility (const LongDeps &, float &measuredDepth, float &expectedDepth, float &expectedSigma)
 
void computePCA (float radius, bool withHalo=true)
 
void computeShowerWidth (float radius, bool withHalo=true)
 
 EGammaPCAHelper ()
 
const TVectorD & eigenValues () const
 
LongDeps energyPerLayer (float radius, bool withHalo=true)
 
void fillHitMap (const HGCRecHitCollection &HGCEERecHits, const HGCRecHitCollection &HGCFHRecHits, const HGCRecHitCollection &HGCBHRecHits)
 to compute from inside - once per event More...
 
std::map< DetId, const HGCRecHit * > * getHitMap ()
 
void pcaInitialComputation ()
 
const TPrincipal & pcaResult ()
 
void printHits (float radius) const
 
void setdEdXWeights (const std::vector< double > &dEdX)
 
void setHitMap (std::map< DetId, const HGCRecHit * > *hitMap)
 to set from outside - once per event More...
 
void setRecHitTools (const hgcal::RecHitTools *recHitTools)
 
double sigmaEE () const
 
double sigmaPP () const
 
const TVectorD & sigmas () const
 
double sigmaUU () const
 
double sigmaVV () const
 
void storeRecHits (const reco::CaloCluster &theCluster)
 
void storeRecHits (const reco::HGCalMultiCluster &cluster)
 
 ~EGammaPCAHelper ()
 

Private Member Functions

bool checkIteration () const
 
float findZFirstLayer (const LongDeps &) const
 
void storeRecHits (const std::vector< std::pair< DetId, float >> &hf)
 

Private Attributes

math::XYZVector axis_
 
math::XYZPoint barycenter_
 
bool debug_
 
std::vector< double > dEdXWeights_
 
std::map< DetId, const HGCRecHit * > * hitMap_
 
int hitMapOrigin_
 
std::vector< double > invThicknessCorrection_
 
std::unique_ptr< TPrincipal > pca_
 
int pcaIteration_
 
bool recHitsStored_
 
const hgcal::RecHitToolsrecHitTools_
 
ShowerDepth showerDepth_
 
double sige_
 
double sigp_
 
double sigu_
 
double sigv_
 
const reco::CaloClustertheCluster_
 
std::vector< SpottheSpots_
 
Transform3D trans_
 

Detailed Description

Definition at line 33 of file EgammaPCAHelper.h.

Member Typedef Documentation

typedef ROOT::Math::Transform3D::Point hgcal::EGammaPCAHelper::Point

Definition at line 37 of file EgammaPCAHelper.h.

typedef ROOT::Math::Transform3D hgcal::EGammaPCAHelper::Transform3D

Definition at line 36 of file EgammaPCAHelper.h.

Constructor & Destructor Documentation

EGammaPCAHelper::EGammaPCAHelper ( )

Definition at line 20 of file EgammaPCAHelper.cc.

20  :
21  // Thickness correction to dEdx weights
22  // (100um, 200um, 300um silicon)
23  // See RecoLocalCalo.HGCalRecProducers.HGCalRecHit_cfi
24  invThicknessCorrection_({1. / 1.132, 1. / 1.092, 1. / 1.084}),
25  pca_(new TPrincipal(3, "D")) {
26  hitMapOrigin_ = 0;
27  hitMap_ = new std::map<DetId, const HGCRecHit *>();
28  debug_ = false;
29 }
std::map< DetId, const HGCRecHit * > * hitMap_
std::unique_ptr< TPrincipal > pca_
std::vector< double > invThicknessCorrection_
EGammaPCAHelper::~EGammaPCAHelper ( )

Definition at line 31 of file EgammaPCAHelper.cc.

References hitMap_, and hitMapOrigin_.

31  {
32  if (hitMapOrigin_ == 2) delete hitMap_;
33 }
std::map< DetId, const HGCRecHit * > * hitMap_

Member Function Documentation

const math::XYZVector& hgcal::EGammaPCAHelper::axis ( ) const
inline

Definition at line 66 of file EgammaPCAHelper.h.

References axis_, and computeShowerWidth().

Referenced by HGCalEgammaIDHelper::axis().

66 {return axis_;}
const math::XYZPoint& hgcal::EGammaPCAHelper::barycenter ( ) const
inline

Definition at line 65 of file EgammaPCAHelper.h.

References barycenter_.

Referenced by HGCalEgammaIDHelper::barycenter().

65 {return barycenter_;}
math::XYZPoint barycenter_
bool EGammaPCAHelper::checkIteration ( ) const
private

Definition at line 233 of file EgammaPCAHelper.cc.

References gather_cfg::cout, debug_, and pcaIteration_.

Referenced by clusterDepthCompatibility(), energyPerLayer(), sigmaEE(), sigmaPP(), sigmas(), sigmaUU(), and sigmaVV().

233  {
234  if (pcaIteration_ == 0) {
235  if(debug_)
236  std::cout << " The PCA has not been run yet " << std::endl;
237  return false;
238  } else if (pcaIteration_ == 1) {
239  if (debug_)
240  std::cout << " The PCA has been run only once - careful " << std::endl;
241  return false;
242  } else if (pcaIteration_ == -1){
243  if (debug_)
244  std::cout << " Not enough layers to perform PCA " << std::endl;
245  return false;
246  }
247  return true;
248 }
void EGammaPCAHelper::clear ( void  )

Definition at line 250 of file EgammaPCAHelper.cc.

References pcaIteration_, sige_, sigp_, sigu_, sigv_, and theSpots_.

Referenced by HGCalEgammaIDHelper::computeHGCAL(), and sigmas().

250  {
251  theSpots_.clear();
252  pcaIteration_ = 0;
253  sigu_ = 0.;
254  sigv_ = 0.;
255  sigp_ = 0.;
256  sige_ = 0.;
257 }
std::vector< Spot > theSpots_
float EGammaPCAHelper::clusterDepthCompatibility ( const LongDeps ld,
float &  measuredDepth,
float &  expectedDepth,
float &  expectedSigma 
)

Definition at line 315 of file EgammaPCAHelper.cc.

References funct::abs(), axis_, barycenter_, checkIteration(), dir, hgcal::LongDeps::energyEE(), findZFirstLayer(), hgcal::ShowerDepth::getClusterDepthCompatibility(), and showerDepth_.

Referenced by HGCalEgammaIDHelper::clusterDepthCompatibility(), and sigmas().

315  {
316  expectedDepth = -999.;
317  expectedSigma = -999.;
318  measuredDepth = -999.;
319  if (!checkIteration()) return -999.;
320 
321  float z = findZFirstLayer(ld);
322  math::XYZVector dir=axis_.unit();
323  measuredDepth = std::abs((z-std::abs(barycenter_.z()))/dir.z());
324  return showerDepth_.getClusterDepthCompatibility(measuredDepth,ld.energyEE(), expectedDepth,expectedSigma);
325 }
math::XYZPoint barycenter_
float energyEE() const
Definition: LongDeps.h:20
float findZFirstLayer(const LongDeps &) const
float getClusterDepthCompatibility(float measuredDepth, float emEnergy, float &expectedDepth, float &expectedSigma) const
Definition: ShowerDepth.cc:6
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
bool checkIteration() const
dbl *** dir
Definition: mlp_gen.cc:35
void EGammaPCAHelper::computePCA ( float  radius,
bool  withHalo = true 
)

Definition at line 133 of file EgammaPCAHelper.cc.

References axis_, barycenter_, gather_cfg::cout, debug_, mps_fire::i, hgcal::RecHitTools::lastLayerEE(), LayerTriplets::layers(), DTRecHitClients_cfi::local, pca_, pcaIteration_, TCMET_cfi::radius, recHitTools_, theSpots_, and trans_.

Referenced by HGCalEgammaIDHelper::computeHGCAL(), and pcaInitialComputation().

133  {
134  // very important - to reset
135  pca_.reset(new TPrincipal(3, "D"));
136  bool initialCalculation = radius < 0;
137  if (debug_)
138  std::cout << " Initial calculation " << initialCalculation << std::endl;
139  if (initialCalculation && withHalo) {
140  edm::LogWarning("EGammaPCAHelper") << "Warning - in the first iteration, the halo hits are excluded " << std::endl;
141  withHalo=false;
142  }
143 
144  float radius2 = radius*radius;
145  if (! initialCalculation) {
146  math::XYZVector mainAxis(axis_);
147  mainAxis.unit();
148  math::XYZVector phiAxis(barycenter_.x(), barycenter_.y(), 0);
149  math::XYZVector udir(mainAxis.Cross(phiAxis));
150  udir = udir.unit();
152  Point(0., 0., 1.), Point(1., 0., 0.));
153  }
154 
155  std::set<int> layers;
156  for (const auto& spot : theSpots_) {
157  if (spot.layer() > recHitTools_->lastLayerEE()) continue;
158  if (!withHalo && (! spot.isCore() ))
159  continue;
160  if (initialCalculation) {
161  // initial calculation, take only core hits
162  if ( ! spot.isCore() ) continue;
163  layers.insert(spot.layer());
164  for (int i = 0; i < spot.multiplicity(); ++i)
165  pca_->AddRow(spot.row());
166  }
167  else {
168  // use a cylinder, include all hits
169  math::XYZPoint local = trans_(Point( spot.row()[0],spot.row()[1],spot.row()[2]));
170  if (local.Perp2() > radius2) continue;
171  layers.insert(spot.layer());
172  for (int i = 0; i < spot.multiplicity(); ++i)
173  pca_->AddRow(spot.row());
174  }
175  }
176  if (debug_)
177  std::cout << " Nlayers " << layers.size() << std::endl;
178  if (layers.size() < 3) {
179  pcaIteration_ = -1;
180  return;
181  }
182  pca_->MakePrincipals();
183  ++pcaIteration_;
184  const TVectorD& means = *(pca_->GetMeanValues());
185  const TMatrixD& eigens = *(pca_->GetEigenVectors());
186 
187  barycenter_ = math::XYZPoint(means[0], means[1], means[2]);
188  axis_ = math::XYZVector(eigens(0, 0), eigens(1, 0), eigens(2, 0));
189  if (axis_.z() * barycenter_.z() < 0.0) {
190  axis_ = -1. * axis_;
191  }
192 }
math::XYZPoint barycenter_
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
std::unique_ptr< TPrincipal > pca_
std::vector< Spot > theSpots_
unsigned int lastLayerEE() const
Definition: RecHitTools.h:58
ROOT::Math::Transform3D Transform3D
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
ROOT::Math::Transform3D::Point Point
const hgcal::RecHitTools * recHitTools_
void EGammaPCAHelper::computeShowerWidth ( float  radius,
bool  withHalo = true 
)

Definition at line 194 of file EgammaPCAHelper.cc.

References hiPixelPairStep_cff::deltaPhi, reco::CaloCluster::eta(), DTRecHitClients_cfi::local, reco::CaloCluster::phi(), TCMET_cfi::radius, sige_, sigp_, sigu_, sigv_, mathSSE::sqrt(), theCluster_, theSpots_, and trans_.

Referenced by axis(), and HGCalEgammaIDHelper::computeHGCAL().

194  {
195  sigu_ = 0.;
196  sigv_ = 0.;
197  sigp_ = 0.;
198  sige_ = 0.;
199  double cyl_ene = 0.;
200 
201  float radius2 = radius * radius;
202  for (const auto& spot : theSpots_) {
203  Point globalPoint(spot.row()[0],spot.row()[1],spot.row()[2]);
204  math::XYZPoint local = trans_(globalPoint);
205  if (local.Perp2() > radius2) continue;
206 
207  // Select halo hits or not
208  if (withHalo && spot.fraction() < 0) continue;
209  if (!withHalo && !(spot.isCore())) continue;
210 
211  sige_ += (globalPoint.eta() - theCluster_->eta()) * (globalPoint.eta() - theCluster_->eta()) * spot.energy();
212  sigp_ += deltaPhi(globalPoint.phi(), theCluster_->phi()) * deltaPhi(globalPoint.phi(), theCluster_->phi()) *
213  spot.energy();
214 
215  sigu_ += local.x() * local.x() * spot.energy();
216  sigv_ += local.y() * local.y() * spot.energy();
217  cyl_ene += spot.energy();
218  }
219 
220  if (cyl_ene > 0.) {
221  const double inv_cyl_ene = 1. / cyl_ene;
222  sigu_ = sigu_ * inv_cyl_ene;
223  sigv_ = sigv_ * inv_cyl_ene;
224  sigp_ = sigp_ * inv_cyl_ene;
225  sige_ = sige_ * inv_cyl_ene;
226  }
227  sigu_ = std::sqrt(sigu_);
228  sigv_ = std::sqrt(sigv_);
229  sigp_ = std::sqrt(sigp_);
230  sige_ = std::sqrt(sige_);
231 }
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:168
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< Spot > theSpots_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
const reco::CaloCluster * theCluster_
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:171
const TVectorD& hgcal::EGammaPCAHelper::eigenValues ( ) const
inline

Definition at line 75 of file EgammaPCAHelper.h.

References pca_.

Referenced by HGCalEgammaIDHelper::eigenValues().

75 {return *pca_->GetEigenValues();}
std::unique_ptr< TPrincipal > pca_
LongDeps EGammaPCAHelper::energyPerLayer ( float  radius,
bool  withHalo = true 
)

Definition at line 259 of file EgammaPCAHelper.cc.

References axis_, barycenter_, checkIteration(), debug_, f, DetId::HGCalEE, DetId::HGCalHSc, DetId::HGCalHSi, HGCEE, HGCHEB, HGCHEF, LayerTriplets::layers(), DTRecHitClients_cfi::local, HGCalClusteringAlgoBase::maxlayer, or, TCMET_cfi::radius, theSpots_, and trans_.

Referenced by HGCalEgammaIDHelper::energyPerLayer(), and sigmas().

259  {
260  if (debug_) checkIteration();
261  std::set<int> layers;
262  float radius2 = radius*radius;
263  std::vector<float> energyPerLayer(HGCalClusteringAlgoBase::maxlayer+1, 0.f);
264  math::XYZVector mainAxis(axis_);
265  mainAxis.unit();
266  math::XYZVector phiAxis(barycenter_.x(), barycenter_.y(), 0);
267  math::XYZVector udir(mainAxis.Cross(phiAxis));
268  udir = udir.unit();
270  Point(0., 0., 1.), Point(1., 0., 0.));
271  float energyEE = 0.;
272  float energyFH = 0.;
273  float energyBH = 0.;
274 
275  for (const auto& spot : theSpots_) {
276  if (!withHalo && ! spot.isCore())
277  continue;
278  math::XYZPoint local = trans_(Point( spot.row()[0],spot.row()[1],spot.row()[2]));
279  if (local.Perp2() > radius2) continue;
280  energyPerLayer[spot.layer()] += spot.energy();
281  layers.insert(spot.layer());
282  if (spot.detId().det() == DetId::HGCalEE or spot.subdet() == HGCEE) { energyEE += spot.energy();}
283  else if (spot.detId().det() == DetId::HGCalHSi or spot.subdet() == HGCHEF) { energyFH += spot.energy();}
284  else if (spot.detId().det() == DetId::HGCalHSc or spot.subdet() == HGCHEB) { energyBH += spot.energy();}
285 
286  }
287  return LongDeps(radius,energyPerLayer,energyEE,energyFH,energyBH,layers);
288 }
math::XYZPoint barycenter_
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
LongDeps energyPerLayer(float radius, bool withHalo=true)
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
double f[11][100]
std::vector< Spot > theSpots_
ROOT::Math::Transform3D Transform3D
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static const unsigned int maxlayer
bool checkIteration() const
ROOT::Math::Transform3D::Point Point
void EGammaPCAHelper::fillHitMap ( const HGCRecHitCollection HGCEERecHits,
const HGCRecHitCollection HGCFHRecHits,
const HGCRecHitCollection HGCBHRecHits 
)

to compute from inside - once per event

Definition at line 45 of file EgammaPCAHelper.cc.

References hitMap_, hitMapOrigin_, and pcaIteration_.

Referenced by HGCalEgammaIDHelper::eventInit().

47  {
48  hitMap_->clear();
49  for (const auto& hit : rechitsEE) {
50  hitMap_->emplace_hint(hitMap_->end(), hit.detid(), &hit);
51  }
52 
53  for (const auto& hit : rechitsFH) {
54  hitMap_->emplace_hint(hitMap_->end(), hit.detid(), &hit);
55  }
56 
57  for (const auto& hit : rechitsBH) {
58  hitMap_->emplace_hint(hitMap_->end(), hit.detid(), &hit);
59  }
60 
61  pcaIteration_ = 0;
62  hitMapOrigin_ = 2;
63 }
std::map< DetId, const HGCRecHit * > * hitMap_
float EGammaPCAHelper::findZFirstLayer ( const LongDeps ld) const
private

Definition at line 303 of file EgammaPCAHelper.cc.

References hgcal::LongDeps::energyPerLayer(), hgcal::RecHitTools::getPositionLayer(), HGCalClusteringAlgoBase::maxlayer, recHitTools_, and PV3DBase< T, PVType, FrameType >::z().

Referenced by clusterDepthCompatibility(), and sigmas().

303  {
304  unsigned int firstLayer = 0;
305  for(unsigned il=1;il<=HGCalClusteringAlgoBase::maxlayer;++il) {
306  if (ld.energyPerLayer()[il] > 0.) {
307  firstLayer = il;
308  break;
309  }
310  }
311  // Make dummy DetId to get abs(z) for layer
312  return recHitTools_->getPositionLayer(firstLayer).z();
313 }
GlobalPoint getPositionLayer(int layer) const
Definition: RecHitTools.cc:124
T z() const
Definition: PV3DBase.h:64
const std::vector< float > & energyPerLayer() const
Definition: LongDeps.h:23
static const unsigned int maxlayer
const hgcal::RecHitTools * recHitTools_
std::map<DetId,const HGCRecHit *>* hgcal::EGammaPCAHelper::getHitMap ( )
inline

Definition at line 54 of file EgammaPCAHelper.h.

References hitMap_, and setRecHitTools().

54 {return hitMap_;}
std::map< DetId, const HGCRecHit * > * hitMap_
void hgcal::EGammaPCAHelper::pcaInitialComputation ( )
inline

Definition at line 60 of file EgammaPCAHelper.h.

References computePCA(), and TCMET_cfi::radius.

Referenced by HGCalEgammaIDHelper::computeHGCAL().

60  {
61  computePCA(-1.,false);
62  }
void computePCA(float radius, bool withHalo=true)
const TPrincipal& hgcal::EGammaPCAHelper::pcaResult ( )
void EGammaPCAHelper::printHits ( float  radius) const

Definition at line 290 of file EgammaPCAHelper.cc.

References gather_cfg::cout, mps_fire::i, DTRecHitClients_cfi::local, TCMET_cfi::radius, hgcal::Spot::row(), mathSSE::sqrt(), theSpots_, and trans_.

Referenced by ntuplePrintersDiff.TrackingParticlePrinter::__call__(), ntuplePrintersDiff.SeedPrinter::diff(), ntuplePrintersDiff.TrackPrinter::diff(), ntuplePrintersDiff.TrackingParticlePrinter::diff(), HGCalEgammaIDHelper::printHits(), ntuplePrintersDiff.SeedPrinter::printSeed(), ntuplePrintersDiff.TrackPrinter::printTrack(), and sigmas().

290  {
291  unsigned nSpots = theSpots_.size();
292  float radius2=radius*radius;
293  for ( unsigned i =0; i< nSpots ; ++i) {
294  Spot spot(theSpots_[i]);
295  math::XYZPoint local = trans_(Point( spot.row()[0],spot.row()[1],spot.row()[2]));
296  if (local.Perp2() < radius2 ) {
297  std::cout << i << " " << theSpots_[i].detId().rawId() << " " << theSpots_[i].layer() << " " << theSpots_[i].energy() << " " <<theSpots_[i].isCore() ;
298  std::cout << " " << std::sqrt(local.Perp2()) << std::endl;
299  }
300  }
301 }
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< Spot > theSpots_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
ROOT::Math::Transform3D::Point Point
void hgcal::EGammaPCAHelper::setdEdXWeights ( const std::vector< double > &  dEdX)
inline

Definition at line 58 of file EgammaPCAHelper.h.

References HGCalRecHit_cfi::dEdX, and dEdXWeights_.

Referenced by HGCalEgammaIDHelper::HGCalEgammaIDHelper().

58 { dEdXWeights_ = dEdX;}
std::vector< double > dEdXWeights_
void EGammaPCAHelper::setHitMap ( std::map< DetId, const HGCRecHit * > *  hitMap)

to set from outside - once per event

Definition at line 35 of file EgammaPCAHelper.cc.

References hitMap_, hitMapOrigin_, and pcaIteration_.

35  {
36  hitMapOrigin_ = 1;
37  hitMap_ = hitMap ;
38  pcaIteration_ = 0;
39 }
std::map< DetId, const HGCRecHit * > * hitMap_
void EGammaPCAHelper::setRecHitTools ( const hgcal::RecHitTools recHitTools)

Definition at line 41 of file EgammaPCAHelper.cc.

References recHitTools_.

Referenced by HGCalEgammaIDHelper::eventInit(), and getHitMap().

41  {
42  recHitTools_ = recHitTools;
43 }
const hgcal::RecHitTools * recHitTools_
double hgcal::EGammaPCAHelper::sigmaEE ( ) const
inline

Definition at line 72 of file EgammaPCAHelper.h.

References checkIteration(), and sige_.

Referenced by HGCalEgammaIDHelper::sigmaEE().

72 { return checkIteration()? sige_ : -1. ;}
bool checkIteration() const
double hgcal::EGammaPCAHelper::sigmaPP ( ) const
inline

Definition at line 73 of file EgammaPCAHelper.h.

References checkIteration(), and sigp_.

Referenced by HGCalEgammaIDHelper::sigmaPP().

73 { return checkIteration()? sigp_ : -1. ;}
bool checkIteration() const
const TVectorD& hgcal::EGammaPCAHelper::sigmas ( ) const
inline
double hgcal::EGammaPCAHelper::sigmaUU ( ) const
inline

Definition at line 70 of file EgammaPCAHelper.h.

References checkIteration(), and sigu_.

Referenced by HGCalEgammaIDHelper::sigmaUU().

70 { return checkIteration()? sigu_ : -1. ;}
bool checkIteration() const
double hgcal::EGammaPCAHelper::sigmaVV ( ) const
inline

Definition at line 71 of file EgammaPCAHelper.h.

References checkIteration(), and sigv_.

Referenced by HGCalEgammaIDHelper::sigmaVV().

71 { return checkIteration()? sigv_ : -1. ;}
bool checkIteration() const
void EGammaPCAHelper::storeRecHits ( const reco::CaloCluster theCluster)

Definition at line 76 of file EgammaPCAHelper.cc.

References reco::CaloCluster::hitsAndFractions(), and theCluster_.

Referenced by HGCalEgammaIDHelper::computeHGCAL(), sigmas(), and storeRecHits().

76  {
77  theCluster_ = &cluster;
78  storeRecHits(cluster.hitsAndFractions());
79 }
void storeRecHits(const reco::CaloCluster &theCluster)
const reco::CaloCluster * theCluster_
void EGammaPCAHelper::storeRecHits ( const reco::HGCalMultiCluster cluster)

Definition at line 65 of file EgammaPCAHelper.cc.

References reco::HGCalMultiCluster::begin(), reco::HGCalMultiCluster::end(), photonIsolationHIProducer_cfi::hf, mps_fire::result, storeRecHits(), and theCluster_.

65  {
66  theCluster_ = &cluster;
67  std::vector<std::pair<DetId, float>> result;
68  for (reco::HGCalMultiCluster::component_iterator it = cluster.begin(); it != cluster.end();
69  it++) {
70  const std::vector<std::pair<DetId, float>> &hf = (*it)->hitsAndFractions();
71  result.insert(result.end(),hf.begin(),hf.end());
72  }
73  storeRecHits(result);
74 }
component_iterator begin() const
void storeRecHits(const reco::CaloCluster &theCluster)
component_iterator end() const
const reco::CaloCluster * theCluster_
void EGammaPCAHelper::storeRecHits ( const std::vector< std::pair< DetId, float >> &  hf)
private

Definition at line 81 of file EgammaPCAHelper.cc.

References gather_cfg::cout, debug_, dEdXWeights_, DetId::det(), plotBeamSpotDB::first, dedxEstimators_cff::fraction, hgcal::RecHitTools::getLayerWithOffset(), hgcal::RecHitTools::getPosition(), hgcal::RecHitTools::getSiThickIndex(), photonIsolationHIProducer_cfi::hf, hitMap_, invThicknessCorrection_, pcaIteration_, DetId::rawId(), recHitTools_, sige_, sigp_, sigu_, sigv_, theSpots_, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

81  {
82  std::vector<double> pcavars;
83  pcavars.resize(3,0.);
84  theSpots_.clear();
85  pcaIteration_ = 0;
86 
87  sigu_ = 0.;
88  sigv_ = 0.;
89  sigp_ = 0.;
90  sige_ = 0.;
91 
92  unsigned hfsize = hf.size();
93  if (debug_)
94  std::cout << "The seed cluster constains " << hfsize << " hits " << std::endl;
95 
96  if (hfsize == 0) return;
97 
98 
99  for (unsigned int j = 0; j < hfsize; j++) {
100  unsigned int layer = recHitTools_->getLayerWithOffset(hf[j].first);
101 
102  const DetId rh_detid = hf[j].first;
103  std::map<DetId,const HGCRecHit *>::const_iterator itcheck= hitMap_->find(rh_detid);
104  if (itcheck == hitMap_->end()) {
105  edm::LogWarning("EgammaPCAHelper") << " Big problem, unable to find a hit " << rh_detid.rawId() << " " << rh_detid.det() << " " << HGCalDetId(rh_detid) << std::endl;
106  continue;
107  }
108  if (debug_) {
109  std::cout << "DetId " << rh_detid.rawId() << " " << layer << " " << itcheck->second->energy() <<std::endl;
110  std::cout << " Hit " << itcheck->second << " " << itcheck->second->energy() << std::endl;
111  }
112  float fraction = hf[j].second;
113 
114  int thickIndex = recHitTools_->getSiThickIndex(rh_detid);
115  double mip = dEdXWeights_[layer] * 0.001; // convert in GeV
116  if(thickIndex>-1 and thickIndex<3) mip *= invThicknessCorrection_[thickIndex];
117 
118  pcavars[0] = recHitTools_->getPosition(rh_detid).x();
119  pcavars[1] = recHitTools_->getPosition(rh_detid).y();
120  pcavars[2] = recHitTools_->getPosition(rh_detid).z();
121  if (pcavars[2] == 0.)
122  edm::LogWarning("EgammaPCAHelper") << " Problem, hit with z =0 ";
123  else {
124  Spot mySpot(rh_detid,itcheck->second->energy(),pcavars,layer,fraction,mip);
125  theSpots_.push_back(mySpot);
126  }
127  }
128  if (debug_) {
129  std::cout << " Stored " << theSpots_.size() << " hits " << std::endl;
130  }
131 }
std::map< DetId, const HGCRecHit * > * hitMap_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
T y() const
Definition: PV3DBase.h:63
T z() const
Definition: PV3DBase.h:64
std::vector< Spot > theSpots_
Definition: DetId.h:18
unsigned int getLayerWithOffset(const DetId &) const
Definition: RecHitTools.cc:283
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:168
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:112
std::vector< double > invThicknessCorrection_
std::vector< double > dEdXWeights_
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

Member Data Documentation

math::XYZVector hgcal::EGammaPCAHelper::axis_
private

Definition at line 104 of file EgammaPCAHelper.h.

Referenced by axis(), clusterDepthCompatibility(), computePCA(), and energyPerLayer().

math::XYZPoint hgcal::EGammaPCAHelper::barycenter_
private
bool hgcal::EGammaPCAHelper::debug_
private

Definition at line 90 of file EgammaPCAHelper.h.

Referenced by checkIteration(), computePCA(), energyPerLayer(), and storeRecHits().

std::vector<double> hgcal::EGammaPCAHelper::dEdXWeights_
private

Definition at line 93 of file EgammaPCAHelper.h.

Referenced by setdEdXWeights(), and storeRecHits().

std::map<DetId, const HGCRecHit *>* hgcal::EGammaPCAHelper::hitMap_
private

Definition at line 98 of file EgammaPCAHelper.h.

Referenced by fillHitMap(), getHitMap(), setHitMap(), storeRecHits(), and ~EGammaPCAHelper().

int hgcal::EGammaPCAHelper::hitMapOrigin_
private

Definition at line 96 of file EgammaPCAHelper.h.

Referenced by fillHitMap(), setHitMap(), and ~EGammaPCAHelper().

std::vector<double> hgcal::EGammaPCAHelper::invThicknessCorrection_
private

Definition at line 94 of file EgammaPCAHelper.h.

Referenced by storeRecHits().

std::unique_ptr<TPrincipal> hgcal::EGammaPCAHelper::pca_
private

Definition at line 110 of file EgammaPCAHelper.h.

Referenced by computePCA(), eigenValues(), and sigmas().

int hgcal::EGammaPCAHelper::pcaIteration_
private
bool hgcal::EGammaPCAHelper::recHitsStored_
private

Definition at line 89 of file EgammaPCAHelper.h.

const hgcal::RecHitTools* hgcal::EGammaPCAHelper::recHitTools_
private

Definition at line 111 of file EgammaPCAHelper.h.

Referenced by computePCA(), findZFirstLayer(), setRecHitTools(), and storeRecHits().

ShowerDepth hgcal::EGammaPCAHelper::showerDepth_
private

Definition at line 112 of file EgammaPCAHelper.h.

Referenced by clusterDepthCompatibility().

double hgcal::EGammaPCAHelper::sige_
private

Definition at line 107 of file EgammaPCAHelper.h.

Referenced by clear(), computeShowerWidth(), sigmaEE(), and storeRecHits().

double hgcal::EGammaPCAHelper::sigp_
private

Definition at line 107 of file EgammaPCAHelper.h.

Referenced by clear(), computeShowerWidth(), sigmaPP(), and storeRecHits().

double hgcal::EGammaPCAHelper::sigu_
private

Definition at line 107 of file EgammaPCAHelper.h.

Referenced by clear(), computeShowerWidth(), sigmaUU(), and storeRecHits().

double hgcal::EGammaPCAHelper::sigv_
private

Definition at line 107 of file EgammaPCAHelper.h.

Referenced by clear(), computeShowerWidth(), sigmaVV(), and storeRecHits().

const reco::CaloCluster* hgcal::EGammaPCAHelper::theCluster_
private

Definition at line 97 of file EgammaPCAHelper.h.

Referenced by computeShowerWidth(), and storeRecHits().

std::vector<Spot> hgcal::EGammaPCAHelper::theSpots_
private
Transform3D hgcal::EGammaPCAHelper::trans_
private

Definition at line 106 of file EgammaPCAHelper.h.

Referenced by computePCA(), computeShowerWidth(), energyPerLayer(), and printHits().