15 using namespace hgcal;
21 invThicknessCorrection_({1. / 1.132, 1. / 1.092, 1. / 1.084}),
22 pca_(
new TPrincipal(3,
"D")) {
24 hitMap_ =
new std::map<DetId, const HGCRecHit *>();
46 for (
const auto&
hit : rechitsEE) {
50 for (
const auto&
hit : rechitsFH) {
54 for (
const auto&
hit : rechitsBH) {
64 std::vector<std::pair<DetId, float>>
result;
67 const std::vector<std::pair<DetId, float>> &
hf = (*it)->hitsAndFractions();
68 result.insert(result.end(),hf.begin(),hf.end());
79 std::vector<double> pcavars;
89 unsigned hfsize =
hf.size();
91 std::cout <<
"The seed cluster constains " << hfsize <<
" hits " << std::endl;
93 if (hfsize == 0)
return;
96 for (
unsigned int j = 0; j < hfsize; j++) {
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;
106 std::cout <<
"DetId " << rh_detid.
rawId() <<
" " << layer <<
" " << itcheck->second->energy() <<std::endl;
107 std::cout <<
" Hit " << itcheck->second <<
" " << itcheck->second->energy() << std::endl;
113 if (thickness > 99. && thickness < 101)
115 else if (thickness > 199 && thickness < 201)
117 else if (thickness > 299 && thickness < 301)
123 if (pcavars[2] == 0.)
126 Spot mySpot(rh_detid,itcheck->second->energy(),pcavars,layer,
fraction,mip);
137 pca_.reset(
new TPrincipal(3,
"D"));
138 bool initialCalculation = radius < 0;
140 std::cout <<
" Initial calculation " << initialCalculation << std::endl;
141 if (initialCalculation && withHalo) {
142 edm::LogWarning(
"EGammaPCAHelper") <<
"Warning - in the first iteration, the halo hits are excluded " << std::endl;
146 float radius2 = radius*
radius;
147 if (! initialCalculation) {
160 if (!withHalo && (! spot.isCore() ))
162 if (initialCalculation) {
164 if ( ! spot.isCore() )
continue;
165 layers.insert(spot.layer());
166 for (
int i = 0;
i < spot.multiplicity(); ++
i)
167 pca_->AddRow(spot.row());
172 if (local.Perp2() > radius2)
continue;
173 layers.insert(spot.layer());
174 for (
int i = 0;
i < spot.multiplicity(); ++
i)
175 pca_->AddRow(spot.row());
179 std::cout <<
" Nlayers " << layers.size() << std::endl;
180 if (layers.size() < 3) {
184 pca_->MakePrincipals();
186 const TVectorD& means = *(
pca_->GetMeanValues());
187 const TMatrixD& eigens = *(
pca_->GetEigenVectors());
203 float radius2 = radius *
radius;
205 Point globalPoint(spot.row()[0],spot.row()[1],spot.row()[2]);
207 if (local.Perp2() > radius2)
continue;
210 if (withHalo && spot.fraction() < 0)
continue;
211 if (!withHalo && !(spot.isCore()))
continue;
217 sigu_ += local.x() * local.x() * spot.energy();
218 sigv_ += local.y() * local.y() * spot.energy();
219 cyl_ene += spot.energy();
223 const double inv_cyl_ene = 1. / cyl_ene;
238 std::cout <<
" The PCA has not been run yet " << std::endl;
242 std::cout <<
" The PCA has been run only once - careful " << std::endl;
246 std::cout <<
" Not enough layers to perform PCA " << std::endl;
264 float radius2 = radius*
radius;
278 if (!withHalo && ! spot.isCore())
281 if (local.Perp2() > radius2)
continue;
282 energyPerLayer[spot.layer()] += spot.energy();
283 layers.insert(spot.layer());
284 if (spot.subdet() ==
HGCEE) { energyEE += spot.energy();}
285 else if (spot.subdet() ==
HGCHEF) { energyFH += spot.energy();}
286 else if (spot.subdet() ==
HGCHEB) { energyBH += spot.energy();}
289 return LongDeps(radius,energyPerLayer,energyEE,energyFH,energyBH,layers);
294 float radius2=radius*
radius;
295 for (
unsigned i =0;
i< nSpots ; ++
i) {
298 if (local.Perp2() < radius2 ) {
306 unsigned int firstLayer = 0;
322 expectedDepth = -999.;
323 expectedSigma = -999.;
324 measuredDepth = -999.;
void computeShowerWidth(float radius, bool withHalo=true)
math::XYZPoint barycenter_
static const unsigned int maxlayer
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
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
ROOT::Math::Transform3DPJ Transform3D
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
float getClusterDepthCompatibility(float measuredDepth, float emEnergy, float &expectedDepth, float &expectedSigma) const
double eta() const
pseudorapidity of cluster centroid
uint32_t rawId() const
get the raw id
void storeRecHits(const reco::CaloCluster &theCluster)
ROOT::Math::Transform3DPJ::Point Point
void setHitMap(std::map< DetId, const HGCRecHit * > *hitMap)
to set from outside - once per event
Abs< T >::type abs(const T &t)
std::vector< Spot > theSpots_
component_iterator end() const
void setRecHitTools(const hgcal::RecHitTools *recHitTools)
float clusterDepthCompatibility(const LongDeps &, float &measuredDepth, float &expectedDepth, float &expectedSigma)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
XYZPointD XYZPoint
point in space with cartesian internal representation
const std::vector< float > & energyPerLayer() const
void printHits(float radius) const
const double * row() const
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
double phi() const
azimuthal angle of cluster centroid
Detector det() const
get the detector field from this detid
const hgcal::RecHitTools * recHitTools_