CMS 3D CMS Logo

EcalRegressionData.cc
Go to the documentation of this file.
2 
9 
12 
15 
17 {
18  float eLeftRightSum = eLeft()+eRight();
19  float eLeftRightDiff = eLeft()-eRight();
20  return eLeftRightSum !=0 ? eLeftRightDiff/eLeftRightSum : 0.;
21 }
22 
24 {
25  float eTopBottomSum = eTop()+eBottom();
26  float eTopBottomDiff = eTop()-eBottom();
27  return eTopBottomSum !=0 ? eTopBottomDiff/eTopBottomSum : 0.;
28 }
29 
30 float EcalRegressionData::subClusRawEnergy(size_t clusNr)const
31 {
32  if(clusNr<subClusRawEnergy_.size()) return subClusRawEnergy_[clusNr];
33  else return 0.;
34 }
35 
36 float EcalRegressionData::subClusDEta(size_t clusNr)const
37 {
38  if(clusNr<subClusDEta_.size()) return subClusDEta_[clusNr];
39  else return 0.;
40 }
41 
42 float EcalRegressionData::subClusDPhi(size_t clusNr)const
43 {
44  if(clusNr<subClusDPhi_.size()) return subClusDPhi_[clusNr];
45  else return 0.;
46 }
47 
49  const EcalRecHitCollection* ebRecHits,const EcalRecHitCollection* eeRecHits,
50  const CaloGeometry* geom,const CaloTopology* topology,
51  int nrVertices)
52 {
53  clear();
54 
55  const DetId& seedid = superClus.seed()->hitsAndFractions().at(0).first;
56  isEB_ = ( seedid.subdetId()==EcalBarrel );
57 
58  // skip HGCal
59  if( EcalTools::isHGCalDet(seedid.det()) ) return;
60 
61  const EcalRecHitCollection* recHits = isEB_ ? ebRecHits : eeRecHits;
62 
63  scRawEnergy_ = superClus.rawEnergy();
64  scCalibEnergy_ = superClus.correctedEnergy();
65  scPreShowerEnergy_ = superClus.preshowerEnergy();
66  scEta_ = superClus.eta();
67  scPhi_ = superClus.phi();
68  scEtaWidth_ = superClus.etaWidth();
69  scPhiWidth_ = superClus.phiWidth();
70  scNrAdditionalClusters_ = superClus.clustersSize()-1;
71 
72  seedClusEnergy_ = superClus.seed()->energy();
73  eMax_ = EcalClusterTools::eMax(*superClus.seed(),recHits);
74  e2nd_ = EcalClusterTools::e2nd(*superClus.seed(),recHits);
75  e3x3_ = EcalClusterTools::e3x3(*superClus.seed(),recHits,topology);
76  eTop_ = EcalClusterTools::eTop(*superClus.seed(),recHits,topology);
77  eBottom_ = EcalClusterTools::eBottom(*superClus.seed(),recHits,topology);
78  eLeft_ = EcalClusterTools::eLeft(*superClus.seed(),recHits,topology);
79  eRight_ = EcalClusterTools::eRight(*superClus.seed(),recHits,topology);
80  std::vector<float> localCovs = EcalClusterTools::localCovariances(*superClus.seed(),recHits,topology);
81  sigmaIEtaIEta_ = edm::isNotFinite(localCovs[0]) ? 0. : std::sqrt(localCovs[0]);
82  sigmaIPhiIPhi_ = edm::isNotFinite(localCovs[2]) ? 0. : std::sqrt(localCovs[2]);
83 
84  if(sigmaIEtaIEta_*sigmaIPhiIPhi_>0) sigmaIEtaIPhi_ = localCovs[1]/(sigmaIEtaIEta_*sigmaIPhiIPhi_);
85  else if(localCovs[1]>0) sigmaIEtaIPhi_ = 1.;
86  else sigmaIEtaIPhi_ = -1.;
87 
89 
90  float dummy;
91  localCoordsFunc( *superClus.seed(),*geom, seedCrysEtaOrX_,seedCrysPhiOrY_,
92  seedCrysIEtaOrIX_,seedCrysIPhiOrIY_, dummy, dummy );
93 
94  for( auto clus = superClus.clustersBegin()+1;clus != superClus.clustersEnd(); ++clus ) {
95  const float dEta = (*clus)->eta() - superClus.seed()->eta();
96  const float dPhi = reco::deltaPhi((*clus)->phi(),superClus.seed()->phi());
97  const float dR2 = dEta*dEta+dPhi*dPhi;
98  if(dR2 > maxSubClusDR2_ || maxSubClusDR2_ == 998001.) {
99  maxSubClusDR2_ = dR2;
102  maxSubClusDRRawEnergy_ = (*clus)->energy();
103  }
104  subClusRawEnergy_.push_back((*clus)->energy());
105  subClusDEta_.push_back(dEta);
106  subClusDPhi_.push_back(dPhi);
107 
108  }
109 
110  nrVtx_ = nrVertices;
111 
112 }
113 
115 {
116  isEB_=false;
117  scRawEnergy_=0.;
118  scCalibEnergy_=0.;
120  scEta_=0.;
121  scPhi_=0.;
122  scEtaWidth_=0.;
123  scPhiWidth_=0.;
125 
126  seedClusEnergy_=0.;
127  eMax_=0.;
128  e2nd_=0.;
129  e3x3_=0.;
130  eTop_=0.;
131  eBottom_=0.;
132  eLeft_=0.;
133  eRight_=0.;
134  sigmaIEtaIEta_=0.;
135  sigmaIEtaIPhi_=0.;
136  sigmaIPhiIPhi_=0.;
137 
138  seedCrysPhiOrY_=0.;
139  seedCrysEtaOrX_=0.;
142 
143  maxSubClusDR2_=998001.;
144  maxSubClusDRDPhi_=999.;
145  maxSubClusDRDEta_=999;
147 
148  subClusRawEnergy_.clear();
149  subClusDPhi_.clear();
150  subClusDEta_.clear();
151 
152  nrVtx_=0;
153 
154 }
155 
156 void EcalRegressionData::fillVec(std::vector<float>& inputVec)const
157 {
158  if(isEB()) fillVecEB_(inputVec);
159  else fillVecEE_(inputVec);
160 }
161 
162 void EcalRegressionData::fillVecEB_(std::vector<float>& inputVec)const
163 {
164  inputVec.clear();
165  inputVec.resize(33);
166  inputVec[0] = nrVtx(); //nVtx
167  inputVec[1] = scEta(); //scEta
168  inputVec[2] = scPhi(); //scPhi
169  inputVec[3] = scEtaWidth(); //scEtaWidth
170  inputVec[4] = scPhiWidth(); //scPhiWidth
171  inputVec[5] = scSeedR9(); //scSeedR9
172  inputVec[6] = seedClusEnergyOverSCRawEnergy(); //scSeedRawEnergy/scRawEnergy
173  inputVec[7] = eMaxOverSCRawEnergy(); //scSeedEmax/scRawEnergy
174  inputVec[8] = e2ndOverSCRawEnergy(); //scSeedE2nd/scRawEnergy
175  inputVec[9] = seedLeftRightAsym();//scSeedLeftRightAsym
176  inputVec[10] = seedTopBottomAsym();//scSeedTopBottomAsym
177  inputVec[11] = sigmaIEtaIEta(); //scSeedSigmaIetaIeta
178  inputVec[12] = sigmaIEtaIPhi(); //scSeedSigmaIetaIphi
179  inputVec[13] = sigmaIPhiIPhi(); //scSeedSigmaIphiIphi
180  inputVec[14] = scNrAdditionalClusters(); //N_ECALClusters
181  inputVec[15] = maxSubClusDR(); //clusterMaxDR
182  inputVec[16] = maxSubClusDRDPhi(); //clusterMaxDRDPhi
183  inputVec[17] = maxSubClusDRDEta(); //clusterMaxDRDEta
184  inputVec[18] = maxSubClusDRRawEnergyOverSCRawEnergy(); //clusMaxDRRawEnergy/scRawEnergy
185  inputVec[19] = subClusRawEnergyOverSCRawEnergy(SubClusNr::C1); //clusterRawEnergy[0]/scRawEnergy
186  inputVec[20] = subClusRawEnergyOverSCRawEnergy(SubClusNr::C2); //clusterRawEnergy[1]/scRawEnergy float scPreShowerEnergy()const{return scPreShowerEnergy_;}
187 
188  inputVec[21] = subClusRawEnergyOverSCRawEnergy(SubClusNr::C3); //clusterRawEnergy[2]/scRawEnergy
189  inputVec[22] = subClusDPhi(SubClusNr::C1); //clusterDPhiToSeed[0]
190  inputVec[23] = subClusDPhi(SubClusNr::C2); //clusterDPhiToSeed[1]
191  inputVec[24] = subClusDPhi(SubClusNr::C3); //clusterDPhiToSeed[2]
192  inputVec[25] = subClusDEta(SubClusNr::C1); //clusterDEtaToSeed[0]
193  inputVec[26] = subClusDEta(SubClusNr::C2); //clusterDEtaToSeed[1]
194  inputVec[27] = subClusDEta(SubClusNr::C3); //clusterDEtaToSeed[2]
195  inputVec[28] = seedCrysEtaOrX(); //scSeedCryEta
196  inputVec[29] = seedCrysPhiOrY(); //scSeedCryPhi
197  inputVec[30] = seedCrysIEtaOrIX(); //scSeedCryIeta
198  inputVec[31] = seedCrysIPhiOrIY(); //scSeedCryIphi
199  inputVec[32] = scCalibEnergy(); //scCalibratedEnergy
200 }
201 
202 void EcalRegressionData::fillVecEE_(std::vector<float>& inputVec)const
203 {
204  inputVec.clear();
205  inputVec.resize(33); //this emulates the old behaviour of RegressionHelper, even if past 29 we dont use elements
206  inputVec[0] = nrVtx(); //nVtx
207  inputVec[1] = scEta(); //scEta
208  inputVec[2] = scPhi(); //scPhi
209  inputVec[3] = scEtaWidth(); //scEtaWidth
210  inputVec[4] = scPhiWidth(); //scPhiWidth
211  inputVec[5] = scSeedR9(); //scSeedR9
212  inputVec[6] = seedClusEnergyOverSCRawEnergy(); //scSeedRawEnergy/scRawEnergy
213  inputVec[7] = eMaxOverSCRawEnergy(); //scSeedEmax/scRawEnergy
214  inputVec[8] = e2ndOverSCRawEnergy(); //scSeedE2nd/scRawEnergy
215  inputVec[9] = seedLeftRightAsym();//scSeedLeftRightAsym
216  inputVec[10] = seedTopBottomAsym();//scSeedTopBottomAsym
217  inputVec[11] = sigmaIEtaIEta(); //scSeedSigmaIetaIeta
218  inputVec[12] = sigmaIEtaIPhi(); //scSeedSigmaIetaIphi
219  inputVec[13] = sigmaIPhiIPhi(); //scSeedSigmaIphiIphi
220  inputVec[14] = scNrAdditionalClusters(); //N_ECALClusters
221  inputVec[15] = maxSubClusDR(); //clusterMaxDR
222  inputVec[16] = maxSubClusDRDPhi(); //clusterMaxDRDPhi
223  inputVec[17] = maxSubClusDRDEta(); //clusterMaxDRDEta
224  inputVec[18] = maxSubClusDRRawEnergyOverSCRawEnergy(); //clusMaxDRRawEnergy/scRawEnergy
225  inputVec[19] = subClusRawEnergyOverSCRawEnergy(SubClusNr::C1); //clusterRawEnergy[0]/scRawEnergy
226  inputVec[20] = subClusRawEnergyOverSCRawEnergy(SubClusNr::C2); //clusterRawEnergy[1]/scRawEnergy
227  inputVec[21] = subClusRawEnergyOverSCRawEnergy(SubClusNr::C3); //clusterRawEnergy[2]/scRawEnergy
228  inputVec[22] = subClusDPhi(SubClusNr::C1); //clusterDPhiToSeed[0]
229  inputVec[23] = subClusDPhi(SubClusNr::C2); //clusterDPhiToSeed[1]
230  inputVec[24] = subClusDPhi(SubClusNr::C3); //clusterDPhiToSeed[2]
231  inputVec[25] = subClusDEta(SubClusNr::C1); //clusterDEtaToSeed[0]
232  inputVec[26] = subClusDEta(SubClusNr::C2); //clusterDEtaToSeed[1]
233  inputVec[27] = subClusDEta(SubClusNr::C3); //clusterDEtaToSeed[2]
234  inputVec[28] = scPreShowerEnergyOverSCRawEnergy();
235  inputVec[29] = scCalibEnergy(); //scCalibratedEnergy
236 }
float sigmaIPhiIPhi() const
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
static bool isHGCalDet(DetId::Detector thedet)
identify HGCal cells
Definition: EcalTools.h:54
const std::vector< float > & subClusRawEnergy() const
float e2ndOverSCRawEnergy() const
void localEcalClusterCoordsEB(const reco::CaloCluster &bclus, const CaloGeometry &geom, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt)
float maxSubClusDRDEta() const
double correctedEnergy() const
Definition: CaloCluster.h:127
float scEtaWidth() const
CaloTopology const * topology(0)
double phiWidth() const
obtain phi and eta width of the Super Cluster
Definition: SuperCluster.h:55
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
void fill(const reco::SuperCluster &superClus, const EcalRecHitCollection *ebRecHits, const EcalRecHitCollection *eeRecHits, const CaloGeometry *geom, const CaloTopology *topology, const reco::VertexCollection *vertices)
int scNrAdditionalClusters() const
float eLeft() const
float seedCrysPhiOrY() const
float sigmaIEtaIEta() const
float sigmaIEtaIPhi() const
float scSeedR9() const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:168
float subClusRawEnergyOverSCRawEnergy(size_t clusNr) const
float seedClusEnergyOverSCRawEnergy() const
std::vector< float > subClusDPhi_
double etaWidth() const
Definition: SuperCluster.h:56
float scPhi() const
float seedCrysIPhiOrIY() const
T sqrt(T t)
Definition: SSEVec.h:18
float maxSubClusDR() const
float seedCrysIEtaOrIX() const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::vector< float > subClusRawEnergy_
float maxSubClusDRDPhi() const
void fillVec(std::vector< float > &inputVec) const
void localEcalClusterCoordsEE(const reco::CaloCluster &bclus, const CaloGeometry &geom, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt)
float scPhiWidth() const
float scPreShowerEnergyOverSCRawEnergy() const
float seedCrysEtaOrX() const
float seedLeftRightAsym() const
Definition: DetId.h:18
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:47
std::vector< float > subClusDEta_
size_t clustersSize() const
number of BasicCluster constituents
Definition: SuperCluster.h:87
float eMaxOverSCRawEnergy() const
float maxSubClusDRRawEnergyOverSCRawEnergy() const
const std::vector< float > & subClusDEta() const
float scEta() const
const std::vector< float > & subClusDPhi() const
float seedTopBottomAsym() const
float eRight() const
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:75
void fillVecEB_(std::vector< float > &inputVec) const
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:66
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:50
float scCalibEnergy() const
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:171
void fillVecEE_(std::vector< float > &inputVec) const
float eBottom() const
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:78
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39