CMS 3D CMS Logo

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