CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EgammaSCEnergyCorrectionAlgo Class Reference

#include <EgammaSCEnergyCorrectionAlgo.h>

Public Member Functions

reco::SuperCluster applyCorrection (const reco::SuperCluster &cl, const EcalRecHitCollection &rhc, reco::CaloCluster::AlgoId theAlgo, const CaloSubdetectorGeometry *geometry, EcalClusterFunctionBaseClass *energyCorrectionFunction, std::string energyCorrectorName_, int modeEB_, int modeEE_)
 
reco::SuperCluster applyCrackCorrection (const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
 
reco::SuperCluster applyLocalContCorrection (const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *localContCorrectionFunction)
 
 EgammaSCEnergyCorrectionAlgo (float noise, reco::CaloCluster::AlgoId theAlgo, const edm::ParameterSet &pset)
 
 ~EgammaSCEnergyCorrectionAlgo ()
 

Private Member Functions

float fNCrystals (int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
 
int nCrystalsGT2Sigma (reco::BasicCluster const &seed, EcalRecHitCollection const &rhc) const
 

Private Attributes

float sigmaElectronicNoise_
 
reco::CaloCluster::AlgoId theAlgo_
 

Detailed Description

Definition at line 20 of file EgammaSCEnergyCorrectionAlgo.h.

Constructor & Destructor Documentation

EgammaSCEnergyCorrectionAlgo::EgammaSCEnergyCorrectionAlgo ( float  noise,
reco::CaloCluster::AlgoId  theAlgo,
const edm::ParameterSet pset 
)
EgammaSCEnergyCorrectionAlgo::~EgammaSCEnergyCorrectionAlgo ( )
inline

Member Function Documentation

reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyCorrection ( const reco::SuperCluster cl,
const EcalRecHitCollection rhc,
reco::CaloCluster::AlgoId  theAlgo,
const CaloSubdetectorGeometry geometry,
EcalClusterFunctionBaseClass energyCorrectionFunction,
std::string  energyCorrectorName_,
int  modeEB_,
int  modeEE_ 
)

Definition at line 19 of file EgammaSCEnergyCorrectionAlgo.cc.

References SuperClusterShapeAlgo::Calculate_Covariances(), GetRecoTauVFromDQM_MC_cff::cl, reco::CaloCluster::dynamicHybrid, reco::CaloCluster::energy(), SuperClusterShapeAlgo::etaWidth(), fNCrystals(), EcalClusterFunctionBaseClass::getValue(), reco::CaloCluster::hybrid, LogTrace, reco::CaloCluster::multi5x5, nCrystalsGT2Sigma(), SuperClusterShapeAlgo::phiWidth(), reco::SuperCluster::preshowerEnergy(), reco::SuperCluster::rawEnergy(), reco::SuperCluster::seed(), reco::CaloCluster::setEnergy(), reco::SuperCluster::setEtaWidth(), reco::SuperCluster::setPhiWidth(), and createJobs::tmp.

Referenced by ~EgammaSCEnergyCorrectionAlgo().

26  {
27  // A little bit of trivial info to be sure all is well
28 
29  {
30  LogTrace("EgammaSCEnergyCorrectionAlgo") << "::applyCorrection" << std::endl;
31  LogTrace() << " SC has energy " << cl.energy() << std::endl;
32  LogTrace() << " Will correct now.... " << std::endl;
33  }
34 
35  // Get the seed cluster
36  const reco::CaloClusterPtr& seedC = cl.seed();
37 
38  LogTrace("EgammaSCEnergyCorrectionAlgo") << " Seed cluster energy... " << seedC->energy() << std::endl;
39 
40  // Find the algorithm used to construct the basic clusters making up the supercluster
41  LogTrace("EgammaSCEnergyCorrectionAlgo") << " The seed cluster used algo " << theAlgo;
42 
43  // Find the detector region of the supercluster
44  // where is the seed cluster?
45  std::vector<std::pair<DetId, float> > seedHits = seedC->hitsAndFractions();
46  EcalSubdetector theBase = EcalSubdetector(seedHits.at(0).first.subdetId());
47 
48  LogTrace("EgammaSCEnergyCorrectionAlgo") << " seed cluster location == " << theBase << std::endl;
49 
50  // Get number of crystals 2sigma above noise in seed basiccluster
51  int nCryGT2Sigma = nCrystalsGT2Sigma(*seedC, rhc);
52 
53  LogTrace("EgammaSCEnergyCorrectionAlgo") << " nCryGT2Sigma " << nCryGT2Sigma << std::endl;
54 
55  // Supercluster enery - seed basiccluster energy
56  float bremsEnergy = cl.energy() - seedC->energy();
57 
58  LogTrace("EgammaSCEnergyCorrectionAlgo") << " bremsEnergy " << bremsEnergy << std::endl;
59 
60  //Create the pointer ot class SuperClusterShapeAlgo
61  //which calculates phiWidth and etaWidth
62  SuperClusterShapeAlgo SCShape(&rhc, geometry);
63 
64  double phiWidth = 0.;
65  double etaWidth = 0.;
66  //Calculate phiWidth & etaWidth for SuperClusters
67  SCShape.Calculate_Covariances(cl);
68  phiWidth = SCShape.phiWidth();
69  etaWidth = SCShape.etaWidth();
70 
71  // Calculate the new supercluster energy
72  //as a function of number of crystals in the seed basiccluster for Endcap
73  //or apply new Enegry SCale correction
74  float newEnergy = 0;
75 
77  tmp.setPhiWidth(phiWidth);
78  tmp.setEtaWidth(etaWidth);
79 
80  if (theAlgo == reco::CaloCluster::hybrid || theAlgo == reco::CaloCluster::dynamicHybrid) {
81  if (energyCorrectorName_ == "EcalClusterEnergyCorrection")
82  newEnergy = tmp.rawEnergy() + energyCorrectionFunction->getValue(tmp, modeEB_);
83  if (energyCorrectorName_ == "EcalClusterEnergyCorrectionObjectSpecific") {
84  //std::cout << "newEnergy="<<newEnergy<<std::endl;
85  newEnergy = energyCorrectionFunction->getValue(tmp, modeEB_);
86  }
87 
88  } else if (theAlgo == reco::CaloCluster::multi5x5) {
89  if (energyCorrectorName_ == "EcalClusterEnergyCorrection")
90  newEnergy = tmp.rawEnergy() + tmp.preshowerEnergy() + energyCorrectionFunction->getValue(tmp, modeEE_);
91  if (energyCorrectorName_ == "EcalClusterEnergyCorrectionObjectSpecific")
92  newEnergy = energyCorrectionFunction->getValue(tmp, modeEE_);
93 
94  } else {
95  //Apply f(nCry) correction on island algo and fixedMatrix algo
96  newEnergy = seedC->energy() / fNCrystals(nCryGT2Sigma, theAlgo, theBase) + bremsEnergy;
97  }
98 
99  // Create a new supercluster with the corrected energy
100 
101  LogTrace("EgammaSCEnergyCorrectionAlgo") << " UNCORRECTED SC has energy... " << cl.energy() << std::endl;
102  LogTrace("EgammaSCEnergyCorrectionAlgo") << " CORRECTED SC has energy... " << newEnergy << std::endl;
103 
104  reco::SuperCluster corrCl = cl;
105 
106  corrCl.setEnergy(newEnergy);
107  corrCl.setPhiWidth(phiWidth);
108  corrCl.setEtaWidth(etaWidth);
109 
110  return corrCl;
111 }
float fNCrystals(int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
void setEnergy(double energy)
Definition: CaloCluster.h:135
void setPhiWidth(double pw)
Definition: SuperCluster.h:73
void setEtaWidth(double ew)
Definition: SuperCluster.h:74
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double energy() const
cluster energy
Definition: CaloCluster.h:148
#define LogTrace(id)
int nCrystalsGT2Sigma(reco::BasicCluster const &seed, EcalRecHitCollection const &rhc) const
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:58
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:77
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:61
tmp
align.sh
Definition: createJobs.py:716
EcalSubdetector
reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyCrackCorrection ( const reco::SuperCluster cl,
EcalClusterFunctionBaseClass crackCorrectionFunction 
)

Definition at line 219 of file EgammaSCEnergyCorrectionAlgo.cc.

References GetRecoTauVFromDQM_MC_cff::cl, reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), reco::CaloCluster::energy(), EcalClusterFunctionBaseClass::getValue(), reco::SuperCluster::rawEnergy(), and reco::CaloCluster::setEnergy().

Referenced by ~EgammaSCEnergyCorrectionAlgo().

220  {
221  double crackcor = 1.;
222 
223  for (reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
224  const reco::CaloClusterPtr cc = *cIt;
225  crackcor *= ((cl.rawEnergy() + cc->energy() * (crackCorrectionFunction->getValue(*cc) - 1.)) / cl.rawEnergy());
226  } // loop on BCs
227 
228  reco::SuperCluster corrCl = cl;
229  corrCl.setEnergy(cl.energy() * crackcor);
230 
231  return corrCl;
232 }
void setEnergy(double energy)
Definition: CaloCluster.h:135
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double energy() const
cluster energy
Definition: CaloCluster.h:148
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:58
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:86
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:89
reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyLocalContCorrection ( const reco::SuperCluster cl,
EcalClusterFunctionBaseClass localContCorrectionFunction 
)

Definition at line 237 of file EgammaSCEnergyCorrectionAlgo.cc.

References GetRecoTauVFromDQM_MC_cff::cl, reco::CaloCluster::energy(), EcalClusterFunctionBaseClass::getValue(), reco::SuperCluster::seed(), and reco::CaloCluster::setEnergy().

Referenced by ~EgammaSCEnergyCorrectionAlgo().

238  {
240 
241  const reco::CaloClusterPtr& seedBC = cl.seed();
242  float seedBCene = seedBC->energy();
243  float correctedSeedBCene = localContCorrectionFunction->getValue(*seedBC, dummy) * seedBCene;
244 
245  reco::SuperCluster correctedSC = cl;
246  correctedSC.setEnergy(cl.energy() - seedBCene + correctedSeedBCene);
247 
248  return correctedSC;
249 }
void setEnergy(double energy)
Definition: CaloCluster.h:135
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double energy() const
cluster energy
Definition: CaloCluster.h:148
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:77
float EgammaSCEnergyCorrectionAlgo::fNCrystals ( int  nCry,
reco::CaloCluster::AlgoId  theAlgo,
EcalSubdetector  theBase 
) const
private

Definition at line 113 of file EgammaSCEnergyCorrectionAlgo.cc.

References EcalBarrel, EcalEndcap, f, reco::CaloCluster::hybrid, reco::CaloCluster::island, LogTrace, p1, p2, p3, p4, mps_fire::result, and x.

Referenced by applyCorrection(), and ~EgammaSCEnergyCorrectionAlgo().

115  {
116  float p0 = 0, p1 = 0, p2 = 0, p3 = 0, p4 = 0;
117  float x = nCry;
118  float result = 1.f;
119 
120  if ((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::hybrid)) {
121  if (nCry <= 10) {
122  p0 = 6.32879e-01f;
123  p1 = 1.14893e-01f;
124  p2 = -2.45705e-02f;
125  p3 = 2.53074e-03f;
126  p4 = -9.29654e-05f;
127  } else if (nCry > 10 && nCry <= 30) {
128  p0 = 6.93196e-01f;
129  p1 = 4.44034e-02f;
130  p2 = -2.82229e-03f;
131  p3 = 8.19495e-05f;
132  p4 = -8.96645e-07f;
133  } else {
134  p0 = 5.65474e+00f;
135  p1 = -6.31640e-01f;
136  p2 = 3.14218e-02f;
137  p3 = -6.84256e-04f;
138  p4 = 5.50659e-06f;
139  }
140  if (x > 40.f)
141  x = 40.f;
142  }
143 
144  else if ((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::hybrid)) {
145  LogTrace("EgammaSCEnergyCorrectionAlgo") << "ERROR! HybridEFRYsc called" << std::endl;
146 
147  return 1.f;
148  }
149 
150  else if ((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
151  p0 = 4.69976e-01f; // extracted from fit to all endcap classes with Ebrem = 0.
152  p1 = 1.45900e-01f;
153  p2 = -1.61359e-02f;
154  p3 = 7.99423e-04f;
155  p4 = -1.47873e-05f;
156  if (x > 16.f)
157  x = 16.f;
158  }
159 
160  else if ((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
161  p0 = 4.69976e-01f; // extracted from fit to all endcap classes with Ebrem = 0.
162  p1 = 1.45900e-01f;
163  p2 = -1.61359e-02f;
164  p3 = 7.99423e-04f;
165  p4 = -1.47873e-05f;
166  if (x > 16.f)
167  x = 16.f;
168  }
169 
170  else {
171  LogTrace("EgammaSCEnergyCorrectionAlgo") << "trying to correct unknown cluster!!!" << std::endl;
172  return 1.f;
173  }
174  result = p0 + x * (p1 + x * (p2 + x * (p3 + x * p4)));
175 
176  //Rescale energy scale correction to take into account change in calibrated
177  //RecHit definition introduced in CMSSW_1_5_0
178  float const ebfact = 1.f / 0.965f;
179  float const eefact = 1.f / 0.975f;
180 
181  if (theBase == EcalBarrel) {
182  result *= ebfact;
183  } else {
184  result *= eefact;
185  }
186 
187  return result;
188 }
double p4[4]
Definition: TauolaWrapper.h:92
double f[11][100]
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
double p1[4]
Definition: TauolaWrapper.h:89
double p3[4]
Definition: TauolaWrapper.h:91
int EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma ( reco::BasicCluster const &  seed,
EcalRecHitCollection const &  rhc 
) const
private

Definition at line 190 of file EgammaSCEnergyCorrectionAlgo.cc.

References edm::SortedCollection< T, SORT >::find(), hfClusterShapes_cfi::hits, LogTrace, and sigmaElectronicNoise_.

Referenced by applyCorrection(), and ~EgammaSCEnergyCorrectionAlgo().

191  {
192  // return number of crystals 2Sigma above
193  // electronic noise
194 
195  std::vector<std::pair<DetId, float> > const& hits = seed.hitsAndFractions();
196 
197  LogTrace("EgammaSCEnergyCorrectionAlgo") << " EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma" << std::endl;
198  LogTrace("EgammaSCEnergyCorrectionAlgo") << " Will calculate number of crystals above 2sigma noise" << std::endl;
199  LogTrace("EgammaSCEnergyCorrectionAlgo") << " sigmaElectronicNoise = " << sigmaElectronicNoise_ << std::endl;
200  LogTrace("EgammaSCEnergyCorrectionAlgo") << " There are " << hits.size() << " recHits" << std::endl;
201 
202  int nCry = 0;
203  std::vector<std::pair<DetId, float> >::const_iterator hit;
205  for (hit = hits.begin(); hit != hits.end(); hit++) {
206  // need to get hit by DetID in order to get energy
207  aHit = rhc.find((*hit).first);
208  if ((*aHit).energy() > 2.f * sigmaElectronicNoise_)
209  nCry++;
210  }
211 
212  LogTrace("EgammaSCEnergyCorrectionAlgo") << " " << nCry << " of these above 2sigma noise" << std::endl;
213 
214  return nCry;
215 }
std::vector< EcalRecHit >::const_iterator const_iterator
#define LogTrace(id)

Member Data Documentation

float EgammaSCEnergyCorrectionAlgo::sigmaElectronicNoise_
private
reco::CaloCluster::AlgoId EgammaSCEnergyCorrectionAlgo::theAlgo_
private

Definition at line 56 of file EgammaSCEnergyCorrectionAlgo.h.