CMS 3D CMS Logo

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