CMS 3D CMS Logo

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

#include <EgammaSCEnergyCorrectionAlgo.h>

Public Types

using BasicClusterFunction = std::function< float(reco::BasicCluster const &, EcalRecHitCollection const &)>
 

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_)
 
 EgammaSCEnergyCorrectionAlgo (float noise)
 

Static Public Member Functions

static reco::SuperCluster applyCrackCorrection (const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
 
static reco::SuperCluster applyLocalContCorrection (const reco::SuperCluster &cl, BasicClusterFunction localContCorrectionFunction)
 

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_
 

Detailed Description

Definition at line 18 of file EgammaSCEnergyCorrectionAlgo.h.

Member Typedef Documentation

◆ BasicClusterFunction

Definition at line 20 of file EgammaSCEnergyCorrectionAlgo.h.

Constructor & Destructor Documentation

◆ EgammaSCEnergyCorrectionAlgo()

EgammaSCEnergyCorrectionAlgo::EgammaSCEnergyCorrectionAlgo ( float  noise)

Member Function Documentation

◆ applyCorrection()

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 15 of file EgammaSCEnergyCorrectionAlgo.cc.

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

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

◆ applyCrackCorrection()

reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyCrackCorrection ( const reco::SuperCluster cl,
EcalClusterFunctionBaseClass crackCorrectionFunction 
)
static

Definition at line 215 of file EgammaSCEnergyCorrectionAlgo.cc.

References gpuPixelDoublets::cc, haddnano::cl, EcalClusterFunctionBaseClass::getValue(), and reco::CaloCluster::setEnergy().

Referenced by EgammaSCCorrectionMaker::produce().

216  {
217  double crackcor = 1.;
218 
219  for (reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
220  const reco::CaloClusterPtr cc = *cIt;
221  crackcor *= ((cl.rawEnergy() + cc->energy() * (crackCorrectionFunction->getValue(*cc) - 1.)) / cl.rawEnergy());
222  } // loop on BCs
223 
224  reco::SuperCluster corrCl = cl;
225  corrCl.setEnergy(cl.energy() * crackcor);
226 
227  return corrCl;
228 }
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
void setEnergy(double energy)
Definition: CaloCluster.h:135
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0

◆ applyLocalContCorrection()

reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyLocalContCorrection ( const reco::SuperCluster cl,
BasicClusterFunction  localContCorrectionFunction 
)
static

Definition at line 233 of file EgammaSCEnergyCorrectionAlgo.cc.

References haddnano::cl, and reco::CaloCluster::setEnergy().

Referenced by EgammaSCCorrectionMaker::produce().

234  {
236 
237  const reco::CaloClusterPtr& seedBC = cl.seed();
238  float seedBCene = seedBC->energy();
239  float correctedSeedBCene = localContCorrectionFunction(*seedBC, dummy) * seedBCene;
240 
241  reco::SuperCluster correctedSC = cl;
242  correctedSC.setEnergy(cl.energy() - seedBCene + correctedSeedBCene);
243 
244  return correctedSC;
245 }
void setEnergy(double energy)
Definition: CaloCluster.h:135

◆ fNCrystals()

float EgammaSCEnergyCorrectionAlgo::fNCrystals ( int  nCry,
reco::CaloCluster::AlgoId  theAlgo,
EcalSubdetector  theBase 
) const
private

Definition at line 109 of file EgammaSCEnergyCorrectionAlgo.cc.

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

Referenced by applyCorrection().

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

◆ nCrystalsGT2Sigma()

int EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma ( reco::BasicCluster const &  seed,
EcalRecHitCollection const &  rhc 
) const
private

Definition at line 186 of file EgammaSCEnergyCorrectionAlgo.cc.

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

Referenced by applyCorrection().

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

Member Data Documentation

◆ sigmaElectronicNoise_

float EgammaSCEnergyCorrectionAlgo::sigmaElectronicNoise_
private

Definition at line 53 of file EgammaSCEnergyCorrectionAlgo.h.

Referenced by nCrystalsGT2Sigma().