CMS 3D CMS Logo

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